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&,
61  const GeometricField<Type, PatchField, GeoMesh>&
62 );
63 
64 template<class Type, template<class> class PatchField, class GeoMesh>
65 Ostream& operator<<
66 (
67  Ostream&,
68  const tmp<GeometricField<Type, PatchField, GeoMesh>>&
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 
134 public:
135 
136  //- Runtime type information
137  TypeName("GeometricField");
138 
139 
140  // Public Typedefs
142  typedef typename Field<Type>::cmptType cmptType;
143 
144  // Static Member Functions
145 
146  //- Return a null geometric field
147  inline static const GeometricField<Type, PatchField, GeoMesh>& null();
148 
149 
150  // Constructors
151 
152  //- Constructor given IOobject, mesh, dimensions and patch field type.
153  // This allocates storage for the field but not values.
154  // Used only within this class to create TEMPORARY variables
156  (
157  const IOobject&,
158  const Mesh&,
159  const dimensionSet&,
160  const word& patchFieldType=PatchField<Type>::calculatedType()
161  );
162 
163  //- Constructor given IOobject, mesh, dimensions and patch field types.
164  // This allocates storage for the field but not values.
165  // Used only within this class to create TEMPORARY variables
167  (
168  const IOobject&,
169  const Mesh&,
170  const dimensionSet&,
171  const wordList& wantedPatchTypes,
172  const wordList& actualPatchTypes = wordList()
173  );
174 
175  //- Constructor given IOobject, mesh, dimensioned<Type>
176  // and patch field type.
178  (
179  const IOobject&,
180  const Mesh&,
181  const dimensioned<Type>&,
182  const word& patchFieldType=PatchField<Type>::calculatedType()
183  );
184 
185  //- Constructor given IOobject, mesh, dimensioned<Type>
186  // and patch field types.
188  (
189  const IOobject&,
190  const Mesh&,
191  const dimensioned<Type>&,
192  const wordList& wantedPatchTypes,
193  const wordList& actualPatchTypes = wordList()
194  );
195 
196  //- Constructor from components
198  (
199  const IOobject&,
200  const Internal&,
201  const PtrList<PatchField<Type>>&
202  );
203 
204  //- Constructor from components
206  (
207  const IOobject&,
208  const Mesh&,
209  const dimensionSet&,
210  const Field<Type>&,
211  const PtrList<PatchField<Type>>&
212  );
213 
214  //- Construct and read given IOobject
216  (
217  const IOobject&,
218  const Mesh&
219  );
220 
221  //- Construct from dictionary
223  (
224  const IOobject&,
225  const Mesh&,
226  const dictionary&
227  );
228 
229  //- Copy constructor
231 
232  //- Move constructor
234 
235  //- Construct as copy of tmp<GeometricField> deleting argument
237 
238  //- Construct as copy resetting IO parameters
240  (
241  const IOobject&,
243  );
244 
245  //- Construct as copy of tmp<GeometricField> resetting IO parameters
247  (
248  const IOobject&,
250  );
251 
252  //- Construct as copy resetting name
254  (
255  const word& newName,
257  );
258 
259  //- Construct as copy resetting name
261  (
262  const word& newName,
264  );
265 
266  //- Construct as copy resetting IO parameters and patch field type
268  (
269  const IOobject&,
271  const word& patchFieldType
272  );
273 
274  //- Construct as copy resetting IO parameters and boundary types
276  (
277  const IOobject&,
279  const wordList& patchFieldTypes,
280  const wordList& actualPatchTypes = wordList()
281  );
282 
283  //- Construct as copy resetting IO parameters and boundary types
285  (
286  const IOobject&,
288  const wordList& patchFieldTypes,
289  const wordList& actualPatchTypes = wordList()
290  );
291 
292  //- Clone
294 
295  //- Clone un-sliced
297 
298  //- Return a temporary field constructed from name,
299  // internal field and list of patch fields
301  (
302  const word& name,
303  const Internal&,
304  const PtrList<PatchField<Type>>&
305  );
306 
307  //- Return a temporary field constructed from name, mesh, dimensionSet
308  // and patch field type.
310  (
311  const word& name,
312  const Mesh&,
313  const dimensionSet&,
314  const word& patchFieldType=PatchField<Type>::calculatedType()
315  );
316 
317  //- Return a temporary field constructed from mesh, dimensioned<Type>
318  // and patch field type.
320  (
321  const word& name,
322  const Mesh&,
323  const dimensioned<Type>&,
324  const word& patchFieldType=PatchField<Type>::calculatedType()
325  );
326 
327  //- Return a temporary field constructed from mesh, dimensioned<Type>
328  // and patch field types.
330  (
331  const word& name,
332  const Mesh&,
333  const dimensioned<Type>&,
334  const wordList& patchFieldTypes,
335  const wordList& actualPatchTypes = wordList()
336  );
337 
338  //- Rename temporary field and return
340  (
341  const word& newName,
343  );
344 
345  //- Rename temporary field, reset patch field type and return
347  (
348  const word& newName,
350  const word&
351  );
352 
353  //- Rename and reset patch fields types of temporary field and return
355  (
356  const word& newName,
358  const wordList& patchFieldTypes,
359  const wordList& actualPatchTypes = wordList()
360  );
361 
362 
363  //- Destructor
364  virtual ~GeometricField();
365 
366 
367  // Member Functions
368 
369  //- Return a reference to the dimensioned internal field
370  // Note: this increments the event counter and checks the
371  // old-time fields; avoid in loops.
372  Internal& ref();
373 
374  //- Return a const-reference to the dimensioned internal field
375  inline const Internal& internalField() const;
376 
377  //- Return a const-reference to the dimensioned internal field
378  // of a "vol" field. Useful in the formulation of source-terms
379  // for FV equations
380  inline const Internal& v() const;
381 
382  //- Return a reference to the internal field
383  // Note: this increments the event counter and checks the
384  // old-time fields; avoid in loops.
386 
387  //- Return a const-reference to the internal field
388  inline const typename Internal::FieldType& primitiveField() const;
389 
390  //- Return a reference to the boundary field
391  // Note: this increments the event counter and checks the
392  // old-time fields; avoid in loops.
393  Boundary& boundaryFieldRef();
394 
395  //- Return const-reference to the boundary field
396  inline const Boundary& boundaryField() const;
397 
398  //- Return the time index of the field
399  inline label timeIndex() const;
400 
401  //- Return the time index of the field
402  inline label& timeIndex();
403 
404  //- Return whether or not this is an old-time field
405  bool isOldTime() const;
406 
407  //- Store the old-time fields
408  void storeOldTimes() const;
409 
410  //- Store the old-time field
411  void storeOldTime() const;
412 
413  //- Return the number of old time fields stored
414  label nOldTimes() const;
415 
416  //- Return old time field
418 
419  //- Return non-const old time field
420  // (Not a good idea but it is used for sub-cycling)
422 
423  //- Return the n-th old time field
425  (
426  const label n
427  ) const;
428 
429  //- Return the n-th non-const old time field
430  // (Not a good idea but it is used for sub-cycling)
432 
433  //- Store the field as the previous iteration value
434  void storePrevIter() const;
435 
436  //- Return previous iteration field
438 
439  //- Delete old time and previous iteration fields
440  void clearOldTimes();
441 
442  //- Correct boundary field
444 
445  //- Reset the field contents to the given field
446  // Used for mesh to mesh mapping
448 
449  //- Does the field need a reference level for solution
450  bool needReference() const;
451 
452  //- Return a component of the field
454  (
455  const direction
456  ) const;
457 
458  //- WriteData member function required by regIOobject
459  bool writeData(Ostream&) const;
460 
461  //- Return transpose (only if it is a tensor field)
463 
464  //- Relax field (for steady-state solution).
465  // alpha = 1 : no relaxation
466  // alpha < 1 : relaxation
467  // alpha = 0 : do nothing
468  void relax(const scalar alpha);
469 
470  //- Relax field (for steady-state solution).
471  // alpha is read from controlDict
472  void relax();
473 
474  //- Select the final iteration parameters if `final' is true
475  // by returning the field name + "Final"
476  // otherwise the standard parameters by returning the field name
477  word select(bool final) const;
478 
479  //- Helper function to write the min and max to an Ostream
480  void writeMinMax(Ostream& os) const;
481 
482 
483  // Member function *this operators
484 
485  void negate();
486 
487  void replace
488  (
489  const direction,
491  );
492 
493  void replace
494  (
495  const direction,
496  const dimensioned<cmptType>&
497  );
498 
499  void max(const dimensioned<Type>&);
500 
501  void min(const dimensioned<Type>&);
502 
503  void maxMin
504  (
505  const dimensioned<Type>& minDt,
506  const dimensioned<Type>& maxDt
507  );
508 
509 
510  // Member Operators
511 
512  //- Return a const-reference to the dimensioned internal field
513  // Useful in the formulation of source-terms for FV equations
514  inline const Internal& operator()() const;
515 
519  void operator=(const dimensioned<Type>&);
520  void operator=(const zero&);
521 
523  void operator==(const dimensioned<Type>&);
524  void operator==(const zero&);
525 
528 
531 
534 
537 
538  void operator+=(const dimensioned<Type>&);
539  void operator-=(const dimensioned<Type>&);
540 
541  void operator*=(const dimensioned<scalar>&);
542  void operator/=(const dimensioned<scalar>&);
543 
544 
545  // Ostream operators
546 
547  friend Ostream& operator<< <Type, PatchField, GeoMesh>
548  (
549  Ostream&,
551  );
552 
553  friend Ostream& operator<< <Type, PatchField, GeoMesh>
554  (
555  Ostream&,
557  );
558 };
559 
560 
561 template<class Type, template<class> class PatchField, class GeoMesh>
562 Ostream& operator<<
563 (
564  Ostream&,
566 );
567 
568 
569 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
570 
571 } // End namespace Foam
572 
573 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
574 
575 #include "GeometricFieldI.H"
576 
577 #ifdef NoRepository
578  #include "GeometricField.C"
579 #endif
580 
581 #include "GeometricFieldFunctions.H"
582 
583 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
584 
585 #endif
586 
587 // ************************************************************************* //
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
const Internal & operator()() const
Return a const-reference to the dimensioned internal field.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
void operator=(const GeometricField< Type, PatchField, GeoMesh > &)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
void clearOldTimes()
Delete old time and previous iteration fields.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
Field< Type >::cmptType cmptType
uint8_t direction
Definition: direction.H:45
pTraits< Type >::cmptType cmptType
Component type.
Definition: Field.H:97
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
void operator/=(const GeometricField< scalar, PatchField, GeoMesh > &)
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
GeoMesh::Mesh Mesh
Type of mesh on which this GeometricField is instantiated.
TypeName("GeometricField")
Runtime type information.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Generic GeometricField class.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
Generic dimensioned Type class.
label nOldTimes() const
Return the number of old time fields stored.
void maxMin(const dimensioned< Type > &minDt, const dimensioned< Type > &maxDt)
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
word select(bool final) const
Select the final iteration parameters if `final&#39; is true.
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
bool isOldTime() const
Return whether or not this is an old-time field.
bool needReference() const
Does the field need a reference level for solution.
Dimension set for the base types.
Definition: dimensionSet.H:121
const Internal & v() const
Return a const-reference to the dimensioned internal field.
Definition: volFieldsI.H:31
void storeOldTimes() const
Store the old-time fields.
void reset(const tmp< GeometricField< Type, PatchField, GeoMesh >> &)
Reset the field contents to the given field.
Pre-declare SubField and related Field type.
Definition: Field.H:56
A class for handling words, derived from string.
Definition: word.H:59
void storeOldTime() const
Store the old-time field.
void operator*=(const GeometricField< scalar, PatchField, GeoMesh > &)
void operator+=(const GeometricField< Type, PatchField, GeoMesh > &)
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
void min(const dimensioned< Type > &)
GeometricField(const IOobject &, const Mesh &, const dimensionSet &, const word &patchFieldType=PatchField< Type >::calculatedType())
Constructor given IOobject, mesh, dimensions and patch field type.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
void operator-=(const GeometricField< Type, PatchField, GeoMesh > &)
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
MESH Mesh
Definition: GeoMesh.H:61
virtual ~GeometricField()
Destructor.
Generic GeometricBoundaryField class.
Definition: fvMesh.H:80
label timeIndex() const
Return the time index of the field.
Internal & ref()
Return a reference to the dimensioned internal field.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
List< word > wordList
A List of words.
Definition: fileName.H:54
void operator==(const tmp< GeometricField< Type, PatchField, GeoMesh >> &)
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
PatchField< Type > Patch
Type of the patch field of which the Boundary is composed.
tmp< GeometricField< Type, PatchField, GeoMesh > > cloneUnSliced() const
Clone un-sliced.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
void max(const dimensioned< Type > &)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void correctBoundaryConditions()
Correct boundary field.
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:49
label n
void storePrevIter() const
Store the field as the previous iteration value.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
A class for managing temporary objects.
Definition: PtrList.H:53
word patchFieldType(const PatchField &pf)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
List of coupled interface fields to be used in coupling.
void relax()
Relax field (for steady-state solution).
Namespace for OpenFOAM.