clouds.H
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) 2021-2024 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 Class
25  Foam::fv::clouds
26 
27 Description
28  This fvModel adds any number of Lagrangian clouds to any single-phase
29  solver. The particles are tracked through, and exchange sources with, the
30  Eulerian flow field.
31 
32  As well as the fvModel controls, properties must be specified for each
33  cloud. For a single cloud, these should be provided in the
34  constant/cloudProperties file. For multiple clouds, the list of cloud names
35  must first be provided in the constant/clouds file. Then, each named cloud
36  has its properties specified in its constant/<cloudName>Properties file.
37 
38  The application of sources to the Eulerian fields is controlled by the
39  solution/coupled switch in each cloud's properties file. If set to "true"
40  then the Eulerian phase will have forces, and heat and mass sources applied
41  to it by the Lagrangian phase. If set to "false" then these will be omitted,
42  and the Lagrangian phase will not affect the Eulerian phase.
43 
44  If this model is used with an incompressible solver, then the density of
45  the Eulerian phase must be specified in the constant/physicalProperties
46  dictionary.
47 
48  Gravity will be read from the constant/g file if present, otherwise it will
49  default to zero.
50 
51 Usage
52  Example usage:
53  \verbatim
54  clouds
55  {
56  libs ("liblagrangianParcel.so");
57  type clouds;
58  }
59  \endverbatim
60 
61  Example usage, for multiple clouds:
62  \verbatim
63  clouds
64  {
65  libs ("liblagrangianParcel.so");
66  type clouds;
67  clouds (cloud1 cloud2 cloud3);
68  }
69  \endverbatim
70 
71  \table
72  Property | Description | Required | Default value
73  type | Type name: clouds | yes |
74  rho | Name of the density field | no | rho
75  U | Name of the velocity field | no | U
76  clouds | Names of the clouds | no | (cloud)
77  \endtable
78 
79 SourceFiles
80  clouds.C
81 
82 \*---------------------------------------------------------------------------*/
83 
84 #ifndef clouds_H
85 #define clouds_H
86 
87 #include "fvModel.H"
88 #include "fluidThermo.H"
89 #include "viscosityModel.H"
91 #include "parcelCloudList.H"
92 
93 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
94 
95 namespace Foam
96 {
97 namespace fv
98 {
99 
100 /*---------------------------------------------------------------------------*\
101  Class clouds Declaration
102 \*---------------------------------------------------------------------------*/
103 
104 class clouds
105 :
106  public fvModel
107 {
108  // Private Data
109 
110  //- Optional acceleration due to gravity
112 
113  //- Flag to indicate whether a thermo model is available for the
114  // carrier
115  const bool carrierHasThermo_;
116 
117  //- Reference to the carrier phase thermo. Valid only if the carrier
118  // has thermo.
119  const tmpNrc<fluidThermo> tCarrierThermo_;
120 
121  //- Reference to the carrier viscosity model. Valid only if the carrier
122  // does not have thermo.
123  const tmpNrc<viscosityModel> tCarrierViscosity_;
124 
125  //- Density field. Valid only if the carrier does not have thermo.
126  const tmp<volScalarField> tRho_;
127 
128  //- Viscosity field. Valid only if the carrier does not have thermo.
129  tmp<volScalarField> tMu_;
130 
131  //- Names of the clouds
132  const wordList cloudNames_;
133 
134  //- Name of the density field
135  const word rhoName_;
136 
137  //- Name of the velocity field
138  const word UName_;
139 
140  //- The Lagrangian cloud list
141  mutable autoPtr<parcelCloudList> cloudsPtr_;
142 
143  //- Current time index (used for updating)
144  mutable label curTimeIndex_;
145 
146 
147 public:
148 
149  //- Runtime type information
150  TypeName("clouds");
151 
152 
153  // Constructors
154 
155  //- Construct from explicit source name and mesh
156  clouds
157  (
158  const word& sourceName,
159  const word& modelType,
160  const fvMesh& mesh,
161  const dictionary& dict
162  );
163 
164  //- Disallow default bitwise copy construction
165  clouds
166  (
167  const clouds&
168  ) = delete;
169 
170 
171  // Member Functions
172 
173  // Checks
174 
175  //- Return the list of fields for which the option adds source term
176  // to the transport equation
177  virtual wordList addSupFields() const;
178 
179 
180  // Correct
181 
182  //- Solve the Lagrangian clouds and update the sources
183  virtual void correct();
184 
185 
186  // Add explicit and implicit contributions to compressible equation
187 
188  //- Add source to continuity equation
189  virtual void addSup
190  (
191  const volScalarField& rho,
192  fvMatrix<scalar>& eqn
193  ) const;
194 
195  //- Add source to enthalpy or species equation
196  virtual void addSup
197  (
198  const volScalarField& rho,
199  const volScalarField& field,
200  fvMatrix<scalar>& eqn
201  ) const;
202 
203  //- Add source to incompressible momentum equation
204  virtual void addSup
205  (
206  const volVectorField& U,
207  fvMatrix<vector>& eqn
208  ) const;
209 
210  //- Add source to compressible momentum equation
211  virtual void addSup
212  (
213  const volScalarField& rho,
214  const volVectorField& U,
215  fvMatrix<vector>& eqn
216  ) const;
217 
218 
219  // Mesh changes
220 
221  //- Prepare for mesh update
222  virtual void preUpdateMesh();
223 
224  //- Update topology using the given map
225  virtual void topoChange(const polyTopoChangeMap&);
226 
227  //- Update from another mesh using the given map
228  virtual void mapMesh(const polyMeshMap&);
229 
230  //- Redistribute or update using the given distribution map
231  virtual void distribute(const polyDistributionMap&);
232 
233  //- Update for mesh motion
234  virtual bool movePoints();
235 
236 
237  // Member Operators
238 
239  //- Disallow default bitwise assignment
240  void operator=(const clouds&) = delete;
241 };
242 
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 } // End namespace fv
247 } // End namespace Foam
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 #endif
252 
253 // ************************************************************************* //
Generic GeometricField class.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
This fvModel adds any number of Lagrangian clouds to any single-phase solver. The particles are track...
Definition: clouds.H:131
virtual bool movePoints()
Update for mesh motion.
Definition: clouds.C:382
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
Definition: clouds.C:158
void operator=(const clouds &)=delete
Disallow default bitwise assignment.
virtual void correct()
Solve the Lagrangian clouds and update the sources.
Definition: clouds.C:189
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn) const
Add source to continuity equation.
Definition: clouds.C:208
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: clouds.C:364
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: clouds.C:376
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: clouds.C:357
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: clouds.C:370
TypeName("clouds")
Runtime type information.
clouds(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: clouds.C:52
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
UniformDimensionedField< vector > uniformDimensionedVectorField
labelList fv(nPoints)
dictionary dict