CloudDerivedField.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 "CloudDerivedField.H"
28 #include "Time.H"
29 
30 /*---------------------------------------------------------------------------*\
31  Class CloudDerivedField::Functor Declaration
32 \*---------------------------------------------------------------------------*/
33 
34 template<class Type>
36 {
37 public:
38 
39  // Constructors
40 
41  //- Construct null
43  {}
44 
45 
46  //- Destructor
47  virtual ~Functor()
48  {}
49 
50 
51  // Member Operators
52 
53  //- Evaluate
54  virtual tmp<LagrangianSubField<Type>> operator()
55  (
56  const LagrangianModelRef& model,
57  const LagrangianSubMesh& subMesh
58  ) const = 0;
59 };
60 
61 
62 /*---------------------------------------------------------------------------*\
63  Class CloudDerivedField::Function Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 template<class Type>
67 template<class F>
69 :
70  public Functor
71 {
72  // Private Member Data
73 
74  //- The function
75  F f_;
76 
77 
78 public:
79 
80  // Constructors
81 
82  //- Construct from a function
83  Function(const F& f)
84  :
85  Functor(),
86  f_(f)
87  {}
88 
89 
90  //- Destructor
91  virtual ~Function()
92  {}
93 
94 
95  // Member Operators
96 
97  //- Evaluate the field
98  virtual tmp<LagrangianSubField<Type>> operator()
99  (
100  const LagrangianModelRef& model,
101  const LagrangianSubMesh& subMesh
102  ) const
103  {
104  return f_(model, subMesh);
105  }
106 };
107 
108 
109 /*---------------------------------------------------------------------------*\
110  Class CloudDerivedField::Method Declaration
111 \*---------------------------------------------------------------------------*/
112 
113 template<class Type>
114 template<class C>
116 :
117  public Functor
118 {
119  // Private Member Data
120 
121  //- The class
122  const C& c_;
123 
124  //- The method
126  (
127  const LagrangianModelRef&,
128  const LagrangianSubMesh&
129  ) const;
130 
131 
132 public:
133 
134  // Constructors
135 
136  //- Construct from a class and a method
138  (
139  const C& c,
141  (
142  const LagrangianModelRef&,
143  const LagrangianSubMesh&
144  ) const
145  )
146  :
147  Functor(),
148  c_(c),
149  m_(m)
150  {}
151 
152 
153  //- Destructor
154  virtual ~Method()
155  {}
156 
157 
158  // Member Operators
159 
160  //- Evaluate the field
161  virtual tmp<LagrangianSubField<Type>> operator()
162  (
163  const LagrangianModelRef& model,
164  const LagrangianSubMesh& subMesh
165  ) const
166  {
167  return (c_.*m_)(model, subMesh);
168  }
169 };
170 
171 
172 /*---------------------------------------------------------------------------*\
173  Class CloudDerivedField::AllFieldToField Declaration
174 \*---------------------------------------------------------------------------*/
175 
176 template<class Type>
178 :
179  public LagrangianInternalField<Type>
180 {
181 public:
182 
183  // Constructors
184 
185  //- Construct from a sub-all-field reference
187  :
189  (
190  IOobject
191  (
192  allField.name(),
193  allField.mesh().mesh().time().name(),
194  allField.mesh().mesh(),
195  IOobject::NO_READ,
196  IOobject::NO_WRITE,
197  false
198  ),
199  allField.mesh().mesh(),
200  allField.dimensions(),
201  NullObjectRef<Field<Type>>()
202  )
203  {
204  this->UList<Type>::shallowCopy(allField);
205  }
206 
207 
208  //- Destructor
210  {
211  this->UList<Type>::shallowCopy(UList<Type>(nullptr, 0));
212  }
213 };
214 
215 
216 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
217 
218 template<class Type>
219 template<class F>
221 :
222  name_(name),
223  functorPtr_(new Function<F>(f))
224 {}
225 
226 
227 template<class Type>
228 template<class F>
230 :
231  name_(NullObjectRef<word>()),
232  functorPtr_(new Function<F>(f))
233 {}
234 
235 
236 template<class Type>
237 template<class C>
239 (
240  const word& name,
241  const C& c,
243  (
244  const LagrangianModelRef&,
245  const LagrangianSubMesh&
246  ) const
247 )
248 :
249  name_(name),
250  functorPtr_(new Method<C>(c, m))
251 {}
252 
253 
254 template<class Type>
255 template<class C>
257 (
258  const C& c,
260  (
261  const LagrangianModelRef&,
262  const LagrangianSubMesh&
263  ) const
264 )
265 :
266  name_(NullObjectRef<word>()),
267  functorPtr_(new Method<C>(c, m))
268 {}
269 
270 
271 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
272 
273 template<class Type>
276 {
277  const LagrangianSubField<Type>& allField =
278  operator()(LagrangianModelRef(), mesh.subAll());
279 
282  (
283  allField.name(),
284  mesh,
285  allField.dimensions()
286  );
287 
288  tResult.ref().primitiveFieldRef() = allField.primitiveField();
289 
290  return tResult;
291 }
292 
293 
294 template<class Type>
296 (
297  const LagrangianSubMesh& subMesh
298 ) const
299 {
300  // Evaluate and store if it doesn't already exist for the sub-mesh
301  if (!psiSubSubPtr_.valid() || psiSubSubMeshIndex_ != subMesh.index())
302  {
303  psiSubSubPtr_.reset
304  (
305  new LagrangianSubSubField<Type>(operator()(subMesh))
306  );
307 
308  psiSubSubUpToDate_ = true;
309  psiSubSubMeshIndex_ = subMesh.index();
310  }
311 
312  // Update the field in-place if the up-to-date flag is not set
313  if (!psiSubSubUpToDate_)
314  {
315  psiSubSubPtr_->UList<Type>::shallowCopy(operator()(subMesh));
316 
317  psiSubSubUpToDate_ = true;
318  }
319 
320  return psiSubSubPtr_();
321 }
322 
323 
324 template<class Type>
326 {
327  psiAllPtr_.clear();
328 
329  // If this is not the final iteration, then retain the fields to be
330  // modified in-place. This is more efficient, and it ensures that old-time
331  // fields are maintained.
332  if (!final)
333  {
334  psiSubUpToDate_ = false;
335  psiSubSubUpToDate_ = false;
336  }
337  else
338  {
339  psiSubPtr_.clear();
340  psiSubSubPtr_.clear();
341  }
342 }
343 
344 
345 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
346 
347 template<class Type>
350 {
351  return
353  (
354  new AllFieldToField
355  (
356  operator()(LagrangianModelRef(), mesh.subAll())
357  ),
358  true
359  );
360 }
361 
362 
363 template<class Type>
366 (
367  const LagrangianSubMesh& subMesh
368 ) const
369 {
370  return operator()(LagrangianModelRef(), subMesh);
371 }
372 
373 
374 template<class Type>
377 (
378  const LagrangianModelRef& model,
379  const LagrangianSubMesh& subMesh
380 ) const
381 {
382  // Evaluate and store if it doesn't already exist for the all-mesh
383  if (&subMesh == &subMesh.mesh().subAll())
384  {
385  if (!psiAllPtr_.valid())
386  {
387  if (notNull(name_))
388  {
389  psiAllPtr_.reset
390  (
392  (
393  name_,
394  functorPtr_()(model, subMesh)
395  )
396  );
397  }
398  else
399  {
400  psiAllPtr_.reset(functorPtr_()(model, subMesh).ptr());
401  }
402  }
403 
404  return psiAllPtr_();
405  }
406 
407  // Evaluate and store if it doesn't already exist for the sub-mesh
408  if (!psiSubPtr_.valid() || psiSubMeshIndex_ != subMesh.index())
409  {
410  if (notNull(name_))
411  {
412  psiSubPtr_.reset
413  (
415  (
416  name_ + ':' + name(subMesh.group()),
417  functorPtr_()(model, subMesh)
418  )
419  );
420  }
421  else
422  {
423  psiSubPtr_.reset(functorPtr_()(model, subMesh).ptr());
424  }
425 
426  psiSubUpToDate_ = true;
427  psiSubMeshIndex_ = subMesh.index();
428  }
429 
430  // Update the field in-place if the up-to-date flag is not set
431  if (!psiSubUpToDate_)
432  {
433  psiSubPtr_() = functorPtr_()(model, subMesh);
434 
435  psiSubUpToDate_ = true;
436  }
437 
438  return psiSubPtr_();
439 }
440 
441 
442 // ************************************************************************* //
Graphite solid properties.
Definition: C.H:51
AllFieldToField(const LagrangianSubField< Type > &allField)
Construct from a sub-all-field reference.
Class to store an evaluation function.
Function(const F &f)
Construct from a function.
Class to store an evaluation method.
Method(const C &c, tmp< LagrangianSubField< Type >>(C::*m)(const LagrangianModelRef &, const LagrangianSubMesh &) const)
Construct from a class and a method.
A field derived from other state fields of the cloud. Stores and virtualises a function or a method w...
tmp< LagrangianInternalField< Type > > operator()(const LagrangianMesh &mesh) const
Compute and return the entire field. This will be a slice of the.
CloudDerivedField(const word &name, const F &f)
Construct from a name and a function.
void clear(const bool final)
Clear.
LagrangianSubSubField< Type > & ref(const LagrangianSubMesh &) const
Access a part of the field.
tmp< LagrangianInternalField< Type > > field(const LagrangianMesh &mesh) const
Compute and return an independent copy of the entire field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
const PrimitiveField< Type > & primitiveField() const
Return a const-reference to the primitive field.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Pre-declare SubField and related Field type.
Definition: Field.H:83
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:307
Class containing Lagrangian geometry and topology.
Simple wrapper to provide an optional reference to a Lagrangian model.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
label index() const
Return the index.
void shallowCopy(const UList< T > &)
Copy the pointer held by the given UList.
Definition: UListI.H:156
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
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 dimensionedScalar F
Faraday constant: default SI units: [C/mol].
const dimensionedScalar c
Speed of light in a vacuum.
static const coefficient C("C", dimTemperature, 234.5)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:64
const T & NullObjectRef()
Return const reference to the nullObject of type T.
Definition: nullObjectI.H:27
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
labelList f(nPoints)