Loading...
Searching...
No Matches
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-2025 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>
35Type Foam::PatchFunction1Types::Sampled<Type>::getAverage
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 const polyPatch& pp
81)
82:
83 PatchFunction1<Type>(rhs, pp),
85 fieldName_(rhs.fieldName_),
94(
95 const Sampled<Type>& rhs
96)
97:
98 Sampled<Type>(rhs, rhs.patch())
99{}
100
101
102// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103
104template<class Type>
107{
109
110 if (this->sameRegion())
111 {
112 const polyMesh& thisMesh =
113 this->mappedPatchBase::patch_.boundaryMesh().mesh();
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 =
132 this->mappedPatchBase::patch_.boundaryMesh().mesh();
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 const int oldTag = UPstream::incrMsgType();
155
156 const fvMesh& thisMesh = refCast<const fvMesh>
157 (
159 );
160 const fvMesh& nbrMesh = refCast<const fvMesh>(this->sampleMesh());
161
162
163 // Result of obtaining remote values
164 auto tnewValues = tmp<Field<Type>>::New();
165 auto& newValues = tnewValues.ref();
166
167 if (!haveSampleField())
168 {
169 UPstream::msgType(oldTag); // Restore tag
170 newValues.resize_nocopy(this->mappedPatchBase::patch_.size());
171 newValues = Zero;
172 return this->transform(tnewValues);
173 }
174
175 switch (this->mode())
176 {
178 {
179 const mapDistribute& distMap = this->map();
180
181 if (interpolationScheme_ != interpolationCell<Type>::typeName)
182 {
183 // Send back sample points to the processor that holds the cell
184 vectorField samples(this->samplePoints());
185 distMap.reverseDistribute
186 (
187 (
188 this->sameRegion()
189 ? thisMesh.nCells()
190 : nbrMesh.nCells()
191 ),
193 samples
194 );
195
196 auto interpolator =
198 (
199 interpolationScheme_,
200 sampleField()
201 );
202
203 const auto& interp = *interpolator;
204
205 newValues.setSize(samples.size(), pTraits<Type>::max);
206 forAll(samples, celli)
207 {
208 if (samples[celli] != point::max)
209 {
210 newValues[celli] = interp.interpolate
211 (
212 samples[celli],
213 celli
214 );
215 }
216 }
217 }
218 else
219 {
220 newValues = sampleField();
221 }
222 distMap.distribute(newValues);
223
224 break;
225 }
228 {
229 const label nbrPatchID =
230 nbrMesh.boundaryMesh().findPatchID(this->samplePatch());
231
232 if (nbrPatchID < 0)
233 {
235 << "Unable to find sample patch " << this->samplePatch()
236 << " in region " << this->sampleRegion()
237 << " for patch " << this->mappedPatchBase::patch_.name() << nl
238 << abort(FatalError);
239 }
240
241 const fieldType& nbrField = sampleField();
242
243 newValues = nbrField.boundaryField()[nbrPatchID];
244 this->distribute(newValues);
245
246 break;
247 }
249 {
250 Field<Type> allValues(nbrMesh.nFaces(), Zero);
251
252 const fieldType& nbrField = sampleField();
253
254 for (const fvPatchField<Type>& pf : nbrField.boundaryField())
255 {
256 label faceStart = pf.patch().start();
257
258 forAll(pf, facei)
259 {
260 allValues[faceStart++] = pf[facei];
261 }
262 }
263
264 this->distribute(allValues);
265 newValues.transfer(allValues);
266
267 break;
268 }
269 default:
270 {
272 << "Unknown sampling mode: " << this->mode() << nl
273 << abort(FatalError);
274 }
275 }
276
277 // Enforce average. Either by scaling (if scaling factor > 0.5) or by
278 // offsetting.
279 if (setAverage_ && returnReduceOr(newValues.size()))
280 {
281 Type averagePsi;
282 if (this->faceValues())
283 {
284 const scalarField magSf
285 (
286 mag(this->mappedPatchBase::patch_.faceAreas())
287 );
288 averagePsi = gWeightedAverage(magSf, newValues);
289 }
290 else
291 {
292 averagePsi = gAverage(newValues);
293 }
294
295 if (mag(averagePsi) > 0.5*mag(average_))
296 {
297 newValues *= mag(average_)/mag(averagePsi);
298 }
299 else
300 {
301 newValues += (average_ - averagePsi);
302 }
303 }
304
305 UPstream::msgType(oldTag); // Restore tag
307 return this->transform(tnewValues);
308}
309
310
311template<class Type>
314(
315 const scalar x1,
316 const scalar x2
317) const
320 return nullptr;
321}
322
323
324template<class Type>
326(
327 Ostream& os
328) const
329{
331
332 os.writeEntry(this->name(), type());
333
335
336 os.writeEntry("field", fieldName_);
337 if (setAverage_)
338 {
339 os.writeEntry("setAverage", "true");
340 os.writeEntry("average", average_);
341 }
342
343 os.writeEntry("interpolationScheme", interpolationScheme_);
344}
345
346
347// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
PatchFunction1 to sample an existing field.
Definition Sampled.H:142
virtual tmp< Field< Type > > integrate(const scalar x1, const scalar x2) const
Integrate between two values.
Definition Sampled.C:307
virtual void writeData(Ostream &os) const
Write in dictionary format.
Definition Sampled.C:319
const Type average_
Average value the mapped field is adjusted to maintain if setAverage_ is set true.
Definition Sampled.H:169
bool haveSampleField() const
Field to sample. Either on my or nbr mesh.
Definition Sampled.C:118
virtual tmp< Field< Type > > value(const scalar x) const
Return sampled value.
Definition Sampled.C:139
Sampled(const polyPatch &pp, const word &redirectType, const word &entryName, const dictionary &dict, const bool faceValues=true)
Construct from entry name and dictionary.
Definition Sampled.C:47
word fieldName_
Name of the field.
Definition Sampled.H:158
const GeometricField< Type, fvPatchField, volMesh > & sampleField() const
Field to sample. Either on my or nbr mesh.
Definition Sampled.C:99
word interpolationScheme_
Interpolation scheme to use for nearestcell mode.
Definition Sampled.H:174
const bool setAverage_
If true adjust the mapped field to maintain average value average_.
Definition Sampled.H:163
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
const polyPatch const word const word & entryName
virtual void writeData(Ostream &os) const
Write in dictionary format.
const polyPatch const word const word const dictionary & dict
const polyPatch & pp
const polyPatch const word const word const dictionary const bool faceValues
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
static const Form max
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Uses the cell value for any location within the cell.
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
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 List data using default commsType, default flip/negate operator.
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
tmp< pointField > samplePoints(const pointField &) const
Get the sample points given the face points.
const mapDistribute & map() const
Return reference to the parallel distribution map.
bool sameRegion() const noexcept
Cached sampleRegion != mesh.name().
const polyPatch & patch_
Patch to sample.
@ NEARESTCELL
nearest cell containing sample
@ NEARESTPATCHFACE
nearest face on selected patch
@ NEARESTPATCHFACEAMI
nearest patch face + AMI interpolation
const polyMesh & sampleMesh() const
Get the region mesh.
virtual void write(Ostream &os) const
Write as a dictionary.
mappedPatchBase(const polyPatch &)
Construct from patch.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE).
sampleMode mode() const noexcept
What to sample.
const word & sampleRegion() const
Region to sample.
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
const polyPatch & patch() const noexcept
Reference to the patch.
label size() const
Number of faces or points on the patch.
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.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
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:75
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition POSIX.C:775
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
scalarField samples(nIntervals, Zero)