2.2.2.1. NonLinearAdvectionDiffusionFoam

\[\frac{\partial U}{\partial t}+ U \bullet \nabla U = \nu \nabla^2 U\]

Where U is Vector velocity field in 3D i.e. (u v w) .

2.2.2.1.2. Family of equations.

Name

Equation

Name

Equation

1D First-order Non Linear Advection

\[\frac{\partial u}{\partial t} + u \frac{\partial u}{\partial x} = 0\]

1D Second-order Non Linear Advection Diffusion

\[\frac{\partial u}{\partial t} + u \frac{\partial u}{\partial x} = \nu \frac{\partial^2 u}{\partial x^2}\]

2D First-order Non Linear Advection

\[ \begin{align}\begin{aligned}\frac{\partial u}{\partial t}+ u \frac{\partial u}{\partial x} + v \frac{\partial u}{\partial y} = 0\\\frac{\partial v}{\partial t}+ u \frac{\partial v}{\partial x} + v \frac{\partial v}{\partial y} = 0\end{aligned}\end{align} \]

2D Second-order Non Linear Advection Diffusion

\[ \begin{align}\begin{aligned}\frac{\partial u}{\partial t} + u \frac{\partial u}{\partial x} + v \frac{\partial u}{\partial y} = \nu ( \frac{\partial^2 u}{\partial x^2} + \frac{\partial^2 u}{\partial y^2} )\\\frac{\partial v}{\partial t} + u \frac{\partial v}{\partial x} + v \frac{\partial v}{\partial y} = \nu ( \frac{\partial^2 v}{\partial x^2} + \frac{\partial^2 v}{\partial y^2} )\end{aligned}\end{align} \]

2.2.2.1.3. Step 1: Understanding modified advection term.

\[\frac{\partial U}{\partial t}+ 0.5 * \nabla \bullet (U U) = \nabla \bullet (\nu \nabla U) = \nu \nabla^2 U\]

In OpenFOAM advection is represented as div(phi,U) where phi = U * face area of cell = flux

2.2.2.1.4. Step 2: Create a base directory for solver.

foamNewApp is script to create template for a new application in easiest way which Creates a directory containing source .C file and Make directory rather than copying/modifying exiting solver code.

For more details

mkdir -p $FOAM_RUN                             // Creates directory "/home/user/OpenFOAM/user-foamVerion/run"
cd $FOAM_RUN                                   // Enter into directory
foamNewApp nonLinearAdvectionDiffusionFoam     // Creates base directory for solver
cd nonLinearAdvectionDiffusionFoam             // Enter into directory
ls */*                                         // List the contents

In nonLinearAdvectionDiffusionFoam directory there are files i.e. Make/files, Make/options, nonLinearAdvectionDiffusionFoam.C

2.2.2.1.5. Step 3: Contents of nonLinearAdvectionDiffusionFoam.C .

Open file with any text editor or use nano in terminal. Basic header declarations are already present. Add following lines and create a createFields.H file (mind the extension .H).

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
 #include "fvCFD.H"
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 int main(int argc, char *argv[])
 {
     #include "setRootCase.H"
     #include "createTime.H"
     #include "createMesh.H"    //<----

     #include "createFields.H"  //<----

     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

     Info<< nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
         << "  ClockTime = " << runTime.elapsedClockTime() << " s"
         << nl << endl;

     Info<< "End\n" << endl;

     return 0;
 }

2.2.2.1.6. Step 4: Create field and variables (U,nu) in createFields.H

Insert the following code in createFields.H file. ‘c’ is defined vector to generalized code.

Note the last line createPhi.H is include which is currently not in base directory.

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
 Info<< "Reading transportProperties\n" << endl;
      // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 IOdictionary transportProperties
 (
     IOobject
     (
         "transportProperties",
         runTime.constant(),
         mesh,
         IOobject::MUST_READ_IF_MODIFIED,
         IOobject::NO_WRITE
     )
 );

 dimensionedScalar nu
 (

     transportProperties.lookup("nu")        //---------- Contents is same
 );                                          //---------- icoFoam solver

 Info<< "Reading field U\n" << endl;

 volVectorField U
 (
     IOobject
     (
         "U",
         runTime.timeName(),
         mesh,
         IOobject::MUST_READ,
         IOobject::AUTO_WRITE
     ),
     mesh
 );


 #include "createPhi.H"               //<------ Note here

2.2.2.1.7. Step 5: About createPhi.H

Although this step is not mandatory but for this tutorial createPhi.H file is copied to current directory without any changes.

cp /opt/openfoam5/src/finiteVolume/cfdTools/incompressible/createPhi.H .
 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
 Info<< "Reading/calculating face flux field phi\n" << endl;

 surfaceScalarField phi
 (
     IOobject
     (
         "phi",
         runTime.timeName(),
         mesh,
         IOobject::READ_IF_PRESENT,
         IOobject::AUTO_WRITE
     ),
     fvc::flux(U)              //<--- set to default and preferred in recent
     //U & mesh.Sf()           //<--- versions, although this works too
 );

2.2.2.1.8. Step 6: Adding main code.

Again open nonLinearAdvectionDiffusionFoam.C file and add/edit following code snippet.

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
 Info<< "\nStarting time loop\n" << endl;

 while (runTime.loop())
 {
     Info<< "Time = " << runTime.timeName() << nl << endl;

     #include "CourantNo.H"

     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

     fvVectorMatrix UEqn
     (
         fvm::ddt(U)
       + 0.5*fvm::div(phi, U)
       - fvm::laplacian(nu, U)
     );

     UEqn.solve();

      runTime.write();

      Info<< nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
          << "  ClockTime = " << runTime.elapsedClockTime() << " s"
          << nl << endl;

 }

2.2.2.1.9. Step 7: Compiling code.

OpenFOAM provides wmake system to build a project or application to avoid complex compilation process. Run following command in terminal.

wmake

2.2.2.1.10. Step 8: Testing

Finally solver is ready and time to build a cases for validation. See Test cases