File: createBlockMesh.H

package info (click to toggle)
openfoam 1812%2Bdfsg1-2
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 220,284 kB
  • sloc: cpp: 1,038,902; sh: 14,536; ansic: 8,240; lex: 657; xml: 387; python: 300; awk: 212; makefile: 94; sed: 88; csh: 3
file content (104 lines) | stat: -rw-r--r-- 2,425 bytes parent folder | download
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
const cellModel& hex = cellModel::ref(cellModel::HEX);

cellShapeList cellShapes;
faceListList boundary;
pointField points;
{
    Info<< "Creating block" << endl;

    block b
    (
        cellShape(hex, identity(8), false),
        pointField
        (
            {
                point(0, 0, 0),
                point(L.x(), 0, 0),
                point(L.x(), L.y(), 0),
                point(0, L.y(), 0),
                point(0, 0, L.z()),
                point(L.x(), 0, L.z()),
                point(L.x(), L.y(), L.z()),
                point(0, L.y(), L.z())
            }
        ),
        blockEdgeList(),
        blockFaceList(),
        N,
        List<gradingDescriptors>(12)
    );

    Info<< "Creating cells" << endl;

    List<FixedList<label, 8>> bCells(b.cells());
    cellShapes.setSize(bCells.size());
    forAll(cellShapes, celli)
    {
        cellShapes[celli] =
            cellShape(hex, labelList(bCells[celli]), false);
    }

    Info<< "Creating boundary faces" << endl;

    boundary.setSize(b.boundaryPatches().size());
    forAll(boundary, patchi)
    {
        faceList faces(b.boundaryPatches()[patchi].size());
        forAll(faces, facei)
        {
            faces[facei] = face(b.boundaryPatches()[patchi][facei]);
        }
        boundary[patchi].transfer(faces);
    }

    points.transfer(const_cast<pointField&>(b.points()));
}

Info<< "Creating patch dictionaries" << endl;
wordList patchNames(boundary.size());
forAll(patchNames, patchi)
{
    patchNames[patchi] = "patch" + Foam::name(patchi);
}

PtrList<dictionary> boundaryDicts(boundary.size());
forAll(boundaryDicts, patchi)
{
    boundaryDicts.set(patchi, new dictionary());
    dictionary& patchDict = boundaryDicts[patchi];
    word nbrPatchName;
    if (patchi % 2 == 0)
    {
        nbrPatchName = "patch" + Foam::name(patchi + 1);
    }
    else
    {
        nbrPatchName = "patch" + Foam::name(patchi - 1); 
    }

    patchDict.add("type", cyclicPolyPatch::typeName);
    patchDict.add("neighbourPatch", nbrPatchName);
}

Info<< "Creating polyMesh" << endl;
polyMesh mesh
(
    IOobject
    (
        polyMesh::defaultRegion,
        runTime.constant(),
        runTime,
        IOobject::NO_READ
    ),
    std::move(points),
    cellShapes,
    boundary,
    patchNames,
    boundaryDicts,
    "defaultFaces",
    cyclicPolyPatch::typeName,
    false
);

Info<< "Writing polyMesh" << endl;
mesh.write();