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-2022 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, diffusion=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  diffusion
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 functionObjects_scalarTransport_H
94 #define functionObjects_scalarTransport_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 \*---------------------------------------------------------------------------*/
110 class scalarTransport
111 :
112  public fvMeshFunctionObject
113 {
114 public:
115 
116  //- Enumeration defining the type of the diffusion
117  enum class diffusionType
118  {
119  none,
120  constant,
121  viscosity
122  };
123 
124 
125 private:
126 
127  // Private Data
128 
129  //- Name of field to process
130  word fieldName_;
131 
132  //- Name of flux field (optional)
133  word phiName_;
134 
135  //- Name of density field (optional)
136  word rhoName_;
137 
138  //- Diffusion type names
139  static const NamedEnum<diffusionType, 3> diffusionTypeNames_;
140 
141  //- The type of diffusion
142  diffusionType diffusion_;
143 
144  //- Constant diffusion coefficient (optional)
145  scalar D_;
146 
147  //- Laminar diffusion coefficient (optional)
148  scalar alphal_;
149 
150  //- Turbulent diffusion 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  //- Scalar volumetric flux
171 
172  //- MULES Correction
173  tmp<surfaceScalarField> tsPhiCorr0_;
174 
175  //- Switch to indicate that s has been restarted for Crank-Nicolson
176  bool sRestart_;
177 
178 
179  // Private Member Functions
180 
181  //- Return the diffusivity field
182  tmp<volScalarField> D() const;
183 
184  void subCycleMULES();
185  void solveMULES();
186 
187 
188 public:
189 
190  //- Runtime type information
191  TypeName("scalarTransport");
192 
193 
194  // Constructors
195 
196  //- Construct from Time and dictionary
198  (
199  const word& name,
200  const Time& runTime,
201  const dictionary& dict
202  );
203 
204  //- Disallow default bitwise copy construction
205  scalarTransport(const scalarTransport&) = delete;
206 
207 
208  //- Destructor
209  virtual ~scalarTransport();
210 
211 
212  // Member Functions
213 
214  //- Read the scalarTransport data
215  virtual bool read(const dictionary&);
216 
217  //- Return the list of fields required
218  virtual wordList fields() const;
219 
220  //- Calculate the scalarTransport
221  virtual bool execute();
222 
223  //- Do nothing.
224  // The volScalarField is registered and written automatically
225  virtual bool write();
226 
227 
228  // Member Operators
229 
230  //- Disallow default bitwise assignment
231  void operator=(const scalarTransport&) = delete;
232 };
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 } // End namespace functionObjects
238 } // End namespace Foam
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 #endif
243 
244 // ************************************************************************* //
dictionary dict
const word & name() const
Return the name of this functionObject.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
diffusionType
Enumeration defining the type of the diffusion.
scalarTransport(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
virtual bool execute()
Calculate the scalarTransport.
Evolves a passive scalar transport equation.
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read(const dictionary &)
Read the scalarTransport data.
Abstract base class for all fluid physical properties.
Definition: viscosity.H:49
void operator=(const scalarTransport &)=delete
Disallow default bitwise assignment.
TypeName("scalarTransport")
Runtime type information.
virtual bool write()
Do nothing.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
A class for managing temporary objects.
Definition: PtrList.H:53
virtual wordList fields() const
Return the list of fields required.
Namespace for OpenFOAM.