dynamicMeshPointInterpolator.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) 2018-2019 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"
28 #include "interpolationWeights.H"
29 #include "uniformInterpolate.H"
30 #include "ReadFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 (
36  const polyMesh& mesh,
37  const dictionary& dict
38 )
39 :
40  mesh_(mesh),
41  fieldName_(dict.lookup("field")),
42  interpolationScheme_(dict.lookup("interpolationScheme"))
43 {
44  const pointMesh& pMesh = pointMesh::New(mesh_);
45 
46  // Read time values
47  instantList allTimes = Time::findTimes(pMesh().time().path());
48 
49  // Only keep those that contain the field
50  DynamicList<word> names(allTimes.size());
51  DynamicList<scalar> values(allTimes.size());
52 
53  forAll(allTimes, i)
54  {
55  IOobject io
56  (
57  fieldName_,
58  allTimes[i].name(),
59  pMesh(),
60  IOobject::MUST_READ,
61  IOobject::NO_WRITE,
62  false
63  );
64  if (io.typeHeaderOk<pointVectorField>(false))
65  {
66  names.append(allTimes[i].name());
67  values.append(allTimes[i].value());
68  }
69  }
70  timeNames_.transfer(names);
71  timeVals_.transfer(values);
72 
73  Info<< mesh_.type() << " : found " << fieldName_ << " for times "
74  << timeNames_ << endl;
75 
76  if (timeNames_.size() < 1)
77  {
79  << "Did not find any times with " << fieldName_
80  << exit(FatalError);
81  }
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
86 
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
95 {
96  if (!interpolatorPtr_.valid())
97  {
98  interpolatorPtr_ = interpolationWeights::New
99  (
100  interpolationScheme_,
101  timeVals_
102  );
103  }
104 
105  const pointMesh& pMesh = pointMesh::New(mesh_);
106 
107  const Time& t = pMesh().time();
108 
109  // Update indices of times and weights
110  bool timesChanged
111  (
112  interpolatorPtr_->valueWeights
113  (
114  t.timeOutputValue(),
115  currentIndices_,
116  currentWeights_
117  )
118  );
119 
120  const wordList currentTimeNames
121  (
122  UIndirectList<word>(timeNames_, currentIndices_)
123  );
124 
125 
126  // Load if necessary fields for this interpolation
127  if (timesChanged)
128  {
129  objectRegistry& fieldsCache = const_cast<objectRegistry&>
130  (
131  pMesh.thisDb().subRegistry("fieldsCache", true)
132  );
133 
134  // Save old times so we now which ones have been loaded and need
135  // 'correctBoundaryConditions'. Bit messy.
136  HashSet<word> oldTimes(fieldsCache.toc());
137 
138  ReadFields<pointVectorField>
139  (
140  fieldName_,
141  pMesh,
142  currentTimeNames
143  );
144 
145  forAllConstIter(objectRegistry, fieldsCache, fieldsCacheIter)
146  {
147  if (!oldTimes.found(fieldsCacheIter.key()))
148  {
149  // Newly loaded fields. Make sure the internal
150  // values are consistent with the boundary conditions.
151  // This is quite often not the case since these
152  // fields typically are constructed 'by hand'
153 
154  const objectRegistry& timeCache = dynamic_cast
155  <
156  const objectRegistry&
157  >(*fieldsCacheIter());
158 
159  timeCache.lookupObjectRef<pointVectorField>(fieldName_)
161  }
162  }
163  }
164 
165 
166  // Interpolate the point field
167  return
168  (
169  uniformInterpolate<pointVectorField>
170  (
171  IOobject
172  (
173  word("uniformInterpolate(") + fieldName_ + ')',
174  pMesh.time().timeName(),
175  pMesh.thisDb(),
178  ),
179  fieldName_,
180  currentTimeNames,
181  currentWeights_
182  )
183  );
184 }
185 
186 
188 {
189  writeEntry(os, "field", fieldName_);
190  writeEntry(os, "interpolationScheme", interpolationScheme_);
191 }
192 
193 
194 // ************************************************************************* //
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
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool typeHeaderOk(const bool checkType=true)
Read header (uses typeFilePath to find file) and check header.
rDeltaT correctBoundaryConditions()
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:158
#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:164
Type & lookupObjectRef(const word &name) const
Lookup and return the object reference of the given Type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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:622
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
dynamicMeshPointInterpolator(const polyMesh &mesh, const dictionary &dict)
Construct from mesh and dictionary.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:115
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
scalar timeOutputValue() const
Return current time value.
Definition: TimeStateI.H:29
const Time & time() const
Return time.
Definition: IOobject.C:328
A List with indirect addressing.
Definition: fvMatrix.H:106
messageStream Info
const objectRegistry & subRegistry(const word &name, const bool forceCreate=false) const
Lookup and return a const sub-objectRegistry. Optionally create.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:202
tmp< pointVectorField > curPointField() const
Return interpolated pointField for the currentTime.
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812
fileName path(UMean.rootPath()/UMean.caseName()/functionObjects::writeFile::outputPrefix/"graphs"/UMean.instance())