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 -------------------------------------------------------------------------------
11 License
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"
34 #include "interpolationWeights.H"
35 #include "uniformInterpolate.H"
36 #include "ReadFields.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
47 (
48  const pointPatch& p,
50 )
51 :
53 {}
54 
55 
58 (
59  const pointPatch& p,
61  const dictionary& dict
62 )
63 :
65  fieldName_(dict.lookup("field")),
66  interpolationScheme_(dict.lookup("interpolationScheme"))
67 {
68  if (!dict.found("value"))
69  {
71  }
72 }
73 
74 
77 (
78  const uniformInterpolatedDisplacementPointPatchVectorField& ptf,
79  const pointPatch& p,
80  const DimensionedField<vector, pointMesh>& iF,
81  const pointPatchFieldMapper& mapper
82 )
83 :
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 
195  ReadFields<pointVectorField>
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 
217  pointVectorField& d =
218  timeCache.lookupObjectRef<pointVectorField>(fieldName_);
219 
221  }
222  }
223  }
224 
225 
226  // Interpolate the whole field
227  pointVectorField result
228  (
229  uniformInterpolate<pointVectorField>
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  (
251  this->patchInternalField(result())
252  );
253 
255 }
256 
257 
259 const
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 // ************************************************************************* //
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
List of names generated by calling name() for each list item and filtered for matches.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
A FixedValue boundary condition for pointField.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
virtual void write(Ostream &os) const
Write.
Field reading functions for post-processing utilities.
Ignore writing from objectRegistry::writeObject()
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
Macros for easy insertion into run-time selection tables.
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
uniformInterpolatedDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Vector< scalar > vector
Definition: vector.H:57
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:73
bool updated() const noexcept
True if the boundary condition has already been updated.
tmp< Field< vector > > patchInternalField() const
Return field created from appropriate internal field values.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
pointPatchField< vector > pointPatchVectorField
static autoPtr< interpolationWeights > New(const word &type, const scalarField &samples)
Return a reference to the selected interpolationWeights.
static instantList findTimes(const fileName &directory, const word &constantDirName="constant")
Search a given directory for valid time directories.
Definition: TimePaths.C:109
List< word > wordList
List of word.
Definition: fileName.H:59
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:61
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
Nothing to be read.
Automatically write from objectRegistry::writeObject()
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
List< instant > instantList
List of instants.
Definition: instantList.H:41
Do not request registration (bool: false)
const DimensionedField< vector, pointMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28