PatchInteractionFields.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) 2020-2023 OpenCFD Ltd.
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 "PatchInteractionFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 template<class CloudType>
35 ({
36  { resetMode::none, "none" },
37  { resetMode::timeStep, "timeStep" },
38  { resetMode::writeTime, "writeTime" },
39 });
40 
41 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
42 
43 template<class CloudType>
45 (
46  autoPtr<volScalarField>& fieldPtr,
47  const word& fieldName,
48  const dimensionSet& dims
49 ) const
50 {
51  if (fieldPtr)
52  {
53  fieldPtr->primitiveFieldRef() = scalar(0);
54  }
55  else
56  {
57  const fvMesh& mesh = this->owner().mesh();
58 
59  fieldPtr.reset
60  (
61  new volScalarField
62  (
63  IOobject
64  (
65  IOobject::scopedName
66  (
67  this->owner().name(),
68  this->modelName(),
69  fieldName
70  ),
71  mesh.time().timeName(),
72  mesh,
73  IOobject::READ_IF_PRESENT,
74  IOobject::NO_WRITE,
75  IOobject::NO_REGISTER
76  ),
77  mesh,
78  dimensionedScalar(dims, Zero)
79  )
80  );
81  }
82 }
83 
84 
85 template<class CloudType>
87 {
88  clearOrReset(massPtr_, "mass", dimMass);
89  clearOrReset(countPtr_, "count", dimless);
90 }
91 
92 
93 template<class CloudType>
95 {
96  if (massPtr_)
97  {
98  massPtr_->write();
99  }
100  else
101  {
103  << "massPtr not valid" << abort(FatalError);
104  }
105 
106  if (countPtr_)
107  {
108  countPtr_->write();
109  }
110  else
111  {
113  << "countPtr not valid" << abort(FatalError);
114  }
115 
116  if (resetMode_ == resetMode::writeTime)
117  {
118  reset();
119  }
120 }
121 
122 
123 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
124 
125 template<class CloudType>
127 (
128  const dictionary& dict,
129  CloudType& owner,
130  const word& modelName
131 )
132 :
133  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
134  massPtr_(nullptr),
135  countPtr_(nullptr),
136  resetMode_
137  (
138  resetModeNames_.getOrDefault
139  (
140  "resetMode",
141  this->coeffDict(),
142  resetMode::none
143  )
144  )
145 {
146  reset();
147 }
148 
149 
150 template<class CloudType>
152 (
154 )
155 :
157  massPtr_(nullptr),
158  countPtr_(nullptr),
159  resetMode_(pii.resetMode_)
160 {}
161 
162 
163 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164 
165 template<class CloudType>
167 (
168  const typename parcelType::trackingData&
169 )
170 {
171  if (resetMode_ == resetMode::timeStep)
172  {
173  reset();
174  }
175 }
176 
177 
178 template<class CloudType>
180 (
181  const parcelType& p,
182  const polyPatch& pp,
183  const typename parcelType::trackingData& td
184 )
185 {
186  const label patchi = pp.index();
187  const label facei = pp.whichFace(p.face());
188 
189  massPtr_->boundaryFieldRef()[patchi][facei] += p.nParticle()*p.mass();
190  countPtr_->boundaryFieldRef()[patchi][facei] += 1;
191 
192  return true;
193 }
194 
195 
196 // ************************************************************************* //
dictionary dict
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
virtual bool postPatch(const parcelType &p, const polyPatch &pp, const typename parcelType::trackingData &td)
Post-patch hook.
const dimensionSet dimless
Dimensionless.
void clearOrReset(autoPtr< volScalarField > &fieldPtr, const word &fieldName, const dimensionSet &dims) const
Helper function to clear or reset fields.
static const Enum< resetMode > resetModeNames_
PatchInteractionFields(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Creates volume fields whose boundaries are used to store patch interaction statistics.
A class for handling words, derived from Foam::string.
Definition: word.H:63
errorManip< error > abort(error &err)
Definition: errorManip.H:139
virtual void write()
Write post-processing info.
virtual void preEvolve(const typename parcelType::trackingData &td)
Pre-evolve hook.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void reset()
Create|read|reset the fields.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
Templated cloud function object base class.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127