Skip to content
🎉 Welcome to Aerosand!.

15_mesh

Important

Visit https://aerosand.cc for the latest updates.

0. Preface

In addition to setRootCase.H and createTime.H, there is another essential header file that inevitably appears: createMesh.H.

This section primarily discusses:

  • Understanding mesh-related classes
  • Understanding the createMesh.H header file
  • Practicing mesh-related methods
  • Compiling and running a mesh project

1. Project Preparation

Run the following commands in the terminal to create the project:

terminal
ofsp
foamNewApp ofsp_15_mesh
cd ofsp_15_mesh
cp -r $FOAM_TUTORIALS/incompressible/icoFoam/cavity/cavity debug_case
code .

Test the initial solver, and provide scripts and documentation.

2. mesh

Let us first not look at the source code of createMesh.H, but start from the needs of the application. That is, we need a mesh for computation based on a specific case. This mesh naturally has its own data (including points, faces, etc.) and its own methods (including returning point lists, returning face centers, etc.).

2.1. primitiveMesh

API page: https://api.openfoam.com/2506/classFoam_1_1primitiveMesh.html

OpenFOAM provides the primitiveMesh class based on geometric elements.

Run the following command in the terminal to view the declaration of the primitiveMesh class:

terminal
find $FOAM_SRC -iname primitiveMesh.H

It can be seen that the primitiveMesh class does not inherit from other classes.

Let us briefly explore the primitiveMesh class by looking at a few parts of the code.

primitiveMesh.H
 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
class primitiveMesh
{
...
    // Constructors

        //- Construct from components
        primitiveMesh
        (
            const label nPoints,
            const label nInternalFaces,
            const label nFaces,
            const label nCells
        );
    ...
    
    // Access
    
	    // Primitive mesh data

