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 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  const pointMesh& pMesh = this->internalField().mesh();
69 
70  // Read time values
71  instantList allTimes = Time::findTimes(pMesh().time().path());
72 
73  // Only keep those that contain the field
74  DynamicList<word> names(allTimes.size());
75  DynamicList<scalar> values(allTimes.size());
76 
77  for (const instant& inst : allTimes)
78  {
79  IOobject io
80  (
81  fieldName_,
82  inst.name(),
83  pMesh(),
86  false
87  );
88  if (io.typeHeaderOk<pointVectorField>(false))
89  {
90  names.append(inst.name());
91  values.append(inst.value());
92  }
93  }
94  timeNames_.transfer(names);
95  timeVals_.transfer(values);
96 
97  Info<< type() << " : found " << fieldName_ << " for times "
98  << timeNames_ << endl;
99 
100  if (timeNames_.size() < 1)
101  {
103  << "Did not find any times with " << fieldName_
104  << exit(FatalError);
105  }
106 
107 
108  if (!dict.found("value"))
109  {
110  updateCoeffs();
111  }
112 }
113 
114 
117 (
119  const pointPatch& p,
121  const pointPatchFieldMapper& mapper
122 )
123 :
124  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
125  fieldName_(ptf.fieldName_),
126  interpolationScheme_(ptf.interpolationScheme_),
127  timeNames_(ptf.timeNames_),
128  timeVals_(ptf.timeVals_),
129  interpolatorPtr_(ptf.interpolatorPtr_)
130 {}
131 
132 
135 (
138 )
139 :
141  fieldName_(ptf.fieldName_),
142  interpolationScheme_(ptf.interpolationScheme_),
143  timeNames_(ptf.timeNames_),
144  timeVals_(ptf.timeVals_),
145  interpolatorPtr_(ptf.interpolatorPtr_)
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
152 {
153  if (this->updated())
154  {
155  return;
156  }
157 
158  if (!interpolatorPtr_)
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  writeEntry("value", os);
265 }
266 
267 
268 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269 
271 (
274 );
275 
276 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277 
278 } // End namespace Foam
279 
280 // ************************************************************************* //
Foam::TimePaths::findTimes
static instantList findTimes(const fileName &directory, const word &constantName="constant")
Search a given directory for valid time directories.
Definition: TimePaths.C:140
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::TimeState::timeOutputValue
scalar timeOutputValue() const
Return current time value.
Definition: TimeStateI.H:31
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::HashTable::toc
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:121
Foam::DynamicList< word >
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::pointMesh::time
const Time & time() const
Return Time from polyMesh.
Definition: pointMesh.H:127
Foam::pointPatch
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:58
polyMesh.H
Foam::HashSet< word, Hash< word > >
Foam::pointPatchField
Abstract base class for point-mesh patch fields.
Definition: pointMVCWeight.H:60
Foam::pointPatchFieldMapper
Foam::pointPatchFieldMapper.
Definition: pointPatchFieldMapper.H:48
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::fixedValuePointPatchField< vector >
Foam::makePointPatchTypeField
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
Foam::uniformInterpolatedDisplacementPointPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.C:258
Foam::valuePointPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: valuePointPatchField.C:135
uniformInterpolate.H
Foam::pointPatchField< vector >::internalField
const DimensionedField< vector, pointMesh > & internalField() const
Return dimensioned internal field reference.
Definition: pointPatchField.H:275
Foam::MeshObject::mesh
const Mesh & mesh() const
Definition: MeshObject.H:122
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::pointPatchField::write
virtual void write(Ostream &) const
Write.
Definition: pointPatchField.C:117
interpolationWeights.H
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:386
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
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:47
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::pointPatchField< vector >::updated
bool updated() const
Return true if the boundary condition has already been updated.
Definition: pointPatchField.H:311
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::Field< vector >::writeEntry
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:608
Foam::uniformInterpolatedDisplacementPointPatchVectorField
Interpolates pre-specified motion.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.H:72
Foam::objectRegistry::lookupObjectRef
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:478
Time.H
uniformInterpolatedDisplacementPointPatchVectorField.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
ReadFields.H
Field reading functions for post-processing utilities.
Foam::interpolationWeights::New
static autoPtr< interpolationWeights > New(const word &type, const scalarField &samples)
Return a reference to the selected interpolationWeights.
Definition: interpolationWeights.C:53
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::objectRegistry::subRegistry
const objectRegistry & subRegistry(const word &name, const bool forceCreate=false, const bool recursive=false) const
Lookup and return a const sub-objectRegistry.
Definition: objectRegistry.C:174
Foam::List< instant >
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::uniformInterpolatedDisplacementPointPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: uniformInterpolatedDisplacementPointPatchVectorField.C:151
Foam::pointMesh::thisDb
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:121
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::pointPatchField< vector >::patchInternalField
tmp< Field< vector > > patchInternalField() const
Return field created from appropriate internal field values.
Definition: pointPatchField.C:130
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::instant
An instant of time. Contains the time value and name.
Definition: instant.H:52
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::GeometricField< vector, pointPatchField, pointMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::PtrListOps::names
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
pointFields.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::IOobject::MUST_READ
Definition: IOobject.H:185