LagrangianModels.H
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 Class
25  Foam::LagrangianModels
26 
27 Description
28  List of Lagrangian models, constructed as a (Lagrangian) mesh object.
29  Provides similar functions to the models themselves and forwards them to
30  each model in turn. This is the high level model interface used by clouds
31  when constructing their injections and transport equations.
32 
33 SourceFiles
34  LagrangianModels.C
35  LagrangianModelsTemplates.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef LagrangianModels_H
40 #define LagrangianModels_H
41 
42 #include "LagrangianModel.H"
43 #include "PtrListDictionary.H"
44 #include "DemandDrivenMeshObject.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class LagrangianModels Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 class LagrangianModels
56 :
58  <
59  LagrangianMesh,
60  TopoChangeableMeshObject,
61  LagrangianModels
62  >,
63  public dictionary,
64  public PtrListDictionary<LagrangianModel>
65 {
66  // Private Member Data
67 
68  //- Time index to check that all defined sources have been applied
69  mutable label checkTimeIndex_;
70 
71  //- Sets of the fields that have had sources added by the
72  // LagrangianModels
73  mutable PtrList<wordHashSet> addSupFields_;
74 
75 
76  // Private Member Functions
77 
78  //- Create IO object for an optional LagrangianModels file
79  IOobject io(const LagrangianMesh& mesh) const;
80 
81  //- Helper for modelTypeFieldSourceTypes. Inserts the field source type
82  // into to the table if the model is of the required type.
83  template<class ... ModelAndFieldSourceTypes>
84  struct modelTypeFieldSourceType;
85 
86  //- Check that all LagrangianModels have been applied
87  void checkApplied() const;
88 
89  //- Return a source for an equation
90  template
91  <
92  class Type,
93  template<class> class PrimitiveField,
94  class ... AlphaRhoFieldTypes
95  >
96  tmp<LagrangianEqn<Type>> sourceTerm
97  (
99  const LagrangianSubScalarField& deltaT,
100  const AlphaRhoFieldTypes& ... alphaRhoFields
101  ) const;
102 
103 
104 protected:
105 
106  //- Let the mesh object use the protected constructors
107  friend class DemandDrivenMeshObject
108  <
112  >;
113 
114 
115  // Protected Constructors
116 
117  //- Construct for a mesh
118  explicit LagrangianModels(const LagrangianMesh& mesh);
119 
120 
121 public:
122 
123  //- Runtime type information
124  TypeName("LagrangianModels");
125 
126 
127  // Constructors
128 
129  //- Disallow default bitwise copy construction
130  LagrangianModels(const LagrangianModels&) = delete;
131 
132  //- Inherit the base New method
134  <
138  >::New;
139 
140 
141  //- Destructor
142  virtual ~LagrangianModels();
143 
144 
145  // Member Functions
146 
147  //- Declare LagrangianModels to be a global dictionary
148  virtual bool global() const
149  {
150  return true;
151  }
152 
153 
154  // Checks
155 
156  //- Return true if the LagrangianModels adds a source term to the
157  // given field's transport equation
158  bool addsSupToField(const word&) const;
159 
160  //- Return true if the LagrangianModels adds a source term to the
161  // given field's transport equation
162  template<class Type, template<class> class PrimitiveField>
163  bool addsSupToField
164  (
166  ) const;
167 
168 
169  //- Correct the LagrangianModels
170  void correct();
171 
172  //- Identify elements in the Lagrangian mesh which are to be
173  // instantaneously modified and put them in contiguous groups
175 
176  //- Instantaneously modify and/or create and remove elements in the
177  // Lagrangian mesh
179  (
181  const LagrangianSubMesh& modifiedMesh
182  ) const;
183 
184  //- Solve equations and/or update continually changing properties
185  void calculate
186  (
187  const LagrangianSubScalarField& deltaT,
188  const bool final
189  );
190 
191  //- Return a table of field source types that are chosen to match given
192  // model types. So, e.g., a zero field source for every injection
193  // model. Useful for programmatically constructing fields. Template
194  // arguments should alternate between model and field types; i.e.,
195  // model-type-1, field-source-type-1, model-type-2,
196  // field-source-type-2, model-type-3, ...
197  template<class ... ModelAndFieldSourceTypes>
199 
200 
201  // Sources
202 
203  //- Return the fractional source
205  (
206  const LagrangianSubScalarField& deltaT
207  ) const;
208 
209  //- Return source for an equation
210  template<class Type, template<class> class PrimitiveField>
212  (
213  const LagrangianSubScalarField& deltaT,
215  ) const;
216 
217  //- Return source for an equation
218  template
219  <
220  class Type,
221  template<class> class PrimitiveField,
222  template<class> class PrimitiveEqnField
223  >
225  (
226  const LagrangianSubScalarField& deltaT,
229  ) const;
230 
231  //- Return source for a mass-weighted equation
232  template<class Type, template<class> class PrimitiveField>
234  (
235  const LagrangianSubScalarField& deltaT,
238  ) const;
239 
240  //- Return source for a mass-weighted equation
241  template
242  <
243  class Type,
244  template<class> class PrimitiveField,
245  template<class> class PrimitiveEqnField
246  >
248  (
249  const LagrangianSubScalarField& deltaT,
253  ) const;
254 
255 
256  // Mesh changes
257 
258  //- Update for mesh motion. Only for mesh object. Does nothing.
259  // Lagrangian evolves continuously across the range of mesh
260  // motion, so no instantaneous update is needed.
261  virtual bool movePoints();
262 
263  //- Update topology using the given map
264  virtual void topoChange(const polyTopoChangeMap&);
265 
266  //- Update from another mesh using the given map
267  virtual void mapMesh(const polyMeshMap&);
268 
269  //- Redistribute or update using the given distribution map
270  virtual void distribute(const polyDistributionMap&);
271 
272 
273  // IO
274 
275  //- Read the LagrangianModels dictionary if it has changed
276  // and update the models
277  virtual bool read();
278 
279  //- ReadData function which reads the LagrangianModels dictionary
280  // required for regIOobject read operation
281  virtual bool readData(Istream&);
282 
283  //- writeData function required by regIOobject but not used
284  // for this class, writeObject is used instead
285  virtual bool writeData(Ostream& os) const;
286 
287  //- Write the LagrangianModels
288  virtual bool writeObject
289  (
293  const bool write
294  ) const;
295 
296 
297  // Member Operators
298 
299  //- Inherit the PtrListDictionary index operators
301 
302  //- Inherit the PtrListDictionary size function
304 
305  //- Disallow default bitwise assignment
306  void operator=(const LagrangianModels&) = delete;
307 };
308 
309 
310 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311 
312 //- Trait for obtaining global status
313 template<>
315 {
316  static const bool global = true;
317 };
318 
319 
320 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
321 
322 #ifdef NoRepository
323  #include "LagrangianModelsTemplates.C"
324 #endif
325 
326 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 
328 } // End namespace Foam
329 
330 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
331 
332 #endif
333 
334 // ************************************************************************* //
Templated abstract base-class for demand-driven mesh objects used to automate their allocation to the...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An STL-conforming hash table.
Definition: HashTable.H:127
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Version number type.
Definition: IOstream.H:97
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:87
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:194
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Class containing Lagrangian geometry and topology.
List of Lagrangian models, constructed as a (Lagrangian) mesh object. Provides similar functions to t...
virtual bool movePoints()
Update for mesh motion. Only for mesh object. Does nothing.
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write) const
Write the LagrangianModels.
void operator=(const LagrangianModels &)=delete
Disallow default bitwise assignment.
tmp< LagrangianSubScalarField > source(const LagrangianSubScalarField &deltaT) const
Return the fractional source.
void correct()
Correct the LagrangianModels.
bool addsSupToField(const word &) const
Return true if the LagrangianModels adds a source term to the.
virtual bool writeData(Ostream &os) const
writeData function required by regIOobject but not used
HashTable< word > modelTypeFieldSourceTypes() const
Return a table of field source types that are chosen to match given.
virtual bool global() const
Declare LagrangianModels to be a global dictionary.
virtual ~LagrangianModels()
Destructor.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
TypeName("LagrangianModels")
Runtime type information.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
LagrangianSubMesh preModify(LagrangianMesh &mesh) const
Identify elements in the Lagrangian mesh which are to be.
tmp< LagrangianEqn< Type > > sourceProxy(const LagrangianSubScalarField &deltaT, const LagrangianSubField< Type, PrimitiveField > &field, const LagrangianSubField< Type, PrimitiveEqnField > &eqnField) const
Return source for an equation.
virtual bool readData(Istream &)
ReadData function which reads the LagrangianModels dictionary.
LagrangianSubMesh modify(LagrangianMesh &mesh, const LagrangianSubMesh &modifiedMesh) const
Instantaneously modify and/or create and remove elements in the.
void calculate(const LagrangianSubScalarField &deltaT, const bool final)
Solve equations and/or update continually changing properties.
LagrangianModels(const LagrangianMesh &mesh)
Construct for a mesh.
virtual bool read()
Read the LagrangianModels dictionary if it has changed.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Template dictionary class which manages the storage associated with it.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
void write(Ostream &, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:226
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
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
Trait for obtaining global status.
Definition: IOobject.H:503
static const bool global
Definition: IOobject.H:504