phaseScalarTransport.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) 2019-2021 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::phaseScalarTransport
26 
27 Description
28  Evolves a passive scalar transport equation within one phase of a
29  multiphase simulation. The scalar is considered to be a phase-intensive
30  property; i.e., its value represents an amount per-unit of the phase. In
31  addition to the scalar, the function also writes out the product of the
32  volume fraction and the scalar, as this provides a phase-extensive field
33  which is often more convenient to post-process.
34 
35  Most entries are the same as for the \c scalarTransport function. Refer to
36  its documentation for details. Entries specific to this function are
37  detailed below. Note that the phase-name will be determined by stripping
38  the extension from the supplied field name.
39 
40  If \c alphaPhi is specified and found in the database then it will be used
41  to transport the field. This is the likely mode of operation if this
42  function is used with Euler-Euler solvers, in which \c alphaPhi -like
43  fluxes are available. If \c alphaPhi is not found, then a pressure-like
44  equation will be solved in order to construct it so that it exactly matches
45  the time-derivative of the phase volume or mass. This is likely to be
46  necessary in volume-of-fluid solvers where \c alphaPhi is not part of the
47  solution procedure. The pressure field name will be required in this case.
48 
49 Usage
50  \table
51  Property | Description | Req'd? | Default
52  alpha | Name of the volume-fraction field | no\\
53  | alpha.<phase-name>
54  alphaPhi | Name of the phase-flux field | no\\
55  | alphaPhi.<phase-name>
56  p | Name of the pressure field | no | p
57  residualAlpha | Small volume fraction used to stabilise the solution\\
58  | no | rootSmall
59  writeAlphaField | Also write out alpha multiplied by the field\\
60  | no | true
61  \endtable
62 
63  Example specification for interFoam:
64  \verbatim
65  phaseScalarTransport1
66  {
67  type phaseScalarTransport;
68  libs ("libsolverFunctionObjects.so");
69 
70  field s.water;
71  p p_rgh;
72  }
73  \endverbatim
74 
75  Example specification for multiphaseEulerFoam:
76  \verbatim
77  phaseScalarTransport1
78  {
79  type phaseScalarTransport;
80  libs ("libsolverFunctionObjects.so");
81 
82  field s.water;
83  alphaPhi alphaRhoPhi.water;
84  rho thermo:rho.water;
85  }
86  \endverbatim
87 
88 See also
89  Foam::functionObjects::fvMeshFunctionObject
90  Foam::functionObjects::scalarTransport
91 
92 SourceFiles
93  phaseScalarTransport.C
94 
95 \*---------------------------------------------------------------------------*/
96 
97 #ifndef functionObjects_phaseScalarTransport_H
98 #define functionObjects_phaseScalarTransport_H
99 
100 #include "fvMeshFunctionObject.H"
101 #include "volFields.H"
102 
103 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104 
105 namespace Foam
106 {
107 namespace functionObjects
108 {
109 
110 /*---------------------------------------------------------------------------*\
111  Class phaseScalarTransport Declaration
112 \*---------------------------------------------------------------------------*/
113 
114 class phaseScalarTransport
115 :
116  public fvMeshFunctionObject
117 {
118  // Private Data
119 
120  //- Name of field to process
121  const word fieldName_;
122 
123  //- Name of the phase in which to solve
124  const word phaseName_;
125 
126  //- Name of phase volume-fraction field (optional)
127  word alphaName_;
128 
129  //- Name of the phase flux field (optional)
130  word alphaPhiName_;
131 
132  //- Name of the mixture flux field (optional)
133  word phiName_;
134 
135  //- Name of density field (optional)
136  word rhoName_;
137 
138  //- Name of the pressure field (optional)
139  word pName_;
140 
141  //- Diffusion coefficient (optional)
142  scalar D_;
143 
144  //- Flag to indicate whether a constant, uniform D_ is specified
145  bool constantD_;
146 
147  //- Laminar diffusion coefficient (optional)
148  scalar alphaD_;
149 
150  //- Turbulent diffusion coefficient (optional)
151  scalar alphaDt_;
152 
153  //- Number of corrector iterations (optional)
154  int nCorr_;
155 
156  //- Residual volume-fraction
157  scalar residualAlpha_;
158 
159  //- Name of field whose schemes are used (optional)
160  word schemesField_;
161 
162  //- Flag to indicate whether to write the field multiplied by the phase
163  // fraction
164  bool writeAlphaField_;
165 
166  //- The field
167  volScalarField s_;
168 
169  //- The field multiplied by the phase fraction
170  autoPtr<volScalarField> alphaSPtr_;
171 
172  //- Potential field used to generate the phase flux
173  autoPtr<volScalarField> PhiPtr_;
174 
175 
176  // Private Member Functions
177 
178  //- Return the potential field used to generate the phase flux.
179  // Constructed on demand.
180  volScalarField& Phi();
181 
182  //- Return the phase flux. Tries to look it up, and generates it if the
183  // lookup fails. The generation creates a crude guess for alphaPhi
184  // then solves a Laplacian to correct the flux to match the time
185  // derivative of alpha.
186  tmp<surfaceScalarField> alphaPhi();
187 
188  //- Return the diffusivity field
189  tmp<volScalarField> D(const surfaceScalarField& alphaPhi) const;
190 
191 
192 public:
193 
194  //- Runtime type information
195  TypeName("phaseScalarTransport");
196 
197 
198  // Constructors
199 
200  //- Construct from Time and dictionary
202  (
203  const word& name,
204  const Time& runTime,
205  const dictionary& dict
206  );
207 
208 
209  //- Destructor
210  virtual ~phaseScalarTransport();
211 
212 
213  // Member Functions
214 
215  //- Read the settings from the given dictionary
216  virtual bool read(const dictionary&);
217 
218  //- Return the list of fields required
219  virtual wordList fields() const;
220 
221  //- Solve for the evolution of the field
222  virtual bool execute();
223 
224  //- Do nothing. The field is registered and written automatically.
225  virtual bool write();
226 };
227 
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 } // End namespace functionObjects
232 } // End namespace Foam
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 #endif
237 
238 // ************************************************************************* //
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
virtual bool read(const dictionary &)
Read the settings from the given dictionary.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
virtual bool write()
Do nothing. The field is registered and written automatically.
A class for handling words, derived from string.
Definition: word.H:59
virtual wordList fields() const
Return the list of fields required.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
virtual bool execute()
Solve for the evolution of the field.
TypeName("phaseScalarTransport")
Runtime type information.
phaseScalarTransport(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Namespace for OpenFOAM.