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-2021 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  {
56  (
57  fieldName_,
58  allTimes[i].name(),
59  pMesh(),
60  IOobject::MUST_READ,
61  IOobject::NO_WRITE,
62  false
63  );
64 
65  if (io.headerOk())
66  {
67  names.append(allTimes[i].name());
68  values.append(allTimes[i].value());
69  }
70  }
71  timeNames_.transfer(names);
72  timeVals_.transfer(values);
73 
74  Info<< mesh_.type() << " : found " << fieldName_ << " for times "
75  << timeNames_ << endl;
76 
77  if (timeNames_.size() < 1)
78  {
80  << "Did not find any times with " << fieldName_
81  << exit(FatalError);
82  }
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
87 
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
96 {
97  if (!interpolatorPtr_.valid())
98  {
99  interpolatorPtr_ = interpolationWeights::New
100  (
101  interpolationScheme_,
102  timeVals_
103  );
104  }
105 
106  const pointMesh& pMesh = pointMesh::New(mesh_);
107 
108  const Time& t = pMesh().time();
109 
110  // Update indices of times and weights
111  bool timesChanged
112  (
113  interpolatorPtr_->valueWeights
114  (
115  t.userTimeValue(),
116  currentIndices_,
117  currentWeights_
118  )
119  );
120 
121  const wordList currentTimeNames
122  (
123  UIndirectList<word>(timeNames_, currentIndices_)
124  );
125 
126 
127  // Load if necessary fields for this interpolation
128  if (timesChanged)
129  {
130  objectRegistry& fieldsCache = const_cast<objectRegistry&>
131  (
132  pMesh.thisDb().subRegistry("fieldsCache", true)
133  );
134 
135  // Save old times so we now which ones have been loaded and need
136  // 'correctBoundaryConditions'. Bit messy.
137  HashSet<word> oldTimes(fieldsCache.toc());
138 
139  ReadFields<pointVectorField>
140  (
141  fieldName_,
142  pMesh,
143  currentTimeNames
144  );
145 
146  forAllConstIter(objectRegistry, fieldsCache, fieldsCacheIter)
147  {
148  if (!oldTimes.found(fieldsCacheIter.key()))
149  {
150  // Newly loaded fields. Make sure the internal
151  // values are consistent with the boundary conditions.
152  // This is quite often not the case since these
153  // fields typically are constructed 'by hand'
154 
155  const objectRegistry& timeCache = dynamic_cast
156  <
157  const objectRegistry&
158  >(*fieldsCacheIter());
159 
160  timeCache.lookupObjectRef<pointVectorField>(fieldName_)
162  }
163  }
164  }
165 
166 
167  // Interpolate the point field
168  return
169  (
170  uniformInterpolate<pointVectorField>
171  (
172  IOobject
173  (
174  word("uniformInterpolate(") + fieldName_ + ')',
175  pMesh.time().timeName(),
176  pMesh.thisDb(),
179  ),
180  fieldName_,
181  currentTimeNames,
182  currentWeights_
183  )
184  );
185 }
186 
187 
189 {
190  writeEntry(os, "field", fieldName_);
191  writeEntry(os, "interpolationScheme", interpolationScheme_);
192 }
193 
194 
195 // ************************************************************************* //
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
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
rDeltaT correctBoundaryConditions()
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Templated form of IOobject providing type information for file reading and header type checking...
Definition: IOobject.H:537
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.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
scalar userTimeValue() const
Return current user time value.
Definition: Time.C:845
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.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
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:109
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
const Time & time() const
Return time.
Definition: IOobject.C:318
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:76
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:98
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864