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 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 
32 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
33 
34 template<class CloudType>
36 {
37  const scalarField z(this->owner().mesh().nCells(), 0);
38 
40  (
41  IOobject
42  (
43  this->owner().name() + ":collisionDensity",
44  this->owner().mesh().time().timeName(),
45  this->owner().mesh()
46  ),
47  this->owner().mesh(),
49  z,
50  collisionDensity_
51  )
52  .write();
53 
55  (
56  IOobject
57  (
58  this->owner().name() + ":collisionDensityRate",
59  this->owner().mesh().time().timeName(),
60  this->owner().mesh()
61  ),
62  this->owner().mesh(),
64  z,
65  (collisionDensity_ - collisionDensity0_)
66  /(this->owner().mesh().time().value() - time0_)
67  )
68  .write();
69 
70  collisionDensity0_ == collisionDensity_;
71  time0_ = this->owner().mesh().time().value();
72 }
73 
74 
75 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
76 
77 template<class CloudType>
79 (
80  const dictionary& dict,
81  CloudType& owner,
82  const word& modelName
83 )
84 :
85  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
86  minSpeed_(dict.lookupOrDefault<scalar>("minSpeed", -1)),
87  collisionDensity_
88  (
89  this->owner().mesh().boundary(),
90  volScalarField::Internal::null(),
92  ),
93  collisionDensity0_
94  (
95  this->owner().mesh().boundary(),
96  volScalarField::Internal::null(),
98  ),
99  time0_(this->owner().mesh().time().value())
100 {
101  collisionDensity_ == 0;
102  collisionDensity0_ == 0;
103 
104  IOobject io
105  (
106  this->owner().name() + ":collisionDensity",
107  this->owner().mesh().time().timeName(),
108  this->owner().mesh(),
109  IOobject::MUST_READ,
110  IOobject::NO_WRITE
111  );
112 
113  if (io.typeHeaderOk<volScalarField>())
114  {
115  const volScalarField collisionDensity(io, this->owner().mesh());
116  collisionDensity_ == collisionDensity.boundaryField();
117  collisionDensity0_ == collisionDensity.boundaryField();
118  }
119 }
120 
121 
122 template<class CloudType>
124 (
126 )
127 :
129  minSpeed_(ppm.minSpeed_),
130  collisionDensity_(ppm.collisionDensity_),
131  collisionDensity0_(ppm.collisionDensity0_),
132  time0_(ppm.time0_)
133 {}
134 
135 
136 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
137 
138 template<class CloudType>
140 {}
141 
142 
143 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144 
145 template<class CloudType>
147 (
148  const parcelType& p,
149  const polyPatch& pp,
150  bool&
151 )
152 {
153  const label patchi = pp.index();
154  const label patchFacei = p.face() - pp.start();
155 
156  vector nw, Up;
157  this->owner().patchData(p, pp, nw, Up);
158 
159  const scalar speed = (p.U() - Up) & nw;
160  if (speed > minSpeed_)
161  {
162  collisionDensity_[patchi][patchFacei] +=
163  1/this->owner().mesh().magSf().boundaryField()[patchi][patchFacei];
164  }
165 }
166 
167 
168 // ************************************************************************* //
dictionary dict
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
bool typeHeaderOk(const bool checkType=true)
Read header (uses typeFilePath to find file) and check header.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Operations on lists of strings.
void write()
Write post-processing info.
Various functions to operate on Lists.
virtual ~PatchCollisionDensity()
Destructor.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
word timeName
Definition: getTimeIndex.H:3
virtual void postPatch(const parcelType &p, const polyPatch &pp, bool &keepParticle)
Post-patch hook.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Function object which generates fields of the number and rate of collisions per unit area on all patc...
label index() const
Return the index of this patch in the boundaryMesh.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
PatchCollisionDensity(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual bool write(const bool valid=true) const
Write using setting from DB.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
Templated cloud function object base class.
label nw
Definition: createFields.H:25
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57