GeometricField.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) 2011-2022 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::GeometricField
26 
27 Description
28  Generic GeometricField class.
29 
30 SourceFiles
31  GeometricFieldI.H
32  GeometricField.C
33  GeometricFieldFunctions.H
34  GeometricFieldFunctions.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef GeometricField_H
39 #define GeometricField_H
40 
43 #include "GeometricBoundaryField.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 class dictionary;
51 
52 // Forward declaration of friend functions and operators
53 
54 template<class Type, template<class> class PatchField, class GeoMesh>
55 class GeometricField;
56 
57 template<class Type, template<class> class PatchField, class GeoMesh>
58 Ostream& operator<<
59 (
60  Ostream&,
62 );
63 
64 template<class Type, template<class> class PatchField, class GeoMesh>
65 Ostream& operator<<
66 (
67  Ostream&,
69 );
70 
71 
72 /*---------------------------------------------------------------------------*\
73  Class GeometricField Declaration
74 \*---------------------------------------------------------------------------*/
75 
76 template<class Type, template<class> class PatchField, class GeoMesh>
77 class GeometricField
78 :
79  public DimensionedField<Type, GeoMesh>
80 {
81  // Private Member Functions
82 
83  //- Read from file if it is present
84  bool readIfPresent();
85 
86  //- Read old time field from file if it is present
87  bool readOldTimeIfPresent();
88 
89 
90 public:
91 
92  // Public Typedefs
93 
94  //- Type of mesh on which this GeometricField is instantiated
95  typedef typename GeoMesh::Mesh Mesh;
96 
97  //- Type of the internal field from which this GeometricField is derived
99 
100  //- Type of the patch field of which the Boundary is composed
101  typedef PatchField<Type> Patch;
102 
103  //- Type of the boundary field
105 
106 
107 private:
108 
109  // Private Data
110 
111  //- Current time index.
112  // Used to trigger the storing of the old-time value
113  mutable label timeIndex_;
114 
115  //- Pointer to old time field
116  mutable GeometricField<Type, PatchField, GeoMesh>* field0Ptr_;
117 
118  //- Pointer to previous iteration (used for under-relaxation)
119  mutable GeometricField<Type, PatchField, GeoMesh>* fieldPrevIterPtr_;
120 
121  //- Boundary Type field containing boundary field values
122  Boundary boundaryField_;
123 
124 
125  // Private Member Functions
126 
127  //- Read the field from the dictionary
128  void readFields(const dictionary&);
129 
130  //- Read the field - create the field dictionary on-the-fly
131  void readFields();
132 
133  //- Store the old-time field
134  void storeOldTime() const;
135 
136  //- Find oldest time and set to nullObject
137  void nullOldestTimeFound();
138 
139 
140 public:
141 
142  //- Runtime type information
143  TypeName("GeometricField");
144 
145 
146  // Public Typedefs
147 
148  typedef typename Field<Type>::cmptType cmptType;
149 
150  // Static Member Functions
151 
152  //- Return a null geometric field
153  inline static const GeometricField<Type, PatchField, GeoMesh>& null();
154 
155 
156  // Constructors
157 
158  //- Constructor given IOobject, mesh, dimensions and patch field type.
159  // This allocates storage for the field but not values.
160  // Used only within this class to create TEMPORARY variables
162  (
163  const IOobject&,
164  const Mesh&,
165  const dimensionSet&,
166  const word& patchFieldType=PatchField<Type>::calculatedType()
167  );
168 
169  //- Constructor given IOobject, mesh, dimensions and patch field types.
170  // This allocates storage for the field but not values.
171  // Used only within this class to create TEMPORARY variables
173  (
174  const IOobject&,
175  const Mesh&,
176  const dimensionSet&,
177  const wordList& wantedPatchTypes,
178  const wordList& actualPatchTypes = wordList()
179  );
180 
181  //- Constructor given IOobject, mesh, dimensioned<Type>
182  // and patch field type.
184  (
185  const IOobject&,
186  const Mesh&,
187  const dimensioned<Type>&,
188  const word& patchFieldType=PatchField<Type>::calculatedType()
189  );
190 
191  //- Constructor given IOobject, mesh, dimensioned<Type>
192  // and patch field types.
194  (
195  const IOobject&,
196  const Mesh&,
197  const dimensioned<Type>&,
198  const wordList& wantedPatchTypes,
199  const wordList& actualPatchTypes = wordList()
200  );
201 
202  //- Constructor from components
204  (
205  const IOobject&,
206  const Internal&,
207  const PtrList<PatchField<Type>>&
208  );
209 
210  //- Constructor from components
212  (
213  const IOobject&,
214  const Mesh&,
215  const dimensionSet&,
216  const Field<Type>&,
217  const PtrList<PatchField<Type>>&
218  );
219 
220  //- Construct and read given IOobject
222  (
223  const IOobject&,
224  const Mesh&
225  );
226 
227  //- Construct from dictionary
229  (
230  const IOobject&,
231  const Mesh&,
232  const dictionary&
233  );
234 
235  //- Copy constructor
237 
238  //- Move constructor
240 
241  //- Construct as copy of tmp<GeometricField> deleting argument
243 
244  //- Construct as copy resetting IO parameters
246  (
247  const IOobject&,
249  );
250 
251  //- Construct as copy of tmp<GeometricField> resetting IO parameters
253  (
254  const IOobject&,
256  );
257 
258  //- Construct as copy resetting name
260  (
261  const word& newName,
263  );
264 
265  //- Construct as copy resetting name
267  (
268  const word& newName,
270  );
271 
272  //- Construct as copy resetting IO parameters and patch field type
274  (
275  const IOobject&,
277  const word& patchFieldType
278  );
279 
280  //- Construct as copy resetting IO parameters and boundary types
282  (
283  const IOobject&,
285  const wordList& patchFieldTypes,
286  const wordList& actualPatchTypes = wordList()
287  );
288 
289  //- Construct as copy resetting IO parameters and boundary types
291  (
292  const IOobject&,
294  const wordList& patchFieldTypes,
295  const wordList& actualPatchTypes = wordList()
296  );
297 
298  //- Clone
300 
301  //- Clone un-sliced
303 
304  //- Return a temporary field constructed from name,
305  // internal field and list of patch fields
307  (
308  const word& name,
309  const Internal&,
310  const PtrList<PatchField<Type>>&
311  );
312 
313  //- Return a temporary field constructed from name, mesh, dimensionSet
314  // and patch field type.
316  (
317  const word& name,
318  const Mesh&,
319  const dimensionSet&,
320  const word& patchFieldType=PatchField<Type>::calculatedType()
321  );
322 
323  //- Return a temporary field constructed from mesh, dimensioned<Type>
324  // and patch field type.
326  (
327  const word& name,
328  const Mesh&,
329  const dimensioned<Type>&,
330  const word& patchFieldType=PatchField<Type>::calculatedType()
331  );
332 
333  //- Return a temporary field constructed from mesh, dimensioned<Type>
334  // and patch field types.
336  (
337  const word& name,
338  const Mesh&,
339  const dimensioned<Type>&,
340  const wordList& patchFieldTypes,
341  const wordList& actualPatchTypes = wordList()
342  );
343 
344  //- Rename temporary field and return
346  (
347  const word& newName,
349  );
350 
351  //- Rename temporary field, reset patch field type and return
353  (
354  const word& newName,
356  const word&
357  );
358 
359  //- Rename and reset patch fields types of temporary field and return
361  (
362  const word& newName,
364  const wordList& patchFieldTypes,
365  const wordList& actualPatchTypes = wordList()
366  );
367 
368 
369  //- Destructor
370  virtual ~GeometricField();
371 
372 
373  // Member Functions
374 
375  //- Return a reference to the dimensioned internal field
376  // Note: this increments the event counter and checks the
377  // old-time fields; avoid in loops.
378  Internal& ref();
379 
380  //- Return a const-reference to the dimensioned internal field
381  inline const Internal& internalField() const;
382 
383  //- Return a const-reference to the dimensioned internal field
384  // of a "vol" field. Useful in the formulation of source-terms
385  // for FV equations
386  inline const Internal& v() const;
387 
388  //- Return a reference to the internal field
389  // Note: this increments the event counter and checks the
390  // old-time fields; avoid in loops.
392 
393  //- Return a const-reference to the internal field
394  inline const typename Internal::FieldType& primitiveField() const;
395 
396  //- Return a reference to the boundary field
397  // Note: this increments the event counter and checks the
398  // old-time fields; avoid in loops.
400 
401  //- Return const-reference to the boundary field
402  inline const Boundary& boundaryField() const;
403 
404  //- Return the time index of the field
405  inline label timeIndex() const;
406 
407  //- Return the time index of the field
408  inline label& timeIndex();
409 
410  //- Return whether or not this is an old-time field
411  bool isOldTime() const;
412 
413  //- Store the old-time fields
414  void storeOldTimes() const;
415 
416  //- Return the number of old time fields stored
417  label nOldTimes() const;
418 
419  //- Return old time field
421 
422  //- Return non-const old time field
423  // (Not a good idea but it is used for sub-cycling)
425 
426  //- Return the n-th old time field
428  (
429  const label n
430  ) const;
431 
432  //- Return the n-th non-const old time field
433  // (Not a good idea but it is used for sub-cycling)
435 
436  //- Store the field as the previous iteration value
437  void storePrevIter() const;
438 
439  //- Return previous iteration field
441 
442  //- Delete old time and previous iteration fields
443  void clearOldTimes();
444 
445  //- Set oldest time field pointer to nullObjectPtr
446  // to remove the oldest time field while maintaining a tag
447  // so that it is reinstated on the next call to oldTime().
448  //
449  // This is used to temporarily remove the oldest-time field
450  // prior to mesh topology change.
451  void nullOldestTime();
452 
453  //- Correct boundary field
455 
456  //- Reset the field contents to the given field
457  // Used for mesh to mesh mapping
459 
460  //- Does the field need a reference level for solution
461  bool needReference() const;
462 
463  //- Return a component of the field
465  (
466  const direction
467  ) const;
468 
469  //- WriteData member function required by regIOobject
470  bool writeData(Ostream&) const;
471 
472  //- Return transpose (only if it is a tensor field)
474 
475  //- Relax field (for steady-state solution).
476  // alpha >= 1 : no relaxation
477  // alpha < 1 : relaxation
478  void relax(const scalar alpha);
479 
480  //- Return the field relaxation factor read from fvSolution
481  // or 1 if not specified
482  scalar relaxationFactor() const;
483 
484  //- Relax current field with respect to the cached previous iteration.
485  // Relaxation factor is read from fvSolution
486  void relax();
487 
488  //- Relax given field with respect to the current field
489  // and reset the field to the result
490  void relax
491  (
493  const scalar alpha
494  );
495 
496  //- Relax given field with respect to the current field
497  // and reset the field to the result
498  // Relaxation factor is read from fvSolution
500 
501  //- Select the final iteration parameters if `final' is true
502  // by returning the field name + "Final"
503  // otherwise the standard parameters by returning the field name
504  word select(bool final) const;
505 
506  //- Helper function to write the min and max to an Ostream
507  void writeMinMax(Ostream& os) const;
508 
509 
510  // Member function *this operators
511 
512  void negate();
513 
514  void replace
515  (
516  const direction,
518  );
519 
520  void replace
521  (
522  const direction,
523  const dimensioned<cmptType>&
524  );
525 
526  void max(const dimensioned<Type>&);
527 
528  void min(const dimensioned<Type>&);
529 
530  void maxMin
531  (
532  const dimensioned<Type>& minDt,
533  const dimensioned<Type>& maxDt
534  );
535 
536 
537  // Member Operators
538 
539  //- Return a const-reference to the dimensioned internal field
540  // Useful in the formulation of source-terms for FV equations
541  inline const Internal& operator()() const;
542 
546  void operator=(const dimensioned<Type>&);
547  void operator=(const zero&);
548 
550  void operator==(const dimensioned<Type>&);
551  void operator==(const zero&);
552 
555 
558 
561 
564 
565  void operator+=(const dimensioned<Type>&);
566  void operator-=(const dimensioned<Type>&);
567 
568  void operator*=(const dimensioned<scalar>&);
569  void operator/=(const dimensioned<scalar>&);
570 
571 
572  // Ostream operators
573 
574  friend Ostream& operator<< <Type, PatchField, GeoMesh>
575  (
576  Ostream&,
578  );
579 
580  friend Ostream& operator<< <Type, PatchField, GeoMesh>
581  (
582  Ostream&,
584  );
585 };
586 
587 
588 template<class Type, template<class> class PatchField, class GeoMesh>
589 Ostream& operator<<
590 (
591  Ostream&,
593 );
594 
595 
596 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
597 
598 } // End namespace Foam
599 
600 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
601 
602 #include "GeometricFieldI.H"
603 
604 #ifdef NoRepository
605  #include "GeometricField.C"
606 #endif
607 
608 #include "GeometricFieldFunctions.H"
609 
610 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
611 
612 #endif
613 
614 // ************************************************************************* //
List of coupled interface fields to be used in coupling.
label n
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
GeoMesh::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
Pre-declare SubField and related Field type.
Definition: Field.H:82
pTraits< Type >::cmptType cmptType
Component type.
Definition: Field.H:97
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
MESH Mesh
Definition: GeoMesh.H:61
Generic GeometricBoundaryField class.
Generic GeometricField class.
void max(const dimensioned< Type > &)
Internal & ref()
Return a reference to the dimensioned internal field.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
void relax()
Relax current field with respect to the cached previous iteration.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
TypeName("GeometricField")
Runtime type information.
void operator/=(const GeometricField< scalar, PatchField, GeoMesh > &)
label timeIndex() const
Return the time index of the field.
void maxMin(const dimensioned< Type > &minDt, const dimensioned< Type > &maxDt)
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
GeoMesh::Mesh Mesh
Type of mesh on which this GeometricField is instantiated.
void operator==(const tmp< GeometricField< Type, PatchField, GeoMesh >> &)
void operator-=(const GeometricField< Type, PatchField, GeoMesh > &)
tmp< GeometricField< Type, PatchField, GeoMesh > > cloneUnSliced() const
Clone un-sliced.
void operator*=(const GeometricField< scalar, PatchField, GeoMesh > &)
void min(const dimensioned< Type > &)
void operator+=(const GeometricField< Type, PatchField, GeoMesh > &)
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
const Internal & v() const
Return a const-reference to the dimensioned internal field.
Definition: volFieldsI.H:31
PatchField< Type > Patch
Type of the patch field of which the Boundary is composed.
Field< Type >::cmptType cmptType
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
GeometricField(const IOobject &, const Mesh &, const dimensionSet &, const word &patchFieldType=PatchField< Type >::calculatedType())
Constructor given IOobject, mesh, dimensions and patch field type.
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
bool needReference() const
Does the field need a reference level for solution.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
virtual ~GeometricField()
Destructor.
scalar relaxationFactor() const
Return the field relaxation factor read from fvSolution.
bool isOldTime() const
Return whether or not this is an old-time field.
void storePrevIter() const
Store the field as the previous iteration value.
label nOldTimes() const
Return the number of old time fields stored.
word select(bool final) const
Select the final iteration parameters if `final' is true.
void storeOldTimes() const
Store the old-time fields.
void correctBoundaryConditions()
Correct boundary field.
void reset(const tmp< GeometricField< Type, PatchField, GeoMesh >> &)
Reset the field contents to the given field.
const Internal & operator()() const
Return a const-reference to the dimensioned internal field.
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
void nullOldestTime()
Set oldest time field pointer to nullObjectPtr.
void operator=(const GeometricField< Type, PatchField, GeoMesh > &)
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
void clearOldTimes()
Delete old time and previous iteration fields.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
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 keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Dimension set for the base types.
Definition: dimensionSet.H:122
Generic dimensioned Type class.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:50
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
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
word patchFieldType(const PatchField &pf)
uint8_t direction
Definition: direction.H:45