mapGeometricFields.C
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) 2011-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 \*---------------------------------------------------------------------------*/
25 
26 #include "mapGeometricFields.H"
27 #include "fvMeshToFvMesh.H"
28 #include "surfaceMesh.H"
29 #include "pointMesh.H"
30 #include "IOobjectList.H"
31 #include "OSspecific.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 template<class Type>
39 void mapVolTypeFields
40 (
41  const fvMeshToFvMesh& interp,
42  const wordReList& cuttingPatches,
43  const HashSet<word>& selectedFields,
44  const IOobjectList& objects
45 )
46 {
47  const fvMesh& srcMesh = static_cast<const fvMesh&>(interp.srcMesh());
48  const fvMesh& tgtMesh = static_cast<const fvMesh&>(interp.tgtMesh());
49 
50  IOobjectList fields = objects.lookupClass(VolField<Type>::typeName);
51 
52  forAllIter(IOobjectList, fields, fieldIter)
53  {
54  const word& fieldName = fieldIter()->name();
55 
56  if (!selectedFields.empty() && !selectedFields.found(fieldName))
57  {
58  continue;
59  }
60 
61  const VolField<Type> fieldSource(*fieldIter(), srcMesh);
62 
63  typeIOobject<VolField<Type>> targetIO
64  (
65  fieldName,
66  tgtMesh.time().name(),
67  tgtMesh,
69  );
70 
71  // Warnings about inconsistent execution
72  if (targetIO.headerOk() && interp.consistent())
73  {
75  << "Mapping of field " << fieldName << " will not utilise "
76  << "the corresponding field in the target case, as the map is "
77  << "consistent (i.e., all patches are mapped)" << endl;
78  }
79  if (!targetIO.headerOk() && !interp.consistent())
80  {
82  << "Cannot map field " << fieldName << " because the "
83  << "map is not consistent (i.e., not all patches are "
84  << "mapped), and there is no corresponding field in "
85  << "the target case" << endl;
86  continue;
87  }
88  if (!targetIO.headerOk() && !cuttingPatches.empty())
89  {
91  << "Cutting patches will not be used for field " << fieldName
92  << " because no there is no corresponding field in the target "
93  << "case" << endl;
94  }
95 
96  if (targetIO.headerOk())
97  {
98  Info<< " mapping into existing field " << fieldName << endl;
99 
100  VolField<Type> fieldTarget(targetIO, tgtMesh);
101 
102  fieldTarget.reset
103  (
104  interp.srcToTgt(fieldSource, fieldTarget, cuttingPatches)
105  );
106 
107  fieldTarget.write();
108  }
109  else
110  {
111  Info<< " creating new field " << fieldName << endl;
112 
113  VolField<Type>(targetIO, interp.srcToTgt(fieldSource)).write();
114  }
115  }
116 }
117 
118 }
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
123 (
124  const fvMeshToFvMesh& interp,
125  const wordReList& cuttingPatches,
126  const HashSet<word>& selectedFields,
127  const bool noLagrangian
128 )
129 {
130  // Search for list of source objects for this time
131  const polyMesh& srcMesh = interp.srcMesh();
132  IOobjectList objects(srcMesh, srcMesh.time().name());
133 
134  // Map the fields
135  #define MapVolTypeFields(Type, nullArg) \
136  mapVolTypeFields<Type> \
137  ( \
138  interp, \
139  cuttingPatches, \
140  selectedFields, \
141  objects \
142  );
143  FOR_ALL_FIELD_TYPES(MapVolTypeFields);
144  #undef MapVolTypeFields
145 }
146 
147 
148 // ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
static const char *const typeName
Definition: Field.H:105
virtual Ostream & write(const char)=0
Write character.
#define FOR_ALL_FIELD_TYPES(Macro,...)
Definition: fieldTypes.H:51
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
Maps geometric fields.
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
void mapGeometricFields(const fvMeshToFvMesh &interp, const wordReList &cuttingPatches, const HashSet< word > &selectedFields, const bool noLagrangian)
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
objects