mappedPatchFieldBase.C
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) 2013-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 #include "mappedPatchFieldBase.H"
27 #include "mappedPatchBase.H"
28 #include "interpolationCell.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
37 template<class Type>
39 (
40  const mappedPatchBase& mapper,
41  const fvPatchField<Type>& patchField,
42  const word& fieldName,
43  const bool setAverage,
44  const Type average,
45  const word& interpolationScheme
46 )
47 :
48  mapper_(mapper),
49  patchField_(patchField),
50  fieldName_(fieldName),
51  setAverage_(setAverage),
52  average_(average),
53  interpolationScheme_(interpolationScheme)
54 {}
55 
56 
57 template<class Type>
59 (
60  const mappedPatchBase& mapper,
61  const fvPatchField<Type>& patchField,
62  const dictionary& dict
63 )
64 :
65  mapper_(mapper),
66  patchField_(patchField),
67  fieldName_
68  (
69  dict.template lookupOrDefault<word>
70  (
71  "field",
72  patchField_.internalField().name()
73  )
74  ),
75  setAverage_(readBool(dict.lookup("setAverage"))),
76  average_(pTraits<Type>(dict.lookup("average"))),
77  interpolationScheme_(interpolationCell<Type>::typeName)
78 {
79  if (mapper_.mode() == mappedPatchBase::NEARESTCELL)
80  {
81  dict.lookup("interpolationScheme") >> interpolationScheme_;
82  }
83 }
84 
85 
86 template<class Type>
88 (
89  const mappedPatchBase& mapper,
90  const fvPatchField<Type>& patchField
91 )
92 :
93  mapper_(mapper),
94  patchField_(patchField),
95  fieldName_(patchField_.internalField().name()),
96  setAverage_(false),
97  average_(Zero),
98  interpolationScheme_(interpolationCell<Type>::typeName)
99 {}
100 
101 
102 template<class Type>
104 (
105  const mappedPatchFieldBase<Type>& mapper
106 )
107 :
108  mapper_(mapper.mapper_),
109  patchField_(mapper.patchField_),
110  fieldName_(mapper.fieldName_),
111  setAverage_(mapper.setAverage_),
112  average_(mapper.average_),
113  interpolationScheme_(mapper.interpolationScheme_)
114 {}
115 
116 
117 template<class Type>
119 (
120  const mappedPatchBase& mapper,
121  const fvPatchField<Type>& patchField,
122  const mappedPatchFieldBase<Type>& base
123 )
124 :
125  mapper_(mapper),
126  patchField_(patchField),
127  fieldName_(base.fieldName_),
128  setAverage_(base.setAverage_),
129  average_(base.average_),
130  interpolationScheme_(base.interpolationScheme_)
131 {}
132 
133 
134 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135 
136 template<class Type>
139 {
141 
142  const fvMesh& nbrMesh = refCast<const fvMesh>(mapper_.sampleMesh());
143 
144  if (mapper_.sameRegion())
145  {
146  if (fieldName_ == patchField_.internalField().name())
147  {
148  // Optimisation: bypass field lookup
149  return
150  dynamic_cast<const fieldType&>
151  (
152  patchField_.internalField()
153  );
154  }
155  else
156  {
157  const fvMesh& thisMesh = patchField_.patch().boundaryMesh().mesh();
158  return thisMesh.template lookupObject<fieldType>(fieldName_);
159  }
160  }
161  else
162  {
163  return nbrMesh.template lookupObject<fieldType>(fieldName_);
164  }
165 }
166 
167 
168 template<class Type>
170 {
172 
173  // Since we're inside initEvaluate/evaluate there might be processor
174  // comms underway. Change the tag we use.
175  int oldTag = UPstream::msgType();
176  UPstream::msgType() = oldTag + 1;
177 
178  const fvMesh& thisMesh = patchField_.patch().boundaryMesh().mesh();
179  const fvMesh& nbrMesh = refCast<const fvMesh>(mapper_.sampleMesh());
180 
181  // Result of obtaining remote values
182  tmp<Field<Type>> tnewValues(new Field<Type>(0));
183  Field<Type>& newValues = tnewValues.ref();
184 
185  switch (mapper_.mode())
186  {
188  {
189  const distributionMap& distMap = mapper_.map();
190 
191  if (interpolationScheme_ != interpolationCell<Type>::typeName)
192  {
193  // Send back sample points to the processor that holds the cell
194  vectorField samples(mapper_.samplePoints());
195  distMap.reverseDistribute
196  (
197  (
198  mapper_.sameRegion()
199  ? thisMesh.nCells()
200  : nbrMesh.nCells()
201  ),
202  point::max,
203  samples
204  );
205 
206  autoPtr<interpolation<Type>> interpolator
207  (
209  (
210  interpolationScheme_,
211  sampleField()
212  )
213  );
214  const interpolation<Type>& interp = interpolator();
215 
216  newValues.setSize(samples.size(), pTraits<Type>::max);
217  forAll(samples, celli)
218  {
219  if (samples[celli] != point::max)
220  {
221  newValues[celli] = interp.interpolate
222  (
223  samples[celli],
224  celli
225  );
226  }
227  }
228  }
229  else
230  {
231  newValues = sampleField();
232  }
233 
234  distMap.distribute(newValues);
235 
236  break;
237  }
240  {
241  const label nbrPatchID =
242  nbrMesh.boundaryMesh().findPatchID(mapper_.samplePatch());
243 
244  if (nbrPatchID < 0)
245  {
247  << "Unable to find sample patch " << mapper_.samplePatch()
248  << " in region " << mapper_.sampleRegion()
249  << " for patch " << patchField_.patch().name() << nl
250  << abort(FatalError);
251  }
252 
253  const fieldType& nbrField = sampleField();
254 
255  newValues = nbrField.boundaryField()[nbrPatchID];
256  mapper_.distribute(newValues);
257 
258  break;
259  }
261  {
262  Field<Type> allValues(nbrMesh.nFaces(), Zero);
263 
264  const fieldType& nbrField = sampleField();
265 
266  forAll(nbrField.boundaryField(), patchi)
267  {
268  const fvPatchField<Type>& pf =
269  nbrField.boundaryField()[patchi];
270  label faceStart = pf.patch().start();
271 
272  forAll(pf, facei)
273  {
274  allValues[faceStart++] = pf[facei];
275  }
276  }
277 
278  mapper_.distribute(allValues);
279  newValues.transfer(allValues);
280 
281  break;
282  }
283  default:
284  {
286  << "Unknown sampling mode: " << mapper_.mode()
287  << nl << abort(FatalError);
288  }
289  }
290 
291  if (setAverage_)
292  {
293  Type averagePsi =
294  gSum(patchField_.patch().magSf()*newValues)
295  /gSum(patchField_.patch().magSf());
296 
297  if (mag(averagePsi)/mag(average_) > 0.5)
298  {
299  newValues *= mag(average_)/mag(averagePsi);
300  }
301  else
302  {
303  newValues += (average_ - averagePsi);
304  }
305  }
306 
307  // Restore tag
308  UPstream::msgType() = oldTag;
309 
310  return tnewValues;
311 }
312 
313 
314 template<class Type>
316 {
317  writeEntry(os, "field", fieldName_);
318  writeEntry(os, "setAverage", setAverage_);
319  writeEntry(os, "average", average_);
320  writeEntry(os, "interpolationScheme", interpolationScheme_);
321 }
322 
323 
324 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
325 
326 } // End namespace Foam
327 
328 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void reverseDistribute(const label constructSize, List< T > &, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
word interpolationScheme_
Interpolation scheme to use for nearestcell mode.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:151
virtual void write(Ostream &) const
Write.
label nFaces() const
const mappedPatchBase & mapper_
Mapping engine.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
label nCells() const
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
label findPatchID(const word &patchName) const
Find patch index given a name.
Generic GeometricField class.
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
bool readBool(Istream &)
Definition: boolIO.C:60
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
virtual Type interpolate(const vector &position, const label celli, const label facei=-1) const =0
Interpolate field to the given point in the given cell.
scalarField samples(nIntervals, 0)
Type gSum(const FieldField< Field, Type > &f)
Pre-declare SubField and related Field type.
Definition: Field.H:56
const bool setAverage_
If true adjust the mapped field to maintain average value average_.
A class for handling words, derived from string.
Definition: word.H:59
const GeometricField< Type, fvPatchField, volMesh > & sampleField() const
Field to sample. Either on my or nbr mesh.
const polyMesh & mesh() const
Return the mesh reference.
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
static const zero Zero
Definition: zero.H:97
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:351
const Type average_
Average value the mapped field is adjusted to maintain if.
const fvPatchField< Type > & patchField_
Underlying patch field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
mappedPatchFieldBase(const mappedPatchBase &mapper, const fvPatchField< Type > &patchField, const word &fieldName, const bool setAverage, const Type average, const word &interpolationScheme)
Construct from components.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Functionality for sampling fields using mappedPatchBase. Every call to mappedField() returns a sample...
void setSize(const label)
Reset size of List.
Definition: List.C:281
Class containing processor-to-processor mapping information.
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
virtual tmp< Field< Type > > mappedField() const
Map sampleField onto *this patch.
Abstract base class for interpolation.
word fieldName_
Name of field to sample.
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
void transfer(List< Type > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
Namespace for OpenFOAM.
Uses the cell value for any point in the cell.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864