MapVolFields.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) 2011-2022 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 #ifndef MapVolFields_H
27 #define MapVolFields_H
28 
29 #include "GeometricField.H"
30 #include "meshToMesh.H"
31 #include "IOobjectList.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 template<class Type>
39 void evaluateConstraintTypes(GeometricField<Type, fvPatchField, volMesh>& fld)
40 {
41  typename GeometricField<Type, fvPatchField, volMesh>::
42  Boundary& fldBf = fld.boundaryFieldRef();
43 
44  if
45  (
48  )
49  {
50  label nReq = Pstream::nRequests();
51 
52  forAll(fldBf, patchi)
53  {
54  fvPatchField<Type>& tgtField = fldBf[patchi];
55 
56  if
57  (
58  tgtField.type() == tgtField.patch().patch().type()
59  && polyPatch::constraintType(tgtField.patch().patch().type())
60  )
61  {
62  tgtField.initEvaluate(Pstream::defaultCommsType);
63  }
64  }
65 
66  // Block for any outstanding requests
67  if
68  (
71  )
72  {
74  }
75 
76  forAll(fldBf, patchi)
77  {
78  fvPatchField<Type>& tgtField = fldBf[patchi];
79 
80  if
81  (
82  tgtField.type() == tgtField.patch().patch().type()
83  && polyPatch::constraintType(tgtField.patch().patch().type())
84  )
85  {
86  tgtField.evaluate(Pstream::defaultCommsType);
87  }
88  }
89  }
91  {
92  const lduSchedule& patchSchedule =
93  fld.mesh().globalData().patchSchedule();
94 
95  forAll(patchSchedule, patchEvali)
96  {
97  label patchi = patchSchedule[patchEvali].patch;
98  fvPatchField<Type>& tgtField = fldBf[patchi];
99 
100  if
101  (
102  tgtField.type() == tgtField.patch().patch().type()
103  && polyPatch::constraintType(tgtField.patch().patch().type())
104  )
105  {
106  if (patchSchedule[patchEvali].init)
107  {
108  tgtField.initEvaluate(Pstream::commsTypes::scheduled);
109  }
110  else
111  {
112  tgtField.evaluate(Pstream::commsTypes::scheduled);
113  }
114  }
115  }
116  }
117 }
118 
119 
120 template<class Type>
121 void MapVolFields
122 (
123  const IOobjectList& objects,
124  const HashSet<word>& selectedFields,
125  const meshToMesh& interp
126 )
127 {
129 
130  const fvMesh& meshSource = static_cast<const fvMesh&>(interp.srcRegion());
131  const fvMesh& meshTarget = static_cast<const fvMesh&>(interp.tgtRegion());
132 
133  IOobjectList fields = objects.lookupClass(fieldType::typeName);
134 
135  forAllIter(IOobjectList, fields, fieldIter)
136  {
137  const word& fieldName = fieldIter()->name();
138 
139  if (selectedFields.empty() || selectedFields.found(fieldName))
140  {
141  const fieldType fieldSource(*fieldIter(), meshSource);
142 
143  typeIOobject<fieldType> targetIO
144  (
145  fieldName,
146  meshTarget.time().timeName(),
147  meshTarget,
149  );
150 
151  if (targetIO.headerOk())
152  {
153  Info<< " interpolating onto existing field "
154  << fieldName << endl;
155  fieldType fieldTarget(targetIO, meshTarget);
156 
157  interp.mapSrcToTgt(fieldSource, fieldTarget);
158 
159  evaluateConstraintTypes(fieldTarget);
160 
161  fieldTarget.write();
162  }
163  else
164  {
165  Info<< " creating new field "
166  << fieldName << endl;
167 
168  targetIO.readOpt() = IOobject::NO_READ;
169 
170  tmp<fieldType> tfieldTarget(interp.mapSrcToTgt(fieldSource));
171 
172  fieldType fieldTarget(targetIO, tfieldTarget);
173 
174  evaluateConstraintTypes(fieldTarget);
175 
176  fieldTarget.write();
177  }
178  }
179  }
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 } // End namespace Foam
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 #endif
190 
191 // ************************************************************************* //
void MapVolFields(const IOobjectList &objects, const meshToMesh0 &meshToMesh0Interp, const meshToMesh0::order &mapOrder)
Definition: MapVolFields.H:40
A HashTable with keys but without contents.
Definition: HashSet.H:59
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
const polyMesh & srcRegion() const
Return const access to the source mesh.
Definition: meshToMeshI.H:30
Templated form of IOobject providing type information for file reading and header type checking...
Definition: IOobject.H:537
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
bool empty() const
Return true if the hash table is empty.
Definition: HashTableI.H:72
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:137
List< lduScheduleEntry > lduSchedule
Definition: lduSchedule.H:74
Generic GeometricField class.
Class to calculate the cell-addressing between two overlapping meshes.
Definition: meshToMesh.H:60
void mapSrcToTgt(const UList< Type > &srcFld, List< Type > &result) const
Map field from src to tgt mesh with defined operation.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
void evaluateConstraintTypes(GeometricField< Type, fvPatchField, volMesh > &fld)
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:265
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
A class for handling words, derived from string.
Definition: word.H:59
objects
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:147
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:190
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, p_rgh, pimple.dict(), thermo.incompressible());mesh.schemes().setFluxRequired(p_rgh.name());hydrostaticInitialisation(p_rgh, p, rho, U, gh, ghf, pRef, thermo, pimple.dict());Info<< "Creating field dpdt\"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar(p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\"<< endl;volScalarField K("K", 0.5 *magSqr(U));dimensionedScalar initialMass=fvc::domainIntegrate(rho);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:131
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
messageStream Info
A class for managing temporary objects.
Definition: PtrList.H:53
const polyMesh & tgtRegion() const
Return const access to the target mesh.
Definition: meshToMeshI.H:36
Namespace for OpenFOAM.