PatchCollisionDensity.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2018-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "PatchCollisionDensity.H"
27 #include "Pstream.H"
28 #include "stringListOps.H"
29 #include "ListOps.H"
30 #include "ListListOps.H"
31 #include "nonConformalFvPatch.H"
32 #include "fvPatchFieldMapper.H"
34 
35 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
36 
37 template<class CloudType>
39 {
40  const scalarField z(this->owner().mesh().nCells(), 0);
41 
43  (
44  IOobject
45  (
46  this->owner().name() + ":numberCollisionDensity",
47  this->owner().mesh().time().name(),
48  this->owner().mesh()
49  ),
50  this->owner().mesh(),
52  z,
53  numberCollisionDensity_
54  )
55  .write();
56 
58  (
59  IOobject
60  (
61  this->owner().name() + ":numberCollisionDensityRate",
62  this->owner().mesh().time().name(),
63  this->owner().mesh()
64  ),
65  this->owner().mesh(),
67  z,
68  (numberCollisionDensity_ - numberCollisionDensity0_)
69  /(this->owner().mesh().time().value() - time0_)
70  )
71  .write();
72 
74  (
75  IOobject
76  (
77  this->owner().name() + ":massCollisionDensity",
78  this->owner().mesh().time().name(),
79  this->owner().mesh()
80  ),
81  this->owner().mesh(),
83  z,
84  massCollisionDensity_
85  )
86  .write();
87 
89  (
90  IOobject
91  (
92  this->owner().name() + ":massCollisionDensityRate",
93  this->owner().mesh().time().name(),
94  this->owner().mesh()
95  ),
96  this->owner().mesh(),
98  z,
99  (massCollisionDensity_ - massCollisionDensity0_)
100  /(this->owner().mesh().time().value() - time0_)
101  )
102  .write();
103 
104  numberCollisionDensity0_ == numberCollisionDensity_;
105  massCollisionDensity0_ == massCollisionDensity_;
106  time0_ = this->owner().mesh().time().value();
107 }
108 
109 
110 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
111 
112 template<class CloudType>
114 (
115  const dictionary& dict,
116  CloudType& owner,
117  const word& modelName
118 )
119 :
120  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
121  minSpeed_(dict.lookupOrDefault<scalar>("minSpeed", -1)),
122  numberCollisionDensity_
123  (
124  this->owner().mesh().boundary(),
125  volScalarField::Internal::null(),
126  calculatedFvPatchField<scalar>::typeName
127  ),
128  numberCollisionDensity0_
129  (
130  this->owner().mesh().boundary(),
131  volScalarField::Internal::null(),
132  calculatedFvPatchField<scalar>::typeName
133  ),
134  massCollisionDensity_
135  (
136  this->owner().mesh().boundary(),
137  volScalarField::Internal::null(),
138  calculatedFvPatchField<scalar>::typeName
139  ),
140  massCollisionDensity0_
141  (
142  this->owner().mesh().boundary(),
143  volScalarField::Internal::null(),
144  calculatedFvPatchField<scalar>::typeName
145  ),
146  time0_(this->owner().mesh().time().value())
147 {
148  numberCollisionDensity_ == 0;
149  numberCollisionDensity0_ == 0;
150  massCollisionDensity_ == 0;
151  massCollisionDensity0_ == 0;
152 
154  (
155  this->owner().name() + ":numberCollisionDensity",
156  this->owner().mesh().time().name(),
157  this->owner().mesh(),
160  );
161 
162  if (numberIo.headerOk())
163  {
164  const volScalarField numberCollisionDensity
165  (
166  numberIo,
167  this->owner().mesh()
168  );
169  numberCollisionDensity_ == numberCollisionDensity.boundaryField();
170  numberCollisionDensity0_ == numberCollisionDensity.boundaryField();
171  }
172 
174  (
175  this->owner().name() + ":massCollisionDensity",
176  this->owner().mesh().time().name(),
177  this->owner().mesh(),
180  );
181 
182  if (massIo.headerOk())
183  {
184  const volScalarField massCollisionDensity
185  (
186  massIo,
187  this->owner().mesh()
188  );
189  massCollisionDensity_ == massCollisionDensity.boundaryField();
190  massCollisionDensity0_ == massCollisionDensity.boundaryField();
191  }
192 }
193 
194 
195 template<class CloudType>
197 (
199 )
200 :
202  minSpeed_(ppm.minSpeed_),
203  numberCollisionDensity_
204  (
205  volScalarField::Internal::null(),
206  ppm.numberCollisionDensity_
207  ),
208  numberCollisionDensity0_
209  (
210  volScalarField::Internal::null(),
211  ppm.numberCollisionDensity0_
212  ),
213  massCollisionDensity_
214  (
215  volScalarField::Internal::null(),
216  ppm.massCollisionDensity_
217  ),
218  massCollisionDensity0_
219  (
220  volScalarField::Internal::null(),
221  ppm.massCollisionDensity0_
222  ),
223  time0_(ppm.time0_)
224 {}
225 
226 
227 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
228 
229 template<class CloudType>
231 {}
232 
233 
234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
235 
236 template<class CloudType>
238 {
239  const fvMesh& mesh = this->owner().mesh();
240 
241  if (!mesh.conformal())
242  {
243  forAll(mesh.boundary(), patchi)
244  {
245  const fvPatch& fvp = mesh.boundary()[patchi];
246 
247  if (isA<nonConformalFvPatch>(fvp))
248  {
249  const setSizeFvPatchFieldMapper mapper(fvp.size());
250 
251  numberCollisionDensity_[patchi].map
252  (
253  numberCollisionDensity_[patchi],
254  mapper
255  );
256  numberCollisionDensity0_[patchi].map
257  (
258  numberCollisionDensity0_[patchi],
259  mapper
260  );
261  massCollisionDensity_[patchi].map
262  (
263  massCollisionDensity_[patchi],
264  mapper
265  );
266  massCollisionDensity0_[patchi].map
267  (
268  massCollisionDensity0_[patchi],
269  mapper
270  );
271 
272  numberCollisionDensity_[patchi] == 0;
273  numberCollisionDensity0_[patchi] == 0;
274  massCollisionDensity_[patchi] == 0;
275  massCollisionDensity0_[patchi] == 0;
276  }
277  }
278  }
279 }
280 
281 
282 template<class CloudType>
284 (
285  const parcelType& p,
286  const polyPatch& pp
287 )
288 {
289  if (pp.coupled()) return;
290 
291  const label patchi = pp.index();
292  const label patchFacei = p.face() - pp.start();
293 
294  vector nw, Up;
295  this->owner().patchData(p, pp, nw, Up);
296 
297  const scalar speed = (p.U() - Up) & nw;
298 
299  if (speed > minSpeed_)
300  {
301  const scalar magSf =
302  this->owner().mesh().magSf().boundaryField()[patchi][patchFacei];
303 
304  numberCollisionDensity_[patchi][patchFacei] +=
305  p.nParticle()/magSf;
306  massCollisionDensity_[patchi][patchFacei] +=
307  p.nParticle()*p.mass()/magSf;
308  }
309 }
310 
311 
312 // ************************************************************************* //
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Templated cloud function object base class.
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Function object which generates fields of the number and mass and rates thereof of collisions per uni...
PatchCollisionDensity(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual void postPatch(const parcelType &p, const polyPatch &pp)
Post-patch hook.
virtual ~PatchCollisionDensity()
Destructor.
void write()
Write post-processing info.
virtual void preEvolve()
Pre-evolve hook.
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
bool conformal() const
Return whether the fvMesh is conformal with the polyMesh.
Definition: fvMesh.C:910
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual label size() const
Return size.
Definition: fvPatch.H:157
label index() const
Return the index of this patch in the boundaryMesh.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:290
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
A class for handling words, derived from string.
Definition: word.H:62
label patchi
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimless
const dimensionSet dimTime
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
const dimensionSet dimMass
const dimensionSet dimArea
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
faceListList boundary(nPatches)
dictionary dict
volScalarField & p
Operations on lists of strings.