fvModels.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) 2021-2023 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 "fvModels.H"
27 #include "fvConstraint.H"
28 #include "fvMesh.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
35 }
36 
37 
38 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
39 
40 Foam::IOobject Foam::fvModels::createIOobject
41 (
42  const fvMesh& mesh
43 ) const
44 {
45  typeIOobject<IOdictionary> io
46  (
47  typeName,
48  mesh.time().constant(),
49  mesh,
52  );
53 
54  if (io.headerOk())
55  {
56  Info<< "Creating fvModels from "
57  << io.instance()/io.name() << nl
58  << endl;
59 
60  io.readOpt() = IOobject::MUST_READ_IF_MODIFIED;
61  return io;
62  }
63  else
64  {
65  // For backward-compatibility
66  // check if the fvOptions file is in constant
67  io.rename("fvOptions");
68 
69  if (io.headerOk())
70  {
71  Warning
72  << "Creating fvModels from "
73  << io.instance()/io.name() << nl
74  << endl;
75 
76  io.readOpt() = IOobject::MUST_READ_IF_MODIFIED;
77  return io;
78  }
79  else
80  {
81  // For backward-compatibility
82  // check if the fvOptions file is in system
83  io.instance() = mesh.time().system();
84 
85  if (io.headerOk())
86  {
87  Warning
88  << "Creating fvModels from "
89  << io.instance()/io.name()
90  << " rather than constant/fvModels"
91  << endl;
92 
93  io.readOpt() = IOobject::MUST_READ_IF_MODIFIED;
94  return io;
95  }
96  else
97  {
98  io.rename(typeName);
99  io.readOpt() = IOobject::NO_READ;
100  return io;
101  }
102  }
103  }
104 }
105 
106 
107 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
108 
109 void Foam::fvModels::checkApplied() const
110 {
111  const label timeIndex =
112  mesh().time().subCycling()
113  ? mesh().time().prevTimeState().timeIndex()
114  : mesh().time().timeIndex();
115 
116  if (timeIndex > checkTimeIndex_)
117  {
118  const PtrListDictionary<fvModel>& modelList(*this);
119 
120  forAll(modelList, i)
121  {
122  const fvModel& model = modelList[i];
123 
124  wordHashSet notAddSupFields(model.addSupFields());
125  notAddSupFields -= addSupFields_[i];
126 
127  forAllConstIter(wordHashSet, notAddSupFields, iter)
128  {
130  << "Model " << model.name()
131  << " defined for field " << iter.key()
132  << " but never used" << endl;
133  }
134  }
135 
136  checkTimeIndex_ = timeIndex;
137  }
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142 
144 (
145  const fvMesh& mesh
146 )
147 :
149  (
150  mesh,
151  createIOobject(mesh)
152  ),
154  checkTimeIndex_(mesh.time().timeIndex() + 1),
155  addSupFields_()
156 {
157  readHeaderOk(IOstream::ASCII, typeName);
158 
159  const bool readFromFvModels(IOobject::name() == typeName);
160 
161  const dictionary& dict(*this);
162 
163  // Count number of active fvModels
164  label count = 0;
166  {
167  if (iter().isDict())
168  {
169  count++;
170  }
171  }
172 
174 
175  addSupFields_.setSize(count);
176 
177  label i = 0;
179  {
180  if (iter().isDict())
181  {
182  const word& name = iter().keyword();
183  const dictionary& modelDict = iter().dict();
184 
185  const word modelType(modelDict.lookup("type"));
186 
187  if
188  (
189  readFromFvModels
190  || !fvConstraint::dictionaryConstructorTablePtr_->found
191  (
192  modelType
193  )
194  )
195  {
197  (
198  i,
199  name,
200  fvModel::New(name, mesh, modelDict).ptr()
201  );
202 
203  addSupFields_.set(i, new wordHashSet());
204 
205  i++;
206  }
207  }
208  }
209 
211  addSupFields_.setSize(i);
212 
213  if (readFromFvModels)
214  {
215  // Add file watch on the fvModels dictionary for
216  // MUST_READ_IF_MODIFIED
217  addWatch();
218  }
219  else
220  {
221  // If the fvOptions file was read re-name to fvModels
222  rename(typeName);
223  }
224 }
225 
226 
227 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
228 
229 bool Foam::fvModels::addsSupToField(const word& fieldName) const
230 {
231  const PtrListDictionary<fvModel>& modelList(*this);
232 
233  forAll(modelList, i)
234  {
235  if (modelList[i].addsSupToField(fieldName))
236  {
237  return true;
238  }
239  }
240 
241  return false;
242 }
243 
244 
245 Foam::scalar Foam::fvModels::maxDeltaT() const
246 {
247  const PtrListDictionary<fvModel>& modelList(*this);
248 
249  scalar maxDeltaT = vGreat;
250 
251  forAll(modelList, i)
252  {
253  maxDeltaT = min(maxDeltaT, modelList[i].maxDeltaT());
254  }
255 
256  return maxDeltaT;
257 }
258 
259 
261 {
262  PtrListDictionary<fvModel>& modelList(*this);
263 
264  forAll(modelList, i)
265  {
266  modelList[i].preUpdateMesh();
267  }
268 }
269 
270 
272 {
273  bool allOk = true;
274 
275  PtrListDictionary<fvModel>& modelList(*this);
276 
277  forAll(modelList, i)
278  {
279  allOk = allOk && modelList[i].movePoints();
280  }
281 
282  return allOk;
283 }
284 
285 
287 {
288  PtrListDictionary<fvModel>& modelList(*this);
289 
290  forAll(modelList, i)
291  {
292  modelList[i].topoChange(map);
293  }
294 }
295 
296 
298 {
299  PtrListDictionary<fvModel>& modelList(*this);
300 
301  forAll(modelList, i)
302  {
303  modelList[i].mapMesh(map);
304  }
305 }
306 
307 
309 {
310  PtrListDictionary<fvModel>& modelList(*this);
311 
312  forAll(modelList, i)
313  {
314  modelList[i].distribute(map);
315  }
316 }
317 
318 
320 {
321  is >> *this;
322  return !is.bad();
323 }
324 
325 
327 {
328  dictionary::write(os, false);
329  return os.good();
330 }
331 
332 
334 {
336  {
337  checkTimeIndex_ = mesh().time().timeIndex() + 1;
338 
339  bool allOk = true;
340 
341  PtrListDictionary<fvModel>& modelList(*this);
342 
343  forAll(modelList, i)
344  {
345  fvModel& model = modelList[i];
346  bool ok = model.read(subDict(model.name()));
347  allOk = (allOk && ok);
348  }
349  return allOk;
350  }
351  else
352  {
353  return false;
354  }
355 }
356 
357 
359 {
360  PtrListDictionary<fvModel>& modelList(*this);
361 
362  forAll(modelList, i)
363  {
364  modelList[i].correct();
365  }
366 }
367 
368 
369 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
370 
371 Foam::Ostream& Foam::operator<<
372 (
373  Ostream& os,
374  const fvModels& models
375 )
376 {
377  models.writeData(os);
378  return os;
379 }
380 
381 
382 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Templated abstract base-class for demand-driven mesh objects used to automate their allocation to the...
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
const word & name() const
Return name.
Definition: IOobject.H:310
bool bad() const
Return true if stream is corrupted.
Definition: IOstream.H:348
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:330
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
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
static const word & system()
Return system name.
Definition: TimePaths.H:112
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:109
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
void write(Ostream &, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:204
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:952
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
Finite volume model abstract base class.
Definition: fvModel.H:59
static autoPtr< fvModel > New(const word &name, const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected fvModel.
Definition: fvModel.C:94
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:187
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:28
Finite volume models.
Definition: fvModels.H:65
virtual bool movePoints()
Update for mesh motion.
Definition: fvModels.C:271
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:358
virtual bool writeData(Ostream &os) const
Write data to Ostream.
Definition: fvModels.C:326
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: fvModels.C:286
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: fvModels.C:308
virtual scalar maxDeltaT() const
Return the maximum time-step for stable operation.
Definition: fvModels.C:245
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:260
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: fvModels.C:297
virtual bool readData(Istream &)
ReadData function which reads the fvModels dictionary.
Definition: fvModels.C:319
virtual bool addsSupToField(const word &fieldName) const
Return true if an fvModel adds a source term to the given.
Definition: fvModels.C:229
fvModels(const fvMesh &mesh)
Construct from components with list of field names.
Definition: fvModels.C:144
virtual bool read()
Read the fvModels dictionary if it has changed.
Definition: fvModels.C:333
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 void rename(const word &newName)
Rename.
Definition: regIOobject.C:417
A class for handling words, derived from string.
Definition: word.H:62
scalar maxDeltaT
#define WarningInFunction
Report a warning using Foam::Warning.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
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:251
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
defineTypeNameAndDebug(combustionModel, 0)
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:208
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
static const char nl
Definition: Ostream.H:260
messageStream Warning
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict