Sampled.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) 2018-2020 OpenCFD Ltd.
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
28#include "fvMesh.H"
29#include "volFields.H"
30#include "interpolationCell.H"
31
32// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33
34template<class Type>
36(
37 const dictionary& dict,
38 const bool mandatory
39)
40{
41 if (mandatory)
42 {
43 return dict.get<Type>("average");
44 }
45
46 return Zero;
47}
48
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
52template<class Type>
54(
55 const polyPatch& pp,
56 const word& redirectType,
57 const word& entryName,
58 const dictionary& dict,
59 const bool faceValues
60)
61:
62 PatchFunction1<Type>(pp, entryName, dict, faceValues),
64 fieldName_(dict.get<word>("field")),
65 setAverage_(dict.getOrDefault("setAverage", false)),
66 average_(getAverage(dict, setAverage_)),
67 interpolationScheme_(interpolationCell<Type>::typeName)
68{
70 {
71 dict.readEntry("interpolationScheme", interpolationScheme_);
72 }
73}
74
75
76template<class Type>
78(
79 const Sampled<Type>& rhs
80)
81:
82 Sampled<Type>(rhs, rhs.patch())
83{}
84
85
86template<class Type>
88(
89 const Sampled<Type>& rhs,
90 const polyPatch& pp
91)
92:
93 PatchFunction1<Type>(rhs, pp),
94 mappedPatchBase(pp, rhs),
95 fieldName_(rhs.fieldName_),
96 setAverage_(rhs.setAverage_),
97 average_(rhs.average_),
98 interpolationScheme_(rhs.interpolationScheme_)
99{}
100
101
102// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103
104template<class Type>
107{
109
110 if (this->sameRegion())
111 {
112 const polyMesh& thisMesh =
114 return thisMesh.template lookupObject<fieldType>(fieldName_);
115 }
116 else
117 {
118 const fvMesh& nbrMesh = refCast<const fvMesh>(this->sampleMesh());
119 return nbrMesh.template lookupObject<fieldType>(fieldName_);
120 }
121}
122
123
124template<class Type>
126{
128
129 if (this->sameRegion())
130 {
131 const polyMesh& thisMesh =
133 return thisMesh.template foundObject<fieldType>(fieldName_);
134 }
135 else
136 {
137 const fvMesh& nbrMesh = refCast<const fvMesh>(this->sampleMesh());
138 return nbrMesh.template foundObject<fieldType>(fieldName_);
139 }
140}
141
142
143template<class Type>
146(
147 const scalar x
148) const
149{
151
152 // Since we're inside initEvaluate/evaluate there might be processor
153 // comms underway. Change the tag we use.
154 int oldTag = UPstream::msgType();
155 UPstream::msgType() = oldTag + 1;
156
157 const fvMesh& thisMesh = refCast<const fvMesh>
158 (
160 );
161 const fvMesh& nbrMesh = refCast<const fvMesh>(this->sampleMesh());
162
163
164 // Result of obtaining remote values
165 auto tnewValues = tmp<Field<Type>>::New();
166 auto& newValues = tnewValues.ref();
167
168 if (!haveSampleField())
169 {
170 // Restore tag
171 UPstream::msgType() = oldTag;
172 newValues.setSize(this->mappedPatchBase::patch_.size());
173 newValues = Zero;
174 return this->transform(tnewValues);
175 }
176
177 switch (this->mode())
178 {
180 {
181 const mapDistribute& distMap = this->map();
182
183 if (interpolationScheme_ != interpolationCell<Type>::typeName)
184 {
185 // Send back sample points to the processor that holds the cell
186 vectorField samples(this->samplePoints());
187 distMap.reverseDistribute
188 (
189 (
190 this->sameRegion()
191 ? thisMesh.nCells()
192 : nbrMesh.nCells()
193 ),
195 samples
196 );
197
198 auto interpolator =
200 (
201 interpolationScheme_,
202 sampleField()
203 );
204
205 const auto& interp = *interpolator;
206
207 newValues.setSize(samples.size(), pTraits<Type>::max);
208 forAll(samples, celli)
209 {
210 if (samples[celli] != point::max)
211 {
212 newValues[celli] = interp.interpolate
213 (
214 samples[celli],
215 celli
216 );
217 }
218 }
219 }
220 else
221 {
222 newValues = sampleField();
223 }
224 distMap.distribute(newValues);
225
226 break;
227 }
230 {
231 const label nbrPatchID =
232 nbrMesh.boundaryMesh().findPatchID(this->samplePatch());
233
234 if (nbrPatchID < 0)
235 {
237 << "Unable to find sample patch " << this->samplePatch()
238 << " in region " << this->sampleRegion()
239 << " for patch " << this->mappedPatchBase::patch_.name() << nl
240 << abort(FatalError);
241 }
242
243 const fieldType& nbrField = sampleField();
244
245 newValues = nbrField.boundaryField()[nbrPatchID];
246 this->distribute(newValues);
247
248 break;
249 }
251 {
252 Field<Type> allValues(nbrMesh.nFaces(), Zero);
253
254 const fieldType& nbrField = sampleField();
255
256 for (const fvPatchField<Type>& pf : nbrField.boundaryField())
257 {
258 label faceStart = pf.patch().start();
259
260 forAll(pf, facei)
261 {
262 allValues[faceStart++] = pf[facei];
263 }
264 }
265
266 this->distribute(allValues);
267 newValues.transfer(allValues);
268
269 break;
270 }
271 default:
272 {
274 << "Unknown sampling mode: " << this->mode() << nl
275 << abort(FatalError);
276 }
277 }
278
279 // Enforce average. Either by scaling (if scaling factor > 0.5) or by
280 // offsetting.
281 if (setAverage_ && returnReduce(newValues.size(), sumOp<label>()))
282 {
283 Type averagePsi;
284 if (this->faceValues())
285 {
286 const scalarField magSf
287 (
288 mag(this->mappedPatchBase::patch_.faceAreas())
289 );
290 averagePsi = gSum(magSf*newValues)/gSum(magSf);
291 }
292 else
293 {
294 averagePsi = gAverage(newValues);
295 }
296
297 if (mag(averagePsi) > 0.5*mag(average_))
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 this->transform(tnewValues);
311}
312
313
314template<class Type>
317(
318 const scalar x1,
319 const scalar x2
320) const
321{
323 return tmp<Field<Type>>(nullptr);
324}
325
326
327template<class Type>
329(
330 Ostream& os
331) const
332{
334
335 os.writeEntry(this->name(), type());
336
338
339 os.writeEntry("field", fieldName_);
340 if (setAverage_)
341 {
342 os.writeEntry("setAverage", "true");
343 os.writeEntry("average", average_);
344 }
345
346 os.writeEntry("interpolationScheme", interpolationScheme_);
347}
348
349
350// ************************************************************************* //
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
PatchFunction1 to sample an existing field.
Definition: Sampled.H:146
virtual tmp< Field< Type > > integrate(const scalar x1, const scalar x2) const
Integrate between two values.
Definition: Sampled.C:317
bool haveSampleField() const
Field to sample. Either on my or nbr mesh.
Definition: Sampled.C:125
const GeometricField< Type, fvPatchField, volMesh > & sampleField() const
Field to sample. Either on my or nbr mesh.
Definition: Sampled.C:106
word interpolationScheme_
Interpolation scheme to use for nearestcell mode.
Definition: Sampled.H:168
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:556
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
virtual bool write()
Write the output fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
Uses the cell value for any location within the cell.
Class containing processor-to-processor mapping information.
void reverseDistribute(const label constructSize, List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
const polyPatch & patch_
Patch to sample.
@ NEARESTCELL
nearest cell containing sample
@ NEARESTPATCHFACE
nearest face on selected patch
@ NEARESTPATCHFACEAMI
nearest patch face + AMI interpolation
@ NEARESTFACE
nearest face
sampleMode mode() const noexcept
What to sample.
static const complex max
complex (VGREAT,VGREAT)
Definition: complex.H:293
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
const word & name() const noexcept
The patch name.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
const polyMesh & mesh() const noexcept
Return the mesh reference.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:315
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
Type gSum(const FieldField< Field, Type > &f)
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition: MSwindows.C:572
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
scalarField samples(nIntervals, Zero)