coupled.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 "autoPtr.H"
27 #include "coupled.H"
28 #include "fvcCurl.H"
29 #include "fvcDdt.H"
30 #include "LagrangianmDdt.H"
31 #include "volFieldsFwd.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace clouds
38 {
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 #define ACCESS_CARRIER_FIELDS(Type, nullArg) \
47 namespace Foam \
48 { \
49  namespace clouds \
50  { \
51  template<> \
52  PtrDictionary<CarrierField<Type>>& coupled::carrierFields() const \
53  { \
54  return CAT3(carrier, CAPITALIZE(Type), Fields_); \
55  } \
56  } \
57 }
59 #undef ACCESS_CARRIER_FIELDS
60 
61 
62 void Foam::clouds::coupled::clearCarrierFields()
63 {
64  #define CLEAR_TYPE_CARRIER_FIELDS(Type, nullArg) \
65  forAllIter \
66  ( \
67  PtrDictionary<CarrierField<Type>>, \
68  carrierFields<Type>(), \
69  iter \
70  ) \
71  { \
72  iter().clear(true); \
73  }
74  FOR_ALL_FIELD_TYPES(CLEAR_TYPE_CARRIER_FIELDS);
75  #undef CLEAR_TYPE_CARRIER_FIELDS
76 }
77 
78 
79 void Foam::clouds::coupled::resetCarrierFields(const bool predict)
80 {
81  #define RESET_TYPE_CARRIER_FIELDS(Type, nullArg) \
82  forAllIter \
83  ( \
84  PtrDictionary<CarrierField<Type>>, \
85  carrierFields<Type>(), \
86  iter \
87  ) \
88  { \
89  iter().reset(predict); \
90  }
92  #undef RESET_TYPE_CARRIER_FIELDS
93 }
94 
95 
96 #define ACCESS_CARRIER_EQNS(Type, nullArg) \
97 namespace Foam \
98 { \
99  namespace clouds \
100  { \
101  template<> \
102  PtrDictionary<CarrierEqn<Type>>& coupled::carrierEqns() const \
103  { \
104  return CAT3(carrier, CAPITALIZE(Type), Eqns_); \
105  } \
106  } \
107 }
109 #undef ACCESS_CARRIER_EQNS
110 
111 
112 Foam::autoPtr<Foam::volVectorField> Foam::clouds::coupled::readDUdtc
113 (
114  const cloud& c
115 ) const
116 {
117  const volVectorField& U =
118  c.mesh().mesh().lookupObject<volVectorField>("U");
119 
120  typeIOobject<volVectorField> io
121  (
122  "ddt(" + U.name() + ")",
123  U.mesh().time().name(),
124  U.mesh(),
127  );
128 
129  return
130  autoPtr<volVectorField>
131  (
132  io.headerOk()
133  ? new volVectorField(io, U.mesh())
134  : nullptr
135  );
136 }
137 
138 
139 const Foam::volVectorField& Foam::clouds::coupled::dUdtc() const
140 {
141  if (Uc.psi().hasStoredOldTimes())
142  {
143  dUdtcPtr_.reset(fvc::ddt(Uc.psi()).ptr());
144  }
145  else if (!dUdtcPtr_.valid())
146  {
147  dUdtcPtr_.set
148  (
149  new volVectorField
150  (
151  IOobject
152  (
153  "ddt(" + Uc.psi().name() + ")",
154  Uc.psi().mesh().time().name(),
155  Uc.psi().mesh(),
158  ),
159  Uc.psi().mesh(),
161  )
162  );
163  }
164 
165  return dUdtcPtr_();
166 }
167 
168 
169 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
170 
172 {
173  #define CLEAR_TYPE_CARRIER_EQNS(Type, nullArg) \
174  forAllIter \
175  ( \
176  PtrDictionary<CarrierEqn<Type>>, \
177  carrierEqns<Type>(), \
178  iter \
179  ) \
180  { \
181  iter().clear(); \
182  }
184  #undef CLEAR_TYPE_CARRIER_EQNS
185 }
186 
187 
188 void Foam::clouds::coupled::initialise(const bool predict)
189 {
190  resetCarrierFields(predict);
191 }
192 
193 
195 {
196  clearCarrierFields();
197 }
198 
199 
200 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
201 
203 :
204  dUdtcPtr_(readDUdtc(c)),
205  Uc
206  (
207  carrierField<vector>
208  (
209  c.mesh().mesh().lookupObject<volVectorField>("U")
210  )
211  ),
212  curlUc
213  (
214  carrierField<vector>
215  (
216  "curlUc",
217  [&]()
218  {
219  return fvc::curl(Uc.psi());
220  }
221  )
222  ),
223  DUDtc
224  (
225  carrierField<vector>
226  (
227  "DUDtc",
228  [&]()
229  {
230  return dUdtc() + (Uc.psi() & fvc::grad(Uc.psi()));
231  }
232  )
233  ),
234  nuc(c.derivedField<scalar>(*this, &coupled::calcNuc))
235 {}
236 
237 
238 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
239 
241 {}
242 
243 
244 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
245 
247 {
248  return
250  (
251  IOobject::member(name) + 'c',
253  );
254 }
255 
256 
257 // ************************************************************************* //
Functions for calculating the time derivative for a Lagrangian equation.
Generic GeometricField class.
void reset(const GeometricField< Type, GeoMesh, PrimitiveField2 > &)
Reset the field contents to the given field.
word group() const
Return group (extension part of name)
Definition: IOobject.C:321
word member() const
Return member (name without the extension)
Definition: IOobject.C:327
static word groupName(Name name, const word &group)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
Definition: cloud.H:63
Base class for clouds which are coupled to a fluid.
Definition: coupled.H:60
void initialise(const bool predict)
Initialisation hook.
Definition: coupled.C:188
static word carrierName(const word &name)
Modify a name to disambiguate it as relating to the carrier.
Definition: coupled.C:246
coupled(const cloud &c)
Construct from a reference to the cloud.
Definition: coupled.C:202
void clearCarrierEqns()
Clear the carrier equations.
Definition: coupled.C:171
const CarrierField< vector > & Uc
Carrier velocity.
Definition: coupled.H:138
void partition()
Partition hook.
Definition: coupled.C:194
virtual ~coupled()
Destructor.
Definition: coupled.C:240
A class for handling words, derived from string.
Definition: word.H:62
#define RESET_TYPE_CARRIER_FIELDS(Type, nullArg)
#define ACCESS_CARRIER_EQNS(Type, nullArg)
Definition: coupled.C:96
#define CLEAR_TYPE_CARRIER_EQNS(Type, nullArg)
#define CLEAR_TYPE_CARRIER_FIELDS(Type, nullArg)
#define ACCESS_CARRIER_FIELDS(Type, nullArg)
Definition: coupled.C:46
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FOR_ALL_FIELD_TYPES(Macro,...)
Definition: fieldTypes.H:51
Calculate the curl of the given volField by constructing the Hodge-dual of the symmetric part of the ...
Calculate the first temporal derivative.
U
Definition: pEqn.H:72
defineTypeNameAndDebug(coupled, 0)
const dimensionedScalar c
Speed of light in a vacuum.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > curl(const VolField< Type > &vf)
Definition: fvcCurl.C:45
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
VolField< vector > volVectorField
Definition: volFieldsFwd.H:63
const dimensionSet dimTime
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimVelocity
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)