fvMeshTools.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvMeshTools.H"
27 #include "processorPolyPatch.H"
28 #include "pointFields.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 // Adds patch if not yet there. Returns patchID.
34 (
35  fvMesh& mesh,
36  const polyPatch& patch,
37  const dictionary& patchFieldDict,
38  const word& defaultPatchFieldType,
39  const bool validBoundary
40 )
41 {
42  polyBoundaryMesh& polyPatches =
43  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
44 
45  label patchi = polyPatches.findPatchID(patch.name());
46  if (patchi != -1)
47  {
48  // Already there
49  return patchi;
50  }
51 
52 
53  // Append at end unless there are processor patches
54  label insertPatchi = polyPatches.size();
55 
56  if (!isA<processorPolyPatch>(patch))
57  {
58  forAll(polyPatches, patchi)
59  {
60  const polyPatch& pp = polyPatches[patchi];
61 
62  if (isA<processorPolyPatch>(pp))
63  {
64  insertPatchi = patchi;
65  break;
66  }
67  }
68  }
69 
70  mesh.addPatch
71  (
72  insertPatchi,
73  patch,
74  patchFieldDict,
75  defaultPatchFieldType,
76  validBoundary
77  );
78 
79  return insertPatchi;
80 }
81 
82 
83 void Foam::fvMeshTools::setPatchFields
84 (
85  fvMesh& mesh,
86  const label patchi,
87  const dictionary& patchFieldDict
88 )
89 {
90  setPatchFields<volScalarField>(mesh, patchi, patchFieldDict);
91  setPatchFields<volVectorField>(mesh, patchi, patchFieldDict);
92  setPatchFields<volSphericalTensorField>(mesh, patchi, patchFieldDict);
93  setPatchFields<volSymmTensorField>(mesh, patchi, patchFieldDict);
94  setPatchFields<volTensorField>(mesh, patchi, patchFieldDict);
95 
96  setPatchFields<surfaceScalarField>(mesh, patchi, patchFieldDict);
97  setPatchFields<surfaceVectorField>(mesh, patchi, patchFieldDict);
98  setPatchFields<surfaceSphericalTensorField>(mesh, patchi, patchFieldDict);
99  setPatchFields<surfaceSymmTensorField>(mesh, patchi, patchFieldDict);
100  setPatchFields<surfaceTensorField>(mesh, patchi, patchFieldDict);
101 
102  if (mesh.foundObject<pointMesh>(pointMesh::typeName))
103  {
104  pointMesh& pm = const_cast<pointMesh&>(pointMesh::New(mesh));
105  setPatchFields<pointScalarField>(pm, patchi, patchFieldDict);
106  setPatchFields<pointVectorField>(pm, patchi, patchFieldDict);
107  setPatchFields<pointSphericalTensorField>(pm, patchi, patchFieldDict);
108  setPatchFields<pointSymmTensorField>(pm, patchi, patchFieldDict);
109  setPatchFields<pointTensorField>(pm, patchi, patchFieldDict);
110  }
111 }
112 
113 
115 {
116  setPatchFields<volScalarField>(mesh, patchi, Zero);
117  setPatchFields<volVectorField>(mesh, patchi, Zero);
118  setPatchFields<volSphericalTensorField>(mesh, patchi, Zero);
119  setPatchFields<volSymmTensorField>(mesh, patchi, Zero);
120  setPatchFields<volTensorField>(mesh, patchi, Zero);
121 
122  setPatchFields<surfaceScalarField>(mesh, patchi, Zero);
123  setPatchFields<surfaceVectorField>(mesh, patchi, Zero);
124  setPatchFields<surfaceSphericalTensorField>(mesh, patchi, Zero);
125  setPatchFields<surfaceSymmTensorField>(mesh, patchi, Zero);
126  setPatchFields<surfaceTensorField>(mesh, patchi, Zero);
127 
128  if (mesh.foundObject<pointMesh>(pointMesh::typeName))
129  {
130  pointMesh& pm = const_cast<pointMesh&>(pointMesh::New(mesh));
131  setPatchFields<pointScalarField>(pm, patchi, Zero);
132  setPatchFields<pointVectorField>(pm, patchi, Zero);
133  setPatchFields<pointSphericalTensorField>(pm, patchi, Zero);
134  setPatchFields<pointSymmTensorField>(pm, patchi, Zero);
135  setPatchFields<pointTensorField>(pm, patchi, Zero);
136  }
137 }
138 
139 
141 (
142  fvMesh& mesh,
143  const labelList& oldToNew,
144  const label nNewPatches,
145  const bool validBoundary
146 )
147 {
148  // Note: oldToNew might have entries beyond nNewPatches so
149  // cannot use invert
150  //const labelList newToOld(invert(nNewPatches, oldToNew));
151  labelList newToOld(nNewPatches, -1);
152  forAll(oldToNew, i)
153  {
154  label newi = oldToNew[i];
155 
156  if (newi >= 0 && newi < nNewPatches)
157  {
158  newToOld[newi] = i;
159  }
160  }
161  mesh.reorderPatches(newToOld, validBoundary);
162 }
163 
164 
165 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches. If validBoundary call is parallel.
Definition: fvMeshTools.C:141
static label addPatch(fvMesh &mesh, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add patch. Inserts patch before all processor patches.
Definition: fvMeshTools.C:34
static void zeroPatchFields(fvMesh &mesh, const label patchi)
Change patchField to zero on registered fields.
Definition: fvMeshTools.C:114
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
virtual void addPatch(const label insertPatchi, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add/insert single patch. If validBoundary the new situation.
Definition: fvMesh.C:1521
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: fvMesh.C:1653
bool foundObject(const word &name) const
Is the named Type in registry.
const word & name() const
Return name.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:52
Foam::polyBoundaryMesh.
label findPatchID(const word &patchName) const
Find patch index given a name.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
A class for handling words, derived from string.
Definition: word.H:62
label patchi
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59