MeshToMeshMapGeometricFields.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) 2022-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 Class
25  Foam::MeshToMeshMapGeometricFields
26 
27 Description
28  Generic internal field mapper. For "real" mapping, add template
29  specialisations for mapping of internal fields depending on mesh
30  type.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef MeshToMeshMapGeometricFields_H
35 #define MeshToMeshMapGeometricFields_H
36 
37 #include "polyMesh.H"
38 #include "fvMeshToFvMesh.H"
39 #include "fvPatchFieldMapper.H"
40 #include "pointPatchFieldMapper.H"
41 #include "setSizeFieldMapper.H"
42 #include "processorPolyPatch.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 template<class Type>
51 (
52  const fvMesh& mesh,
53  const fvMeshToFvMesh& mapper
54 )
55 {
57  (
58  mesh.objectRegistry::template lookupClass<VolField<Type>>()
59  );
60 
61  // Deleted old time fields
62  for
63  (
64  typename HashTable<const VolField<Type>*>::
65  iterator fieldIter = fields.begin();
66  fieldIter != fields.end();
67  ++fieldIter
68  )
69  {
70  VolField<Type>& field =
71  const_cast<VolField<Type>&>
72  (*fieldIter());
73 
74  field.clearOldTimes();
75  }
76 
77  fields =
78  (
79  mesh.objectRegistry::template lookupClass<VolField<Type>>()
80  );
81 
82  for
83  (
84  typename HashTable<const VolField<Type>*>::
85  iterator fieldIter = fields.begin();
86  fieldIter != fields.end();
87  ++fieldIter
88  )
89  {
90  VolField<Type>& field =
91  const_cast<VolField<Type>&>(*fieldIter());
92 
93  if (fvMeshToFvMesh::debug)
94  {
95  Info<< "Mapping " << field.typeName << ' ' << field.name()
96  << endl;
97  }
98 
99  field.reset(mapper.srcToTgt(field));
100 
101  field.instance() = field.time().name();
102  }
103 }
104 
105 
106 template<class Type>
108 (
109  const fvMesh& mesh,
110  const fvMeshToFvMesh& mapper
111 )
112 {
113  // Note: use strict flag on lookupClass to avoid picking up volFields
115  (
116  mesh.objectRegistry::template lookupClass<VolInternalField<Type>>(true)
117  );
118 
119  for
120  (
121  typename HashTable<const VolInternalField<Type>*>::
122  iterator fieldIter = fields.begin();
123  fieldIter != fields.end();
124  ++fieldIter
125  )
126  {
127  VolInternalField<Type>& field =
128  const_cast<VolInternalField<Type>&>(*fieldIter());
129 
130  if (fvMeshToFvMesh::debug)
131  {
132  Info<< "Mapping " << field.typeName << ' ' << field.name()
133  << endl;
134  }
135 
136  field.reset(mapper.srcToTgt<Type>(field));
137 
138  field.instance() = field.time().name();
139  }
140 }
141 
142 
143 template
144 <
145  class Type,
146  template<class> class PatchField,
147  class GeoMesh,
148  class SetSizePatchFieldMapper
149 >
151 (
152  const fvMesh& mesh,
153  const fvMeshToFvMesh& mapper
154 )
155 {
157 
159  (
160  mesh.objectRegistry::template lookupClass<GField>()
161  );
162 
163  // Deleted old time fields
164  for
165  (
166  typename HashTable<const GField*>::
167  iterator fieldIter = fields.begin();
168  fieldIter != fields.end();
169  ++fieldIter
170  )
171  {
172  GField& field = const_cast<GField&>(*fieldIter());
173  field.clearOldTimes();
174  }
175 
176  fields = mesh.objectRegistry::template lookupClass<GField>();
177 
178  Type NaN;
179 
180  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
181  {
182  setComponent(NaN, cmpt) = std::numeric_limits<scalar>::signaling_NaN();
183  }
184 
185  for
186  (
187  typename HashTable<const GField*>::iterator fieldIter = fields.begin();
188  fieldIter != fields.end();
189  ++fieldIter
190  )
191  {
192  GField& field = const_cast<GField&>(*fieldIter());
193 
194  if (fvMeshToFvMesh::debug)
195  {
196  Info<< "Setting to NaN " << field.typeName << ' ' << field.name()
197  << endl;
198  }
199 
200  const typename GField::Mesh& mesh = field.mesh();
201 
202  field.primitiveFieldRef().setSize(GeoMesh::size(mesh));
203  field.primitiveFieldRef() = NaN;
204 
205  field.boundaryFieldRef().setSize(mesh.boundary().size());
206 
207  forAll(mesh.boundary(), patchi)
208  {
209  if (isA<processorPolyPatch>(mesh().boundaryMesh()[patchi]))
210  {
211  field.boundaryFieldRef().set
212  (
213  patchi,
215  (
217  mesh.boundary()[patchi],
218  field
219  )
220  );
221  }
222  else
223  {
224  typename GField::Patch& pf = field.boundaryFieldRef()[patchi];
225  pf.map(pf, SetSizePatchFieldMapper(pf.patch().size()));
226  }
227 
228  field.boundaryFieldRef()[patchi] == NaN;
229  }
230 
231  field.instance() = field.time().name();
232  }
233 }
234 
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 } // End namespace Foam
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 #endif
243 
244 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static const char *const typeName
Definition: Field.H:105
Generic GeometricField class.
void reset(const tmp< GeometricField< Type, PatchField, GeoMesh >> &)
Reset the field contents to the given field.
void clearOldTimes()
Delete old time and previous iteration fields.
An STL-conforming iterator.
Definition: HashTable.H:429
An STL-conforming hash table.
Definition: HashTable.H:127
const Time & time() const
Return time.
Definition: IOobject.C:318
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:355
const word & name() const
Return name.
Definition: IOobject.H:310
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:294
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
const word & name() const
Return const reference to name.
tmp< VolField< Type > > srcToTgt(const VolField< Type > &srcFld) const
Interpolate a source vol field to the target with no left.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
label patchi
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
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static Type NaN()
Return a primitive with all components set to NaN.
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
label & setComponent(label &l, const direction)
Definition: label.H:86
messageStream Info
void MeshToMeshMapVolFields(const fvMesh &mesh, const fvMeshToFvMesh &mapper)
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:58
void MeshToMeshMapVolInternalFields(const fvMesh &mesh, const fvMeshToFvMesh &mapper)
uint8_t direction
Definition: direction.H:45
void NaNGeometricFields(const fvMesh &mesh, const fvMeshToFvMesh &mapper)