ScalarTransportFoam
Contents
1 Introduction and applications
The scalarTransportFoam is a basic solver which resolves a transport equation for a passive scalar, using a user-specified stationary velocity field.
Typical applications are:
- Solution of a scalar convection-diffusion problem on a given velocity field.
2 scalarTransportFoam features, capabilities and limitations
3 scalarTransportFoam capabilities
The scalarTransportFoam solver implements and solves a convection-diffusion scalar transport equation without source terms.
The main features of the solver are:
- Solution of a convection-diffusion equation with user-specified boundary conditions
- Arbitrary velocity field provided by the user and read at runtime
4 scalarTransportFoam limitations
The main limitations of the solver are:
- The diffusion coefficient is assumed to be a constant scalar
- An option to solve for the flow coupled with the scalar transport is not available
5 scalarTransportFoam theory
The scalarTransportFoam solver uses a complete convection-diffusion equation, in the incompressible form (the equation is divided by the density)
where is the transported scalar, is the fluid velocity, and is the diffusion coefficient divided by the fluid density, both supposed to be constant.
6 scalarTransportFoam implementation
6.1 Numerical methodology and solution procedure
The numerical methodology used to solve the scalar trasnport equation in scalarTransportFoam is based on the finite volume technique. In particular, all the terms of the equations are discretized implicitly, using the fvm family of operators.
6.2 Code representation
The implementation of a scalar transport equation in OpenFOAM(r) is straightforward due to its advanced C++ syntax. A few preliminary steps are necessary however to define the fields involved in the terms of the transport equation.
Since this is a basic solver, it is worth to spend a few words on the implementation details, which will be useful to interpret the code of more complex solvers in other guides. In the main body of the solver reported below
#include "fvCFD.H" int main(int argc, char *argv[]) { // Setting the root directory of the case # include "setRootCase.H" // Creating the time object "runTime" # include "createTime.H" // Creating the mesh object # include "createMesh.H" // Creating fields and reading properties from dictionaries # include "createFields.H" Info<< "\nCalculating scalar transport\n" << endl; // Computing the Courant number # include "CourantNo.H" // Starting the time loop while (runTime.loop()) { // Printing the current time Info<< "Time = " << runTime.timeName() << nl << endl; // Reading control parameters for the SIMPLE algorithm # include "readSIMPLEControls.H" // Performing non-orthogonal correction loop for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++) { // Building and solving the scalar tranport equation solve ( fvm::ddt(T) // Unsteady term + fvm::div(phi, T) // Convective term - fvm::laplacian(DT, T) // Diffusive term ); } // Writing data runTime.write(); } Info<< "End\n" << endl; return 0; }
we notice that one main header, fvCFD.H, is included at the top of the code. This header contains a series of other include directives, whose purpose is to make basic OpenFOAM functionalities available in the code. In the following lines inside the main C++ function, a set of header files are included to perform standard tasks:
- define the root directory of the simulation case
- create the time object, which carries information about the flow time, the time step, and some additional information on the solution loop, like output time and data write time.
- create the mesh object and read the mesh
6.2.1 Fields definition
In order to solve the scalar transport equation, the scalar field for and the vector field for must be defined, as well as the diffusion coefficient has to be declared and its value must be read. This is managed in the code in the createFields.H header file as follows:
Info<< "Reading field T\n" << endl; volScalarField T ( IOobject ( "T", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh ); Info<< "Reading field U\n" << endl; volVectorField U ( IOobject ( "U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh );
It is worth to take a deeper look at how the fields are defined and initialized. Both fields are created using an IOobject and the mesh. The IOobject is the standard way in OpenFOAM to manage the input and the output of fields automatically. The structure of the IOobject used in this case is the following
IOobject ( "U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE )
where
- "U" is the name of the file that contains the field on the disk
- runTime.timeName() defines from which time directory the file has to be read, according to what specified in controlDict by the user
- mesh represents the mesh object
- IOobject::MUST_READ specifies the field has to be read from disk. If not available, an error message will be printed and the solver will terminate.
- IOobject::AUTO_WRITE specified the field has to be written automatically, according to what specified in controlDict by the user
Another important point is the definition of a dictionary, which contains, in this case, the physical properties. The dictionary is declared and read as follows:
Info<< "Reading transportProperties\n" << endl; IOdictionary transportProperties ( IOobject ( "transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE ) );
An IOdictionary object is defined to contain the dictionary, and it is initialized using an IOobject, exactly as done for the fields. The only difference is in the instruction runTime.constant(), which specifies that the dictionary file is physically located in the case directory, inside the constant subdirectory. Since the dictionary is read at the beginning of the simulation and it is not changed by the solver, the IOobject::NO_WRITE option is set, instead than IOobject::AUTO_WRITE, so that the solver won't write the file to disk.
The scalar diffusivity is contained in the transportProperties dictionary, and its dimensioned value is read as follows:
Info<< "Reading diffusivity D\n" << endl; dimensionedScalar DT ( transportProperties.lookup("DT") ); # include "createPhi.H"
The diffusivity is declared as a dimensioned scalar, and its value is recovered from the transportProperties dictionary using the lookup command. Note that the name of the variable between " " is the same name specified in the dictionary, and it is not related to the name assigned to the dimensionedScalar inside the code.
6.3 Communicating with the user
The Info command can be considered as the analogous of cout, with the difference that it produces only one output when the code is run in parallel. Its usage is illustrated, for example, by the line
Info<< "\nCalculating scalar transport\n" << endl;
6.4 Solution of the scalar equation
The actual code that solves the scalar transport equation is represented by the following lines
solve ( fvm::ddt(T) + fvm::div(phi, T) - fvm::laplacian(DT, T) );
where all the terms are discretised implicitly, since methods of the fvm namespace are used. It is worth to notice that this syntax could be generalized as follows
fvScalarMatrix TEqn ( fvm::ddt(T) + fvm::div(phi, T) - fvm::laplacian(DT, T) ); TEqn.solve();
where a fvScalarMatrix object is used to contain the matrix originating from the discretization of the scalar transport equation, which is then solved using the method solve(). This syntax is more general, since it allows the easy addition of terms to the transport equations, according to the user needs. Let us assume, for example, that the equation has an additional term if a condition is satisfied. It is possible to add the specific term simply with the following lines, to be added before the solve() method is called:
if (yourCondition) { TEqn += <your term here> }
The above lines are general, and can be used with both explicitly and implicitly discretized terms. In case of explicitly discretized terms, a faster approach is to put the terms directly in the call to the solve() method
solve ( TEqn == <your term here> );
7 scalarTransportFoam setup and solution strategy
7.1 Case structure and user input
The case structure for scalarTransportFoam is the following
. |-- 0 | |-- T | `-- U |-- constant | |-- polyMesh | | |-- blockMeshDict | | `-- boundary | `-- transportProperties `-- system |-- controlDict |-- fvSchemes `-- fvSolutionThe directory
0
contains two files
-
T
, where the initial and boundary conditions for the scalar transport equation are specified -
U
, where the fixed velocity field is contained
constantdirectory contains the physical model settings and the mesh files, in particular
-
polyMesh/blockMeshDict
is the mesh dictionary for blockMesh -
polyMesh/boundary
is the file containing the boundary specifications to construct the mesh -
transportProperties
contains the value of the scalar diffusivity
systemdirectory contains the solver settings:
-
controlDict
contains information on the solution global iterations, and it allows the flags for data storage to be set up -
fvSchemes
contains the definitions of the numerical schemes used during the simulation -
fvSolution
contains the linear solver settings and tolerances
Since the solver has no specific settings, the reader is invited to refer to the official documentation for further information on how to build computational meshes with blockMesh, and on the numerical schemes for the different terms of the transport equations.
8 References
- OpenCFD. OpenFOAM - The Open Source CFD Toolbox - Programmer’s Guide. OpenCFD Ltd., United Kingdom, 1.6 edition, 2009.
- OpenCFD. OpenFOAM - The Open Source CFD Toolbox - User’s Guide. OpenCFD Ltd., United Kingdom, 1.6 edition, 2009.