mappedPatchFieldBase.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 #include "mappedPatchFieldBase.H"
30 #include "mappedPatchBase.H"
31 #include "interpolationCell.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const dictionary& dict,
39  const bool mandatory
40 )
41 {
42  if (mandatory)
43  {
44  return dict.get<Type>("average");
45  }
46 
47  return Zero;
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
53 template<class Type>
55 (
56  const mappedPatchBase& mapper,
57  const fvPatchField<Type>& patchField,
58  const word& fieldName,
59  const bool setAverage,
60  const Type average,
61  const word& interpolationScheme
62 )
63 :
64  mapper_(mapper),
65  patchField_(patchField),
66  fieldName_(fieldName),
67  setAverage_(setAverage),
68  average_(average),
69  interpolationScheme_(interpolationScheme)
70 {}
71 
72 
73 template<class Type>
75 (
76  const mappedPatchBase& mapper,
77  const fvPatchField<Type>& patchField,
78  const dictionary& dict
79 )
80 :
81  mapper_(mapper),
82  patchField_(patchField),
83  fieldName_
84  (
85  dict.template lookupOrDefault<word>
86  (
87  "field",
88  patchField_.internalField().name()
89  )
90  ),
91  setAverage_(dict.lookupOrDefault("setAverage", false)),
92  average_(getAverage(dict, setAverage_)),
93  interpolationScheme_(interpolationCell<Type>::typeName)
94 {
95  if (mapper_.mode() == mappedPatchBase::NEARESTCELL)
96  {
97  dict.readEntry("interpolationScheme", interpolationScheme_);
98  }
99 }
100 
101 
102 template<class Type>
104 (
105  const mappedPatchBase& mapper,
106  const fvPatchField<Type>& patchField
107 )
108 :
109  mapper_(mapper),
110  patchField_(patchField),
111  fieldName_(patchField_.internalField().name()),
112  setAverage_(false),
113  average_(Zero),
114  interpolationScheme_(interpolationCell<Type>::typeName)
115 {}
116 
117 
118 template<class Type>
120 (
121  const mappedPatchFieldBase<Type>& mapper
122 )
123 :
124  mapper_(mapper.mapper_),
125  patchField_(mapper.patchField_),
126  fieldName_(mapper.fieldName_),
127  setAverage_(mapper.setAverage_),
128  average_(mapper.average_),
129  interpolationScheme_(mapper.interpolationScheme_)
130 {}
131 
132 
133 template<class Type>
135 (
136  const mappedPatchBase& mapper,
137  const fvPatchField<Type>& patchField,
138  const mappedPatchFieldBase<Type>& base
139 )
140 :
141  mapper_(mapper),
142  patchField_(patchField),
143  fieldName_(base.fieldName_),
144  setAverage_(base.setAverage_),
145  average_(base.average_),
146  interpolationScheme_(base.interpolationScheme_)
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
152 template<class Type>
155 {
157 
158  if (mapper_.sameRegion())
159  {
160  if (fieldName_ == patchField_.internalField().name())
161  {
162  // Optimisation: bypass field lookup
163  return
164  dynamic_cast<const fieldType&>
165  (
166  patchField_.internalField()
167  );
168  }
169  else
170  {
171  const fvMesh& thisMesh = patchField_.patch().boundaryMesh().mesh();
172  return thisMesh.template lookupObject<fieldType>(fieldName_);
173  }
174  }
175 
176  const fvMesh& nbrMesh = refCast<const fvMesh>(mapper_.sampleMesh());
177  return nbrMesh.template lookupObject<fieldType>(fieldName_);
178 }
179 
180 
181 template<class Type>
184 {
186 
187  // Since we're inside initEvaluate/evaluate there might be processor
188  // comms underway. Change the tag we use.
189  int oldTag = UPstream::msgType();
190  UPstream::msgType() = oldTag + 1;
191 
192  const fvMesh& thisMesh = patchField_.patch().boundaryMesh().mesh();
193  const fvMesh& nbrMesh = refCast<const fvMesh>(mapper_.sampleMesh());
194 
195  // Result of obtaining remote values
196  auto tnewValues = tmp<Field<Type>>::New();
197  auto& newValues = tnewValues.ref();
198 
199  switch (mapper_.mode())
200  {
201  case mappedPatchBase::NEARESTCELL:
202  {
203  const mapDistribute& distMap = mapper_.map();
204 
205  if (interpolationScheme_ != interpolationCell<Type>::typeName)
206  {
207  // Send back sample points to the processor that holds the cell
208  vectorField samples(mapper_.samplePoints());
209  distMap.reverseDistribute
210  (
211  (
212  mapper_.sameRegion()
213  ? thisMesh.nCells()
214  : nbrMesh.nCells()
215  ),
216  point::max,
217  samples
218  );
219 
220  auto interpolator =
222  (
223  interpolationScheme_,
224  sampleField()
225  );
226 
227  const auto& interp = *interpolator;
228 
229  newValues.setSize(samples.size(), pTraits<Type>::max);
230  forAll(samples, celli)
231  {
232  if (samples[celli] != point::max)
233  {
234  newValues[celli] = interp.interpolate
235  (
236  samples[celli],
237  celli
238  );
239  }
240  }
241  }
242  else
243  {
244  newValues = sampleField();
245  }
246 
247  distMap.distribute(newValues);
248 
249  break;
250  }
251  case mappedPatchBase::NEARESTPATCHFACE:
252  case mappedPatchBase::NEARESTPATCHFACEAMI:
253  {
254  const label nbrPatchID =
255  nbrMesh.boundaryMesh().findPatchID(mapper_.samplePatch());
256 
257  if (nbrPatchID < 0)
258  {
260  << "Unable to find sample patch " << mapper_.samplePatch()
261  << " in region " << mapper_.sampleRegion()
262  << " for patch " << patchField_.patch().name() << nl
263  << abort(FatalError);
264  }
265 
266  const fieldType& nbrField = sampleField();
267 
268  newValues = nbrField.boundaryField()[nbrPatchID];
269  mapper_.distribute(newValues);
270 
271  break;
272  }
273  case mappedPatchBase::NEARESTFACE:
274  {
275  Field<Type> allValues(nbrMesh.nFaces(), Zero);
276 
277  const fieldType& nbrField = sampleField();
278 
279  for (const fvPatchField<Type>& pf : nbrField.boundaryField())
280  {
281  label faceStart = pf.patch().start();
282 
283  forAll(pf, facei)
284  {
285  allValues[faceStart++] = pf[facei];
286  }
287  }
288 
289  mapper_.distribute(allValues);
290  newValues.transfer(allValues);
291 
292  break;
293  }
294  default:
295  {
297  << "Unknown sampling mode: " << mapper_.mode() << nl
298  << abort(FatalError);
299  }
300  }
301 
302  if (setAverage_)
303  {
304  Type averagePsi =
305  gSum(patchField_.patch().magSf()*newValues)
306  /gSum(patchField_.patch().magSf());
307 
308  if (mag(averagePsi)/mag(average_) > 0.5)
309  {
310  newValues *= mag(average_)/mag(averagePsi);
311  }
312  else
313  {
314  newValues += (average_ - averagePsi);
315  }
316  }
317 
318  // Restore tag
319  UPstream::msgType() = oldTag;
320 
321  return tnewValues;
322 }
323 
324 
325 template<class Type>
327 {
328  os.writeEntry("field", fieldName_);
329 
330  if (setAverage_)
331  {
332  os.writeEntry("setAverage", "true");
333  os.writeEntry("average", average_);
334  }
335 
336  os.writeEntry("interpolationScheme", interpolationScheme_);
337 }
338 
339 
340 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
Foam::mappedPatchFieldBase::average_
const Type average_
Definition: mappedPatchFieldBase.H:127
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::mappedPatchFieldBase::mappedPatchFieldBase
mappedPatchFieldBase(const mappedPatchBase &mapper, const fvPatchField< Type > &patchField, const word &fieldName, const bool setAverage, const Type average, const word &interpolationScheme)
Construct from components.
Definition: mappedPatchFieldBase.C:55
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::mappedPatchFieldBase::patchField_
const fvPatchField< Type > & patchField_
Underlying patch field.
Definition: mappedPatchFieldBase.H:117
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
interpolationCell.H
Foam::mappedPatchFieldBase::mappedField
virtual tmp< Field< Type > > mappedField() const
Map sampleField onto *this patch.
Definition: mappedPatchFieldBase.C:183
mappedPatchFieldBase.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:105
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::interpolationCell
Uses the cell value for any point in the cell.
Definition: interpolationCell.H:50
Foam::mappedPatchFieldBase::fieldName_
word fieldName_
Name of field to sample.
Definition: mappedPatchFieldBase.H:120
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::label
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:62
Foam::Field< vector >
Foam::fac::average
tmp< GeometricField< Type, faPatchField, areaMesh > > average(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Area-weighted average a edgeField creating a areaField.
Definition: facAverage.C:47
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
samples
scalarField samples(nIntervals, Zero)
Foam::polyBoundaryMesh::mesh
const polyMesh & mesh() const
Return the mesh reference.
Definition: polyBoundaryMesh.H:144
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Definition: polyBoundaryMesh.C:766
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::mapDistribute::distribute
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Definition: mapDistributeTemplates.C:152
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mappedPatchFieldBase
Functionality for sampling fields using mappedPatchBase. Every call to mappedField() returns a sample...
Definition: mappedPatchFieldBase.H:102
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::mappedPatchFieldBase::interpolationScheme_
word interpolationScheme_
Interpolation scheme to use for nearestcell mode.
Definition: mappedPatchFieldBase.H:130
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
Foam::mappedPatchFieldBase::write
virtual void write(Ostream &os) const
Write.
Definition: mappedPatchFieldBase.C:326
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::mappedPatchFieldBase::mapper_
const mappedPatchBase & mapper_
Mapping engine.
Definition: mappedPatchFieldBase.H:114
Foam::mappedPatchFieldBase::sampleField
const GeometricField< Type, fvPatchField, volMesh > & sampleField() const
Field to sample. Either on my or nbr mesh.
Definition: mappedPatchFieldBase.C:154
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:52
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::mapDistribute::reverseDistribute
void reverseDistribute(const label constructSize, List< T > &, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
Definition: mapDistributeTemplates.C:182
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh >
Foam::mappedPatchFieldBase::setAverage_
const bool setAverage_
If true adjust the mapped field to maintain average value average_.
Definition: mappedPatchFieldBase.H:123
mappedPatchBase.H