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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "PatchCollisionDensity.H"
29 #include "Pstream.H"
30 #include "stringListOps.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 
42  (
43  IOobject
44  (
45  this->owner().name() + ":collisionDensity",
46  this->owner().mesh().time().timeName(),
47  this->owner().mesh()
48  ),
49  this->owner().mesh(),
51  z,
52  collisionDensity_
53  )
54  .write();
55 
57  (
58  IOobject
59  (
60  this->owner().name() + ":collisionDensityRate",
61  this->owner().mesh().time().timeName(),
62  this->owner().mesh()
63  ),
64  this->owner().mesh(),
66  z,
67  (collisionDensity_ - collisionDensity0_)
68  /(this->owner().mesh().time().value() - time0_)
69  )
70  .write();
71 
72  collisionDensity0_ == collisionDensity_;
73  time0_ = this->owner().mesh().time().value();
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78 
79 template<class CloudType>
81 (
82  const dictionary& dict,
83  CloudType& owner,
84  const word& modelName
85 )
86 :
87  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
88  minSpeed_(dict.lookupOrDefault<scalar>("minSpeed", -1)),
89  collisionDensity_
90  (
91  this->owner().mesh().boundary(),
92  volScalarField::Internal::null(),
94  ),
95  collisionDensity0_
96  (
97  this->owner().mesh().boundary(),
98  volScalarField::Internal::null(),
100  ),
101  time0_(this->owner().mesh().time().value())
102 {
103  collisionDensity_ == 0;
104  collisionDensity0_ == 0;
105 
106  IOobject io
107  (
108  this->owner().name() + ":collisionDensity",
109  this->owner().mesh().time().timeName(),
110  this->owner().mesh(),
111  IOobject::MUST_READ,
112  IOobject::NO_WRITE
113  );
114 
115  if (io.typeHeaderOk<volScalarField>())
116  {
117  const volScalarField collisionDensity(io, this->owner().mesh());
118  collisionDensity_ == collisionDensity.boundaryField();
119  collisionDensity0_ == collisionDensity.boundaryField();
120  }
121 }
122 
123 
124 template<class CloudType>
126 (
128 )
129 :
131  minSpeed_(ppm.minSpeed_),
132  collisionDensity_(ppm.collisionDensity_),
133  collisionDensity0_(ppm.collisionDensity0_),
134  time0_(ppm.time0_)
135 {}
136 
137 
138 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
139 
140 template<class CloudType>
142 {}
143 
144 
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
146 
147 template<class CloudType>
149 (
150  const parcelType& p,
151  const polyPatch& pp,
152  bool&
153 )
154 {
155  const label patchi = pp.index();
156  const label patchFacei = p.face() - pp.start();
157 
158  vector nw, Up;
159  this->owner().patchData(p, pp, nw, Up);
160 
161  const scalar speed = (p.U() - Up) & nw;
162  if (speed > minSpeed_)
163  {
164  collisionDensity_[patchi][patchFacei] +=
165  1/this->owner().mesh().magSf().boundaryField()[patchi][patchFacei];
166  }
167 }
168 
169 
170 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::PatchCollisionDensity::write
void write()
Write post-processing info.
Definition: PatchCollisionDensity.C:37
PatchCollisionDensity.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::IOobject::typeHeaderOk
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
Definition: IOobjectTemplates.C:39
ListListOps.H
Foam::PatchCollisionDensity
Function object which generates fields of the number and rate of collisions per unit area on all patc...
Definition: PatchCollisionDensity.H:62
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:404
Foam::PatchCollisionDensity::postPatch
virtual void postPatch(const parcelType &p, const polyPatch &pp, bool &keepParticle)
Post-patch hook.
Definition: PatchCollisionDensity.C:149
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::label
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:62
Foam::Field< scalar >
Foam::fvMesh::magSf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Definition: fvMeshGeometry.C:347
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::PatchCollisionDensity::~PatchCollisionDensity
virtual ~PatchCollisionDensity()
Destructor.
Definition: PatchCollisionDensity.C:141
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Pstream.H
Foam::calculatedFvPatchField
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Definition: calculatedFvPatchField.H:66
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:311
Foam::PatchCollisionDensity::PatchCollisionDensity
PatchCollisionDensity(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: PatchCollisionDensity.C:81
Foam::Vector< scalar >
Foam::CloudFunctionObject
Templated cloud function object base class.
Definition: CloudFunctionObject.H:61
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
stringListOps.H
Operations on lists of strings.
ListOps.H
Various functions to operate on Lists.
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::patchIdentifier::index
label index() const
The index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:133
boundary
faceListList boundary
Definition: createBlockMesh.H:4
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62