cloudBoundaryCollisionFlux.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "cloud.H"
28 #include "CompactListList.H"
30 #include "functionName.H"
31 #include "functionObject.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const word& name,
49  const Time& runTime,
50  const dictionary& dict,
51  const word& phiName,
52  const dimensionSet& phiDims
53 )
54 :
56  faceiPtr_(nullptr),
57  qPtr_(nullptr),
58  phiName_(phiName),
59  phiDims_(phiDims),
60  phib_
61  (
62  mesh().boundary(),
63  surfaceScalarField::Internal::null(),
64  calculatedFvsPatchField<scalar>::typeName
65  )
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
70 
72 {}
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
78 {
79  return wordList::null();
80 }
81 
82 
84 {
85  return false;
86 }
87 
88 
90 {
91  return true;
92 }
93 
94 
96 {
97  phib_ == Zero;
98 }
99 
100 
102 (
104 )
105 {
106  faceiPtr_.set
107  (
109  (
110  IOobject
111  (
112  name() + ":facei",
113  time_.name(),
114  cloud().mesh()
115  ),
116  cloud().mesh(),
118  )
119  );
120 
121  qPtr_.set
122  (
124  (
125  IOobject
126  (
127  name() + ":q",
128  time_.name(),
129  cloud().mesh()
130  ),
131  cloud().mesh(),
132  dimensionedScalar(phiDims_*dimTime, scalar(0))
133  )
134  );
135 }
136 
137 
139 (
140  const LagrangianSubScalarSubField& fraction
141 )
142 {
143  const LagrangianSubMesh& subMesh = fraction.mesh();
144 
145  if
146  (
147  faceiPtr_.empty()
148  || subMesh.group() == LagrangianGroup::none
149  || subMesh.group() == LagrangianGroup::inInternalMesh
150  ) return;
151 
152  const Foam::cloud& c = cloud();
153 
154  SubField<label> facei(subMesh.sub(faceiPtr_->primitiveFieldRef()));
155  facei = subMesh.sub(c.mesh().facei());
156 
157  LagrangianSubScalarSubField q(subMesh.sub(qPtr_()));
158  q += this->q(fraction, +1);
159 }
160 
161 
163 (
164  const LagrangianSubScalarSubField& fraction
165 )
166 {
167  const LagrangianSubMesh& subMesh = fraction.mesh();
168 
169  if
170  (
171  faceiPtr_.empty()
172  || subMesh.group() == LagrangianGroup::none
173  || subMesh.group() == LagrangianGroup::inInternalMesh
174  ) return;
175 
176  const Foam::cloud& c = cloud();
177 
178  const SubField<label> facei0(subMesh.sub(faceiPtr_->primitiveFieldRef()));
179  const SubField<label> facei = subMesh.sub(c.mesh().facei());
180 
181  LagrangianSubScalarSubField q(subMesh.sub(qPtr_()));
182  q += this->q(fraction, -1);
183 
185 
186  forAll(subMesh, subi)
187  {
188  if (facei0[subi] != facei[subi]) continue;
189 
190  const LagrangianState state = c.mesh().state(subi + subMesh.start());
191 
192  if (state != LagrangianState::inCell) continue;
193 
194  const label bFacei = facei[subi] - mesh().nInternalFaces();
195 
196  const labelUList patchis = mesh().polyBFacePatches()[bFacei];
197  const labelUList patchFaceis = mesh().polyBFacePatchFaces()[bFacei];
198 
199  forAll(patchis, i)
200  {
201  phib_[patchis[i]][patchFaceis[i]] +=
202  q[subi]
203  /time_.deltaTValue()
204  *magSfb[patchis[i]][patchFaceis[i]]
205  /mesh().magFaceAreas()[facei[subi]];
206  }
207  }
208 }
209 
210 
212 (
214 )
215 {
216  faceiPtr_.clear();
217  qPtr_.clear();
218 }
219 
220 
222 {
223  return
225  (
226  IOobject
227  (
228  cloud().mesh().name() + ":" + phiName_ + "Coll",
229  time_.name(),
230  mesh(),
233  ),
234  mesh(),
235  phiDims_,
236  scalarField(mesh().nInternalFaces(), scalar(0)),
237  phib_
238  ).write();
239 }
240 
241 
243 {
244  return true;
245 }
246 
247 
248 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
SubList< Type > sub(const List< Type > &list) const
Return a sub-list corresponding to this sub-mesh.
LagrangianGroup group() const
Return the group.
label start() const
Return start.
static const List< word > & null()
Return a null List.
Definition: ListI.H:118
Pre-declare related SubField type.
Definition: SubField.H:63
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
Foam::calculatedFvsPatchField.
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
Definition: cloud.H:63
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
Generic dimensioned Type class.
Base class for functions which generate a boundary collision flux for a cloud.
virtual wordList fields() const
Return the list of fields required.
cloudBoundaryCollisionFlux(const word &name, const Time &runTime, const dictionary &dict, const word &phiName, const dimensionSet &phiDims)
Construct from Time and dictionary.
virtual bool executeAtStart() const
Return false so this function does not execute at the start.
virtual void postCrossFaces(const LagrangianSubScalarSubField &fraction)
Hook following face crossings of a specific sub-mesh.
virtual void preSolve()
Hook before solution steps.
virtual bool execute()
Do nothing. Everything happens in faces crossing hooks.
virtual void preCrossFaces(const LagrangianScalarInternalDynamicField &fraction)
Hook before all face crossings.
Base class for function objects that refer to an fvMesh and a cloud. Used, for example,...
Function object that solves for the evolution of a cloud. Only provides one-way coupling with a finit...
const UCompactListList< label > & polyBFacePatches() const
Return poly-bFace-patch addressing.
Definition: fvMesh.C:1016
const UCompactListList< label > & polyBFacePatchFaces() const
Return poly-bFace-patch-face addressing.
Definition: fvMesh.C:1092
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
label nInternalFaces() const
const scalarField & magFaceAreas() const
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const dimensionedScalar c
Speed of light in a vacuum.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
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:59
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
LagrangianState
Lagrangian state enumeration.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTime
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
faceListList boundary(nPatches)
dictionary dict