MapGeometricFields.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::MapInternalField
29 
30 Description
31  Generic internal field mapper. For "real" mapping, add template
32  specialisations for mapping of internal fields depending on mesh
33  type.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef MapGeometricFields_H
38 #define MapGeometricFields_H
39 
40 #include "polyMesh.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 template<class Type, class MeshMapper, class GeoMesh>
48 class MapInternalField
49 {
50 public:
51 
53  {}
54 
55  void operator()
56  (
58  const MeshMapper& mapper
59  ) const;
60 };
61 
62 
63 //- Generic Geometric field mapper.
64 // For "real" mapping, add template specialisations
65 // for mapping of internal fields depending on mesh type.
66 template
67 <
68  class Type,
69  template<class> class PatchField,
70  class MeshMapper,
71  class GeoMesh
72 >
74 (
75  const MeshMapper& mapper
76 )
77 {
79 
81  (
82  mapper.thisDb().objectRegistry::template
83  lookupClass<FieldType>()
84  );
85 
86  // It is necessary to enforce that all old-time fields are stored
87  // before the mapping is performed. Otherwise, if the
88  // old-time-level field is mapped before the field itself, sizes
89  // will not match.
90 
91  forAllConstIters(fields, fieldIter)
92  {
93  FieldType& field = const_cast<FieldType&>(*fieldIter());
94 
95  //Note: check can be removed once pointFields are actually stored on
96  // the pointMesh instead of now on the polyMesh!
97  if (&field.mesh() == &mapper.mesh())
98  {
99  field.storeOldTimes();
100  }
101  }
102 
103  forAllConstIters(fields, fieldIter)
104  {
105  FieldType& field = const_cast<FieldType&>(*fieldIter());
106 
107  if (&field.mesh() == &mapper.mesh())
108  {
109  if (polyMesh::debug)
110  {
111  Info<< "Mapping " << field.typeName << ' ' << field.name()
112  << endl;
113  }
114 
115  // Map the internal field
117  (
118  field.internalFieldRef(),
119  mapper
120  );
121 
122  // Map the patch fields
123  auto& bfield = field.boundaryFieldRef();
124 
125  forAll(bfield, patchi)
126  {
127  // Cannot check sizes for patch fields because of
128  // empty fields in FV and because point fields get their size
129  // from the patch which has already been resized
130  //
131 
132  bfield[patchi].autoMap(mapper.boundaryMap()[patchi]);
133  }
134 
135  field.instance() = field.time().timeName();
136  }
137  else if (polyMesh::debug)
138  {
139  Info<< "Not mapping " << field.typeName << ' ' << field.name()
140  << " since originating mesh differs from that of mapper."
141  << endl;
142  }
143  }
144 }
145 
146 
147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 
149 } // End namespace Foam
150 
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153 #endif
154 
155 // ************************************************************************* //
rDeltaTY field()
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
void MapGeometricFields(const MeshMapper &mapper)
Generic Geometric field mapper.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
A HashTable similar to std::unordered_map.
Definition: HashTable.H:102
int debug
Static debugging option.
Generic internal field mapper. For "real" mapping, add template specialisations for mapping of intern...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
messageStream Info
Information stream (stdout output on master, null elsewhere)
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:42
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28