GeometricField.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::GeometricField
29 
30 Description
31  Generic GeometricField class.
32 
33 SourceFiles
34  GeometricFieldI.H
35  GeometricField.C
36  GeometricFieldFunctions.H
37  GeometricFieldFunctions.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef Foam_GeometricField_H
42 #define Foam_GeometricField_H
43 
44 #include "regIOobject.H"
45 #include "GeometricBoundaryField.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward Declarations
53 class dictionary;
54 
55 template<class Type, template<class> class PatchField, class GeoMesh>
56 class GeometricField;
57 
58 template<class Type, template<class> class PatchField, class GeoMesh>
59 Ostream& operator<<
60 (
61  Ostream&,
62  const GeometricField<Type, PatchField, GeoMesh>&
63 );
64 
65 template<class Type, template<class> class PatchField, class GeoMesh>
66 Ostream& operator<<
67 (
68  Ostream&,
69  const tmp<GeometricField<Type, PatchField, GeoMesh>>&
70 );
71 
72 
73 /*---------------------------------------------------------------------------*\
74  Class GeometricField Declaration
75 \*---------------------------------------------------------------------------*/
76 
77 template<class Type, template<class> class PatchField, class GeoMesh>
78 class GeometricField
79 :
80  public DimensionedField<Type, GeoMesh>
81 {
82 public:
83 
84  // Public Typedefs
85 
86  //- The mesh type for the GeometricField
87  typedef typename GeoMesh::Mesh Mesh;
88 
89  //- The boundary mesh type for the boundary fields
90  typedef typename GeoMesh::BoundaryMesh BoundaryMesh;
91 
92  //- The internal field type from which this GeometricField is derived
94 
95  //- Type of boundary fields
97 
98  //- The patch field type for the GeometricBoundaryField
99  typedef PatchField<Type> Patch;
100 
101  //- Component type of the field elements
102  typedef typename Field<Type>::cmptType cmptType;
103 
105 private:
106 
107  // Private Data
108 
109  //- Current time index.
110  // Used to trigger the storing of the old-time value
111  mutable label timeIndex_;
112 
113  //- Pointer to old time field
114  mutable GeometricField<Type, PatchField, GeoMesh>* field0Ptr_;
115 
116  //- Pointer to previous iteration (used for under-relaxation)
117  mutable GeometricField<Type, PatchField, GeoMesh>* fieldPrevIterPtr_;
118 
119  //- Boundary field containing boundary field values
120  Boundary boundaryField_;
121 
122 
123  // Private Member Functions
124 
125  //- Read from file if it is present
126  bool readIfPresent();
127 
128  //- Read old time field from file if it is present
129  bool readOldTimeIfPresent();
130 
131  //- Read the field from the dictionary
132  void readFields(const dictionary& dict);
133 
134  //- Read the field - create the field dictionary on-the-fly
135  void readFields();
136 
137 
138 public:
139 
140  //- Runtime type information
141  TypeName("GeometricField");
142 
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  //- Construct given IOobject, mesh, dimensions and patch type.
153  // This allocates storage for the field but does not set values.
154  // Used only within this class to create TEMPORARY variables
156  (
157  const IOobject& io,
158  const Mesh& mesh,
159  const dimensionSet& ds,
160  const word& patchFieldType = PatchField<Type>::calculatedType()
161  );
162 
163  //- Construct given IOobject, mesh, dimensions and patch types.
164  // This allocates storage for the field but does not set values.
165  // Used only within this class to create TEMPORARY variables
167  (
168  const IOobject& io,
169  const Mesh& mesh,
170  const dimensionSet& ds,
171  const wordList& wantedPatchTypes,
172  const wordList& actualPatchTypes = wordList()
173  );
174 
175  //- Construct given IOobject, mesh, dimensioned<Type> and patch type.
176  // This assigns both dimensions and values.
177  // The name of the dimensioned<Type> has no influence.
179  (
180  const IOobject& io,
181  const Mesh& mesh,
182  const dimensioned<Type>& dt,
183  const word& patchFieldType = PatchField<Type>::calculatedType()
184  );
185 
186  //- Construct given IOobject, mesh, dimensioned<Type> and patch types.
187  // This assigns both dimensions and values.
188  // The name of the dimensioned<Type> has no influence.
190  (
191  const IOobject& io,
192  const Mesh& mesh,
193  const dimensioned<Type>& dt,
194  const wordList& wantedPatchTypes,
195  const wordList& actualPatchTypes = wordList()
196  );
197 
198  //- Copy construct from internal field and a patch list to clone
200  (
201  const IOobject& io,
202  const Internal& diField,
203  const PtrList<PatchField<Type>>& ptfl
204  );
205 
206  //- Move construct from internal field and a patch list to clone
208  (
209  const IOobject& io,
210  Internal&& diField,
211  const PtrList<PatchField<Type>>& ptfl
212  );
213 
214  //- Move construct from internal field and a patch list to clone
216  (
217  const IOobject& io,
218  const tmp<Internal>& tfield,
219  const PtrList<PatchField<Type>>& ptfl
220  );
221 
222  //- Copy construct from internal field and a patch list to clone
224  (
225  const Internal& diField,
226  const PtrList<PatchField<Type>>& ptfl
227  );
228 
229  //- Move construct from internal field and a patch list to clone
231  (
232  Internal&& diField,
233  const PtrList<PatchField<Type>>& ptfl
234  );
235 
236  //- Copy construct from internal field, with specified patch type
238  (
239  const IOobject& io,
240  const Mesh& mesh,
241  const dimensionSet& ds,
242  const Field<Type>& iField,
243  const word& patchFieldType = PatchField<Type>::calculatedType()
244  );
245 
246  //- Move construct from internal field, with specified patch type
248  (
249  const IOobject& io,
250  const Mesh& mesh,
251  const dimensionSet& ds,
252  Field<Type>&& iField,
253  const word& patchFieldType = PatchField<Type>::calculatedType()
254  );
255 
256  //- Copy construct from components
258  (
259  const IOobject& io,
260  const Mesh& mesh,
261  const dimensionSet& ds,
262  const Field<Type>& iField,
263  const PtrList<PatchField<Type>>& ptfl
264  );
265 
266  //- Move construct from internal field and a patch list to clone
268  (
269  const IOobject& io,
270  const Mesh& mesh,
271  const dimensionSet& ds,
272  Field<Type>&& iField,
273  const PtrList<PatchField<Type>>& ptfl
274  );
275 
276  //- Copy construct from components
278  (
279  const IOobject& io,
280  const Mesh& mesh,
281  const dimensionSet& ds,
282  const tmp<Field<Type>>& tiField,
283  const PtrList<PatchField<Type>>& ptfl
284  );
285 
286  //- Construct and read given IOobject
288  (
289  const IOobject& io,
290  const Mesh& mesh,
291  const bool readOldTime = true
292  );
293 
294  //- Construct from dictionary
296  (
297  const IOobject& io,
298  const Mesh& mesh,
299  const dictionary& dict
300  );
301 
302  //- Copy construct
304  (
306  );
307 
308  //- Construct from tmp<GeometricField> deleting argument
310  (
312  );
313 
314  //- Construct as copy resetting IO parameters
316  (
317  const IOobject& io,
319  );
320 
321  //- Construct from tmp<GeometricField> resetting IO parameters
323  (
324  const IOobject& io,
326  );
327 
328  //- Copy construct with a new name
330  (
331  const word& newName,
333  );
334 
335  //- Construct with a new name from tmp<GeometricField>
337  (
338  const word& newName,
340  );
341 
342  //- Construct as copy resetting IO parameters and patch type
344  (
345  const IOobject& io,
347  const word& patchFieldType
348  );
349 
350  //- Construct as copy resetting IO parameters and boundary type
351  //- for selected patchIDs
353  (
354  const IOobject& io,
356  const labelList& patchIDs,
357  const word& patchFieldType
358  );
359 
360  //- Construct as copy resetting IO parameters and boundary types
362  (
363  const IOobject& io,
365  const wordList& patchFieldTypes,
366  const wordList& actualPatchTypes = wordList()
367  );
368 
369  //- Construct as copy resetting IO parameters and boundary types
371  (
372  const IOobject& io,
374  const wordList& patchFieldTypes,
375  const wordList& actualPatchTypes = wordList()
376  );
377 
378  //- Clone
380 
381 
382  // Factory Methods
383 
384  //- Return tmp field from name, mesh, dimensions and patch type.
385  // The field is NO_READ, NO_WRITE, unregistered and uses the
386  // current timeName from the mesh registry
388  (
389  const word& name,
390  const Mesh& mesh,
391  const dimensionSet& ds,
392  const word& patchFieldType = PatchField<Type>::calculatedType()
393  );
394 
395  //- Return tmp field from name, mesh, dimensions,
396  //- copy of internal field, with specified patch type.
397  // The field is NO_READ, NO_WRITE, unregistered and uses the
398  // current timeName from the mesh registry
400  (
401  const word& name,
402  const Mesh& mesh,
403  const dimensionSet& ds,
404  const Field<Type>& iField,
405  const word& patchFieldType = PatchField<Type>::calculatedType()
406  );
407 
408  //- Return tmp field from name, mesh, dimensions,
409  //- moved internal field contents, with specified patch type.
410  // The field is NO_READ, NO_WRITE, unregistered and uses the
411  // current timeName from the mesh registry
413  (
414  const word& name,
415  const Mesh& mesh,
416  const dimensionSet& ds,
417  Field<Type>&& iField,
418  const word& patchFieldType = PatchField<Type>::calculatedType()
419  );
420 
421  //- Return tmp field from name, mesh, dimensioned<Type>
422  //- and patch type.
423  // The field is NO_READ, NO_WRITE, unregistered and uses the
424  // current timeName from the mesh registry
426  (
427  const word& name,
428  const Mesh& mesh,
429  const dimensioned<Type>& dt,
430  const word& patchFieldType = PatchField<Type>::calculatedType()
431  );
432 
433  //- Return tmp field from name, mesh, dimensioned<Type>
434  //- and patch types.
435  // The field is NO_READ, NO_WRITE, unregistered and uses the
436  // current timeName from the mesh registry
438  (
439  const word& name,
440  const Mesh& mesh,
441  const dimensioned<Type>& dt,
442  const wordList& patchFieldTypes,
443  const wordList& actualPatchTypes = wordList()
444  );
445 
446  //- Return renamed tmp field
447  // The field is NO_READ, NO_WRITE, unregistered and uses the
448  // current timeName from the mesh registry
450  (
451  const word& newName,
453  );
454 
455  //- Rename tmp field and reset patch field types and return
456  // The field is NO_READ, NO_WRITE, unregistered and uses the
457  // current timeName from the mesh registry
459  (
460  const word& newName,
462  const wordList& patchFieldTypes,
463  const wordList& actualPatchTypes = wordList()
464  );
465 
466  //- Construct tmp field based on mesh/registry information from
467  //- an existing field.
468  // Created NO_READ, NO_WRITE, NO_REGISTER, using the instance
469  // from the field
470  template<class AnyType>
472  (
474  const word& name,
475  const dimensionSet& dims,
476  const word& patchFieldType = PatchField<Type>::calculatedType()
477  );
478 
479  //- Construct tmp field based on mesh/registry information from
480  //- an existing field and initialise with value.
481  // Created NO_READ, NO_WRITE, NO_REGISTER, using the instance
482  // from the field
483  template<class AnyType>
485  (
487  const word& name,
488  const dimensioned<Type>& dt,
489  const word& patchFieldType = PatchField<Type>::calculatedType()
490  );
491 
492 
493  //- Destructor
494  virtual ~GeometricField();
495 
496 
497  // Member Functions
498 
499  //- Return a const-reference to the dimensioned internal field.
500  inline const Internal& internalField() const noexcept;
501 
502  //- Return a reference to the dimensioned internal field.
503  // \param updateAccessTime update event counter and check
504  // old-time fields
505  //
506  // \note Should avoid using updateAccessTime = true within loops.
507  Internal& internalFieldRef(const bool updateAccessTime = true);
508 
509  //- Same as internalFieldRef()
510  Internal& ref(const bool updateAccessTime = true)
511  {
512  return this->internalFieldRef(updateAccessTime);
513  }
514 
515  //- Return a const-reference to the dimensioned internal field
516  //- of a "vol" field.
517  // Useful in the formulation of source-terms for FV equations
518  //
519  // \note definition in finiteVolume/fields/volFields/volFieldsI.H
520  inline const Internal& v() const;
521 
522  //- Return a const-reference to the internal field values.
523  inline const typename Internal::FieldType& primitiveField()
524  const noexcept;
525 
526  //- Return a reference to the internal field values.
527  // \param updateAccessTime update event counter and check
528  // old-time fields
529  //
530  // \note Should avoid using updateAccessTime = true within loops.
532  (
533  const bool updateAccessTime = true
534  );
535 
536  //- Return const-reference to the boundary field
537  inline const Boundary& boundaryField() const noexcept;
538 
539  //- Return a reference to the boundary field
540  // \param updateAccessTime update event counter and check
541  // old-time fields
542  //
543  // \note Should avoid using updateAccessTime = true within loops.
544  Boundary& boundaryFieldRef(const bool updateAccessTime = true);
545 
546  //- Return the time index of the field
547  inline label timeIndex() const noexcept;
548 
549  //- Write-access to the time index of the field
550  inline label& timeIndex() noexcept;
551 
552  //- Store the old-time fields
553  void storeOldTimes() const;
554 
555  //- Store the old-time field
556  void storeOldTime() const;
557 
558  //- Return the number of old time fields stored
559  label nOldTimes() const;
560 
561  //- Return old time field
562  const GeometricField<Type, PatchField, GeoMesh>& oldTime() const;
563 
564  //- Return non-const old time field
565  // (Not a good idea but it is used for sub-cycling)
566  GeometricField<Type, PatchField, GeoMesh>& oldTime();
567 
568  //- Store the field as the previous iteration value
569  void storePrevIter() const;
570 
571  //- Return previous iteration field
572  const GeometricField<Type, PatchField, GeoMesh>& prevIter() const;
573 
574  //- Correct boundary field
576 
577  //- Does the field need a reference level for solution
578  bool needReference() const;
579 
580  //- Return a component of the field
582  (
583  const direction
584  ) const;
585 
586  //- WriteData member function required by regIOobject
587  bool writeData(Ostream&) const;
588 
589  //- Return transpose (only if it is a tensor field)
590  tmp<GeometricField<Type, PatchField, GeoMesh>> T() const;
591 
592  //- Relax field (for steady-state solution).
593  // alpha = 1 : no relaxation
594  // alpha < 1 : relaxation
595  // alpha = 0 : do nothing
596  void relax(const scalar alpha);
597 
598  //- Relax field (for steady-state solution).
599  // alpha is read from controlDict
600  void relax();
601 
602  //- Select the final iteration parameters if `final' is true
603  // by returning the field name + "Final"
604  // otherwise the standard parameters by returning the field name
605  word select(bool final) const;
606 
607  //- Helper function to write the min and max to an Ostream
608  void writeMinMax(Ostream& os) const;
609 
610 
611  // Member Function *this Operators
612 
613  //- Negate the field inplace. See notes in Field
614  void negate();
615 
616  //- Normalise the field inplace. See notes in Field
617  void normalise();
618 
619  //- Replace specified field component with content from another field
620  void replace
621  (
622  const direction d,
623  const GeometricField<cmptType, PatchField, GeoMesh>& gcf
624  );
625 
626  //- Replace specified field component with specified value
627  void replace
628  (
629  const direction d,
630  const dimensioned<cmptType>& ds
631  );
633  //- Use the minimum of the field and specified value
634  // This sets the \em ceiling on the field values
635  void min(const dimensioned<Type>& dt);
636 
637  //- Use the maximum of the field and specified value
638  // This sets the \em floor on the field values
639  void max(const dimensioned<Type>& dt);
640 
641  //- Clip the field to be bounded within the specified range
642  void clip(const dimensioned<MinMax<Type>>& range);
643 
644  //- Clip the field to be bounded within the specified range
645  void clip
646  (
647  const dimensioned<Type>& minVal,
648  const dimensioned<Type>& maxVal
649  );
650 
651  //- Deprecated(2019-01) identical to clip()
652  // \deprecated(2019-01) identical to clip()
653  FOAM_DEPRECATED_FOR(2019-01, "clip() method")
654  void maxMin
655  (
656  const dimensioned<Type>& minVal,
657  const dimensioned<Type>& maxVal
658  );
659 
660 
661  // Member Operators
662 
663  //- Return a const-reference to the dimensioned internal field.
664  //- Same as internalField().
665  // Useful in the formulation of source-terms for FV equations
666  const Internal& operator()() const { return *this; }
667 
670  void operator=(const dimensioned<Type>&);
671 
673  void operator==(const dimensioned<Type>&);
674 
677 
680 
683 
686 
687  void operator+=(const dimensioned<Type>&);
688  void operator-=(const dimensioned<Type>&);
689 
690  void operator*=(const dimensioned<scalar>&);
691  void operator/=(const dimensioned<scalar>&);
692 
693 
694  // Ostream operators
695 
696  friend Ostream& operator<< <Type, PatchField, GeoMesh>
697  (
698  Ostream&,
700  );
701 
702  friend Ostream& operator<< <Type, PatchField, GeoMesh>
703  (
704  Ostream&,
706  );
707 };
708 
709 
710 template<class Type, template<class> class PatchField, class GeoMesh>
711 Ostream& operator<<
712 (
713  Ostream&,
714  const typename GeometricField<Type, PatchField, GeoMesh>::
715  Boundary&
716 );
717 
718 
719 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
720 
721 } // End namespace Foam
722 
723 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
724 
725 #include "GeometricFieldI.H"
726 
727 #ifdef NoRepository
728  #include "GeometricField.C"
729 #endif
730 
731 #include "GeometricFieldFunctions.H"
732 
733 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
734 
735 #endif
736 
737 // ************************************************************************* //
dictionary dict
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
uint8_t direction
Definition: direction.H:46
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
GeoMesh::BoundaryMesh BoundaryMesh
The boundary mesh type for the boundary fields.
void operator=(const GeometricField< Type, PatchField, GeoMesh > &)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
Field< Type >::cmptType cmptType
Component type of the field elements.
MESH::BoundaryMesh BoundaryMesh
Definition: GeoMesh.H:59
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
pTraits< Type >::cmptType cmptType
Component type.
Definition: Field.H:123
A min/max value pair with additional methods. In addition to conveniently storing values...
Definition: HashSet.H:72
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
The mesh type for the GeometricField.
void clip(const dimensioned< MinMax< Type >> &range)
Clip the field to be bounded within the specified range.
TypeName("GeometricField")
Runtime type information.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of boundary fields.
Generic dimensioned Type class.
label nOldTimes() const
Return the number of old time fields stored.
GeometricField(const IOobject &io, const Mesh &mesh, const dimensionSet &ds, const word &patchFieldType=PatchField< Type >::calculatedType())
Construct given IOobject, mesh, dimensions and patch type.
DimensionedField< Type, GeoMesh > Internal
The internal field type from which this GeometricField is derived.
Field< Type > FieldType
Type of the field from which this DimensionedField is derived.
word select(bool final) const
Select the final iteration parameters if `final&#39; is true.
scalar range
class FOAM_DEPRECATED_FOR(2017-05, "Foam::Enum") NamedEnum
Definition: NamedEnum.H:65
bool needReference() const
Does the field need a reference level for solution.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &ds, const word &patchFieldType=PatchField< Type >::calculatedType())
Return tmp field from name, mesh, dimensions and patch type.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
const Internal & v() const
Return a const-reference to the dimensioned internal field of a "vol" field.
Definition: volFieldsI.H:26
void storeOldTimes() const
Store the old-time fields.
Generic templated field type.
Definition: Field.H:61
A class for handling words, derived from Foam::string.
Definition: word.H:63
void min(const dimensioned< Type > &dt)
Use the minimum of the field and specified value.
void storeOldTime() const
Store the old-time field.
void operator*=(const GeometricField< scalar, PatchField, GeoMesh > &)
void maxMin(const dimensioned< Type > &minVal, const dimensioned< Type > &maxVal)
Deprecated(2019-01) identical to clip()
void operator+=(const GeometricField< Type, PatchField, GeoMesh > &)
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
GeoMesh::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
void negate()
Negate the field inplace. See notes in Field.
const Mesh & mesh() const noexcept
Return mesh.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
const direction noexcept
Definition: Scalar.H:258
Internal & internalFieldRef(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
void operator-=(const GeometricField< Type, PatchField, GeoMesh > &)
void normalise()
Normalise the field inplace. See notes in Field.
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
OBJstream os(runTime.globalPath()/outputName)
MESH Mesh
Definition: GeoMesh.H:58
virtual ~GeometricField()
Destructor.
Generic GeometricBoundaryField class.
Definition: areaFieldsFwd.H:46
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
List< word > wordList
A List of words.
Definition: fileName.H:58
void operator==(const tmp< GeometricField< Type, PatchField, GeoMesh >> &)
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
label timeIndex() const noexcept
Return the time index of the field.
PatchField< Type > Patch
The patch field type for the GeometricBoundaryField.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
friend Ostream & operator(Ostream &, const GeometricField< Type, PatchField, GeoMesh > &)
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
void storePrevIter() const
Store the field as the previous iteration value.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:42
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
void relax()
Relax field (for steady-state solution).
Namespace for OpenFOAM.
void replace(const direction d, const GeometricField< cmptType, PatchField, GeoMesh > &gcf)
Replace specified field component with content from another field.