MapVolFields.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2018-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
28
29#ifndef MapConsistentVolFields_H
30#define MapConsistentVolFields_H
31
32#include "GeometricField.H"
33#include "meshToMesh.H"
34#include "IOobjectList.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40
41template<class Type>
43{
44 auto& fldBf = fld.boundaryFieldRef();
45
47
48 if
49 (
52 )
53 {
54 const label startOfRequests = UPstream::nRequests();
55
56 forAll(fldBf, patchi)
57 {
58 fvPatchField<Type>& tgtField = fldBf[patchi];
59
60 if
61 (
62 tgtField.type() == tgtField.patch().patch().type()
63 && polyPatch::constraintType(tgtField.patch().patch().type())
64 )
65 {
66 tgtField.initEvaluate(commsType);
67 }
68 }
69
70 // Wait for outstanding requests
71 if
72 (
75 )
76 {
77 UPstream::waitRequests(startOfRequests);
78 }
79
80 forAll(fldBf, patchi)
81 {
82 fvPatchField<Type>& tgtField = fldBf[patchi];
83
84 if
85 (
86 tgtField.type() == tgtField.patch().patch().type()
87 && polyPatch::constraintType(tgtField.patch().patch().type())
88 )
89 {
90 tgtField.evaluate(commsType);
91 }
92 }
93 }
94 else if (commsType == UPstream::commsTypes::scheduled)
95 {
96 const lduSchedule& patchSchedule =
97 fld.mesh().globalData().patchSchedule();
98
99 for (const auto& schedEval : patchSchedule)
100 {
101 const label patchi = schedEval.patch;
102
103 fvPatchField<Type>& tgtField = fldBf[patchi];
104
105 if
106 (
107 tgtField.type() == tgtField.patch().patch().type()
108 && polyPatch::constraintType(tgtField.patch().patch().type())
109 )
110 {
111 if (schedEval.init)
112 {
113 tgtField.initEvaluate(commsType);
114 }
115 else
116 {
117 tgtField.evaluate(commsType);
118 }
119 }
120 }
121 }
122}
123
124
125template<class Type, class CombineOp>
127(
128 const IOobjectList& objects,
129 const wordRes& selectedFields,
130 const meshToMesh& interp,
131 const CombineOp& cop
132)
133{
135
136 const fvMesh& meshSource = static_cast<const fvMesh&>(interp.srcRegion());
137 const fvMesh& meshTarget = static_cast<const fvMesh&>(interp.tgtRegion());
138
139 // Available fields, sorted order
140 const wordList fieldNames =
141 (
142 selectedFields.empty()
143 ? objects.sortedNames<fieldType>()
144 : objects.sortedNames<fieldType>(selectedFields)
145 );
146
147 for (const word& fieldName : fieldNames)
148 {
149 const fieldType fieldSource(*(objects[fieldName]), meshSource, false);
150
151 IOobject targetIO
152 (
153 fieldName,
154 meshTarget.time().timeName(),
155 meshTarget,
157 );
158
159 if (targetIO.typeHeaderOk<fieldType>(true))
160 {
161 Info<< " interpolating onto existing field "
162 << fieldName << endl;
163
164 fieldType fieldTarget(targetIO, meshTarget, false);
165
166 interp.mapSrcToTgt(fieldSource, cop, fieldTarget);
167
168 evaluateConstraintTypes(fieldTarget);
169
170 fieldTarget.write();
171 }
172 else
173 {
174 Info<< " creating new field "
175 << fieldName << endl;
176
177 targetIO.readOpt(IOobject::NO_READ);
178
179 tmp<fieldType> tfieldTarget
180 (
181 interp.mapSrcToTgt(fieldSource, cop)
182 );
183
184 fieldType fieldTarget(targetIO, tfieldTarget);
185
186 evaluateConstraintTypes(fieldTarget);
187
188 fieldTarget.write();
189 }
190 }
191}
192
193
194// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195
196} // End namespace Foam
197
198// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199
200#endif
201
202// ************************************************************************* //
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic GeometricField class.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
wordList sortedNames() const
The sorted names of the IOobjects.
Definition: IOobjectList.C:383
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
readOption readOpt() const noexcept
The read option.
Definition: IOobjectI.H:164
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
commsTypes
Types of communications.
Definition: UPstream.H:67
@ nonBlocking
"nonBlocking"
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:90
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:100
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:281
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
virtual const word & constraintType() const
Return the constraint type this pointPatch implements.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
virtual void initEvaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Initialise the evaluation of the patch field.
Definition: fvPatchField.H:461
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:345
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:362
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:167
Class to calculate the cell-addressing between two overlapping meshes.
Definition: meshToMesh.H:65
const polyMesh & srcRegion() const
Return const access to the source mesh.
Definition: meshToMeshI.H:33
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:39
A class for managing temporary objects.
Definition: tmp.H:65
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
Namespace for OpenFOAM.
void MapVolFields(const IOobjectList &objects, const meshToMesh0 &meshToMesh0Interp, const meshToMesh0::order &mapOrder, const CombineOp &cop)
Definition: MapVolFields.H:42
messageStream Info
Information stream (stdout output on master, null elsewhere)
void evaluateConstraintTypes(GeometricField< Type, fvPatchField, volMesh > &fld)
Definition: MapVolFields.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333