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-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 \*---------------------------------------------------------------------------*/
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>
40 {
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  {
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  {
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  {
109  }
110  else
111  {
113  }
114  }
115  }
116  }
117 }
118 
119 
120 template<class Type, class CombineOp>
121 void MapVolFields
122 (
123  const IOobjectList& objects,
124  const HashSet<word>& selectedFields,
125  const meshToMesh& interp,
126  const CombineOp& cop
127 )
128 {
130 
131  const fvMesh& meshSource = static_cast<const fvMesh&>(interp.srcRegion());
132  const fvMesh& meshTarget = static_cast<const fvMesh&>(interp.tgtRegion());
133 
134  IOobjectList fields = objects.lookupClass(fieldType::typeName);
135 
136  forAllIter(IOobjectList, fields, fieldIter)
137  {
138  const word& fieldName = fieldIter()->name();
139 
140  if (selectedFields.empty() || selectedFields.found(fieldName))
141  {
142  const fieldType fieldSource(*fieldIter(), meshSource);
143 
144  IOobject targetIO
145  (
146  fieldName,
147  meshTarget.time().timeName(),
148  meshTarget,
150  );
151 
152  if (targetIO.typeHeaderOk<fieldType>(true))
153  {
154  Info<< " interpolating onto existing field "
155  << fieldName << endl;
156  fieldType fieldTarget(targetIO, meshTarget);
157 
158  interp.mapSrcToTgt(fieldSource, cop, fieldTarget);
159 
160  evaluateConstraintTypes(fieldTarget);
161 
162  fieldTarget.write();
163  }
164  else
165  {
166  Info<< " creating new field "
167  << fieldName << endl;
168 
169  targetIO.readOpt() = IOobject::NO_READ;
170 
172  tfieldTarget(interp.mapSrcToTgt(fieldSource, cop));
173 
174  fieldType fieldTarget(targetIO, tfieldTarget);
175 
176  evaluateConstraintTypes(fieldTarget);
177 
178  fieldTarget.write();
179  }
180  }
181  }
182 }
183 
184 
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186 
187 } // End namespace Foam
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 #endif
192 
193 // ************************************************************************* //
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
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
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
#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
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void initEvaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Initialise the evaluation of the patch field.
Definition: fvPatchField.H:419
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:137
Generic GeometricField class.
Class to calculate the cell-addressing between two overlapping meshes.
Definition: meshToMesh.H:60
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
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
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
void evaluateConstraintTypes(GeometricField< Type, fvPatchField, volMesh > &fld)
Definition: MapVolFields.H:39
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:259
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
A class for handling words, derived from string.
Definition: word.H:59
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:138
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:337
const Mesh & mesh() const
Return mesh.
objects
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:147
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:229
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:192
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
messageStream Info
void MapVolFields(const IOobjectList &objects, const meshToMesh0 &meshToMesh0Interp, const meshToMesh0::order &mapOrder, const CombineOp &cop)
Definition: MapVolFields.H:40
A class for managing temporary objects.
Definition: PtrList.H:53
void mapSrcToTgt(const UList< Type > &srcFld, const CombineOp &cop, List< Type > &result) const
Map field from src to tgt mesh with defined operation.
const polyMesh & tgtRegion() const
Return const access to the target mesh.
Definition: meshToMeshI.H:36
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.