CloudStateField.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 "CloudStateField.H"
27 #include "toSubField.H"
28 
29 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const LagrangianSubMesh& subMesh
35 )
36 {
37  // Evaluate and store if it doesn't already exist for the sub-mesh
38  if (!psiSubSubPtr_.valid() || psiSubSubMeshIndex_ != subMesh.index())
39  {
40  psiSubSubPtr_.reset(subMesh.sub(*this).ptr());
41 
42  psiSubSubMeshIndex_ = subMesh.index();
43  }
44 
45  return psiSubSubPtr_();
46 }
47 
48 
49 template<class Type>
51 {
52  psiSubSubPtr_.clear();
53 }
54 
55 
56 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
57 
58 template<class Type>
61 {
62  return *this;
63 }
64 
65 
66 template<class Type>
69 (
70  const LagrangianSubMesh& subMesh
71 ) const
72 {
73  return subMesh.sub(*this);
74 }
75 
76 
77 template<class Type>
80 (
81  const LagrangianModel& model,
82  const LagrangianSubMesh& subMesh
83 ) const
84 {
85  return this->sources()[model.name()].value(model, subMesh);
86 }
87 
88 
89 template<class Type>
92 (
93  const LagrangianModelRef& model,
94  const LagrangianSubMesh& subMesh
95 ) const
96 {
97  if (!model.valid())
98  {
99  return subMesh.sub(*this);
100  }
101  else
102  {
103  return toSubField(operator()(model(), subMesh));
104  }
105 }
106 
107 
108 // ************************************************************************* //
LagrangianSubSubField< Type > & ref(const LagrangianSubMesh &subMesh)
Access a part of the field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
const Internal & operator()() const
Return a const-reference to the dimensioned internal field.
friend Ostream & operator(Ostream &, const GeometricField< Type, GeoMesh, PrimitiveField > &)
Class containing Lagrangian geometry and topology.
Simple wrapper to provide an optional reference to a Lagrangian model.
bool valid() const
Check if the reference is valid.
Base class for Lagrangian models.
const word & name() const
The source 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.
label index() const
Return the index.
A class for managing temporary objects.
Definition: tmp.H:55
tmp< DimensionedField< Type, GeoMesh, SubField > > toSubField(const DimensionedField< Type, GeoMesh, Field > &)
Return a temporary sub-field from a reference to a field.
Functions to cast/convert dimensioned field references and temporaries based on a primitive field to ...