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-2023 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 "ListOps.H"
32 #include "ListListOps.H"
33 
34 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
35 
36 template<class CloudType>
38 {
39  const scalarField z(this->owner().mesh().nCells(), Zero);
40 
41  {
43  (
44  this->owner().mesh().newIOobject
45  (
46  IOobject::scopedName
47  (
48  this->owner().name(),
49  "collisionDensity"
50  )
51  ),
52  this->owner().mesh(),
54  z,
55  collisionDensity_
56  ).write();
57  }
58 
59  {
61  (
62  this->owner().mesh().newIOobject
63  (
64  IOobject::scopedName
65  (
66  this->owner().name(),
67  "collisionDensityRate"
68  )
69  ),
70  this->owner().mesh(),
72  z,
73  (collisionDensity_ - collisionDensity0_)
74  /(this->owner().mesh().time().value() - time0_)
75  ).write();
76  }
77 
78  collisionDensity0_ == collisionDensity_;
79  time0_ = this->owner().mesh().time().value();
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
85 template<class CloudType>
87 (
88  const dictionary& dict,
89  CloudType& owner,
90  const word& modelName
91 )
92 :
93  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
94  minSpeed_(dict.getOrDefault<scalar>("minSpeed", -1)),
95  collisionDensity_
96  (
97  this->owner().mesh().boundary(),
98  volScalarField::Internal::null(),
100  ),
101  collisionDensity0_
102  (
103  this->owner().mesh().boundary(),
104  volScalarField::Internal::null(),
106  ),
107  time0_(this->owner().mesh().time().value())
108 {
109  collisionDensity_ == 0;
110  collisionDensity0_ == 0;
111 
112  IOobject io
113  (
114  this->owner().mesh().newIOobject
115  (
116  IOobject::scopedName(this->owner().name(), "collisionDensity"),
118  )
119  );
120 
121  if (io.typeHeaderOk<volScalarField>())
122  {
123  const volScalarField collisionDensity(io, this->owner().mesh());
124  collisionDensity_ == collisionDensity.boundaryField();
125  collisionDensity0_ == collisionDensity.boundaryField();
126  }
127 }
128 
129 
130 template<class CloudType>
132 (
133  const PatchCollisionDensity<CloudType>& ppm
134 )
135 :
136  CloudFunctionObject<CloudType>(ppm),
137  minSpeed_(ppm.minSpeed_),
138  collisionDensity_(ppm.collisionDensity_),
139  collisionDensity0_(ppm.collisionDensity0_),
140  time0_(ppm.time0_)
141 {}
142 
143 
144 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
145 
146 template<class CloudType>
148 (
149  const parcelType& p,
150  const polyPatch& pp,
151  const typename parcelType::trackingData& td
152 )
153 {
154  const label patchi = pp.index();
155  const label patchFacei = p.face() - pp.start();
156 
157  vector nw, Up;
158  this->owner().patchData(p, pp, nw, Up);
159 
160  const scalar speed = (p.U() - Up) & nw;
161  if (speed > minSpeed_)
162  {
163  collisionDensity_[patchi][patchFacei] +=
164  1/this->owner().mesh().magSf().boundaryField()[patchi][patchFacei];
165  }
166 
167  return true;
168 }
169 
170 
171 // ************************************************************************* //
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
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
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.
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:72
const word calculatedType
A calculated patch field type.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
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:180
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