Loading...
Searching...
No Matches
uniformInterpolatedDisplacementPointPatchVectorField.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-2016 OpenFOAM Foundation
9 Copyright (C) 2019-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 "pointFields.H"
32#include "Time.H"
33#include "polyMesh.H"
35#include "uniformInterpolate.H"
36#include "ReadFields.H"
37
38// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39
40namespace Foam
41{
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
47(
48 const pointPatch& p,
66 interpolationScheme_(dict.lookup("interpolationScheme"))
67{
84 fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
85 fieldName_(ptf.fieldName_),
86 interpolationScheme_(ptf.interpolationScheme_),
87 timeNames_(ptf.timeNames_),
88 timeVals_(ptf.timeVals_),
89 interpolatorPtr_(ptf.interpolatorPtr_)
90{}
91
92
95(
98)
99:
101 fieldName_(ptf.fieldName_),
102 interpolationScheme_(ptf.interpolationScheme_),
103 timeNames_(ptf.timeNames_),
104 timeVals_(ptf.timeVals_),
105 interpolatorPtr_(ptf.interpolatorPtr_)
106{}
107
108
109// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110
112{
113 if (this->updated())
114 {
115 return;
116 }
117
118 if (!interpolatorPtr_)
119 {
120 const pointMesh& pMesh = this->internalField().mesh();
121
122 // Read time values
123 const instantList allTimes = Time::findTimes(pMesh().time().path());
124
125 // Only keep those that contain the field
126 DynamicList<word> names(allTimes.size());
127 DynamicList<scalar> values(allTimes.size());
128
129 for (const instant& inst : allTimes)
130 {
131 IOobject io
132 (
133 fieldName_,
134 inst.name(),
135 pMesh(),
139 );
140 if (io.typeHeaderOk<pointVectorField>(false))
141 {
142 names.append(inst.name());
143 values.append(inst.value());
144 }
145 }
146 timeNames_.transfer(names);
147 timeVals_.transfer(values);
148
149 Info<< type() << " : found " << fieldName_ << " for times "
150 << timeNames_ << endl;
151
152 if (timeNames_.size() < 1)
153 {
155 << "Did not find any times with " << fieldName_
156 << exit(FatalError);
157 }
158
159
160 interpolatorPtr_ = interpolationWeights::New
161 (
162 interpolationScheme_,
163 timeVals_
164 );
165 }
166
167 const pointMesh& pMesh = this->internalField().mesh();
168 const Time& t = pMesh().time();
169
170 // Update indices of times and weights
171 bool timesChanged = interpolatorPtr_->valueWeights
172 (
173 t.timeOutputValue(),
174 currentIndices_,
175 currentWeights_
176 );
177
178 const wordList currentTimeNames
179 (
180 UIndirectList<word>(timeNames_, currentIndices_)
181 );
182
183
184 // Load if necessary fields for this interpolation
185 if (timesChanged)
186 {
187 objectRegistry& fieldsCache = const_cast<objectRegistry&>
188 (
189 pMesh.thisDb().subRegistry("fieldsCache", true)
190 );
191 // Save old times so we know which ones have been loaded and need
192 // 'correctBoundaryConditions'. Bit messy.
193 wordHashSet oldTimes(fieldsCache.toc());
194
196 (
197 fieldName_,
198 pMesh,
199 currentTimeNames
200 );
201
202 forAllConstIters(fieldsCache, fieldsCacheIter)
203 {
204 if (!oldTimes.found(fieldsCacheIter.key()))
205 {
206 // Newly loaded fields. Make sure the internal
207 // values are consistent with the boundary conditions.
208 // This is quite often not the case since these
209 // fields typically are constructed 'by hand'
210
211 const objectRegistry& timeCache = dynamic_cast
212 <
213 const objectRegistry&
214 >(*fieldsCacheIter());
215
216
218 timeCache.lookupObjectRef<pointVectorField>(fieldName_);
219
221 }
222 }
223 }
224
225
226 // Interpolate the whole field
227 pointVectorField result
228 (
230 (
231 IOobject
232 (
233 word("uniformInterpolate(")
234 + this->internalField().name()
235 + ')',
236 pMesh.time().timeName(),
237 pMesh.thisDb(),
240 ),
241 fieldName_,
242 currentTimeNames,
243 currentWeights_
244 )
245 );
246
247
248 // Extract back from the internal field
249 this->operator==
250 (
259const
260{
262 os.writeEntry("field", fieldName_);
263 os.writeEntry("interpolationScheme", interpolationScheme_);
264 this->writeValueEntry(os);
265}
266
267
268// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269
271(
273 uniformInterpolatedDisplacementPointPatchVectorField
274);
275
276// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277
278} // End namespace Foam
279
280// ************************************************************************* //
Field reading functions for post-processing utilities.
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const noexcept
Return const reference to mesh.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void correctBoundaryConditions()
Correct boundary field.
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition HashTable.C:141
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static instantList findTimes(const fileName &directory, const word &constantDirName="constant")
Search a given directory for valid time directories.
Definition TimePaths.C:109
scalar timeOutputValue() const
Return the current user-time value. (ie, after applying any timeToUserTime() conversion).
Definition TimeStateI.H:24
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FixedValue boundary condition for pointField.
fixedValuePointPatchField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name.
Definition instant.H:56
static autoPtr< interpolationWeights > New(const word &type, const scalarField &samples)
Return a reference to the selected interpolationWeights.
Registry of regIOobjects.
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Lookup and return non-const reference to the object of the given Type. Fatal if not found or the wron...
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
const Time & time() const
Return Time from polyMesh.
Definition pointMesh.H:209
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition pointMesh.H:201
bool updated() const noexcept
True if the boundary condition has already been updated.
Foam::pointPatchFieldMapper.
tmp< Field< vector > > patchInternalField() const
virtual void write(Ostream &os) const
Write.
const DimensionedField< vector, pointMesh > & internalField() const noexcept
Basic pointPatch represents a set of points from the mesh.
Definition pointPatch.H:67
Lookup type of boundary radiation properties.
Definition lookup.H:60
uniformInterpolatedDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from Foam::string.
Definition word.H:66
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const auto & io
auto & name
auto & names
Namespace for OpenFOAM.
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
List< word > wordList
List of word.
Definition fileName.H:60
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition HashSet.H:80
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
List< instant > instantList
List of instants.
Definition instantList.H:41
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
tmp< GeoField > uniformInterpolate(const HashPtrTable< GeoField, label, Hash< label > > &fields, const labelList &indices, const scalarField &weights)
Interpolate selected fields (given by indices and corresponding weights).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
pointPatchField< vector > pointPatchVectorField
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete pointPatchField type and add to run-time tables Example, (pointPatchScalarField,...
dictionary dict
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235