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 "fieldMapper.H"
40 #include "setSizeFieldMapper.H"
41 #include "processorPolyPatch.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 template<class Type>
50 (
51  const fvMesh& mesh,
52  const fvMeshToFvMesh& mapper
53 )
54 {
56 
57  forAll(fields, i)
58  {
59  VolField<Type>& field = fields[i];
60 
61  // Delete old time fields
62  field.clearOldTimes();
63 
64  if (fvMeshToFvMesh::debug)
65  {
66  Info<< "Mapping " << field.typeName << ' ' << field.name()
67  << endl;
68  }
69 
70  field.reset(mapper.srcToTgt(field));
71 
72  field.instance() = field.time().name();
73  }
74 }
75 
76 
77 template<class Type>
79 (
80  const fvMesh& mesh,
81  const fvMeshToFvMesh& mapper
82 )
83 {
85  (
86  mesh.fields<VolInternalField<Type>>(true)
87  );
88 
89  forAll(fields, i)
90  {
91  VolInternalField<Type>& field = fields[i];
92 
93  if (fvMeshToFvMesh::debug)
94  {
95  Info<< "Mapping " << field.typeName << ' ' << field.name()
96  << endl;
97  }
98 
99  field.reset(mapper.srcToTgt<Type>(field));
100 
101  field.instance() = field.time().name();
102  }
103 }
104 
105 
106 template<class Type, template<class> class PatchField, class GeoMesh>
108 (
109  const fvMesh& mesh,
110  const fvMeshToFvMesh& mapper
111 )
112 {
114 
115  UPtrList<GField> fields(mesh.curFields<GField>());
116 
117  forAll(fields, i)
118  {
119  GField& field = fields[i];
120 
121  // Delete old time fields
122  field.clearOldTimes();
123 
124  if (fvMeshToFvMesh::debug)
125  {
126  Info<< "Setting to NaN " << field.typeName << ' ' << field.name()
127  << endl;
128  }
129 
130  const typename GField::Mesh& mesh = field.mesh();
131 
132  field.primitiveFieldRef().setSize(GeoMesh::size(mesh));
133  field.primitiveFieldRef() = pTraits<Type>::nan;
134 
135  field.boundaryFieldRef().setSize(mesh.boundary().size());
136 
137  forAll(mesh.boundary(), patchi)
138  {
139  if (isA<processorPolyPatch>(mesh().boundaryMesh()[patchi]))
140  {
141  field.boundaryFieldRef().set
142  (
143  patchi,
145  (
147  mesh.boundary()[patchi],
148  field
149  )
150  );
151  }
152  else
153  {
154  typename GField::Patch& pf = field.boundaryFieldRef()[patchi];
155  pf.map(pf, setSizeFieldMapper(pf.patch().size()));
156  }
157 
158  field.boundaryFieldRef()[patchi] == pTraits<Type>::nan;
159  }
160 
161  field.instance() = field.time().name();
162  }
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 } // End namespace Foam
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 #endif
173 
174 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static const char *const typeName
Definition: Field.H:106
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
Generic GeometricField class.
void reset(const tmp< GeometricField< Type, PatchField, GeoMesh >> &)
Reset the field contents to the given field.
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
void clearOldTimes()
Clear old time fields.
Definition: OldTimeField.C:271
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:99
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.
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:441
Traits class for primitives.
Definition: pTraits.H:53
Mapper which sets the field size. It does not actually map values.
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:228
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
void MeshToMeshMapVolFields(const fvMesh &mesh, const fvMeshToFvMesh &mapper)
void NaNGeometricFields(const fvMesh &mesh, const fvMeshToFvMesh &mapper)
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:61
void MeshToMeshMapVolInternalFields(const fvMesh &mesh, const fvMeshToFvMesh &mapper)