MapDimensionedFields.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) 2012-2018 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 Description
25  Generic internal field mapper for dimensioned fields. For "real" mapping,
26  add template specialisations for mapping of internal fields depending on
27  mesh type.
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #ifndef MapDimensionedFields_H
32 #define MapDimensionedFields_H
33 
34 #include "polyMesh.H"
35 #include "MapFvVolField.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 
42 template<class Type, class MeshMapper, class GeoMesh>
43 void MapDimensionedFields(const MeshMapper& mapper)
44 {
45  typedef DimensionedField<Type, GeoMesh> FieldType;
46  typedef HashTable<const FieldType*> TableType;
47 
48  TableType fields(mapper.thisDb().template lookupClass<FieldType>(true));
49 
50  forAllConstIter(typename TableType, fields, fieldIter)
51  {
52  FieldType& field = const_cast<FieldType&>(*fieldIter());
53 
54  if (&field.mesh() == &mapper.mesh())
55  {
56  if (polyMesh::debug)
57  {
58  Info<< "Mapping " << field.typeName << ' ' << field.name()
59  << endl;
60  }
61 
63 
64  field.instance() = field.time().timeName();
65  }
66  else if (polyMesh::debug)
67  {
68  Info<< "Not mapping " << field.typeName << ' ' << field.name()
69  << " since originating mesh differs from that of mapper."
70  << endl;
71  }
72  }
73 }
74 
75 
76 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77 
78 } // End namespace Foam
79 
80 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81 
82 #endif
83 
84 // ************************************************************************* //
void MapDimensionedFields(const MeshMapper &mapper)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Info<< "Calculating turbulent flame speed field St\"<< endl;volScalarField St(IOobject("St", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:228
An STL-conforming hash table.
Definition: HashTable.H:61
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
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...
messageStream Info
rDeltaTY field()
Map volume internal field on topology change. This is a partial template specialisation, see MapGeometricFields.
Namespace for OpenFOAM.