LagrangianModel.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::LagrangianModel
26 
27 Description
28  Base class for Lagrangian models
29 
30 SourceFiles
31  LagrangianModelI.H
32  LagrangianModel.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef LagrangianModel_H
37 #define LagrangianModel_H
38 
39 #include "LagrangianEqn.H"
40 #include "LagrangianModelM.H"
41 #include "stateModel.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 class LagrangianMesh;
49 class polyTopoChangeMap;
50 class polyMeshMap;
51 class polyDistributionMap;
52 
53 /*---------------------------------------------------------------------------*\
54  Class LagrangianModel Declaration
55 \*---------------------------------------------------------------------------*/
56 
57 class LagrangianModel
58 :
59  public stateModel
60 {
61 public:
62 
63  // Public Enumerations
64 
65  //- Enumeration of the types of instantaneous modification
66  enum class modification : label
67  {
68  change = 0,
69  remove = 1
70  };
71 
72 
73  // Public Type Definitions
74 
75  //- Class containing an element-index and a modification-enumeration
77 
78 
79 private:
80 
81  // Private Member Data
82 
83  //- Model name
84  const word name_;
85 
86  //- Reference to the mesh
87  const LagrangianMesh& mesh_;
88 
89 
90 protected:
91 
92  // Protected Member Functions
93 
94  //- Add a source term to an equation
95  template<class Type>
96  void addSupType
97  (
98  const LagrangianSubScalarField& deltaT,
99  const LagrangianSubSubField<Type>& field,
101  ) const;
102 
103  //- Add a source term to a mass-weighted equation
104  template<class Type>
105  void addSupType
106  (
107  const LagrangianSubScalarField& deltaT,
109  const LagrangianSubSubField<Type>& field,
111  ) const;
112 
113 
114 public:
115 
116  //- Runtime type information
117  TypeName("LagrangianModel");
118 
119 
120  //- Declare run-time constructor selection table
122  (
123  autoPtr,
125  dictionary,
126  (
127  const word& name,
128  const LagrangianMesh& mesh,
129  const dictionary& modelDict,
130  const dictionary& stateDict
131  ),
132  (name, mesh, modelDict, stateDict)
133  );
134 
135 
136  // Static Member Functions
137 
138  //- Return the name of the field associated with a source term
139  template<class AlphaRhoFieldType, class ... AlphaRhoFieldTypes>
140  static word fieldName
141  (
142  const AlphaRhoFieldType& alphaRhoField,
143  const AlphaRhoFieldTypes& ... alphaRhoFields
144  );
145 
146  //- Return the name of the field associated with a source term (base
147  // condition for the above)
148  template<class AlphaRhoFieldType>
149  static word fieldName(const AlphaRhoFieldType& alphaRhoField);
150 
151  //- Return the name of the product of the fields associated with a
152  // source term
153  template<class AlphaRhoFieldType, class ... AlphaRhoFieldTypes>
154  static word fieldsName
155  (
156  const AlphaRhoFieldType& alphaRhoField,
157  const AlphaRhoFieldTypes& ... alphaRhoFields
158  );
159 
160  //- Return the name of the product of the fields associated with a
161  // source term (base condition for the above)
162  template<class AlphaRhoFieldType>
163  static word fieldsName(const AlphaRhoFieldType& alphaRhoField);
164 
165 
166  // Constructors
167 
168  //- Construct from components
170  (
171  const word& name,
172  const LagrangianMesh& mesh
173  );
174 
175  //- Disallow default bitwise copy construction
176  LagrangianModel(const LagrangianModel&) = delete;
177 
178  //- Clone
180  {
182  return autoPtr<LagrangianModel>(nullptr);
183  }
184 
185  //- List construction class
186  class iNew
187  {
188  const word& name_;
189 
190  const LagrangianMesh& mesh_;
191 
192  public:
193 
194  iNew(const word& name, const LagrangianMesh& mesh)
195  :
196  name_(name),
197  mesh_(mesh)
198  {}
199 
201  {
203  (
204  LagrangianModel::New(name_, mesh_, dictionary(is))
205  );
206  }
207  };
208 
209 
210  //- Selector
212  (
213  const word& name,
214  const LagrangianMesh& mesh,
215  const dictionary& modelDict
216  );
217 
218 
219  //- Destructor
220  virtual ~LagrangianModel();
221 
222 
223  // Member Functions
224 
225  // Access
226 
227  //- The source name
228  inline const word& name() const;
229 
230  //- The database
231  inline const objectRegistry& db() const;
232 
233  //- The mesh
234  inline const LagrangianMesh& mesh() const;
235 
236 
237  // Checks
238 
239  //- Return the list of fields for which the LagrangianModel adds
240  // source term to the transport equation
241  virtual wordList addSupFields() const;
242 
243  //- Return true if the LagrangianModel adds a source term to the
244  // given field's transport equation
245  virtual bool addsSupToField(const word&) const;
246 
247  //- Return true if the LagrangianModels adds a source term to the
248  // given field's transport equation
249  template<class Type, template<class> class PrimitiveField>
250  bool addsSupToField
251  (
253  ) const;
254 
255 
256  //- Do post construction steps which require access to other models
257  virtual void postConstruct();
258 
259  //- Correct the LagrangianModel
260  virtual void correct();
261 
262  //- Identify elements in the Lagrangian mesh which are to be
263  // instantaneously modified or removed
264  virtual void preModify
265  (
266  const LagrangianMesh& mesh,
267  DynamicList<elementModification>& elementModifications
268  ) const;
269 
270  //- Instantaneously modify and/or create and remove elements in the
271  // Lagrangian mesh
272  virtual LagrangianSubMesh modify
273  (
275  const LagrangianSubMesh& modifiedMesh
276  ) const;
277 
278  //- Solve equations and/or update continually changing properties
279  virtual void calculate
280  (
281  const LagrangianSubScalarField& deltaT,
282  const bool final
283  );
284 
285 
286  // Sources
287 
288  //- Add a fractional source term
289  virtual void addSup
290  (
291  const LagrangianSubScalarField& deltaT,
293  ) const;
294 
295  //- Add a source term to an equation
297 
298  //- Add a source term to a mass-weighted equation
300 
301 
302  // Mesh changes
303 
304  //- Update topology using the given map
305  virtual void topoChange(const polyTopoChangeMap&);
306 
307  //- Update from another mesh using the given map
308  virtual void mapMesh(const polyMeshMap&);
309 
310  //- Redistribute or update using the given distribution map
311  virtual void distribute(const polyDistributionMap&);
312 
313 
314  // IO
315 
316  //- Read dictionary
317  virtual bool read(const dictionary& modelDict);
318 
319  //- Write data
320  virtual bool write(const bool write) const;
321 
322 
323  // Member Operators
324 
325  //- Disallow default bitwise assignment
326  void operator=(const LagrangianModel&) = delete;
327 };
328 
329 
330 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
331 
332 } // End namespace Foam
333 
334 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
335 
336 #include "LagrangianModelI.H"
337 
338 #ifdef NoRepository
339  #include "LagrangianModelTemplates.C"
340 #endif
341 
342 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
343 
344 #endif
345 
346 // ************************************************************************* //
#define DEFINE_LAGRANGIAN_MODEL_ADD_M_FIELD_SUP(Type, nullArg)
#define DEFINE_LAGRANGIAN_MODEL_ADD_FIELD_SUP(Type, nullArg)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
This class stores the coefficients of a Lagrangian equation, and facilitates solving that equation an...
Definition: LagrangianEqn.H:56
Class to hold element-group indices, and associate the group.
Class containing Lagrangian geometry and topology.
List construction class.
iNew(const word &name, const LagrangianMesh &mesh)
autoPtr< LagrangianModel > operator()(Istream &is) const
Base class for Lagrangian models.
virtual ~LagrangianModel()
Destructor.
virtual wordList addSupFields() const
Return the list of fields for which the LagrangianModel adds.
virtual void addSup(const LagrangianSubScalarField &deltaT, LagrangianSubScalarField &S) const
Add a fractional source term.
virtual void postConstruct()
Do post construction steps which require access to other models.
modification
Enumeration of the types of instantaneous modification.
void addSupType(const LagrangianSubScalarField &deltaT, const LagrangianSubSubField< Type > &field, LagrangianEqn< Type > &eqn) const
Add a source term to an equation.
virtual void preModify(const LagrangianMesh &mesh, DynamicList< elementModification > &elementModifications) const
Identify elements in the Lagrangian mesh which are to be.
virtual void correct()
Correct the LagrangianModel.
virtual bool addsSupToField(const word &) const
Return true if the LagrangianModel adds a source term to the.
declareRunTimeSelectionTable(autoPtr, LagrangianModel, dictionary,(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict, const dictionary &stateDict),(name, mesh, modelDict, stateDict))
Declare run-time constructor selection table.
virtual void topoChange(const polyTopoChangeMap &)
Add a source term to an equation.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
const objectRegistry & db() const
The database.
void operator=(const LagrangianModel &)=delete
Disallow default bitwise assignment.
TypeName("LagrangianModel")
Runtime type information.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
autoPtr< LagrangianModel > clone() const
Clone.
virtual bool read(const dictionary &modelDict)
Read dictionary.
static word fieldsName(const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
Return the name of the product of the fields associated with a.
const LagrangianMesh & mesh() const
The mesh.
virtual bool write(const bool write) const
Write data.
LagrangianModel(const word &name, const LagrangianMesh &mesh)
Construct from components.
virtual LagrangianSubMesh modify(LagrangianMesh &mesh, const LagrangianSubMesh &modifiedMesh) const
Instantaneously modify and/or create and remove elements in the.
static word fieldName(const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
Return the name of the field associated with a source term.
const word & name() const
The source name.
LagrangianMesh::elementGroup< modification > elementModification
Class containing an element-index and a modification-enumeration.
virtual void calculate(const LagrangianSubScalarField &deltaT, const bool final)
Solve equations and/or update continually changing properties.
static autoPtr< LagrangianModel > New(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict)
Selector.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Registry of regIOobjects.
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.
Base class for models with state.
Definition: stateModel.H:49
static dictionary stateDict(const word &name, const objectRegistry &db)
Construct and return the state dictionary for reading.
Definition: stateModel.C:137
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
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
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)