LagrangianFieldSource.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 "LagrangianFieldSource.H"
27 #include "LagrangianMesh.H"
28 #include "LagrangianSubMesh.H"
29 #include "LagrangianSource.H"
30 #include "LagrangianInjection.H"
31 #include "dictionary.H"
32 #include "regIOobject.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class Type>
38 (
39  const regIOobject& iIo
40 )
41 :
42  libs_(fileNameList::null()),
43  internalIo_(iIo),
44  internalField_
45  (
47  ),
48  internalNonDynamicField_
49  (
50  refCastNull<const LagrangianInternalField<Type>>(iIo)
51  )
52 {}
53 
54 
55 template<class Type>
57 (
58  const regIOobject& iIo,
59  const dictionary& dict
60 )
61 :
62  libs_(dict.lookupOrDefault("libs", fileNameList::null())),
63  internalIo_(iIo),
64  internalField_
65  (
67  ),
68  internalNonDynamicField_
69  (
70  refCastNull<const LagrangianInternalField<Type>>(iIo)
71  )
72 {}
73 
74 
75 template<class Type>
77 (
78  const LagrangianFieldSource<Type>& field,
79  const regIOobject& iIo
80 )
81 :
82  libs_(field.libs_),
83  internalIo_(iIo),
84  internalField_
85  (
87  ),
88  internalNonDynamicField_
89  (
90  refCastNull<const LagrangianInternalField<Type>>(iIo)
91  )
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
96 
97 template<class Type>
100 (
101  const word& fieldSourceType,
102  const regIOobject& iIo
103 )
104 {
105  typename nullConstructorTable::iterator cstrIter =
106  nullConstructorTablePtr_->find(fieldSourceType);
107 
108  if (cstrIter == nullConstructorTablePtr_->end())
109  {
111  << "Unknown null-constructable fieldSource type " << fieldSourceType
112  << nl << nl << "Valid fieldSource types are :" << endl
113  << nullConstructorTablePtr_->sortedToc()
114  << exit(FatalError);
115  }
116 
117  return cstrIter()(iIo);
118 }
119 
120 
121 template<class Type>
124 (
125  const regIOobject& iIo,
126  const dictionary& dict
127 )
128 {
129  const word fieldSourceType(dict.lookup("type"));
130 
131  libs.open(dict, "libs", dictionaryConstructorTablePtr_);
132 
133  typename dictionaryConstructorTable::iterator cstrIter =
134  dictionaryConstructorTablePtr_->find(fieldSourceType);
135 
136  if (cstrIter == dictionaryConstructorTablePtr_->end())
137  {
138  if (!disallowGenericLagrangianFieldSource)
139  {
140  cstrIter = dictionaryConstructorTablePtr_->find("generic");
141  }
142 
143  if (cstrIter == dictionaryConstructorTablePtr_->end())
144  {
146  << "Unknown fieldSource type " << fieldSourceType
147  << " for model " << dict.dictName() << nl << nl
148  << "Valid fieldSource types are :" << endl
149  << dictionaryConstructorTablePtr_->sortedToc()
150  << exit(FatalIOError);
151  }
152  }
153 
154  return cstrIter()(iIo, dict);
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
159 
160 template<class Type>
162 {}
163 
164 
165 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
166 
167 template<class Type>
169 {
170  return internalIo_.db();
171 }
172 
173 
174 template<class Type>
175 const Foam::dimensionSet&
177 {
178  if (notNull(internalField_))
179  {
180  return internalField_.dimensions();
181  }
182 
183  if (notNull(internalNonDynamicField_))
184  {
185  return internalNonDynamicField_.dimensions();
186  }
187 
189  << "Dimensions of internal object " << internalIo_.name()
190  << " could not be determined"
191  << exit(FatalError);
192 
193  return NullObjectRef<dimensionSet>();
194 }
195 
196 
197 template<class Type>
200 {
201  if (notNull(internalField_))
202  {
203  return internalField_;
204  }
205 
207  << "Internal field " << internalIo_.name() << " is not of type "
209  << exit(FatalError);
210 
211  return NullObjectRef<LagrangianInternalDynamicField<Type>>();
212 }
213 
214 
215 template<class Type>
218 (
219  const LagrangianSource& source,
220  const LagrangianSubMesh& subMesh
221 ) const
222 {
224  << "The " << typeName << " " << source.name()
225  << " does not define a value for a continuous source "
226  << exit(FatalError);
227 
228  return tmp<LagrangianSubField<Type>>(nullptr);
229 }
230 
231 
232 template<class Type>
235 (
236  const LagrangianSource& source,
237  const LagrangianSubMesh& subMesh
238 ) const
239 {
241  << "The " << typeName << " " << source.name()
242  << " does not define a value for a continuous source"
243  << exit(FatalError);
244 
245  return tmp<LagrangianSubScalarField>(nullptr);
246 }
247 
248 
249 template<class Type>
252 (
253  const LagrangianSource& source,
254  const LagrangianSubMesh& subMesh
255 ) const
256 {
257  return (1 - internalCoeff(source, subMesh))*sourceValue(source, subMesh);
258 }
259 
260 
261 template<class Type>
264 (
265  const LagrangianSource& source,
266  const LagrangianSubMesh& subMesh
267 ) const
268 {
269  return
270  sourceCoeff(source, subMesh)
271  + internalCoeff(source, subMesh)*subMesh.sub(internalField());
272 }
273 
274 
275 template<class Type>
278 (
279  const LagrangianInjection& injection,
280  const LagrangianSubMesh& subMesh
281 ) const
282 {
284  << "The " << typeName << " " << injection.name()
285  << " does not define a value for an instantaneous injection"
286  << exit(FatalError);
287 
288  return tmp<LagrangianSubField<Type>>(nullptr);
289 }
290 
291 
292 template<class Type>
295 (
296  const LagrangianModel& model,
297  const LagrangianSubMesh& subMesh
298 ) const
299 {
300  if (isA<LagrangianSource>(model))
301  {
302  return value(refCast<const LagrangianSource>(model), subMesh);
303  }
304 
305  if (isA<LagrangianInjection>(model))
306  {
307  return value(refCast<const LagrangianInjection>(model), subMesh);
308  }
309 
311  << "The " << LagrangianModel::typeName << " " << model.name()
312  << " is neither a " << LagrangianSource::typeName << " nor a "
313  << LagrangianInjection::typeName << " so source values cannot "
314  << "be generated for it"
315  << exit(FatalError);
316 
317  return tmp<LagrangianSubField<Type>>(nullptr);
318 }
319 
320 
321 template<class Type>
322 template<class OtherType>
325 (
326  const word& name,
327  const LagrangianModel& model
328 ) const
329 {
331  db().template lookupObject<LagrangianDynamicField<OtherType>>(name);
332 
333  return lf.sources()[model.name()];
334 }
335 
336 
337 template<class Type>
338 template<class OtherType, class OtherFieldSourceType>
340 (
341  const word& name,
342  const LagrangianModel& model
343 ) const
344 {
346  fieldSource<OtherType>(name, model);
347 
348  if (!isA<OtherFieldSourceType>(lfs))
349  {
351  << "The '" << type() << "' source of field '"
352  << (db().dbDir()/internalIo_.name()).c_str()
353  << "' for the '" << model.type() << "' Lagrangian model '"
354  << model.name() << "' requires the corresponding source of field '"
355  << (db().dbDir()/name).c_str()
356  << "' to be of type '" << OtherFieldSourceType::typeName
357  << "' (or a derivation thereof), rather than '" << lfs.type()
358  << "'" << exit(FatalError);
359  }
360 
361  return refCast<const OtherFieldSourceType>(lfs);
362 }
363 
364 
365 template<class Type>
366 template<class OtherModelType>
368 (
369  const LagrangianModel& model
370 ) const
371 {
372  if (!isA<OtherModelType>(model))
373  {
375  << "The '" << type() << "' source of field '"
376  << (db().dbDir()/internalIo_.name()).c_str()
377  << "' for the Lagrangian model '" << model.name()
378  << "' requires a model of type '" << OtherModelType::typeName
379  << "' (or a derivation thereof), rather than '" << model.type()
380  << "'" << exit(FatalError);
381  }
382 
383  return refCast<const OtherModelType>(model);
384 }
385 
386 
387 template<class Type>
388 template<class OtherType>
391 (
392  const word& name,
393  const LagrangianSource& source,
394  const LagrangianSubMesh& subMesh
395 ) const
396 {
397  return fieldSource<OtherType>(name, source).value(source, subMesh);
398 }
399 
400 
401 template<class Type>
402 template<class OtherType>
405 (
406  const word& name,
407  const LagrangianInjection& injection,
408  const LagrangianSubMesh& subMesh
409 ) const
410 {
411  return fieldSource<OtherType>(name, injection).value(injection, subMesh);
412 }
413 
414 
415 template<class Type>
417 {
418  writeEntry(os, "type", type());
419 }
420 
421 
422 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
423 
424 template<class Type>
425 Foam::Ostream& Foam::operator<<
426 (
427  Ostream& os,
428  const LagrangianFieldSource<Type>& ptf
429 )
430 {
431  ptf.write(os);
432 
433  os.check
434  (
435  "Ostream& operator<<(Ostream&, const LagrangianFieldSource<Type>&)"
436  );
437 
438  return os;
439 }
440 
441 
442 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
const Sources & sources() const
Return const-reference to the sources.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:309
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
Base class for Lagrangian source conditions.
const OtherFieldSourceType & fieldSourceCast(const word &name, const LagrangianModel &) const
Lookup and return another field source and then cast to a.
const dimensionSet & internalDimensions() const
Return internal dimensions reference.
virtual void write(Ostream &) const
Write.
static autoPtr< LagrangianFieldSource< Type > > New(const word &fieldSourceType, const regIOobject &)
Return a pointer to a new field source.
const objectRegistry & db() const
Return local objectRegistry.
tmp< LagrangianSubField< Type > > value(const LagrangianSource &, const LagrangianSubMesh &) const
Return the value for a continuous source.
const OtherModelType & modelCast(const LagrangianModel &) const
Cast a model to the given type. Handle errors.
tmp< LagrangianSubField< Type > > sourceCoeff(const LagrangianSource &, const LagrangianSubMesh &) const
Return the source coefficient.
virtual tmp< LagrangianSubField< Type > > sourceValue(const LagrangianSource &, const LagrangianSubMesh &) const
Return the source value.
virtual ~LagrangianFieldSource()
Destructor.
virtual tmp< LagrangianSubScalarField > internalCoeff(const LagrangianSource &, const LagrangianSubMesh &) const
Return the internal coefficient.
LagrangianFieldSource(const regIOobject &)
Construct from internal field.
const LagrangianFieldSource< OtherType > & fieldSource(const word &name, const LagrangianModel &) const
Lookup and return another field source.
const LagrangianInternalDynamicField< Type > & internalField() const
Return internal field reference.
Base class for Lagrangian injections. Minimal wrapper over LagrangianSource. Implements some utility ...
Base class for Lagrangian models.
const word & name() const
The source name.
Base class for Lagrangian sources. Minimal wrapper over LagrangianModel that provides an interface to...
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.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual Ostream & write(const char)=0
Write character.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
Dimension set for the base types.
Definition: dimensionSet.H:125
bool open(const fileName &libName, const bool verbose=true)
Open the named library, optionally with warnings if problems occur.
Registry of regIOobjects.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dlLibraryTable libs
Table of loaded dynamic libraries.
To & refCastNull(From &r)
Reference type cast template function,.
Definition: typeInfo.H:155
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:64
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
static const char nl
Definition: Ostream.H:267
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict