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 -------------------------------------------------------------------------------
10 License
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 
28 template<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  int oldTag = UPstream::msgType();
40  UPstream::msgType() = oldTag + 1;
41 
42  const AMIPatchToPatchInterpolation& ami =
43  interRegionAMI(nbrRegion, regionPatchi, nbrPatchi, flip);
44 
45  tmp<Field<Type>> tresult(ami.interpolateToSource(nbrField));
46 
47  UPstream::msgType() = oldTag;
48 
49  return tresult;
50 }
51 
52 
53 template<class Type>
56 (
57  const regionModel& nbrRegion,
58  const word& fieldName,
59  const label regionPatchi,
60  const bool flip
61 ) const
62 {
64 
65  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
66 
67  if (nbrRegionMesh.foundObject<fieldType>(fieldName))
68  {
69  const label nbrPatchi = nbrCoupledPatchID(nbrRegion, regionPatchi);
70 
71  int oldTag = UPstream::msgType();
72  UPstream::msgType() = oldTag + 1;
73 
74  const AMIPatchToPatchInterpolation& ami =
75  interRegionAMI(nbrRegion, regionPatchi, nbrPatchi, flip);
76 
77  const fieldType& nbrField =
78  nbrRegionMesh.lookupObject<fieldType>(fieldName);
79 
80  const Field<Type>& nbrFieldp = nbrField.boundaryField()[nbrPatchi];
81 
82  tmp<Field<Type>> tresult(ami.interpolateToSource(nbrFieldp));
83 
84  UPstream::msgType() = oldTag;
85 
86  return tresult;
87  }
88  else
89  {
90  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
91 
92  return tmp<Field<Type>>::New(p.size(), Zero);
93  }
94 }
95 
96 
97 template<class Type>
100 (
101  const regionModel& nbrRegion,
102  const word& fieldName,
103  const label regionPatchi,
104  const bool flip
105 ) const
106 {
108 
109  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
110 
111  if (nbrRegionMesh.foundObject<fieldType>(fieldName))
112  {
113  const label nbrPatchi = nbrCoupledPatchID(nbrRegion, regionPatchi);
114 
115  int oldTag = UPstream::msgType();
116  UPstream::msgType() = oldTag + 1;
117 
118  const AMIPatchToPatchInterpolation& ami =
119  interRegionAMI(nbrRegion, regionPatchi, nbrPatchi, flip);
120 
121  const fieldType& nbrField =
122  nbrRegionMesh.lookupObject<fieldType>(fieldName);
123 
124  const fvPatchField<Type>& nbrFieldp =
125  nbrField.boundaryField()[nbrPatchi];
126 
127  tmp<Field<Type>> tresult
128  (
129  ami.interpolateToSource(nbrFieldp.patchInternalField())
130  );
131 
132  UPstream::msgType() = oldTag;
133 
134  return tresult;
135  }
136  else
137  {
138  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
139 
140  return tmp<Field<Type>>::New(p.size(), Zero);
141  }
142 }
143 
144 
145 template<class Type>
147 (
148  const label regionPatchi,
149  List<Type>& regionField
150 ) const
151 {
152  forAll(intCoupledPatchIDs_, i)
153  {
154  if (intCoupledPatchIDs_[i] == regionPatchi)
155  {
156  const mappedPatchBase& mpb =
157  refCast<const mappedPatchBase>
158  (
159  regionMesh().boundaryMesh()[regionPatchi]
160  );
161  mpb.reverseDistribute(regionField);
162  return;
163  }
164  }
165 
167  << "Region patch ID " << regionPatchi << " not found in region mesh"
168  << abort(FatalError);
169 }
170 
171 
172 template<class Type>
174 (
175  const label regionPatchi,
176  List<Type>& primaryField
177 ) const
178 {
179  forAll(intCoupledPatchIDs_, i)
180  {
181  if (intCoupledPatchIDs_[i] == regionPatchi)
182  {
183  const mappedPatchBase& mpb =
184  refCast<const mappedPatchBase>
185  (
186  regionMesh().boundaryMesh()[regionPatchi]
187  );
188  mpb.distribute(primaryField);
189  return;
190  }
191  }
192 
194  << "Region patch ID " << regionPatchi << " not found in region mesh"
195  << abort(FatalError);
196 }
197 
198 
199 template<class Type, class CombineOp>
201 (
202  const label regionPatchi,
203  List<Type>& regionField,
204  const CombineOp& cop
205 ) const
206 {
207  forAll(intCoupledPatchIDs_, i)
208  {
209  if (intCoupledPatchIDs_[i] == regionPatchi)
210  {
211  const mappedPatchBase& mpb =
212  refCast<const mappedPatchBase>
213  (
214  regionMesh().boundaryMesh()[regionPatchi]
215  );
216  mpb.reverseDistribute(regionField, cop);
217  return;
218  }
219  }
220 
222  << "Region patch ID " << regionPatchi << " not found in region mesh"
223  << abort(FatalError);
224 }
225 
226 
227 template<class Type, class CombineOp>
229 (
230  const label regionPatchi,
231  List<Type>& primaryField,
232  const CombineOp& cop
233 ) const
234 {
235  forAll(intCoupledPatchIDs_, i)
236  {
237  if (intCoupledPatchIDs_[i] == regionPatchi)
238  {
239  const mappedPatchBase& mpb =
240  refCast<const mappedPatchBase>
241  (
242  regionMesh().boundaryMesh()[regionPatchi]
243  );
244  mpb.distribute(primaryField, cop);
245  return;
246  }
247  }
248 
250  << "Region patch ID " << regionPatchi << " not found in region mesh"
251  << abort(FatalError);
252 }
253 
254 
255 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:112
Foam::objectRegistry::foundObject
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Definition: objectRegistryTemplates.C:379
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::regionModels::regionModel
Base class for region models.
Definition: regionModel.H:60
Foam::regionModels::regionModel::toPrimary
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
Definition: regionModelTemplates.C:147
Foam::AMIInterpolation::interpolateToSource
void interpolateToSource(const UList< Type > &fld, const CombineOp &cop, List< Type > &result, const UList< Type > &defaultValues=UList< Type >::null()) const
Definition: AMIInterpolationTemplates.C:122
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::mappedPatchBase::distribute
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
Definition: mappedPatchBaseTemplates.C:29
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::regionModels::regionModel::mapRegionPatchInternalField
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.
Foam::FatalError
error FatalError
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::regionModels::regionModel::mapRegionPatchField
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.
Foam::List< Type >
Foam::AMIInterpolation
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
Definition: AMIInterpolation.H:79
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
Foam::regionModels::regionModel::toRegion
void toRegion(const label regionPatchi, List< Type > &primaryFieldField) const
Convert a primary region field to the local region.
Definition: regionModelTemplates.C:174
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::mappedPatchBase::reverseDistribute
void reverseDistribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
Definition: mappedPatchBaseTemplates.C:96