CarrierField.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) 2025 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 
26 #include "CarrierField.H"
27 #include "coupled.H"
28 
29 /*---------------------------------------------------------------------------*\
30  Class CloudDerivedField::Functor Declaration
31 \*---------------------------------------------------------------------------*/
32 
33 template<class Type>
35 {
36 public:
37 
38  // Constructors
39 
40  //- Construct null
42  {}
43 
44 
45  //- Destructor
46  virtual ~Functor()
47  {}
48 
49 
50  // Member Operators
51 
52  //- Evaluate the field
53  virtual tmp<VolField<Type>> operator()() const = 0;
54 };
55 
56 
57 /*---------------------------------------------------------------------------*\
58  Class CloudDerivedField::Function Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 template<class Type>
62 template<class F>
64 :
65  public Functor
66 {
67  // Private Member Data
68 
69  //- The name
70  const word& name_;
71 
72  //- The function
73  F f_;
74 
75 
76 public:
77 
78  // Constructors
79 
80  //- Construct from a name and a function
81  Function(const word& name, const F& f)
82  :
83  Functor(),
84  name_(name),
85  f_(f)
86  {}
87 
88 
89  //- Destructor
90  virtual ~Function()
91  {}
92 
93 
94  // Member Operators
95 
96  //- Evaluate the field
97  virtual tmp<VolField<Type>> operator()() const
98  {
99  tmp<VolField<Type>> tpsi = f_();
100  return tpsi.isTmp() ? VolField<Type>::New(name_, tpsi()) : tpsi;
101  }
102 };
103 
104 
105 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
106 
107 template<class Type>
110 (
112  const LagrangianSubMesh& subMesh
113 ) const
114 {
115  const LagrangianMesh& mesh = subMesh.mesh();
116 
117  // Determine whether old-time values are available, and whether or not they
118  // should be used
119  const bool useOldTime =
120  psi().hasStoredOldTimes()
121  && psi().nOldTimes(false)
123  (
124  LagrangianMesh::fractionName
125  );
126 
127  /*
128  Info<< "--> Interpolation of carrier field " << name_ << ' '
129  << (useOldTime ? "*IS*" : "is *NOT*") << " using old-time values"
130  << endl;
131  */
132 
133  // Construct the interpolation on demand
134  if (!psiInterpolationPtr_.valid())
135  {
136  psiInterpolationPtr_ =
138  (
139  word(subMesh.mesh().schemes().interpolation(name())),
140  psi()
141  );
142  }
143 
144  // Construct the old-time interpolation on demand, if necessary
145  if (useOldTime && !psi0InterpolationPtr_.valid())
146  {
147  psi0InterpolationPtr_ =
149  (
150  word(subMesh.mesh().schemes().interpolation(name())),
151  psi().oldTime()
152  );
153  }
154 
155  // Interpolate the values
156  tmp<Field<Type>> tpsic =
157  psiInterpolationPtr_->interpolate
158  (
159  subMesh.sub(mesh.coordinates()),
160  subMesh.sub(mesh.celli()),
161  subMesh.sub(mesh.facei()),
162  subMesh.sub(mesh.faceTrii())
163  );
164 
165  // If using old-time values, interpolate them also, and combine with the
166  // current-time values using the current tracking fraction
167  if (useOldTime)
168  {
169  tmp<Field<Type>> tpsi0c =
170  psi0InterpolationPtr_->interpolate
171  (
172  subMesh.sub(mesh.coordinates()),
173  subMesh.sub(mesh.celli()),
174  subMesh.sub(mesh.facei()),
175  subMesh.sub(mesh.faceTrii())
176  );
177 
178  SubField<scalar> fraction =
179  subMesh.sub
180  (
181  mesh
182  .lookupObject<LagrangianScalarInternalDynamicField>
183  (
184  LagrangianMesh::fractionName
185  ).primitiveField()
186  );
187 
188  tpsic = (1 - fraction)*tpsi0c + fraction*tpsic;
189  }
190 
191  // Build the dimensioned field and return
192  return tmp<LagrangianSubField<Type>>
193  (
194  new LagrangianSubField<Type>
195  (
196  IOobject
197  (
198  name() + ':' + Foam::name(subMesh.group()),
199  mesh.time().name(),
200  mesh,
201  IOobject::NO_READ,
202  IOobject::NO_WRITE,
203  false
204  ),
205  subMesh,
206  psi().dimensions(),
207  tpsic
208  )
209  );
210 }
211 
212 
213 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
214 
215 template<class Type>
217 :
219  name_(clouds::coupled::carrierName(psi.name())),
220  functorPtr_(nullptr),
221  tpsi_(psi),
222  psiInterpolationPtr_(nullptr),
223  psi0InterpolationPtr_(nullptr)
224 {}
225 
226 
227 template<class Type>
228 template<class F>
230 :
232  name_(name),
233  functorPtr_(new Function<F>(name_, f)),
234  tpsi_(nullptr),
235  psiInterpolationPtr_(nullptr),
236  psi0InterpolationPtr_(nullptr)
237 {}
238 
239 
240 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
241 
242 template<class Type>
244 {
245  if (!tpsi_.valid())
246  {
247  tpsi_ = functorPtr_->operator()();
248  }
249 
250  return tpsi_();
251 }
252 
253 
254 template<class Type>
256 {
257  return name_;
258 }
259 
260 
261 template<class Type>
262 void Foam::CarrierField<Type>::reset(const bool predict)
263 {
265 
266  if (tpsi_.valid() && tpsi_.isTmp())
267  {
268  if (predict)
269  {
270  tpsi_.ref().nullOldestTime();
271  }
272  else
273  {
274  tpsi_.ref().oldTime();
275  tpsi_.ref() = functorPtr_->operator()();
276  }
277  }
278 
279  psiInterpolationPtr_.clear();
280  psi0InterpolationPtr_.clear();
281 }
282 
283 
284 // ************************************************************************* //
Class to store an evaluation function.
Definition: CarrierField.C:66
Function(const word &name, const F &f)
Construct from a name and a function.
Definition: CarrierField.C:81
virtual tmp< VolField< Type > > operator()() const
Evaluate the field.
Definition: CarrierField.C:97
virtual ~Function()
Destructor.
Definition: CarrierField.C:90
Functor()
Construct null.
Definition: CarrierField.C:41
virtual ~Functor()
Destructor.
Definition: CarrierField.C:46
virtual tmp< VolField< Type > > operator()() const =0
Evaluate the field.
A field interpolated from the carrier to the cloud. Uses CloudDerivedField to provide flexible access...
Definition: CarrierField.H:57
const VolField< Type > & psi() const
Access the carrier field.
Definition: CarrierField.C:243
void reset(const bool predict)
Reset.
Definition: CarrierField.C:262
const word & name() const
Carrier field name.
Definition: CarrierField.C:255
CarrierField(const VolField< Type > &)
Construct from a reference to a carrier field.
Definition: CarrierField.C:216
A field derived from other state fields of the cloud. Stores and virtualises a function or a method w...
void clear(const bool final)
Clear.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Class containing Lagrangian geometry and topology.
const LagrangianSchemes & schemes() const
Access the schemes.
Simple wrapper to provide an optional reference to a Lagrangian model.
ITstream & interpolation(const word &name) const
Get the interpolation scheme for the given field name.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
SubList< Type > sub(const List< Type > &list) const
Return a sub-list corresponding to this sub-mesh.
LagrangianGroup group() const
Return the group.
const LagrangianMesh & mesh() const
Return the mesh.
const word & name() const
Return const reference to name.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
bool foundObject(const word &name) const
Is the named Type in registry.
A class for managing temporary objects.
Definition: tmp.H:55
bool isTmp() const
Return true if this is really a temporary object.
Definition: tmpI.H:169
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const volScalarField & psi
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
labelList f(nPoints)