The following article is related to the numerical schemes implemented in the open source toolbox OpenFOAM®. The article will describe the implementation of different classes and demonstrate how the code accesses different schemes. The description is kept clear and does not explain everything. In particular, different constructors, steps, classes, etc. are not mentioned. The intention is related to a general understanding of the OpenFOAM® workflow (which is complex if we are not working in the top-level).

The Code.

Understanding OpenFOAM® concerning discretization, the analysis of a convection term is given. OpenFOAM® offers the possibility to discretize each term by its own. The equation of interest:

\[\frac{\partial \rho \phi}{\partial t} + \nabla \bullet (\rho \textbf{U} \phi) = 0\]

This equation is the same as given in the section Equation and Conditions. However, as we can see
\[\phi\]
is the quantity of interest. However, in OpenFOAM
\[\phi\]
is commonly used for the cell face fluxes and represents the product of the density and velocity. To keep clearance, we replace the quantity
\[\phi\]
using S. Thus, the equation of interest is written as:
\[\frac{\partial \rho S}{\partial t} + \nabla \bullet (\rho \textbf{U} S) = 0\]

If we implement this equation in FOAM, it follows:

fvm::ddt(rho,S) + fvm::div(phi,S) == 0

Again, the quantity phi is nothing else than the product of the density and the velocity at the cell faces.
Generally, we are solving a linear system of equations. To be more precise and conform with the programming style in OpenFOAM®, we are building a new matrix that is solved using an appropriate linear solver.

fvScalarMatrix SEqn
(
    fvm::ddt(rho,S)
  + fvm::div(phi,S)
);

As we are constructing a scalar matrix named SEqn, we can say that the single terms represent single matrices such as the timeMatrix and the convectionMatrix, which are summed up. To illustrate this, we can write:

//- Exemplify
fvScalarMatrix SEqn
(
    timeMatrix
  + convectionMatrix
);

It should be clear that we just substituted the code and thus, the following statement holds:

//- Exemplify
timeMatrix = fvm::ddt(rho, S)
convectionMatrix = fvm::div(phi, S)

The key message is, that the fvm::ddt() and fvm::div() functions have to return a matrix. This can be further proven if we check the constructor of the fvScalarMatrix class. The fvScalarMatrix is derived from the template class named fvMatrix. Searching for fvMatrix in Doxygen, we get the following information for the constructors:

//- Doxygen snippet 
fvMatrix (const GeometricField< Type, fvPatchField, volMesh > &, const dimensionSet &)
Construct given a field to solve for. More...

fvMatrix (const fvMatrix< Type > &)
Construct as copy. More...

fvMatrix (const tmp< fvMatrix< Type >> &)
Construct as copy of tmp> deleting argument. More...

fvMatrix (const GeometricField< Type, fvPatchField, volMesh > &, Istream &)
Construct from Istream given field to solve for. More...

As it can be seen, the constructor of the fvMatrix class takes either:

  • a GeometricField and a dimensionSet
  • a fvMatrix (e.g., a fvMatrix equals to fvScalarMatrix or fvMatrix equals to fvVectorMatrix, and so on)
  • a fvMatrix which is a tmp object (a temporary object - special class in OpenFOAM®)
  • a GeometricField and a Istream

For the fvm::ddt() and fvm::div() functions, the second or third constructors will be called. To understand this statement, the fvm::ddt() or fvm::div() functions have to be analyzed. It is worth to mention, that fvm is a namespace in OpenFOAM®.
Let us consider the fvm::div() function in the following:

fvm::div(phi, S)

The first exciting part is the scope operator :: which tells us, that we are calling the function div() from the namespace fvm. As a tip, the phi field is the flux field and represents the fluxes on the surface. Thus, it is a surfaceScalarField. The quantity S is the unknown variable of type volScalarField. The analysis where these two fields are built previously is not further investigated. Knowing the two field, we know that we are calling the function div() which takes two arguments:

//- As illustration
//  The div() is a function that takes two arguments
//  div(anyType argument1, anyType argument2)
div(surfaceScalarField& phi, volScalarField& vf)

