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 | Copyright (C) 2013-2015 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  "fieldName",
72  patchField_.dimensionedInternalField().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_.dimensionedInternalField().name()),
96  setAverage_(false),
97  average_(pTraits<Type>::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_.dimensionedInternalField().name())
147  {
148  // Optimisation: bypass field lookup
149  return
150  dynamic_cast<const fieldType&>
151  (
152  patchField_.dimensionedInternalField()
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();
184 
185  switch (mapper_.mode())
186  {
188  {
189  const mapDistribute& 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(), pTraits<Type>::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  os.writeKeyword("fieldName") << fieldName_ << token::END_STATEMENT << nl;
318  os.writeKeyword("setAverage") << setAverage_ << token::END_STATEMENT << nl;
319  os.writeKeyword("average") << average_ << token::END_STATEMENT << nl;
320  os.writeKeyword("interpolationScheme") << interpolationScheme_
321  << token::END_STATEMENT << nl;
322 }
323 
324 
325 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
326 
327 } // End namespace Foam
328 
329 // ************************************************************************* //
Foam::fvPatchField< Type >
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::mappedPatchFieldBase::average_
const Type average_
Average value the mapped field is adjusted to maintain if.
Definition: mappedPatchFieldBase.H:85
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name.
Definition: polyBoundaryMesh.C:678
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::mappedPatchBase::NEARESTPATCHFACEAMI
@ NEARESTPATCHFACEAMI
Definition: mappedPatchBase.H:113
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:39
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::mappedPatchFieldBase::patchField_
const fvPatchField< Type > & patchField_
Underlying patch field.
Definition: mappedPatchFieldBase.H:75
Foam::mappedPatchBase::NEARESTFACE
@ NEARESTFACE
Definition: mappedPatchBase.H:115
interpolationCell.H
Foam::Vector< scalar >::max
static const Vector max
Definition: Vector.H:82
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::interpolation::interpolate
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.
mappedPatchFieldBase.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::mappedPatchFieldBase::mappedField
virtual tmp< Field< Type > > mappedField() const
Map sampleField onto *this patch.
Definition: mappedPatchFieldBase.C:169
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:101
Foam::interpolationCell
Uses the cell value for any point in the cell.
Definition: interpolationCell.H:48
Foam::mappedPatchFieldBase::fieldName_
word fieldName_
Name of field to sample.
Definition: mappedPatchFieldBase.H:78
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
samples
scalarField samples(nIntervals, 0)
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:59
Foam::Field< Type >
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:152
Foam::polyBoundaryMesh::mesh
const polyMesh & mesh() const
Return the mesh reference.
Definition: polyBoundaryMesh.H:140
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:155
Foam::interpolation
Abstract base class for interpolation.
Definition: mappedPatchFieldBase.H:57
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:64
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::mappedPatchFieldBase::interpolationScheme_
word interpolationScheme_
Interpolation scheme to use for nearestcell mode.
Definition: mappedPatchFieldBase.H:88
Foam::mappedPatchFieldBase::write
virtual void write(Ostream &) const
Write.
Definition: mappedPatchFieldBase.C:315
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:452
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::mappedPatchFieldBase::mapper_
const mappedPatchBase & mapper_
Mapping engine.
Definition: mappedPatchFieldBase.H:72
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::mappedPatchBase::NEARESTPATCHFACE
@ NEARESTPATCHFACE
Definition: mappedPatchBase.H:112
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::fvPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
Foam::mappedPatchBase::NEARESTCELL
@ NEARESTCELL
Definition: mappedPatchBase.H:111
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:187
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:300
Foam::readBool
bool readBool(Istream &)
Definition: boolIO.C:60
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::mappedPatchFieldBase::sampleField
const GeometricField< Type, fvPatchField, volMesh > & sampleField() const
Field to sample. Either on my or nbr mesh.
Definition: mappedPatchFieldBase.C:138
Foam::mappedPatchFieldBase::setAverage_
const bool setAverage_
If true adjust the mapped field to maintain average value average_.
Definition: mappedPatchFieldBase.H:81
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:335
mappedPatchBase.H