PatchCollisionDensity.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2018 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "PatchCollisionDensity.H"
30 #include "Pstream.H"
31 #include "stringListOps.H"
32 #include "ListOps.H"
33 #include "ListListOps.H"
34 
35 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
36 
37 template<class CloudType>
39 {
40  const scalarField z(this->owner().mesh().nCells(), Zero);
41 
43  (
44  IOobject
45  (
46  this->owner().name() + ":collisionDensity",
47  this->owner().mesh().time().timeName(),
48  this->owner().mesh()
49  ),
50  this->owner().mesh(),
52  z,
53  collisionDensity_
54  )
55  .write();
56 
58  (
59  IOobject
60  (
61  this->owner().name() + ":collisionDensityRate",
62  this->owner().mesh().time().timeName(),
63  this->owner().mesh()
64  ),
65  this->owner().mesh(),
67  z,
68  (collisionDensity_ - collisionDensity0_)
69  /(this->owner().mesh().time().value() - time0_)
70  )
71  .write();
72 
73  collisionDensity0_ == collisionDensity_;
74  time0_ = this->owner().mesh().time().value();
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 
80 template<class CloudType>
82 (
83  const dictionary& dict,
84  CloudType& owner,
85  const word& modelName
86 )
87 :
88  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
89  minSpeed_(dict.getOrDefault<scalar>("minSpeed", -1)),
90  collisionDensity_
91  (
92  this->owner().mesh().boundary(),
93  volScalarField::Internal::null(),
95  ),
96  collisionDensity0_
97  (
98  this->owner().mesh().boundary(),
99  volScalarField::Internal::null(),
101  ),
102  time0_(this->owner().mesh().time().value())
103 {
104  collisionDensity_ == 0;
105  collisionDensity0_ == 0;
106 
107  IOobject io
108  (
109  this->owner().name() + ":collisionDensity",
110  this->owner().mesh().time().timeName(),
111  this->owner().mesh(),
114  );
115 
116  if (io.typeHeaderOk<volScalarField>())
117  {
118  const volScalarField collisionDensity(io, this->owner().mesh());
119  collisionDensity_ == collisionDensity.boundaryField();
120  collisionDensity0_ == collisionDensity.boundaryField();
121  }
122 }
123 
124 
125 template<class CloudType>
127 (
128  const PatchCollisionDensity<CloudType>& ppm
129 )
130 :
131  CloudFunctionObject<CloudType>(ppm),
132  minSpeed_(ppm.minSpeed_),
133  collisionDensity_(ppm.collisionDensity_),
134  collisionDensity0_(ppm.collisionDensity0_),
135  time0_(ppm.time0_)
136 {}
137 
138 
139 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140 
141 template<class CloudType>
143 (
144  const parcelType& p,
145  const polyPatch& pp,
146  const typename parcelType::trackingData& td
147 )
148 {
149  const label patchi = pp.index();
150  const label patchFacei = p.face() - pp.start();
151 
152  vector nw, Up;
153  this->owner().patchData(p, pp, nw, Up);
154 
155  const scalar speed = (p.U() - Up) & nw;
156  if (speed > minSpeed_)
157  {
158  collisionDensity_[patchi][patchFacei] +=
159  1/this->owner().mesh().magSf().boundaryField()[patchi][patchFacei];
160  }
161 
162  return true;
163 }
164 
165 
166 // ************************************************************************* //
faceListList boundary
dictionary dict
DSMCCloud< dsmcParcel > CloudType
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Template invariant parts for fvPatchField.
Definition: fvPatchField.H:77
virtual bool postPatch(const parcelType &p, const polyPatch &pp, const typename parcelType::trackingData &td)
Post-patch hook.
Operations on lists of strings.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
void write()
Write post-processing info.
Various functions to operate on Lists.
const CloudType & owner() const
Return const access to the owner cloud.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const word calculatedType
A calculated patch field type.
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:37
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
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:69
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
Templated cloud function object base class.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127