LagrangianMeshTemplates.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 "LagrangianMesh.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Type, template<class> class GeoField>
31 void Foam::LagrangianMesh::appendSpecifiedField
32 (
33  const LagrangianSubMesh& appendMesh,
34  GeoField<Type>& geoField,
35  const Field<Type>& field
36 ) const
37 {
38  geoField.resize(size());
39 
40  appendMesh.sub(geoField).ref().primitiveFieldRef() = field;
41 }
42 
43 
44 template<class Type, template<class> class GeoField>
45 bool Foam::LagrangianMesh::appendSpecifiedField
46 (
47  const LagrangianSubMesh& appendMesh,
48  const word& fieldName,
49  const Field<Type>& field
50 ) const
51 {
52  if (foundObject<GeoField<Type>>(fieldName))
53  {
54  appendSpecifiedField<Type, GeoField>
55  (
56  appendMesh,
57  lookupObjectRef<GeoField<Type>>(fieldName),
58  field
59  );
60 
61  return true;
62  }
63  else
64  {
65  return false;
66  }
67 }
68 
69 
70 template<class Type, class ... FieldNamesAndFields>
71 Foam::wordHashSet Foam::LagrangianMesh::appendSpecifiedFields
72 (
73  const LagrangianSubMesh& appendMesh,
74  const word& fieldName,
75  const Field<Type>& field,
76  const FieldNamesAndFields& ... fieldNamesAndFields
77 ) const
78 {
79  const bool found =
80  appendSpecifiedField<Type, LagrangianField>
81  (
82  appendMesh,
83  fieldName,
84  field
85  )
86  || appendSpecifiedField<Type, LagrangianDynamicField>
87  (
88  appendMesh,
89  fieldName,
90  field
91  )
92  || appendSpecifiedField<Type, LagrangianInternalField>
93  (
94  appendMesh,
95  fieldName,
96  field
97  )
98  || appendSpecifiedField<Type, LagrangianInternalDynamicField>
99  (
100  appendMesh,
101  fieldName,
102  field
103  );
104 
105  if (!found)
106  {
107  static const word& TypeName = pTraits<Type>::typeName;
108 
114  sort(fieldNames);
115 
117  << TypeName.capitalise() << " values were specified for "
118  << fieldName << " but a " << TypeName << " Lagrangian field of "
119  << "that name was not found." << nl << nl
120  << "Available " << TypeName << " Lagrangian fields are:" << nl
121  << fieldNames << exit(FatalError);
122  }
123 
124  wordHashSet result
125  (
126  appendSpecifiedFields(appendMesh, fieldNamesAndFields ...)
127  );
128 
129  result.insert(fieldName);
130 
131  return result;
132 }
133 
134 
135 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
136 
137 template<class GeoField>
140 {
141  HashTable<const GeoField*> fields(lookupClass<GeoField>(strict));
142 
143  forAllIter(typename HashTable<GeoField*>, fields, iter)
144  {
145  if (iter()->isOldTime())
146  {
147  fields.erase(iter);
148  }
149  }
150 
151  return fields;
152 }
153 
154 
155 template<class GeoField>
158 {
159  HashTable<GeoField*> fields(lookupClass<GeoField>(strict));
160 
161  forAllIter(typename HashTable<GeoField*>, fields, iter)
162  {
163  if (iter()->isOldTime())
164  {
165  fields.erase(iter);
166  }
167  }
168 
169  return fields;
170 }
171 
172 
173 template<class ... FieldNamesAndFields>
175 (
177  const labelField& celli,
178  const labelField& facei,
179  const labelField& faceTrii,
180  const FieldNamesAndFields& ... fieldNamesAndFields
181 )
182 {
183  const LagrangianSubMesh injectionMesh =
184  append(coordinates, celli, facei, faceTrii);
185 
186  const wordHashSet specifiedFieldNames =
187  appendSpecifiedFields(injectionMesh, fieldNamesAndFields ...);
188 
189  injectUnspecifiedFields(injectionMesh, specifiedFieldNames);
190 
191  return injectionMesh;
192 }
193 
194 
195 template<class ... FieldNamesAndFields>
197 (
198  const LagrangianInjection& injection,
200  const labelField& celli,
201  const labelField& facei,
202  const labelField& faceTrii,
203  const FieldNamesAndFields& ... fieldNamesAndFields
204 )
205 {
206  const LagrangianSubMesh injectionMesh =
207  append(coordinates, celli, facei, faceTrii);
208 
209  const wordHashSet specifiedFieldNames =
210  appendSpecifiedFields(injectionMesh, fieldNamesAndFields ...);
211 
212  injectUnspecifiedFields(injection, injectionMesh, specifiedFieldNames);
213 
214  return injectionMesh;
215 }
216 
217 
218 template<class ... FieldNamesAndFields>
220 (
221  const labelList& parents,
222  const FieldNamesAndFields& ... fieldNamesAndFields
223 )
224 {
225  const LagrangianSubMesh birthMesh = append(parents);
226 
227  const wordHashSet specifiedFieldNames =
228  appendSpecifiedFields(birthMesh, fieldNamesAndFields ...);
229 
230  birthUnspecifiedFields(parents, birthMesh, specifiedFieldNames);
231 
232  return birthMesh;
233 }
234 
235 
236 template<class Enumeration>
238 (
239  const label nGroups,
240  const UList<elementGroup<Enumeration>>& elementsGroups
241 )
242 {
243  return
244  partition
245  (
246  nGroups,
247  reinterpret_cast<const UList<labelPair>&>(elementsGroups)
248  );
249 }
250 
251 
252 // ************************************************************************* //
bool found
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:458
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pre-declare SubField and related Field type.
Definition: Field.H:83
Generic GeometricField class.
A HashTable with keys but without contents.
Definition: HashSet.H:62
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
An STL-conforming hash table.
Definition: HashTable.H:127
Base class for Lagrangian injections. Minimal wrapper over LagrangianSource. Implements some utility ...
Class to hold element-group indices, and associate the group.
LagrangianSubMesh inject(const barycentricField &coordinates, const labelField &celli, const labelField &facei, const labelField &faceTrii, const FieldNamesAndFields &... fieldNamesAndFields)
Inject specified elements into the mesh. This method does not.
label size() const
Get the number of elements.
LagrangianSubMesh birth(const labelList &parents, const FieldNamesAndFields &... fieldNamesAndFields)
Birth specified elements into the mesh. Fields are expanded.
HashTable< const GeoField * > lookupCurrentFields(const bool strict=false) const
Lookup all current-time fields of the given type.
void partition()
Partition the mesh such that the groups are contiguous in memory.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:74
Traits class for primitives.
Definition: pTraits.H:53
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
static List< word > fieldNames
Definition: globalFoam.H:46
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:234
barycentric coordinates(const polyMesh &mesh, const point &position, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the coordinates given the position and tet topology.
Definition: tracking.C:1259
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
void sort(UList< T > &)
Definition: UList.C:115
error FatalError
static const char nl
Definition: Ostream.H:267
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:74