LagrangianModels.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 
26 #include "LagrangianModels.H"
27 #include "LagrangianMesh.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
34 }
35 
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 Foam::IOobject Foam::LagrangianModels::io(const LagrangianMesh& mesh) const
40 {
41  typeIOobject<IOdictionary> result
42  (
43  typeName,
44  mesh.time().constant(),
45  mesh,
48  );
49 
50  if (!result.headerOk())
51  {
52  result.readOpt() = IOobject::NO_READ;
53  }
54 
55  return result;
56 }
57 
58 
59 void Foam::LagrangianModels::checkApplied() const
60 {
61  const label timeIndex =
62  mesh().time().subCycling()
64  : mesh().time().timeIndex();
65 
66  if (timeIndex > checkTimeIndex_)
67  {
68  const PtrListDictionary<LagrangianModel>& modelList(*this);
69 
70  forAll(modelList, i)
71  {
72  const LagrangianModel& model = modelList[i];
73 
74  wordHashSet notAddSupFields(model.addSupFields());
75  notAddSupFields -= addSupFields_[i];
76 
77  forAllConstIter(wordHashSet, notAddSupFields, iter)
78  {
80  << "Model " << model.name()
81  << " defined for field " << iter.key()
82  << " but never used" << endl;
83  }
84  }
85 
86  checkTimeIndex_ = timeIndex;
87  }
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
94 :
96  <
100  >
101  (
102  io(mesh),
103  mesh
104  ),
105  dictionary(),
107  checkTimeIndex_(mesh.time().timeIndex() + 1),
108  addSupFields_()
109 {
110  readHeaderOk(IOstream::ASCII, typeName);
111 
112  const dictionary& dict(*this);
113 
114  // Iterate through the dictionary to determine the number of models
115  label i = 0;
117  {
118  i += iter().isDict();
119  }
120 
121  // Size the storage
123  addSupFields_.setSize(i);
124 
125  // Iterate through the dictionary to construct the models
126  i = 0;
128  {
129  if (!iter().isDict()) continue;
130 
131  const word& name = iter().keyword();
132  const dictionary& modelDict = iter().dict();
133 
135  (
136  i,
137  name,
138  LagrangianModel::New(name, mesh, modelDict).ptr()
139  );
140 
141  addSupFields_.set(i, new wordHashSet());
142 
143  i ++;
144  }
145 
146  // Do post-construction
147  {
148  PtrListDictionary<LagrangianModel>& modelList(*this);
149 
150  forAll(modelList, i)
151  {
152  modelList[i].postConstruct();
153  }
154  }
155 
156  // Enable re-reading
157  addWatch();
158 }
159 
160 
161 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
162 
164 {}
165 
166 
167 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
168 
169 bool Foam::LagrangianModels::addsSupToField(const word& fieldName) const
170 {
171  const PtrListDictionary<LagrangianModel>& modelList(*this);
172 
173  forAll(modelList, i)
174  {
175  const LagrangianModel& model = modelList[i];
176 
177  if (model.addsSupToField(fieldName))
178  {
179  return true;
180  }
181  }
182 
183  return false;
184 }
185 
186 
188 {
189  PtrListDictionary<LagrangianModel>& modelList(*this);
190 
191  forAll(modelList, i)
192  {
193  modelList[i].correct();
194  }
195 }
196 
197 
199 (
201 ) const
202 {
203  const PtrListDictionary<LagrangianModel>& modelList(*this);
204 
205  // Get the models to produce a list of which elements are to be modified
206  // and which are to be removed
208  forAll(modelList, i)
209  {
210  modelList[i].preModify(mesh, elementModifications);
211  }
212 
213  // Partition the mesh into modified and removed blocks of elements (and
214  // unmodified elements, at the start)
215  const labelList offsets = mesh.partition(2, elementModifications);
216 
217  // Return the sub-mesh associated with the modified block
218  return
220  (
221  mesh,
223  offsets[1] - offsets[0],
224  offsets[0]
225  );
226 
227  // Note: We don't remove the removed elements here. The outer code might
228  // need to do something with them (e.g., remove them from averages, or
229  // transfer their mass to the carrier) before they get deleted. They will
230  // actually get removed in LagrangianModels::modify.
231 }
232 
233 
235 (
237  const LagrangianSubMesh& modifiedMesh
238 ) const
239 {
240  const PtrListDictionary<LagrangianModel>& modelList(*this);
241 
242  // Remove the removed block from the mesh
243  mesh.remove(mesh.size() - modifiedMesh.end());
244 
245  // Get the models to modify and create elements
246  LagrangianSubMesh modifiedAndCreatedMesh(mesh.subNone());
247  modifiedAndCreatedMesh += modifiedMesh;
248  forAll(modelList, i)
249  {
250  modifiedAndCreatedMesh += modelList[i].modify(mesh, modifiedMesh);
251  }
252 
253  // Return a sub-mesh identifying all modified and created elements
254  return modifiedAndCreatedMesh;
255 }
256 
257 
259 (
260  const LagrangianSubScalarField& deltaT,
261  const bool final
262 )
263 {
264  PtrListDictionary<LagrangianModel>& modelList(*this);
265 
266  forAll(modelList, i)
267  {
268  modelList[i].calculate(deltaT, final);
269  }
270 }
271 
272 
275 {
276  checkApplied();
277 
280  (
281  word::null,
282  deltaT.mesh(),
284  );
286 
287  const PtrListDictionary<LagrangianModel>& modelList(*this);
288 
289  forAll(modelList, i)
290  {
291  const LagrangianModel& model = modelList[i];
292 
293  if (model.addsSupToField(word::null))
294  {
295  addSupFields_[i].insert(word::null);
296 
297  model.addSup(deltaT, S);
298  }
299  }
300 
301  return tS;
302 }
303 
304 
306 {
307  return true;
308 }
309 
310 
312 {
313  PtrListDictionary<LagrangianModel>& modelList(*this);
314 
315  forAll(modelList, i)
316  {
317  modelList[i].topoChange(map);
318  }
319 }
320 
321 
323 {
324  PtrListDictionary<LagrangianModel>& modelList(*this);
325 
326  forAll(modelList, i)
327  {
328  modelList[i].mapMesh(map);
329  }
330 }
331 
332 
334 {
335  PtrListDictionary<LagrangianModel>& modelList(*this);
336 
337  forAll(modelList, i)
338  {
339  modelList[i].distribute(map);
340  }
341 }
342 
343 
345 {
346  if (regIOobject::read())
347  {
348  checkTimeIndex_ = mesh().time().timeIndex() + 1;
349 
350  bool allOk = true;
351 
352  PtrListDictionary<LagrangianModel>& modelList(*this);
353 
354  forAll(modelList, i)
355  {
356  const dictionary& modelDict = subDict(modelList[i].name());
357 
358  const bool ok =
359  modelList[i].read
360  (
361  modelDict.optionalSubDict(modelList[i].type() + "Coeffs")
362  );
363  allOk = allOk && ok;
364  }
365 
366  return allOk;
367  }
368  else
369  {
370  return false;
371  }
372 }
373 
374 
376 {
377  is >> *this;
378  return !is.bad();
379 }
380 
381 
383 {
384  dictionary::write(os, false);
385  return os.good();
386 }
387 
388 
390 (
394  const bool write
395 ) const
396 {
397  bool allOk = true;
398 
399  const PtrListDictionary<LagrangianModel>& modelList(*this);
400 
401  forAll(modelList, i)
402  {
403  const bool ok = modelList[i].write(write);
404  allOk = allOk && ok;
405 
406  const bool okState = modelList[i].writeState(write);
407  allOk = allOk && okState;
408  }
409 
410  return allOk;
411 }
412 
413 
414 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Templated abstract base-class for demand-driven mesh objects used to automate their allocation to the...
void insert(const word &, T *)
Add at head of dictionary.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
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
bool bad() const
Return true if stream is corrupted.
Definition: IOstream.H:351
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:194
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
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.
const Time & time() const
Return time.
Base class for Lagrangian models.
virtual void addSup(const LagrangianSubScalarField &deltaT, LagrangianSubScalarField &S) const
Add a fractional source term.
virtual bool addsSupToField(const word &) const
Return true if the LagrangianModel adds a source term to the.
static autoPtr< LagrangianModel > New(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict)
Selector.
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.
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
virtual ~LagrangianModels()
Destructor.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
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.
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...
label end() const
Return end.
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.
autoPtr< T > set(const label, const word &key, T *)
Set element to pointer provided and return old element.
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
const TimeState & prevTimeState() const
Return previous TimeState if time is being sub-cycled.
Definition: Time.H:438
bool subCycling() const
Return true if time currently being sub-cycled, otherwise false.
Definition: Time.H:432
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:111
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
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:927
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:803
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
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.
bool readHeaderOk(const IOstream::streamFormat defaultFormat, const word &typeName)
Read header, check readOpt flags and read data if necessary.
void addWatch()
Add file watch on object (if registered and READ_IF_MODIFIED)
Definition: regIOobject.C:259
virtual bool read()
Read object.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimRate
defineTypeNameAndDebug(combustionModel, 0)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:210
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict