fvMeshTemplates.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) 2015-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 "fvMesh.H"
27 
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 template<class Type>
32 {
33  return pow
34  (
35  this->solutionD(),
36  pTraits
37  <
38  typename powProduct<Vector<label>,
40  >::zero
41  );
42 }
43 
44 
45 template<class GeoField>
47 {
49  (
50  const_cast<fvMesh&>(*this).lookupClass<GeoField>(strict)
51  );
52  UPtrList<GeoField> curFields(fields.size());
53 
54  label i = 0;
55  forAllIter(typename HashTable<GeoField*>, fields, iter)
56  {
57  if (!geometryFields.found(iter()->name()))
58  {
59  curFields.set(i++, iter());
60  }
61  }
62  curFields.setSize(i);
63 
64  return curFields;
65 }
66 
67 
68 template<class GeoField>
70 {
72  (
73  const_cast<fvMesh&>(*this).lookupClass<GeoField>()
74  );
75  UPtrList<GeoField> curFields(fields.size());
76 
77  label i = 0;
78  forAllIter(typename HashTable<GeoField*>, fields, iter)
79  {
80  if (!geometryFields.found(iter()->name()) && !iter()->isOldTime())
81  {
82  curFields.set(i++, iter());
83  }
84  }
85  curFields.setSize(i);
86 
87  return curFields;
88 }
89 
90 
91 template<class Type, template<class> class GeoField>
92 void Foam::fvMesh::storeOldTimeFields()
93 {
94  UPtrList<GeoField<Type>> curFields(this->curFields<GeoField<Type>>());
95 
96  forAll(curFields, i)
97  {
98  curFields[i].storeOldTimes();
99  }
100 }
101 
102 
103 template<template<class> class GeoField>
104 void Foam::fvMesh::storeOldTimeFields()
105 {
106  #define StoreOldTimeFields(Type, nullArg) \
107  storeOldTimeFields<Type, GeoField>();
109  #undef StoreOldTimeFields
110 }
111 
112 
113 template<class Type, template<class> class GeoField>
114 void Foam::fvMesh::nullOldestTimeFields()
115 {
116  UPtrList<GeoField<Type>> curFields(this->curFields<GeoField<Type>>());
117 
118  forAll(curFields, i)
119  {
120  curFields[i].nullOldestTime();
121  }
122 }
123 
124 
125 template<template<class> class GeoField>
126 void Foam::fvMesh::nullOldestTimeFields()
127 {
128  #define nullOldestTimeFields(Type, nullArg) \
129  nullOldestTimeFields<Type, GeoField>();
131  #undef nullOldestTimeFields
132 }
133 
134 
135 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
An STL-conforming hash table.
Definition: HashTable.H:127
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
bool set(const label) const
Is element set.
Definition: UPtrListI.H:78
void setSize(const label)
Reset size of UPtrList. This can only be used to set the size.
Definition: UPtrList.C:54
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
pTraits< Type >::labelType validComponents() const
Return a labelType of valid component indicators.
UPtrList< GeoField > fields(const bool strict=false) const
Return the list of fields of type GeoField.
UPtrList< GeoField > curFields() const
Return the list of current fields of type GeoField.
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
Traits class for primitives.
Definition: pTraits.H:53
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:1060
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:50
#define StoreOldTimeFields(Type, nullArg)
#define nullOldestTimeFields(Type, nullArg)
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
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 name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
FOR_ALL_FIELD_TYPES(DefineContiguousFvWallLocationDataType)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488