			//- Return mesh points		
			virtual const pointField& points() const = 0;
...

From this constructor, it can be seen that the primitiveMesh class constructs an object from the basic geometric elements: points, faces, and cells.

The geometric elements—points, faces, and cells—are generated by blockMesh or by third-party mesh software and then converted to OpenFOAM format. In any case, the points, faces, and cells are generated in the case/constant/polyMesh/ folder. The geometric data here are just data and do not have any methods.

Note that this class is a pure virtual class (abstract class) and cannot be instantiated directly.

2.2. polyMesh

Based on the pure virtual class primitiveMesh, OpenFOAM provides the derived class polyMesh. The polyMesh class is essentially still geometric and topological, but additionally provides some geometric and topological methods.

API page: https://api.openfoam.com/2506/classFoam_1_1polyMesh.html

Run the following command in the terminal to view the declaration of polyMesh:

terminal
find $FOAM_SRC -iname polyMesh.H

Let us briefly explore the polyMesh class by looking at a few parts of the code.

polyMesh.H
 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
...
class polyMesh
:
	public objectRegistry, // Inherits from these two classes
	public primitiveMesh
{
public:
	...
private:
	...
public:
	...
	// Constructors
		//- Read construct from IOobject
		explicit polyMesh(const IOobject& io, const bool doInit = true);
		// Constructs a polyMesh object based on an IOobject
	...
...

Using the polyMesh class, we can construct a mesh object and perform some mesh operations.

Geometric elements can be encapsulated for subsequent use.

OpenFOAM provides the IOobject class to encapsulate access to these data.

API page: https://api.openfoam.com/2506/classFoam_1_1IOobject.html

Run the following command in the terminal to find IOobject.H:

terminal
find $FOAM_SRC -iname IOobject.H

Let us briefly explore the IOobject class by looking at a few parts of the code.

IOobject.H
 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
...
class IOobject
:
	public IOobjectOption
{
public:
	...
private:
	...
protected:
	...
public:
	...
	// Constructors
	inline IOobject
	(
		const word& name, // Name of the object
		const fileName& instance, // File name of the object
		const objectRegistry& registry, // Registry-related, not needed for now
		IOobjectOption::readOption rOpt, // Read option
		IOobjectOption::writeOption wOpt = IOobjectOption::NO_WRITE, // Write option
		bool registerObject = true, // Has a default value, can be ignored
		bool globalObject = false // Has a default value, can be ignored
		// Understand the usage for now
	);
	...
...

Modify the main source code as follows:

ofsp_15_mesh/ofsp_15_mesh.C
 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
38
#include "fvCFD.H"

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

int main(int argc, char *argv[])
{
    #include "setRootCase.H"
    #include "createTime.H"

    Foam::word regionName(Foam::polyMesh::defaultRegion);

    Foam::polyMesh mesh
    (
        IOobject
        (
            // fvMesh::defaultRegion,
            // You can directly use defaultRegion, or create a new regionName
            regionName,
            runTime.timeName(),
            runTime,
            IOobject::MUST_READ
        )
    );

    Info<< "Max cell centre: " << max(mesh.cellCentres()) << endl;
    Info<< "Max cell volumes: " << max(mesh.cellVolumes()) << endl;
    Info<< "Max cell face cetres: " << max(mesh.faceCentres()) << endl;
	// To reduce output, a max function is applied to each

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

    Info<< nl;
    runTime.printExecutionTime(Info);

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

    return 0;
}

Some readers may notice that the construction of IOobject is based on runTime, rather than directly on constant/polyMesh. Recalling the search mechanism of the Time class discussed in the previous section, it can be understood that the Time class has a complex search mechanism.

Tip

This involves implementation details of the source code. At the current stage, it is not necessary or recommended to delve too deeply into it.

Compile and run.

The terminal output is as follows:

terminal
Create time

Max cell centre: (0.0975 0.0975 0.005)
Max cell volumes: 2.5e-07
Max cell face cetres: (0.1 0.1 0.01)

ExecutionTime = 0 s  ClockTime = 0 s

End

2.3. fvMesh

Based on the polyMesh class, OpenFOAM adds finite volume methods to derive the fvMesh class.

API page: https://api.openfoam.com/2506/classFoam_1_1fvMesh.html

Run the following command in the terminal to view the declaration of the fvMesh class:

terminal
find $FOAM_SRC -iname fvMesh.H

Let us briefly explore the fvMesh class by looking at a few parts of the code.

fvMesh.H
 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
...
class fvMesh
:
    public polyMesh, // Inherits from polyMesh
    public lduMesh,
    public fvSchemes,
    public surfaceInterpolation,    // needs input from fvSchemes
    public fvSolution,
    public data
{
protected:
	...
public:
	...
    // Constructors

        //- Construct from IOobject
        explicit fvMesh(const IOobject& io, const bool doInit=true);
        // Constructs an object based on IOobject
	...
...

Using the fvMesh class, we can construct a mesh object and perform various operations, including mesh operations.

Modify the main source code as follows:

ofsp_15_mesh/ofsp_15_mesh.C
 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
#include "fvCFD.H"

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

int main(int argc, char *argv[])
{
    #include "setRootCase.H"
    #include "createTime.H"

    Foam::word regionName(Foam::polyMesh::defaultRegion);

    Foam::fvMesh mesh
    (
        IOobject
        (
            regionName,
            runTime.timeName(),
            runTime,
            IOobject::MUST_READ
        )
    );

    Info<< "Max cell centre: " << max(mesh.C()) << endl;
    Info<< "Max cell volumes: " << max(mesh.V()) << endl;
    Info<< "Max cell face cetres: " << max(mesh.Cf()) << endl;


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

    Info<< nl;
    runTime.printExecutionTime(Info);

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

    return 0;
}

Compile and run.

The terminal output is as follows:

terminal
Create time

Max cell centre: max(C) [0 1 0 0 0 0 0] (0.1 0.1 0.005)
Max cell volumes: max(V) [0 3 0 0 0 0 0] 2.5e-07
Max cell face cetres: max(Cf) [0 1 0 0 0 0 0] (0.1 0.1 0.005)

ExecutionTime = 0 s  ClockTime = 0 s

End

Readers can also remove the max function and recompile to view the results. Based on the discussion of these three subclasses and base classes, we can imagine what the main code statements in createMesh.H are.

3. createMesh.H

API page: https://api.openfoam.com/2506/createMesh_8H.html

Run the following command in the terminal to view this header file:

terminal
find $FOAM_SRC -iname createMesh.H

For ease of understanding, let us examine the code from the OpenFOAM 2.0x version.

The GitHub repository file link is: https://github.com/OpenFOAM/OpenFOAM-2.0.x/blob/master/src/OpenFOAM/include/createMesh.H

The code content is as follows:

createMesh.H
 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
//
// createMesh.H
// ~~~~~~~~~~~~

    Foam::Info
        << "Create mesh for time = "
        << runTime.timeName() << Foam::nl << Foam::endl;

    Foam::fvMesh mesh
    (
        Foam::IOobject
        (
            Foam::fvMesh::defaultRegion,
            runTime.timeName(),
            runTime,
            Foam::IOobject::MUST_READ
        )
    );

It can be seen that this is the same as the code we wrote earlier.

The code in the modern version is as follows:

createMesh.H
 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
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
Foam::autoPtr<Foam::fvMesh> meshPtr(nullptr);
// Create an automatic pointer to manage the fvMesh

// "getRegionOption.H"
Foam::word regionName
(
    args.getOrDefault<word>("region", Foam::polyMesh::defaultRegion)
);
// Get the region name, default to defaultRegion

if (args.dryRun() || args.found("dry-run-write")) // Whether in dry-run test mode
{
    Foam::Info
        << "Operating in 'dry-run' mode: case will run for 1 time step.  "
        << "All checks assumed OK on a clean exit" << Foam::endl;

    Foam::FieldBase::allowConstructFromLargerSize = true;

    // Create a simplified 1D mesh and attempt to re-create boundary conditions
    meshPtr.reset
    (
        new Foam::simplifiedMeshes::columnFvMesh(runTime, regionName)
    );
    // Create a simplified 1D mesh

    // Stop after 1 iteration of the simplified mesh

    if (args.found("dry-run-write")) // Whether to write
    {
        // Using saWriteNow triggers function objects execute(), write()
        runTime.stopAt(Foam::Time::saWriteNow);

        // Make sure mesh gets output to the current time (since instance
        // no longer constant)
        meshPtr().setInstance(runTime.timeName());
    }
    else
    {
        // Using saNoWriteNow triggers function objects execute(),
        // but not write()
        runTime.stopAt(Foam::Time::saNoWriteNow);
    }

    Foam::functionObject::outputPrefix = "postProcessing-dry-run";
}
else // Non-test mode
{
    Foam::Info << "Create mesh";
    if (!Foam::polyMesh::regionName(regionName).empty())
    {
        Foam::Info << ' ' << regionName;
    }
    Foam::Info << " for time = " << runTime.timeName() << Foam::nl;
    // If the region name is non-empty (the default region is usually region0)

	// Create an fvMesh object via the pointer
    meshPtr.reset
    (
        new Foam::fvMesh
        (
            Foam::IOobject
            (
                regionName,
                runTime.timeName(),
                runTime,
                Foam::IOobject::MUST_READ
            ),
            false
        )
    );
    meshPtr().init(true);   // initialise all (lower levels and current)

    Foam::Info << Foam::endl;
}

Foam::fvMesh& mesh = meshPtr();

Comparing the two versions, we can roughly understand the main content of this file and what mechanisms have been added in the modern version.

4. Mesh Information

We can integrate the above discussion and use mesh information in the project.

Modify the main source code as follows:

ofsp_15_mesh/ofsp_15_mesh.C
  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
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
#include "fvCFD.H"

#include "IOmanip.H" // Output format control

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

int main(int argc, char *argv[])
{
    #include "setRootCase.H"
    #include "createTime.H"

    #include "createMesh.H"

    Info<< "Mesh directory: " << mesh.meshDir() << endl; // methods from polyMesh
	// Output the folder path where the mesh is stored

    const pointField& p = mesh.points(); // Points stored in the mesh object
    Info<< "Number of points: " << p.size() << endl;
    for (int i=0; i<3; ++i) // Display the first 3 points
    {
        Info<< "("
            << setf(ios_base::scientific)
            << setw(15)
            << p[i][0]; // Coordinate component 1 of the i-th point
        Info<< ", "
            << setf(ios_base::scientific)
            << setw(15)
            << p[i][1]; // Coordinate component 2 of the i-th point
        Info<< ", "
            << setf(ios_base::scientific)
            << setw(15)
            << p[i][2] // Coordinate component 3 of the i-th point
            << ")" << endl;
    }

    const faceList& f = mesh.faces(); // Faces in the mesh object
    Info<< "Number of faces: " << f.size() << endl;
    forAll(f,i) // Traversal method provided by OpenFOAM
    {
        Info<< "("
            << setw(6) << f[i][0] // Index of the 1st point composing the i-th face
            << setw(6) << f[i][1] // Index of the 2nd point composing the i-th face
            << setw(6) << f[i][2] // Index of the 3rd point composing the i-th face
            << setw(6) << f[i][3] // Index of the 4th point composing the i-th face
            << ")" << endl;
    }
    // Can open debug_case/constant/polyMesh/faces for comparison

    const labelList& fOwner = mesh.faceOwner();
    Info<< "Number of face owner: " << fOwner.size() << endl;
    forAll(fOwner,i)
    {
        Info<< setw(6) << fOwner[i] << endl;
    }
    // Output the cell index of the owner for each face in natural list order
	// For example:
	// 0                    // Face 0 is owned by cell 0
	// 0                    // Face 1 is owned by cell 0
	// 1                    // Face 2 is owned by cell 1
	// 2                    // Face 3 is owned by cell 2
	// ...
	// And so on

    const labelList& fNeigh = mesh.faceNeighbour();
    Info<< "Number of face neighbour: " << fNeigh.size() << endl;
    forAll(fNeigh,i)
    {
        Info<< setw(6) << fNeigh[i] << endl;
    }
    // Output the cell index of the neighbour for each face in natural list order
    // For example:
    // 1                   // Face 0 has neighbour cell 1
    // 2                   // Face 1 has neighbour cell 2
    // 3                   // Face 2 has neighbour cell 3
    // 3                   // Face 3 has neighbour cell 3
    // Other faces are not neighbours of any cell

    const polyBoundaryMesh& bm = mesh.boundaryMesh();
    Info<< "Number of boundary mesh: " << bm.size() << endl;
    forAll(bm,i)
    {
        Info<< "Boundary name: " << bm[i].name()
            << "\tBoundary type: " << bm[i].type()
            << endl;
    }
    // Output the configured boundaries


    Info<< nl << endl;

    Info<< "Bounding box: " << mesh.bounds() << endl;
    Info<< "Mesh volume: " << sum(mesh.V()).value() << endl;

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

    Info<< nl;
    runTime.printExecutionTime(Info);

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

    return 0;
}

To facilitate result display, modify the mesh resolution settings in the test case file debug_case/system/blockMeshDict:

debug_case/system/blockMeshDict
1
2
3
4
5
6
...
blocks
(
    hex (0 1 2 3 4 5 6 7) (2 2 1) simpleGrading (1 1 1)
);
...

Compile and run.

The terminal output is as follows:

terminal
Create time

Create mesh for time = 0

Mesh directory: "polyMesh"
Number of points: 18
(   0.000000e+00,    0.000000e+00,    0.000000e+00)
(   5.000000e-02,    0.000000e+00,    0.000000e+00)
(   1.000000e-01,    0.000000e+00,    0.000000e+00)
Number of faces: 20
(     1     4    13    10)
(     3    12    13     4)
(     4    13    14     5)
(     4     7    16    13)
(     6    15    16     7)
(     7    16    17     8)
(     0     9    12     3)
(     3    12    15     6)
(     2     5    14    11)
(     5     8    17    14)
(     0     1    10     9)
(     1     2    11    10)
(     0     3     4     1)
(     3     6     7     4)
(     1     4     5     2)
(     4     7     8     5)
(     9    10    13    12)
(    12    13    16    15)
(    10    11    14    13)
(    13    14    17    16)
Number of face owner: 20
     0
     0
     1
     2
     2
     3
     0
     2
     1
     3
     0
     1
     0
     2
     1
     3
     0
     2
     1
     3
Number of face neighbour: 4
     1
     2
     3
     3
Number of boundary mesh: 3
Boundary name: movingWall       Boundary type: wall
Boundary name: fixedWalls       Boundary type: wall
Boundary name: frontAndBack     Boundary type: empty


Bounding box: (0.000000e+00 0.000000e+00 0.000000e+00) (1.000000e-01 1.000000e-01 1.000000e-02)
Mesh volume: 1.000000e-04

ExecutionTime = 0.000000e+00 s  ClockTime = 0.000000e+00 s

End

5. Mesh Methods

Let us practice using more mesh methods (member methods).

Modify the main source code as follows:

ofsp_15_mesh/ofsp_15_mesh.C
  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
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
#include "fvCFD.H"

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

int main(int argc, char *argv[])
{
    #include "setRootCase.H"
    #include "createTime.H"

    #include "createMesh.H"

    Info<< "Time                    : " << runTime.timeName() << nl
        << "Number of mesh cells    : " << mesh.C().size() << nl
        << "Number of internal faces: " << mesh.Cf().size() << nl
        << endl;

    for (label cellI=0; cellI<mesh.C().size(); ++cellI) // Native C++ traversal
    {
        if (cellI%100 == 0)
        {
            Info<< "Cell " << cellI
                << " with center at " << mesh.C()[cellI] // Return cell center coordinates
                << endl;
        }
    }
    Info<< nl << endl;


    forAll(mesh.owner(),faceI) // OpenFOAM-style traversal
    {
        if (faceI%200 == 0)
        {
            Info<< "Internal face " << faceI
                << " with center at " << mesh.Cf()[faceI] // Return face center coordinates
                << " with owner cell " << mesh.owner()[faceI] // Return cell index
                << " and neighbour cell " << mesh.neighbour()[faceI] // Return cell index
                << endl;
        }
    }
    Info<< nl << endl;

    forAll(mesh.boundaryMesh(), patchI)
    {
        Info<< "Patch " << patchI
            << " is " << mesh.boundary()[patchI].name() // Return boundary name
            << " with " << mesh.boundary()[patchI].Cf().size() << " faces. "
            // Return the number of boundary faces
            << "Start from face " << mesh.boundary()[patchI].start()
            // Return the starting face index of the boundary
            << endl;
    }
    Info<< nl << endl;

    label nIndex(0); // Custom index
    forAll(mesh.boundaryMesh(), patchI) // Traverse boundaries
    {
        Info<< "Patch " << patchI << nl
            << "\tits face " << nIndex
            << " adjacent to cell " << mesh.boundary()[patchI].patch().faceCells()[nIndex] << nl
            // Return the adjacent cell index
            << "\tits normal vector " << mesh.boundary()[patchI].Sf()[nIndex] << nl
            // Return the face normal vector of the boundary face
            << "\tits surface area " << mag(mesh.boundary()[patchI].Sf()[nIndex]) << nl
            // Return the face area
            << endl;
    }
    Info<< nl << endl;


    const faceList& fcs = mesh.faces();
    const pointField& pts = mesh.points(); // Can also use List<point>& type
    const List<point>& cents = mesh.faceCentres(); // All face centers

    forAll(fcs,faceI) // Traverse all faces
    {
        if (faceI%200 == 0)
        {
            if (faceI < mesh.Cf().size()) // If it is an internal face
            {
                Info<< "Internal face ";
            } else // Otherwise, it is a boundary face
            {
                forAll(mesh.boundary(),patchI)
                {
                    if ( (mesh.boundary()[patchI].start() <= faceI) &&
                        (faceI < mesh.boundary()[patchI].start() +
                        mesh.boundary()[patchI].Cf().size()))
                        // If within this boundary
                    {
                        Info<< "Face on patch " << patchI << ", faceI ";
                        break;
                    }
                }
            }
            Info<< faceI << " with centre at " << cents[faceI] // Return face center coordinates
                << " has " << fcs[faceI].size() << " vertices: "; // Return number of face vertices
            forAll(fcs[faceI],vertexI) // Traverse vertices of this face
            {
                Info<< " " << pts[fcs[faceI][vertexI]]; // Output vertex coordinates
            }
            Info<< endl;
        }
    }
    Info<< nl << endl;

	// Generally, empty boundaries require special attention and sometimes need specific handling
	// Need to identify faces of empty boundaries
    forAll(mesh.boundaryMesh(),patchI) // Traverse all faces composing the boundary
    {
        const polyPatch& pp = mesh.boundaryMesh()[patchI];
        if (isA<emptyPolyPatch>(pp)) // If the face is empty
        {
            Info<< "Patch " << patchI
                << ": " << mesh.boundary()[patchI].name() // Return boundary name
                << " is empty."
                << endl;
        }
    }
    Info<< nl << endl;

    word patchName("movingWall");
    label patchID = mesh.boundaryMesh().findPatchID(patchName);
    // Find the boundary face matching the name
    Info<< "Retrived patch " << patchName
        << " at index " << patchID
        << " using its name only."
        << endl;
    Info<< nl << endl;


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

    Info<< nl;
    runTime.printExecutionTime(Info);

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

    return 0;
}

Compile and run.

The terminal output is as follows:

terminal
Create time

Create mesh for time = 0

Time                    : 0
Number of mesh cells    : 4
Number of internal faces: 4

Cell 0 with center at (0.025 0.025 0.005)


Internal face 0 with center at (0.05 0.025 0.005) with owner cell 0 and neighbour cell 1


Patch 0 is movingWall with 2 faces. Start from face 4
Patch 1 is fixedWalls with 6 faces. Start from face 6
Patch 2 is frontAndBack with 0 faces. Start from face 12


Patch 0
        its face 0 adjacent to cell 2
        its normal vector (0 0.0005 0)
        its surface area 0.0005

Patch 1
        its face 0 adjacent to cell 0
        its normal vector (-0.0005 0 0)
        its surface area 0.0005

Patch 2
        its face 0 adjacent to cell 0
        its normal vector (0 0 -0.0025)
        its surface area 0.0025



Internal face 0 with centre at (0.05 0.025 0.005) has 4 vertices:  (0.05 0 0) (0.05 0.05 0) (0.05 0.05 0.01) (0.05 0 0.01)


Patch 2: frontAndBack is empty.


Retrived patch movingWall at index 0 using its name only.



ExecutionTime = 0.01 s  ClockTime = 0 s

End

6. Summary

Through this discussion, we have gained a brief understanding of mesh-related classes and what createMesh.H is.

At this point, the essential header files—setRootCase.HcreateTime.H, and createMesh.H—have all been discussed.

This section has completed the following discussions:

  • Understanding mesh-related classes
  • Understanding the createMesh.H header file
  • Practicing mesh-related methods
  • Compiling and running a mesh project

Support us

Tip

Hopefully, the sharing here can be helpful to you.

If you find this content helpful, your comments or donations would be greatly appreciated. Your support helps ensure the ongoing updates, corrections, refinements, and improvements to this and future series, ultimately benefiting new readers as well.

The information and message provided during donation will be displayed as an acknowledgment of your support.

Copyright @ 2026 Aerosand

Last updated on • Aerosand