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 | Copyright (C) 2012-2015 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("fieldName")),
63  interpolationScheme_(dict.lookup("interpolationScheme"))
64 {
65  const pointMesh& pMesh = this->dimensionedInternalField().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.headerOk())
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  {
158  (
160  timeVals_
161  );
162  }
163 
164  const pointMesh& pMesh = this->dimensionedInternalField().mesh();
165  const Time& t = pMesh().time();
166 
167  // Update indices of times and weights
168  bool timesChanged = interpolatorPtr_->valueWeights
169  (
170  t.timeOutputValue(),
173  );
174 
175  const wordList currentTimeNames
176  (
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 
214  pointVectorField& d = const_cast<pointVectorField&>
215  (
216  timeCache.lookupObject<pointVectorField>
217  (
218  fieldName_
219  )
220  );
222  }
223  }
224  }
225 
226 
227  // Interpolate the whole field
228  pointVectorField result
229  (
230  uniformInterpolate<pointVectorField>
231  (
232  IOobject
233  (
234  word("uniformInterpolate(")
235  + this->dimensionedInternalField().name()
236  + ')',
237  pMesh.time().timeName(),
238  pMesh.thisDb(),
241  ),
242  fieldName_,
243  currentTimeNames,
245  )
246  );
247 
248 
249  // Extract back from the internal field
250  this->operator==
251  (
253  );
254 
256 }
257 
258 
260 const
261 {
263  os.writeKeyword("fieldName")
265  os.writeKeyword("interpolationScheme")
267  writeEntry("value", os);
268 }
269 
270 
271 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 
274 (
277 );
278 
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280 
281 } // End namespace Foam
282 
283 // ************************************************************************* //
Foam::TimeState::timeOutputValue
scalar timeOutputValue() const
Return current time value.
Definition: TimeState.C:67
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Foam::uniformInterpolatedDisplacementPointPatchVectorField::currentWeights_
scalarField currentWeights_
Cached interpolation weights.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.H:91
Foam::GeometricField::dimensionedInternalField
DimensionedInternalField & dimensionedInternalField()
Return dimensioned internal field.
Definition: GeometricField.C:713
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::interpolationWeights::New
static autoPtr< interpolationWeights > New(const word &type, const scalarField &samples)
Return a reference to the selected interpolationWeights.
Definition: interpolationWeights.C:55
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::HashTable::toc
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:201
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynamicList< word >
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::uniformInterpolatedDisplacementPointPatchVectorField::interpolationScheme_
const word interpolationScheme_
Definition: uniformInterpolatedDisplacementPointPatchVectorField.H:75
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::pointPatchField< vector >::patchInternalField
tmp< Field< Type > > patchInternalField() const
Return field created from appropriate internal field values.
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::pointPatch
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
Foam::uniformInterpolatedDisplacementPointPatchVectorField::currentIndices_
labelList currentIndices_
Cached interpolation times.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.H:88
polyMesh.H
Foam::HashSet
A HashTable with keys but without contents.
Definition: HashSet.H:59
Foam::IOobject::time
const Time & time() const
Return time.
Definition: IOobject.C:245
Foam::pointPatchField
Abstract base class for point-mesh patch fields.
Definition: pointMVCWeight.H:58
Foam::uniformInterpolatedDisplacementPointPatchVectorField::timeNames_
wordList timeNames_
Times with pre-specified displacement.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.H:78
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::pointPatchFieldMapper
Foam::pointPatchFieldMapper.
Definition: pointPatchFieldMapper.H:46
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::IOobject::headerOk
bool headerOk()
Read and check header info.
Definition: IOobject.C:439
Foam::objectRegistry::subRegistry
const objectRegistry & subRegistry(const word &name, const bool forceCreate=false) const
Lookup and return a const sub-objectRegistry. Optionally create.
Definition: objectRegistry.C:156
Foam::fixedValuePointPatchField< vector >
Foam::makePointPatchTypeField
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::uniformInterpolatedDisplacementPointPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.C:259
Foam::valuePointPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: valuePointPatchField.C:150
uniformInterpolate.H
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::pointPatchField::write
virtual void write(Ostream &) const
Write.
Definition: pointPatchField.C:119
interpolationWeights.H
Foam::uniformInterpolatedDisplacementPointPatchVectorField::fieldName_
const word fieldName_
Name of displacement field.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.H:73
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::uniformInterpolatedDisplacementPointPatchVectorField::uniformInterpolatedDisplacementPointPatchVectorField
uniformInterpolatedDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.C:44
Foam::HashTable< nil, word, string::hash >::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::pointPatchField< vector >::updated
bool updated() const
Return true if the boundary condition has already been updated.
Definition: pointPatchField.H:314
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
Foam::Field< vector >::writeEntry
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:700
Foam::uniformInterpolatedDisplacementPointPatchVectorField
Interpolates pre-specified motion.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.H:66
uniformInterpolatedDisplacementPointPatchVectorField.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::uniformInterpolatedDisplacementPointPatchVectorField::interpolatorPtr_
autoPtr< interpolationWeights > interpolatorPtr_
User-specified interpolator.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.H:84
ReadFields.H
Helper routine to read fields.
Foam::uniformInterpolatedDisplacementPointPatchVectorField::timeVals_
scalarField timeVals_
Times with pre-specified displacement.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.H:81
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::Time::findTimes
static instantList findTimes(const fileName &, const word &constantName="constant")
Search a given directory for valid time directories.
Definition: findTimes.C:38
Foam::uniformInterpolatedDisplacementPointPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.C:148
Foam::pointMesh::thisDb
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:118
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::pointPatchField< vector >::dimensionedInternalField
const DimensionedField< Type, pointMesh > & dimensionedInternalField() const
Return dimensioned internal field reference.
Definition: pointPatchField.H:278
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
pointFields.H
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51