uniformInterpolatedDisplacementPointPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "pointFields.H"
29 #include "Time.H"
30 #include "polyMesh.H"
31 #include "interpolationWeights.H"
32 #include "uniformInterpolate.H"
33 #include "ReadFields.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
44 (
45  const pointPatch& p,
47 )
48 :
50 {}
51 
52 
55 (
56  const pointPatch& p,
58  const dictionary& dict
59 )
60 :
62  fieldName_(dict.lookup("field")),
63  interpolationScheme_(dict.lookup("interpolationScheme"))
64 {
65  const pointMesh& pMesh = this->internalField().mesh();
66 
67  // Read time values
68  instantList allTimes = Time::findTimes(pMesh().time().path());
69 
70  // Only keep those that contain the field
71  DynamicList<word> names(allTimes.size());
72  DynamicList<scalar> values(allTimes.size());
73 
74  forAll(allTimes, i)
75  {
76  IOobject io
77  (
78  fieldName_,
79  allTimes[i].name(),
80  pMesh(),
83  false
84  );
85  if (io.typeHeaderOk<pointVectorField>(false))
86  {
87  names.append(allTimes[i].name());
88  values.append(allTimes[i].value());
89  }
90  }
91  timeNames_.transfer(names);
92  timeVals_.transfer(values);
93 
94  Info<< type() << " : found " << fieldName_ << " for times "
95  << timeNames_ << endl;
96 
97  if (timeNames_.size() < 1)
98  {
100  << "Did not find any times with " << fieldName_
101  << exit(FatalError);
102  }
103 
104 
105  if (!dict.found("value"))
106  {
107  updateCoeffs();
108  }
109 }
110 
111 
114 (
116  const pointPatch& p,
118  const pointPatchFieldMapper& mapper
119 )
120 :
121  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
122  fieldName_(ptf.fieldName_),
123  interpolationScheme_(ptf.interpolationScheme_),
124  timeNames_(ptf.timeNames_),
125  timeVals_(ptf.timeVals_),
126  interpolatorPtr_(ptf.interpolatorPtr_)
127 {}
128 
129 
132 (
135 )
136 :
138  fieldName_(ptf.fieldName_),
139  interpolationScheme_(ptf.interpolationScheme_),
140  timeNames_(ptf.timeNames_),
141  timeVals_(ptf.timeVals_),
142  interpolatorPtr_(ptf.interpolatorPtr_)
143 {}
144 
145 
146 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
147 
149 {
150  if (this->updated())
151  {
152  return;
153  }
154 
155  if (!interpolatorPtr_.valid())
156  {
157  interpolatorPtr_ = interpolationWeights::New
158  (
159  interpolationScheme_,
160  timeVals_
161  );
162  }
163 
164  const pointMesh& pMesh = this->internalField().mesh();
165  const Time& t = pMesh().time();
166 
167  // Update indices of times and weights
168  bool timesChanged = interpolatorPtr_->valueWeights
169  (
170  t.timeOutputValue(),
171  currentIndices_,
172  currentWeights_
173  );
174 
175  const wordList currentTimeNames
176  (
177  UIndirectList<word>(timeNames_, currentIndices_)
178  );
179 
180 
181  // Load if necessary fields for this interpolation
182  if (timesChanged)
183  {
184  objectRegistry& fieldsCache = const_cast<objectRegistry&>
185  (
186  pMesh.thisDb().subRegistry("fieldsCache", true)
187  );
188  // Save old times so we now which ones have been loaded and need
189  // 'correctBoundaryConditions'. Bit messy.
190  HashSet<word> oldTimes(fieldsCache.toc());
191 
192  ReadFields<pointVectorField>
193  (
194  fieldName_,
195  pMesh,
196  currentTimeNames
197  );
198 
199  forAllConstIter(objectRegistry, fieldsCache, fieldsCacheIter)
200  {
201  if (!oldTimes.found(fieldsCacheIter.key()))
202  {
203  // Newly loaded fields. Make sure the internal
204  // values are consistent with the boundary conditions.
205  // This is quite often not the case since these
206  // fields typically are constructed 'by hand'
207 
208  const objectRegistry& timeCache = dynamic_cast
209  <
210  const objectRegistry&
211  >(*fieldsCacheIter());
212 
213  timeCache.lookupObjectRef<pointVectorField>(fieldName_)
215  }
216  }
217  }
218 
219 
220  // Interpolate the whole field
221  pointVectorField result
222  (
223  uniformInterpolate<pointVectorField>
224  (
225  IOobject
226  (
227  word("uniformInterpolate(")
228  + this->internalField().name()
229  + ')',
230  pMesh.time().timeName(),
231  pMesh.thisDb(),
234  ),
235  fieldName_,
236  currentTimeNames,
237  currentWeights_
238  )
239  );
240 
241 
242  // Extract back from the internal field
243  this->operator==
244  (
245  this->patchInternalField(result())
246  );
247 
249 }
250 
251 
253 const
254 {
256  os.writeKeyword("field")
257  << fieldName_ << token::END_STATEMENT << nl;
258  os.writeKeyword("interpolationScheme")
259  << interpolationScheme_ << token::END_STATEMENT << nl;
260  writeEntry("value", os);
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 
267 (
270 );
271 
272 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
273 
274 } // End namespace Foam
275 
276 // ************************************************************************* //
static autoPtr< interpolationWeights > New(const word &type, const scalarField &samples)
Return a reference to the selected interpolationWeights.
A HashTable with keys but without contents.
Definition: HashSet.H:59
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
bool typeHeaderOk(const bool checkType=true)
Read header (uses typeFilePath to find file) and check header.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Type & lookupObjectRef(const word &name) const
Lookup and return the object reference of the given Type.
Foam::pointPatchFieldMapper.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:726
Field reading functions for post-processing utilities.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
A class for handling words, derived from string.
Definition: word.H:59
uniformInterpolatedDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:118
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
U correctBoundaryConditions()
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
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:53
static const char nl
Definition: Ostream.H:265
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const DimensionedField< vector, pointMesh > & internalField() const
Return dimensioned internal field reference.
static instantList findTimes(const fileName &, const word &constantName="constant")
Search a given directory for valid time directories.
Definition: findTimes.C:36
virtual void write(Ostream &) const
Write.
const Mesh & mesh() const
Definition: MeshObject.H:148
scalar timeOutputValue() const
Return current time value.
Definition: TimeStateI.H:29
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
const Time & time() const
Return time.
Definition: IOobject.C:367
A List with indirect addressing.
Definition: fvMatrix.H:106
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
bool updated() const
Return true if the boundary condition has already been updated.
messageStream Info
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
const objectRegistry & subRegistry(const word &name, const bool forceCreate=false) const
Lookup and return a const sub-objectRegistry. Optionally create.
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:202
volScalarField & p
Registry of regIOobjects.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576