The ampersand sign is not discussed here. However, in FOAM the abbreviation vf is related to volumeField and phi represents the fluxes on the surfaces in general.
Gathering the information, we know that we have to search in the fvm namespace for the function div() which takes two arguments. Thus, searching the namespace fvm in Doxygen, the following information are presented:

template< class Type >
tmp< fvMatrix< Type > > 	d2dt2 (const GeometricField< Type, fvPatchField, volMesh > &vf)
.
.
. 
template< class Type >
tmp< fvMatrix< Type > > 	div (const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
 
template< class Type >
tmp< fvMatrix< Type > > 	div (const tmp< surfaceScalarField > &tflux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
 
template< class Type >
tmp< fvMatrix< Type > > 	div (const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template< class Type >
tmp< fvMatrix< Type > > 	div (const tmp< surfaceScalarField > &tflux, const GeometricField< Type, fvPatchField, volMesh > &vf)
.
.
. 
template< class Type >
tmp< fvMatrix< Type > > 	SuSp (const tmp< volScalarField > &, const GeometricField< Type, fvPatchField, volMesh > &)

Hence, we can see that we are calling the function with the following declaration:

template< class Type >
tmp< fvMatrix< Type > > 	div (const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf)

Analyzing the function, we can see, that the function returns an object of  fvMatrix.  Due to the fact that we built a fvScalarMatrix (SEqn), the function returns the same Type (The type is a scalar; thus: tmp> which is identical to tmp). If we analyze the fvm::ddt(rho, S) function, we will end up with the same result, and the statement we did before is correct:

//- Exemplify
fvScalarMatrix SEqn
(
    timeMatrix
  + convectionMatrix
);

However, we still don't know how the discretization of the divergence term will be performed. To understand the implementation in this particular manner, we have to dive more into the code. The next step is to check the source code of the fvm:div() function we mentioned before. Doxygen tells us, that we will find the definition of the function in the fvm.C file in line 75. The corresponding code is:

template< class Type>
tmp< fvMatrix< Type>>
div
(
    const surfaceScalarField& flux,
    const GeometricField& vf
)
{
    return fvm::div(flux, vf, "div("+flux.name()+','+vf.name()+')');
}

Now things get more interesting. The function itself calls another function with the same name div(). However, now we provide three arguments. Formatting the return statement of the div() function, it follows:

{
    return fvm::div
    (
        flux,                                 // Argument #1 (the phi field - volSurfaceField)
        vf,                                   // Argument #2 (the S field - volScalarField)
        "div("+flux.name()+','+vf.name()+')'  // Argument #3 (a string)
    );
}

The argument #3 is a string which is constructed from the fields named flux and fv. As we can see, we use these two objects and call the function name(). The string we built is:

//- Constructed string (argument #3)
"div(phi,S)"

//- The field flux is the phi field built in the createPhi.H (see Doxygen)
surfaceScalarField phi
(
    IOobject
    (
        "phi",                   // <-- this is the name and is returned by flux.name()
        runTime.timeName(),
        mesh,
        IOobject::READ_IF_PRESENT,
        IOobject::AUTO_WRITE
    ),
    fvc::flux(U)
);

//- The field vf is the S field built in the createFields.H of the solver
volScalarField S
(
    IOobject
    (
        "S",                     // <-- this is the name and is returned by vf.name()
        runTime.timeName(),
        mesh,
        IOobject::MUST_READ,
        IOobject::AUTO_WRITE
    ),
    mesh
);

The resulting string in the above-mentioned example would be: "div(phi,S)". After we know the string that is constructed for the third argument, we can use Doxygen to analyze the function div(), which has three arguments. The function is defined in line 46 in the fvm.H

  • surfaceScalarField
  • volScalarField
  • word
template< class Type>
tmp< fvMatrix< Type>>
div
(
    const surfaceScalarField& flux,
    const GeometricField& vf,
    const word& name
)
{
    return fv::convectionScheme::New
    (
        vf.mesh(),
        flux,
        vf.mesh().divScheme(name)
    )().fvmDiv(flux, vf);
}

This function calls the New() function of the fv::convectionScheme class (fv is another namespace in FOAM). This function takes three arguments. Additionally, the object (that is returned by the New() function) is used to call its function named fvmDiv(). Based on the fact that things can get complex now, we split it into single steps. First, we have to check what the New() function will return. Therefore, we use Doxygen again:

static tmp< convectionScheme< Type > > 	New (const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
 	Return a pointer to a new convectionScheme created on freestore. More...

The New() function will return an object of the class convectionScheme< Type >. This class further has the function fvmDiv() implemented. As we can see in the listing below, the function is implemented as a virtual function which returns a fvMatrix. A virtual function is used if the function should be defined in other classes. E.g., the discretization of the divergence term can be done differently. If we make different classes that use different discretization techniques such as bounded or limited schemes, we can implement the fvmDiv() function in these classes differently. This allow flexibility, sustainability, and easy maintenance. It should be clear, that the statement given above is just a way of expression and does not replace the real virtual methodology. 

virtual tmp< fvMatrix< Type > > 	fvmDiv (const surfaceScalarField &, const GeometricField< Type, fvPatchField, volMesh > &) const =0

Even though we know that the fvmDiv() function has to be implemented in some other classes, we still do not know which function will be called. Therefore, we have to understand the New() function of the convectionScheme class. In fact, there are two New() functions declared and defined in the convectionScheme. For now, we are focusing on the one we have mentioned above.

template< class Type>
tmp< convectionScheme< Type>> convectionScheme< Type>::New
(
    const fvMesh& mesh,
    const surfaceScalarField& faceFlux,
    Istream& schemeData
)
{
    if (fv::debug)
    {
        InfoInFunction << "Constructing convectionScheme< Type>" << endl;
    }

    if (schemeData.eof())
    {
        FatalIOErrorInFunction
        (
            schemeData
        )   << "Convection scheme not specified" << endl << endl
            << "Valid convection schemes are :" << endl
            << IstreamConstructorTablePtr_->sortedToc()
            << exit(FatalIOError);
    }

    const word schemeName(schemeData);

    typename IstreamConstructorTable::iterator cstrIter =
        IstreamConstructorTablePtr_->find(schemeName);

    if (cstrIter == IstreamConstructorTablePtr_->end())
    {
        FatalIOErrorInFunction
        (
            schemeData
        )   << "Unknown convection scheme " << schemeName << nl << nl
            << "Valid convection schemes are :" << endl
            << IstreamConstructorTablePtr_->sortedToc()
            << exit(FatalIOError);
    }

    return cstrIter()(mesh, faceFlux, schemeData);
}

The function takes three arguments which comes from the fvm::div() function. The fvm::div() function uses the field vf and calls the mesh() function. As we already know, the vf field is derived from the GeometricField class. Thus, the mesh() function has to be implemented in this class or in any public parent class that is used by the GeometricField class. Doxygen will tell us, that the mesh() function is implemented in the public member function of the DimensionedField < Type, GeoMesh >> class. The function returns a reference to the mesh. The second argument is the flux field, and the third one is an Istream object named schemeData. This comes from the following line in the fvm::div() function: 

vf.mesh().divScheme(name)

Again, the vf object is based on the GeometricField class. First, we call the mesh() function that will return the reference to the mesh object. The returned object is used to apply the divScheme() function. The argument »name« was built in the previous function - "div(phi,S)" - and will be used for the divScheme() function in the fvMesh class. Analyzing the fvMesh class, we will see that this function is not defined directly in the fvMesh class. However, the fvMesh class inherits the fvSchemes class which includes the divScheme() function.

Foam::ITstream& Foam::fvSchemes::divScheme(const word& name) const
{
    if (debug)
    {
        Info<< "Lookup divScheme for " << name << endl;
    }

    if (divSchemes_.found(name) || defaultDivScheme_.empty())
    {
        return divSchemes_.lookup(name);
    }
    else
    {
        const_cast(defaultDivScheme_).rewind();
        return const_cast(defaultDivScheme_);
    }
}

The function is searching in the divSchemes_ object for the corresponding word - div(phi,S). The object divSchemes_ is a dictionary declared in the header file of the fvSchemes class and will be initialized while calling the constructor of the class.
To keep it short, the defined divergence schemes in the system/fvSchemes file are inside the d ivSchemes_ object. The code search for the entry - div(phi,S) - and if this is found, the corresponding ITstream will be returned. Otherwise, the default scheme is returned. See the lookup() function in the dictionary class. The conversion of Istream and ITstream is possible based on the class implementations (see source code). Knowing the Istream inside the New() function of the convectionScheme class, we return the corresponding convectionScheme object which is further used to call the function fvmDiv().

    const word schemeName(schemeData);

    typename IstreamConstructorTable::iterator cstrIter =
        IstreamConstructorTablePtr_->find(schemeName);

The convectionScheme can be of type:

  • boundedConvectionScheme
  • gaussConvectionScheme
  • multivariateGaussConvectionScheme

Considering a Gauss scheme in the system/fvSchemes file, we get a gaussConvectionScheme object. This object calls the fvmDiv() function.

template< class Type>
tmp< fvMatrix< Type>>
gaussConvectionScheme::fvmDiv
(
    const surfaceScalarField& faceFlux,
    const GeometricField& vf
) const
{
    tmp< surfaceScalarField> tweights = tinterpScheme_().weights(vf);
    const surfaceScalarField& weights = tweights();

    tmp< fvMatrix< Type>> tfvm
    (
        new fvMatrix< Type>
        (
            vf,
            faceFlux.dimensions()*vf.dimensions()
        )
    );
    fvMatrix< Type>& fvm = tfvm.ref();

    fvm.lower() = -weights.primitiveField()*faceFlux.primitiveField();
    fvm.upper() = fvm.lower() + faceFlux.primitiveField();
    fvm.negSumDiag();

    forAll(vf.boundaryField(), patchi)
    {
        const fvPatchField< Type>& psf = vf.boundaryField()[patchi];
        const fvsPatchScalarField& patchFlux = faceFlux.boundaryField()[patchi];
        const fvsPatchScalarField& pw = weights.boundaryField()[patchi];

        fvm.internalCoeffs()[patchi] = patchFlux*psf.valueInternalCoeffs(pw);
        fvm.boundaryCoeffs()[patchi] = -patchFlux*psf.valueBoundaryCoeffs(pw);
    }

    if (tinterpScheme_().corrected())
    {
        fvm += fvc::surfaceIntegrate(faceFlux*tinterpScheme_().correction(vf));
    }

    return tfvm;
}

The function builts the finite volume matrix and takes the interpolation schemes specified in the fvSchemes file. The scheme in use is called by tinterpScheme_(). This function is defined in the gaussConvectionScheme class.

template< class Type>
const surfaceInterpolationScheme< Type>&
gaussConvectionScheme< Type>::interpScheme() const
{
    return tinterpScheme_();
}

This function returns the surface interpolation scheme we set in the system/fvSchemes file. As we can see, the returned object is related to the class surfaceInterpolationScheme. Doxygen tells us a bunch of different interpolation schemes that OpenFOAM® offers here and we can use. For example, the upwind scheme is inherited by the limitedSurfaceInterpolationScheme. Depending on the scheme in use, different functions that are called by the tinterpScheme_ object are evaluated. Furthermore, we can see above, that the following functions are in use inside the gaussConvectionScheme::div() function.

  • tinterpScheme_().weights()
  • tinterpScheme_().corrected()
  • tinterpScheme_().correction()

This functions are related to the surfaceInterpolationScheme class which are related to the scheme in use. E.g., the weights() function is not implemented in the surfaceInterpolationScheme class. This is even not implemented in the limitedSurfaceInterpolationScheme class but in the upwind class. As one can see, things get complex now. Thus, a fundamental understanding of C++ is required to analyze the code.

As shown in this article, Doxygen is the tool to use to understand the code.