Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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
29
30// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31
32template<class CloudType>
35({
36 { resetMode::none, "none" },
37 { resetMode::timeStep, "timeStep" },
38 { resetMode::writeTime, "writeTime" },
39});
40
41// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
42
43template<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
85template<class CloudType>
88 clearOrReset(massPtr_, "mass", dimMass);
89 clearOrReset(countPtr_, "count", dimless);
90}
91
92
93template<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();
120}
121
122
123// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
124
125template<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 )
145{
146 reset();
147}
148
149
150template<class CloudType>
152(
154)
155:
157 massPtr_(nullptr),
158 countPtr_(nullptr),
159 resetMode_(pii.resetMode_)
160{}
161
162
163// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164
165template<class CloudType>
167(
168 const typename parcelType::trackingData&
169)
170{
171 if (resetMode_ == resetMode::timeStep)
173 reset();
174 }
175}
176
177
178template<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// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Templated cloud function object base class.
CloudFunctionObject(CloudType &owner)
Construct null from owner.
const CloudType & owner() const
Return const access to the owner cloud.
particle::trackingData trackingData
Definition DSMCParcel.H:155
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
@ NO_REGISTER
Do not request registration (bool: false).
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
Creates volume fields whose boundaries are used to store patch interaction statistics.
void clearOrReset(autoPtr< volScalarField > &fieldPtr, const word &fieldName, const dimensionSet &dims) const
Helper function to clear or reset fields.
static const Enum< resetMode > resetModeNames_
virtual void write()
Write post-processing info.
PatchInteractionFields(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual void preEvolve(const typename parcelType::trackingData &td)
Pre-evolve hook.
void reset()
Create|read|reset the fields.
virtual bool postPatch(const parcelType &p, const polyPatch &pp, const typename parcelType::trackingData &td)
Post-patch hook.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition autoPtrI.H:37
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
limits reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL))
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
DSMCCloud< dsmcParcel > CloudType
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
dictionary dict