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-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::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 the solver does not provide an \c alphaPhi flux, or that flux is for
41  some reason unreliable, then the \c solveAlphaPhi switch can be used to
42  make this function solve a pressure-like equation from which \c alphaPhi is
43  recovered.
44 
45 Usage
46  \table
47  Property | Description | Req'd? | Default
48  alpha | Name of the volume-fraction field | no\\
49  | alpha.<phase-name>
50  alphaPhi | Name of the phase-flux field | no\\
51  | alphaPhi.<phase-name>
52  solveAlphaPhi | Solve for the alphaPhi flux, rather than looking it\\
53  up | no | false
54  p | Name of the pressure field | no | p
55  residualAlpha | Small volume fraction used to stabilise the solution\\
56  | no | rootSmall
57  writeAlphaField | Also write out alpha multiplied by the field\\
58  | no | true
59  \endtable
60 
61  Example specification for incompressibleVoF:
62  \verbatim
63  phaseScalarTransport1
64  {
65  type phaseScalarTransport;
66  libs ("libsolverFunctionObjects.so");
67 
68  field s.water;
69  }
70  \endverbatim
71 
72  Example specification for multiphaseEuler:
73  \verbatim
74  phaseScalarTransport1
75  {
76  type phaseScalarTransport;
77  libs ("libsolverFunctionObjects.so");
78 
79  field s.water;
80  alphaPhi alphaRhoPhi.water;
81  rho rho.water;
82  }
83  \endverbatim
84 
85 See also
86  Foam::functionObjects::fvMeshFunctionObject
87  Foam::functionObjects::scalarTransport
88 
89 SourceFiles
90  phaseScalarTransport.C
91 
92 \*---------------------------------------------------------------------------*/
93 
94 #ifndef phaseScalarTransport_functionObject_H
95 #define phaseScalarTransport_functionObject_H
96 
97 #include "scalarTransport.H"
98 
99 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
100 
101 namespace Foam
102 {
103 namespace functionObjects
104 {
105 
106 /*---------------------------------------------------------------------------*\
107  Class phaseScalarTransport Declaration
108 \*---------------------------------------------------------------------------*/
109 
110 class phaseScalarTransport
111 :
112  public fvMeshFunctionObject
113 {
114  // Private Data
115 
116  //- Name of field to process
117  const word fieldName_;
118 
119  //- Name of the phase in which to solve
120  const word phaseName_;
121 
122  //- Solve for the alphaPhi flux, rather than looking it up (optional)
123  bool solveAlphaPhi_;
124 
125  //- Name of phase volume-fraction field (optional)
126  word alphaName_;
127 
128  //- Name of the phase flux field (optional)
129  word alphaPhiName_;
130 
131  //- Name of the mixture flux field (optional)
132  word phiName_;
133 
134  //- Name of density field (optional)
135  word rhoName_;
136 
137  //- Name of the pressure field (optional)
138  word pName_;
139 
140  //- The type of diffusivity
142 
143  //- Constant diffusivity coefficient (optional)
144  scalar D_;
145 
146  //- Laminar diffusivity coefficient (optional)
147  scalar alphal_;
148 
149  //- Turbulent diffusivity coefficient (optional)
150  scalar alphat_;
151 
152  //- Number of corrector iterations (optional)
153  int nCorr_;
154 
155  //- Residual volume-fraction
156  scalar residualAlpha_;
157 
158  //- Name of field whose schemes are used (optional)
159  word schemesField_;
160 
161  //- Flag to indicate whether to write the field multiplied by the phase
162  // fraction
163  bool writeAlphaField_;
164 
165  //- The field
166  volScalarField s_;
167 
168  //- Potential field used to generate the phase flux
169  autoPtr<volScalarField> PhiPtr_;
170 
171 
172  // Private Member Functions
173 
174  //- Return the potential field used to generate the phase flux.
175  // Constructed on demand.
176  volScalarField& Phi();
177 
178  //- Return the phase flux. Tries to look it up, and generates it if the
179  // lookup fails. The generation creates a crude guess for alphaPhi
180  // then solves a Laplacian to correct the flux to match the time
181  // derivative of alpha.
182  tmp<surfaceScalarField> alphaPhi();
183 
184  //- Return the diffusivity field
185  tmp<volScalarField> D(const surfaceScalarField& alphaPhi) const;
186 
187 
188 public:
189 
190  //- Runtime type information
191  TypeName("phaseScalarTransport");
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 
205  //- Destructor
206  virtual ~phaseScalarTransport();
207 
208 
209  // Member Functions
210 
211  //- Read the settings from the given dictionary
212  virtual bool read(const dictionary&);
213 
214  //- Return the list of fields required
215  virtual wordList fields() const;
216 
217  //- Do not execute at the start of the run
218  virtual bool executeAtStart() const
219  {
220  return false;
221  }
222 
223  //- Solve for the evolution of the field
224  virtual bool execute();
225 
226  //- Do nothing. The field is registered and written automatically.
227  virtual bool write();
228 };
229 
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 } // End namespace functionObjects
234 } // End namespace Foam
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 #endif
239 
240 // ************************************************************************* //
Generic GeometricField class.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
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
const word & name() const
Return the name of this functionObject.
TypeName("phaseScalarTransport")
Runtime type information.
virtual wordList fields() const
Return the list of fields required.
virtual bool executeAtStart() const
Do not execute at the start of the run.
virtual bool execute()
Solve for the evolution of the field.
virtual bool write()
Do nothing. The field is registered and written automatically.
virtual bool read(const dictionary &)
Read the settings from the given dictionary.
phaseScalarTransport(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
diffusivityType
Enumeration defining the type of the diffusivity.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
dictionary dict