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  forAll(fields, i)
59  {
60  VolField<Type>& field = fields[i];
61 
62  // Delete old time fields
63  field.clearOldTimes();
64 
65  if (fvMeshToFvMesh::debug)
66  {
67  Info<< "Mapping " << field.typeName << ' ' << field.name()
68  << endl;
69  }
70 
71  field.reset(mapper.srcToTgt(field));
72 
73  field.instance() = field.time().name();
74  }
75 }
76 
77 
78 template<class Type>
80 (
81  const fvMesh& mesh,
82  const fvMeshToFvMesh& mapper
83 )
84 {
86  (
87  mesh.fields<VolInternalField<Type>>(true)
88  );
89 
90  forAll(fields, i)
91  {
92  VolInternalField<Type>& field = fields[i];
93 
94  if (fvMeshToFvMesh::debug)
95  {
96  Info<< "Mapping " << field.typeName << ' ' << field.name()
97  << endl;
98  }
99 
100  field.reset(mapper.srcToTgt<Type>(field));
101 
102  field.instance() = field.time().name();
103  }
104 }
105 
106 
107 template
108 <
109  class Type,
110  template<class> class PatchField,
111  class GeoMesh,
112  class SetSizePatchFieldMapper
113 >
115 (
116  const fvMesh& mesh,
117  const fvMeshToFvMesh& mapper
118 )
119 {
121 
122  Type NaN;
123 
124  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
125  {
126  setComponent(NaN, cmpt) = std::numeric_limits<scalar>::signaling_NaN();
127  }
128 
129  UPtrList<GField> fields(mesh.curFields<GField>());
130 
131  forAll(fields, i)
132  {
133  GField& field = fields[i];
134 
135  // Delete old time fields
136  field.clearOldTimes();
137 
138  if (fvMeshToFvMesh::debug)
139  {
140  Info<< "Setting to NaN " << field.typeName << ' ' << field.name()
141  << endl;
142  }
143 
144  const typename GField::Mesh& mesh = field.mesh();
145 
146  field.primitiveFieldRef().setSize(GeoMesh::size(mesh));
147  field.primitiveFieldRef() = NaN;
148 
149  field.boundaryFieldRef().setSize(mesh.boundary().size());
150 
151  forAll(mesh.boundary(), patchi)
152  {
153  if (isA<processorPolyPatch>(mesh().boundaryMesh()[patchi]))
154  {
155  field.boundaryFieldRef().set
156  (
157  patchi,
159  (
161  mesh.boundary()[patchi],
162  field
163  )
164  );
165  }
166  else
167  {
168  typename GField::Patch& pf = field.boundaryFieldRef()[patchi];
169  pf.map(pf, SetSizePatchFieldMapper(pf.patch().size()));
170  }
171 
172  field.boundaryFieldRef()[patchi] == NaN;
173  }
174 
175  field.instance() = field.time().name();
176  }
177 }
178 
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 } // End namespace Foam
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 #endif
187 
188 // ************************************************************************* //
#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.
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
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
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
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.
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)