scalarTransport.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) 2012-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 Class
25  Foam::functionObjects::scalarTransport
26 
27 Description
28  Evolves a passive scalar transport equation.
29 
30  - To specify the field name set the \c field entry
31  - To employ the same numerical schemes as another field set
32  the \c schemesField entry,
33  - The \c diffusivity entry can be set to \c none, \c constant, \c viscosity
34  - A constant diffusivity is specified with the \c D entry,
35  - If a momentum transport model is available and the \c viscosity
36  diffusivity option specified an effective diffusivity may be constructed
37  from the laminar and turbulent viscosities using the diffusivity
38  coefficients \c alphal and \c alphat:
39  \verbatim
40  D = alphal*nu + alphat*nut
41  \endverbatim
42 
43  Example:
44  \verbatim
45  #includeFunc scalarTransport(T, alphal=1, alphat=1)
46  \endverbatim
47 
48  For incompressible flow the passive scalar may optionally be solved with the
49  MULES limiter and sub-cycling or semi-implicit in order to maintain
50  boundedness, particularly if a compressive, PLIC or MPLIC convection
51  scheme is used.
52 
53  Example:
54  \verbatim
55  #includeFunc scalarTransport(tracer, diffusivity=none)
56 
57  with scheme specification:
58  div(phi,tracer) Gauss interfaceCompression vanLeer 1;
59 
60  and solver specification:
61  tracer
62  {
63  nCorr 1;
64  nSubCycles 3;
65 
66  MULESCorr no;
67  nLimiterIter 5;
68  applyPrevCorr yes;
69 
70  solver smoothSolver;
71  smoother symGaussSeidel;
72  tolerance 1e-8;
73  relTol 0;
74 
75  diffusivity
76  {
77  solver smoothSolver;
78  smoother symGaussSeidel;
79  tolerance 1e-8;
80  relTol 0;
81  }
82  }
83  \endverbatim
84 
85 See also
86  Foam::functionObjects::fvMeshFunctionObject
87 
88 SourceFiles
89  scalarTransport.C
90 
91 \*---------------------------------------------------------------------------*/
92 
93 #ifndef scalarTransport_functionObject_H
94 #define scalarTransport_functionObject_H
95 
96 #include "fvMeshFunctionObject.H"
97 #include "volFields.H"
98 
99 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
100 
101 namespace Foam
102 {
103 namespace functionObjects
104 {
105 
106 /*---------------------------------------------------------------------------*\
107  Class scalarTransport Declaration
108 \*---------------------------------------------------------------------------*/
109 
110 class scalarTransport
111 :
112  public fvMeshFunctionObject
113 {
114 public:
115 
116  //- Enumeration defining the type of the diffusivity
117  enum class diffusivityType
118  {
119  none,
120  constant,
121  viscosity
122  };
123 
124  //- Diffusivity type names
126 
127 
128 private:
129 
130  // Private Data
131 
132  //- Name of field to process
133  word fieldName_;
134 
135  //- Name of flux field (optional)
136  word phiName_;
137 
138  //- Name of density field (optional)
139  word rhoName_;
140 
141  //- The type of diffusivity
142  diffusivityType diffusivity_;
143 
144  //- Constant diffusivity coefficient (optional)
145  scalar D_;
146 
147  //- Laminar diffusivity coefficient (optional)
148  scalar alphal_;
149 
150  //- Turbulent diffusivity coefficient (optional)
151  scalar alphat_;
152 
153  //- Number of corrector iterations (optional)
154  int nCorr_;
155 
156  //- Name of field whose schemes are used (optional)
157  word schemesField_;
158 
159  //- The scalar field
160  volScalarField s_;
161 
162  //- Switch for MULES limited solution
163  bool MULES_;
164 
165  //- Stabilisation for normalisation of the interface normal
166  // needed if a compressive convection scheme is used
167  const dimensionedScalar deltaN_;
168 
169  //- MULES Correction
170  tmp<surfaceScalarField> tsPhiCorr0_;
171 
172 
173  // Private Member Functions
174 
175  //- Return the diffusivity field
176  tmp<volScalarField> D() const;
177 
178  void subCycleMULES();
179  void solveMULES();
180 
181 
182 public:
183 
184  //- Runtime type information
185  TypeName("scalarTransport");
186 
187 
188  // Constructors
189 
190  //- Construct from Time and dictionary
192  (
193  const word& name,
194  const Time& runTime,
195  const dictionary& dict
196  );
197 
198  //- Disallow default bitwise copy construction
199  scalarTransport(const scalarTransport&) = delete;
200 
201 
202  //- Destructor
203  virtual ~scalarTransport();
204 
205 
206  // Member Functions
207 
208  //- Read the scalarTransport data
209  virtual bool read(const dictionary&);
210 
211  //- Return the list of fields required
212  virtual wordList fields() const;
213 
214  //- Do not execute at the start of the run
215  virtual bool executeAtStart() const
216  {
217  return false;
218  }
219 
220  //- Calculate the scalarTransport
221  virtual bool execute();
222 
223  //- Write the updated scalar field
224  virtual bool write();
225 
226 
227  // Member Operators
228 
229  //- Disallow default bitwise assignment
230  void operator=(const scalarTransport&) = delete;
231 };
232 
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 } // End namespace functionObjects
237 } // End namespace Foam
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 #endif
242 
243 // ************************************************************************* //
Generic GeometricField class.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const word & name() const
Return the name of this functionObject.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Evolves a passive scalar transport equation.
TypeName("scalarTransport")
Runtime type information.
static const NamedEnum< diffusivityType, 3 > diffusivityTypeNames_
Diffusivity type names.
virtual wordList fields() const
Return the list of fields required.
diffusivityType
Enumeration defining the type of the diffusivity.
virtual bool executeAtStart() const
Do not execute at the start of the run.
virtual bool execute()
Calculate the scalarTransport.
virtual bool write()
Write the updated scalar field.
void operator=(const scalarTransport &)=delete
Disallow default bitwise assignment.
virtual bool read(const dictionary &)
Read the scalarTransport data.
scalarTransport(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
dictionary dict