Loading...
Searching...
No Matches
regionModelTemplates.C
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-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28template<class Type>
31(
32 const regionModel& nbrRegion,
33 const label regionPatchi,
34 const label nbrPatchi,
35 const Field<Type>& nbrField,
36 const bool flip
37) const
38{
39 const int oldTag = UPstream::incrMsgType();
40
42 interRegionAMI(nbrRegion, regionPatchi, nbrPatchi, flip);
43
44 tmp<Field<Type>> tresult(ami.interpolateToSource(nbrField));
45
46 UPstream::msgType(oldTag); // Restore tag
48 return tresult;
49}
50
51
52template<class Type>
55(
56 const regionModel& nbrRegion,
57 const word& fieldName,
58 const label regionPatchi,
59 const bool flip
60) const
61{
63
64 const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
65
66 if (nbrRegionMesh.foundObject<fieldType>(fieldName))
67 {
68 const label nbrPatchi = nbrCoupledPatchID(nbrRegion, regionPatchi);
69
70 const int oldTag = UPstream::incrMsgType();
71
73 interRegionAMI(nbrRegion, regionPatchi, nbrPatchi, flip);
74
75 const fieldType& nbrField =
76 nbrRegionMesh.lookupObject<fieldType>(fieldName);
77
78 const Field<Type>& nbrFieldp = nbrField.boundaryField()[nbrPatchi];
79
80 tmp<Field<Type>> tresult(ami.interpolateToSource(nbrFieldp));
81
82 UPstream::msgType(oldTag); // Restore tag
83
84 return tresult;
85 }
86 else
87 {
88 const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
89
90 return tmp<Field<Type>>::New(p.size(), Zero);
91 }
92}
93
94
95template<class Type>
98(
99 const regionModel& nbrRegion,
100 const word& fieldName,
101 const label regionPatchi,
102 const bool flip
103) const
104{
106
107 const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
108
109 if (nbrRegionMesh.foundObject<fieldType>(fieldName))
110 {
111 const label nbrPatchi = nbrCoupledPatchID(nbrRegion, regionPatchi);
112
113 const int oldTag = UPstream::incrMsgType();
114
116 interRegionAMI(nbrRegion, regionPatchi, nbrPatchi, flip);
117
118 const fieldType& nbrField =
119 nbrRegionMesh.lookupObject<fieldType>(fieldName);
120
121 const fvPatchField<Type>& nbrFieldp =
122 nbrField.boundaryField()[nbrPatchi];
123
124 tmp<Field<Type>> tresult
125 (
126 ami.interpolateToSource(nbrFieldp.patchInternalField())
127 );
128
129 UPstream::msgType(oldTag); // Restore tag
130
131 return tresult;
132 }
133 else
134 {
135 const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
137 return tmp<Field<Type>>::New(p.size(), Zero);
138 }
139}
140
141
142template<class Type>
144(
145 const label regionPatchi,
146 List<Type>& regionField
147) const
148{
149 forAll(intCoupledPatchIDs_, i)
150 {
151 if (intCoupledPatchIDs_[i] == regionPatchi)
152 {
153 const mappedPatchBase& mpb =
155 (
156 regionMesh().boundaryMesh()[regionPatchi]
157 );
158 mpb.reverseDistribute(regionField);
159 return;
160 }
161 }
162
164 << "Region patch ID " << regionPatchi << " not found in region mesh"
165 << abort(FatalError);
166}
167
168
169template<class Type>
171(
172 const label regionPatchi,
173 List<Type>& primaryField
174) const
175{
176 forAll(intCoupledPatchIDs_, i)
177 {
178 if (intCoupledPatchIDs_[i] == regionPatchi)
179 {
180 const mappedPatchBase& mpb =
182 (
183 regionMesh().boundaryMesh()[regionPatchi]
184 );
185 mpb.distribute(primaryField);
186 return;
187 }
188 }
189
191 << "Region patch ID " << regionPatchi << " not found in region mesh"
192 << abort(FatalError);
193}
194
195
196template<class Type, class CombineOp>
198(
199 const label regionPatchi,
200 List<Type>& regionField,
201 const CombineOp& cop
202) const
203{
204 forAll(intCoupledPatchIDs_, i)
205 {
206 if (intCoupledPatchIDs_[i] == regionPatchi)
207 {
208 const mappedPatchBase& mpb =
210 (
211 regionMesh().boundaryMesh()[regionPatchi]
212 );
213 mpb.reverseDistribute(regionField, cop);
214 return;
215 }
216 }
217
219 << "Region patch ID " << regionPatchi << " not found in region mesh"
220 << abort(FatalError);
221}
222
223
224template<class Type, class CombineOp>
226(
227 const label regionPatchi,
228 List<Type>& primaryField,
229 const CombineOp& cop
230) const
231{
232 forAll(intCoupledPatchIDs_, i)
233 {
234 if (intCoupledPatchIDs_[i] == regionPatchi)
235 {
236 const mappedPatchBase& mpb =
238 (
239 regionMesh().boundaryMesh()[regionPatchi]
240 );
241 mpb.distribute(primaryField, cop);
242 return;
243 }
244 }
245
247 << "Region patch ID " << regionPatchi << " not found in region mesh"
248 << abort(FatalError);
249}
250
251
252// ************************************************************************* //
void interpolateToSource(const UList< Type > &fld, const CombineOp &cop, List< Type > &result, const UList< Type > &defaultValues=UList< Type >::null()) const
Interpolate from target to source with supplied op to combine existing value with remote value and we...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
void reverseDistribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
bool foundObject(const word &name, const bool recursive=false) const
Contains the named Type?
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Base class for region models.
Definition regionModel.H:59
void toRegion(const label regionPatchi, List< Type > &primaryFieldField) const
Convert a primary region field to the local region.
label nbrCoupledPatchID(const regionModel &nbrRegion, const label regionPatchi) const
Return the coupled patch ID paired with coupled patch.
virtual const AMIPatchToPatchInterpolation & interRegionAMI(const regionModel &nbrRegion, const label regionPatchi, const label nbrPatchi, const bool flip) const
Create or return a new inter-region AMI object.
const fvMesh & regionMesh() const
Return the region mesh database.
tmp< Field< Type > > mapRegionPatchInternalField(const regionModel &nbrRegion, const word &fieldName, const label regionPatchi, const bool flip=false) const
Map patch internal field from another region model to local.
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
tmp< Foam::Field< Type > > mapRegionPatchField(const regionModel &nbrRegion, const label regionPatchi, const label nbrPatchi, const Field< Type > &nbrField, const bool flip=false) const
Map patch field from another region model to local patch.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
AMIInterpolation AMIPatchToPatchInterpolation
Patch-to-patch interpolation == Foam::AMIInterpolation.
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299