MapDimensionedFields.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) 2012 OpenFOAM Foundation
9  Copyright (C) 2019-2023 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 Description
28  Generic internal field mapper for dimensioned fields. For "real" mapping,
29  add template specialisations for mapping of internal fields depending on
30  mesh type.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef MapDimensionedFields_H
35 #define MapDimensionedFields_H
36 
37 #include "polyMesh.H"
38 #include "MapFvVolField.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 template<class Type, class MeshMapper, class GeoMesh>
46 void MapDimensionedFields(const MeshMapper& mapper)
47 {
48  typedef DimensionedField<Type, GeoMesh> FieldType;
49 
50  // Note: use strict=true to avoid picking up volFields
52  (
53  mapper.thisDb().template csorted<FieldType, true>()
54  );
55 
56  for (const auto& fld : fields)
57  {
58  FieldType& field = const_cast<FieldType&>(fld);
59 
60  if (&field.mesh() == &mapper.mesh())
61  {
62  if (polyMesh::debug)
63  {
64  Info<< "Mapping "
65  << FieldType::typeName << ' ' << field.name() << endl;
66  }
67 
69 
70  field.instance() = field.time().timeName();
71  }
72  else if (polyMesh::debug)
73  {
74  Info<< "Not mapping "
75  << FieldType::typeName << ' ' << field.name()
76  << " since originating mesh differs from that of mapper."
77  << endl;
78  }
79  }
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85 } // End namespace Foam
86 
87 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
88 
89 #endif
90 
91 // ************************************************************************* //
rDeltaTY field()
void MapDimensionedFields(const MeshMapper &mapper)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
int debug
Static debugging option.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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)
Map volume internal field on topology change. This is a partial template specialisation, see MapGeometricFields.
Namespace for OpenFOAM.