Loading...
Searching...
No Matches
timeVaryingMappedFixedValuePointPatchField.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) 2012-2017 OpenFOAM Foundation
9 Copyright (C) 2020-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
30#include "Time.H"
31#include "rawIOField.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
35template<class Type>
38(
39 const pointPatch& p,
41)
42:
44 setAverage_(false),
45 perturb_(0),
46 fieldTableName_(iF.name()),
47 pointsName_("points"),
48 mapMethod_(),
49 mapperPtr_(nullptr),
50 sampleTimes_(),
51 sampleIndex_(-1, -1),
52 sampleAverage_(Zero, Zero),
53 sampleValues_(),
54 offset_(nullptr)
55{}
56
57
58template<class Type>
61(
62 const pointPatch& p,
64 const dictionary& dict
65)
66:
67 fixedValuePointPatchField<Type>(p, iF, dict, IOobjectOption::NO_READ),
68 setAverage_(dict.getOrDefault("setAverage", false)),
69 perturb_(dict.getOrDefault("perturb", 1e-5)),
70 fieldTableName_(iF.name()),
71 pointsName_(dict.getOrDefault<word>("points", "points")),
72 mapMethod_(),
73 mapperPtr_(nullptr),
74 sampleTimes_(),
75 sampleIndex_(-1, -1),
76 sampleAverage_(Zero, Zero),
77 sampleValues_(),
78 offset_
79 (
80 Function1<Type>::NewIfPresent("offset", dict, word::null, &this->db())
81 )
82{
83 if
84 (
85 dict.readIfPresent("mapMethod", mapMethod_)
86 && !mapMethod_.empty()
87 && mapMethod_ != "nearest"
88 && !mapMethod_.starts_with("planar")
89 )
90 {
91 FatalIOErrorInFunction(dict)
92 << "Unknown mapMethod type " << mapMethod_
93 << "\n\nValid mapMethod types :\n"
94 << "(nearest planar)" << nl
95 << exit(FatalIOError);
96 }
97
98 dict.readIfPresentCompat
99 (
100 "fieldTable", {{"fieldTableName", 2206}},
101 fieldTableName_
102 );
103
104 if (!this->readValueEntry(dict))
105 {
106 // Note: use evaluate to do updateCoeffs followed by a reset
107 // of the pointPatchField::updated_ flag. This is
108 // so if first use is in the next time step it retriggers
109 // a new update.
111 }
112}
113
114
115template<class Type>
118(
119 const timeVaryingMappedFixedValuePointPatchField<Type>& ptf,
120 const pointPatch& p,
121 const DimensionedField<Type, pointMesh>& iF,
122 const pointPatchFieldMapper& mapper
123)
124:
125 fixedValuePointPatchField<Type>(ptf, p, iF, mapper),
126 setAverage_(ptf.setAverage_),
127 perturb_(ptf.perturb_),
128 fieldTableName_(ptf.fieldTableName_),
129 pointsName_(ptf.pointsName_),
130 mapMethod_(ptf.mapMethod_),
131 mapperPtr_(nullptr),
132 sampleTimes_(),
133 sampleIndex_(-1, -1),
134 sampleAverage_(Zero, Zero),
135 sampleValues_(),
136 offset_(ptf.offset_.clone())
137{}
138
139
140template<class Type>
143(
145)
146:
147 fixedValuePointPatchField<Type>(ptf),
148 setAverage_(ptf.setAverage_),
149 perturb_(ptf.perturb_),
150 fieldTableName_(ptf.fieldTableName_),
151 pointsName_(ptf.pointsName_),
152 mapMethod_(ptf.mapMethod_),
153 mapperPtr_(ptf.mapperPtr_),
154 sampleTimes_(ptf.sampleTimes_),
155 sampleIndex_(ptf.sampleIndex_),
156 sampleAverage_(ptf.sampleAverage_),
157 sampleValues_(ptf.sampleValues_),
158 offset_(ptf.offset_.clone())
159{}
160
161
162template<class Type>
165(
168)
169:
170 fixedValuePointPatchField<Type>(ptf, iF),
171 setAverage_(ptf.setAverage_),
172 perturb_(ptf.perturb_),
173 fieldTableName_(ptf.fieldTableName_),
174 pointsName_(ptf.pointsName_),
175 mapMethod_(ptf.mapMethod_),
176 mapperPtr_(ptf.mapperPtr_),
177 sampleTimes_(ptf.sampleTimes_),
178 sampleIndex_(ptf.sampleIndex_),
179 sampleAverage_(ptf.sampleAverage_),
180 sampleValues_(ptf.sampleValues_),
181 offset_(ptf.offset_.clone())
182{}
183
184
185// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186
187template<class Type>
189(
190 const pointPatchFieldMapper& m
191)
192{
194
195 if (sampleValues_.first().size())
196 {
197 sampleValues_.first().autoMap(m);
198 }
199
200 if (sampleValues_.second().size())
201 {
202 sampleValues_.second().autoMap(m);
203 }
204
205 // Clear interpolator
206 mapperPtr_.reset(nullptr);
207 sampleIndex_ = labelPair(-1, -1);
208}
209
210
211template<class Type>
213(
214 const pointPatchField<Type>& ptf,
215 const labelList& addr
216)
217{
219
222
223 sampleValues_.first().rmap(tiptf.sampleValues_.first(), addr);
224 sampleValues_.second().rmap(tiptf.sampleValues_.second(), addr);
225
226 // Clear interpolator
227 mapperPtr_.reset(nullptr);
228 sampleIndex_ = labelPair(-1, -1);
229}
230
231
232template<class Type>
233void Foam::timeVaryingMappedFixedValuePointPatchField<Type>::updateSampledValues
234(
235 const label sampleIndex,
236 Field<Type>& field,
237 Type& avg
238) const
239{
240 tmp<Field<Type>> tvalues;
241
242 // Update sampled data fields
243 {
244 const Time& time = this->db().time();
245
246 if (debug)
247 {
248 Pout<< "checkTable : Reading values from "
249 <<
250 (
251 "boundaryData"
252 / this->patch().name()
253 / sampleTimes_[sampleIndex].name()
254 / fieldTableName_
255 ) << endl;
256 }
257
258 // Reread values and interpolate
259 const fileName valsFile
260 (
261 time.caseConstant()
262 /"boundaryData"
263 /this->patch().name()
264 /sampleTimes_[sampleIndex].name()
265 /fieldTableName_
266 );
267
269 (
270 valsFile, // absolute path
271 time,
275 true // is global object (currently not used)
276 );
277
278 rawIOField<Type> vals(io, setAverage_);
279
280 if (vals.hasAverage())
281 {
282 avg = vals.average();
283 }
284
285 if (vals.size() != mapperPtr_().sourceSize())
286 {
288 << "Number of values (" << vals.size()
289 << ") differs from the number of points ("
290 << mapperPtr_().sourceSize()
291 << ") in file " << valsFile << exit(FatalError);
292 }
293
294 tvalues = tmp<Field<Type>>::New(std::move(vals.field()));
295 }
296
297 // From input values to interpolated (sampled) positions
298 field = mapperPtr_().interpolate(tvalues);
299}
300
301
302template<class Type>
303void Foam::timeVaryingMappedFixedValuePointPatchField<Type>::checkTable
304(
305 const scalar t
306)
307{
308 const Time& time = this->db().time();
309
310 // Initialise
311 if (sampleIndex_.first() == -1 && sampleIndex_.second() == -1)
312 {
313 const polyMesh& pMesh = this->patch().boundaryMesh().mesh()();
314
315 // Read the initial point position
316 pointField meshPts;
317
318 if (pMesh.pointsInstance() == pMesh.facesInstance())
319 {
320 meshPts = pointField(pMesh.points(), this->patch().meshPoints());
321 }
322 else
323 {
324 // Load points from facesInstance
326 << "Reloading points0 from " << pMesh.facesInstance()
327 << endl;
328
330 (
332 (
333 "points",
334 pMesh.facesInstance(),
336 pMesh,
340 )
341 );
342 meshPts = pointField(points0, this->patch().meshPoints());
343 }
344
345 // Reread values and interpolate
346 const fileName samplePointsFile
347 (
348 time.caseConstant()
349 /"boundaryData"
350 /this->patch().name()
351 /pointsName_
352 );
353
355 (
356 samplePointsFile, // absolute path
357 time,
361 true // is global object (currently not used)
362 );
363
364 // Read data (no average value!)
365 const rawIOField<point> samplePoints(io);
366
367 // tbd: run-time selection
368 const bool nearestOnly =
369 (
370 !mapMethod_.empty() && !mapMethod_.starts_with("planar")
371 );
372
373 // Allocate the interpolator
374 mapperPtr_.reset
375 (
377 (
378 samplePoints,
379 meshPts,
380 perturb_,
381 nearestOnly
382 )
383 );
384
385
386 // Read the times for which data is available
387 const fileName samplePointsDir = samplePointsFile.path();
388 sampleTimes_ = Time::findTimes(samplePointsDir);
389
391 << "Found times "
393 << endl;
394 }
395
396 // Find range of current time indices in sampleTimes
397 Pair<label> timeIndices = instant::findRange
398 (
399 sampleTimes_,
400 t, // time.value(),
401 sampleIndex_.first()
402 );
403
404 if (timeIndices.first() < 0)
405 {
407 << "Cannot find starting sampling values for current time "
408 << time.value() << nl
409 << "Have sampling values for times "
411 << "In directory "
412 << time.constant()/"boundaryData"/this->patch().name()
413 << "\n on patch " << this->patch().name()
414 << " of field " << fieldTableName_
415 << exit(FatalError);
416 }
417
418
419 // Update sampled data fields.
420
421 if (sampleIndex_.first() != timeIndices.first())
422 {
423 sampleIndex_.first() = timeIndices.first();
424
425 if (sampleIndex_.first() == sampleIndex_.second())
426 {
427 // Can reuse previous end values
428 sampleValues_.first() = sampleValues_.second();
429 sampleAverage_.first() = sampleAverage_.second();
430 }
431 else
432 {
433 // Update first() values
434 this->updateSampledValues
435 (
436 sampleIndex_.first(),
437 sampleValues_.first(),
438 sampleAverage_.first()
439 );
440 }
441 }
442
443 if (sampleIndex_.second() != timeIndices.second())
444 {
445 sampleIndex_.second() = timeIndices.second();
446
447 if (sampleIndex_.second() == -1)
448 {
449 // Index no longer valid - clear values accordingly
450 sampleValues_.second().clear();
451 }
452 else
453 {
454 // Update second() values
455 this->updateSampledValues
456 (
457 sampleIndex_.second(),
458 sampleValues_.second(),
459 sampleAverage_.second()
460 );
461 }
462 }
463}
464
465
466template<class Type>
468{
469 if (this->updated())
470 {
471 return;
472 }
473
474 // Current time value
475 const scalar x = this->db().time().value();
476
477 checkTable(x);
478
479 // Interpolate between the sampled data
480 auto& fld = static_cast<Field<Type>&>(*this);
481 Type wantedAverage;
482
483 if (sampleIndex_.second() == -1)
484 {
485 // Only start value
486 fld = sampleValues_.first();
487 wantedAverage = sampleAverage_.first();
488 }
489 else
490 {
491 const scalar beg = sampleTimes_[sampleIndex_.first()].value();
492 const scalar end = sampleTimes_[sampleIndex_.second()].value();
493 const scalar s = (x - beg)/(end - beg);
494
495 fld = (1 - s)*sampleValues_.first() + s*sampleValues_.second();
496
497 wantedAverage
498 = (1 - s)*sampleAverage_.first() + s*sampleAverage_.second();
499
501 << "Sampled, interpolated values"
502 << " between time:"
503 << sampleTimes_[sampleIndex_.first()].name()
504 << " and time:" << sampleTimes_[sampleIndex_.second()].name()
505 << " with weight:" << s << endl;
506 }
507
508 // Enforce average. Either by scaling (if scaling factor > 0.5) or by
509 // offsetting.
510 if (setAverage_)
511 {
512 Type averagePsi = gAverage(fld);
513
514 if (debug)
515 {
516 Pout<< "updateCoeffs :"
517 << " actual average:" << averagePsi
518 << " wanted average:" << wantedAverage
519 << endl;
520 }
521
522 if (mag(averagePsi) < VSMALL)
523 {
524 // Field too small to scale. Offset instead.
525 const Type offset = wantedAverage - averagePsi;
526 if (debug)
527 {
528 Pout<< "updateCoeffs :"
529 << " offsetting with:" << offset << endl;
530 }
531 fld += offset;
532 }
533 else
534 {
535 const scalar scale = mag(wantedAverage)/mag(averagePsi);
536
537 if (debug)
538 {
539 Pout<< "updateCoeffs :"
540 << " scaling with:" << scale << endl;
541 }
542 fld *= scale;
543 }
544 }
545
546 // Apply offset to mapped values
547 if (offset_)
548 {
549 const scalar t = this->db().time().timeOutputValue();
550 fld += offset_->value(t);
551 }
552
553 if (debug)
554 {
555 auto limits = gMinMax(*this);
556 auto avg = gAverage(*this);
557
558 Pout<< "updateCoeffs : set fixedValue to min:" << limits.min()
559 << " max:" << limits.max()
560 << " avg:" << avg << endl;
562
564}
565
566
567template<class Type>
569(
570 Ostream& os
571) const
572{
574
575 os.writeEntryIfDifferent
576 (
577 "fieldTable",
578 this->internalField().name(),
579 fieldTableName_
580 );
581
582 if (!pointsName_.empty())
583 {
584 os.writeEntryIfDifferent<word>("points", "points", pointsName_);
585 }
586
587 if (!mapMethod_.empty() && !mapMethod_.starts_with("planar"))
588 {
589 os.writeEntry("mapMethod", mapMethod_);
590 }
591
592 if (setAverage_)
593 {
594 os.writeEntry("setAverage", setAverage_);
595 }
596
597 os.writeEntryIfDifferent<scalar>("perturb", 1e-5, perturb_);
598
599 if (offset_)
600 {
601 offset_->writeData(os);
602 }
603}
604
605
606// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition Field.C:465
static const Field< Type > & null() noexcept
Return a null Field (reference to a nullObject). Behaves like an empty Field.
Definition Field.H:192
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition Field.C:528
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition Ostream.H:331
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition Ostream.H:346
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
static instantList findTimes(const fileName &directory, const word &constantDirName="constant")
Search a given directory for valid time directories.
Definition TimePaths.C:109
fileName caseConstant() const
Return the constant name for the case, which is ../constant() for parallel runs.
Definition TimePathsI.H:143
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
iterator end() noexcept
Definition UListI.H:454
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
Definition UPstream.H:82
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A class for handling file names.
Definition fileName.H:75
static std::string path(const std::string &str)
Return directory path name (part before last /).
Definition fileNameI.H:169
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition fileNameI.H:192
A FixedValue boundary condition for pointField.
fixedValuePointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
static Pair< label > findRange(const UList< instant > &times, const scalar timeVal, const label start=-1)
Find lower/upper indices for given time value in list of instances (linear search) continuing after t...
Definition instant.C:54
const Time & time() const noexcept
Return time registry.
const objectRegistry & db() const
The associated objectRegistry.
bool updated() const noexcept
True if the boundary condition has already been updated.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
const DimensionedField< Type, pointMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field, sets updated() to false.
Basic pointPatch represents a set of points from the mesh.
Definition pointPatch.H:67
Interpolates between two sets of unstructured points using 2D Delaunay triangulation....
static wordList timeNames(const instantList &times)
Helper: extract words of times.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
Like IOField but falls back to raw IFstream if no header found. Optionally reads average value....
Definition rawIOField.H:52
bool starts_with(char c) const
True if string starts with given character (cf. C++20).
Definition string.H:436
A time-varying form of a mapped fixed value boundary condition.
timeVaryingMappedFixedValuePointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual autoPtr< pointPatchField< Type > > clone() const
Return a clone.
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given PointPatchField onto.
A class for managing temporary objects.
Definition tmp.H:75
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
rDeltaTY field()
auto limits
Definition setRDeltaT.H:186
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const auto & io
auto & name
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define DebugInFunction
Report an information message using Foam::Info.
const char * end
Definition SVGTools.H:223
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
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
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
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.
List< label > labelList
A List of labels.
Definition List.H:62
vectorIOField pointIOField
pointIOField is a vectorIOField.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
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...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
volScalarField & e
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))