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-2019 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 "pointFields.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 // Adds patch if not yet there. Returns patchID.
33 (
34  fvMesh& mesh,
35  const polyPatch& patch,
36  const dictionary& patchFieldDict,
37  const word& defaultPatchFieldType,
38  const bool validBoundary
39 )
40 {
41  polyBoundaryMesh& polyPatches =
42  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
43 
44  label patchi = polyPatches.findPatchID(patch.name());
45  if (patchi != -1)
46  {
47  // Already there
48  return patchi;
49  }
50 
51 
52  // Append at end unless there are processor patches
53  label insertPatchi = polyPatches.size();
54 
55  if (!isA<processorPolyPatch>(patch))
56  {
57  forAll(polyPatches, patchi)
58  {
59  const polyPatch& pp = polyPatches[patchi];
60 
61  if (isA<processorPolyPatch>(pp))
62  {
63  insertPatchi = patchi;
64  break;
65  }
66  }
67  }
68 
69  mesh.addPatch
70  (
71  insertPatchi,
72  patch,
73  patchFieldDict,
74  defaultPatchFieldType,
75  validBoundary
76  );
77 
78  return insertPatchi;
79 }
80 
81 
82 void Foam::fvMeshTools::setPatchFields
83 (
84  fvMesh& mesh,
85  const label patchi,
86  const dictionary& patchFieldDict
87 )
88 {
89  setPatchFields<volScalarField>(mesh, patchi, patchFieldDict);
90  setPatchFields<volVectorField>(mesh, patchi, patchFieldDict);
91  setPatchFields<volSphericalTensorField>(mesh, patchi, patchFieldDict);
92  setPatchFields<volSymmTensorField>(mesh, patchi, patchFieldDict);
93  setPatchFields<volTensorField>(mesh, patchi, patchFieldDict);
94 
95  setPatchFields<surfaceScalarField>(mesh, patchi, patchFieldDict);
96  setPatchFields<surfaceVectorField>(mesh, patchi, patchFieldDict);
97  setPatchFields<surfaceSphericalTensorField>(mesh, patchi, patchFieldDict);
98  setPatchFields<surfaceSymmTensorField>(mesh, patchi, patchFieldDict);
99  setPatchFields<surfaceTensorField>(mesh, patchi, patchFieldDict);
100 
101  if (mesh.foundObject<pointMesh>(pointMesh::typeName))
102  {
103  pointMesh& pm = const_cast<pointMesh&>(pointMesh::New(mesh));
104  setPatchFields<pointScalarField>(pm, patchi, patchFieldDict);
105  setPatchFields<pointVectorField>(pm, patchi, patchFieldDict);
106  setPatchFields<pointSphericalTensorField>(pm, patchi, patchFieldDict);
107  setPatchFields<pointSymmTensorField>(pm, patchi, patchFieldDict);
108  setPatchFields<pointTensorField>(pm, patchi, patchFieldDict);
109  }
110 }
111 
112 
114 {
115  setPatchFields<volScalarField>(mesh, patchi, Zero);
116  setPatchFields<volVectorField>(mesh, patchi, Zero);
117  setPatchFields<volSphericalTensorField>(mesh, patchi, Zero);
118  setPatchFields<volSymmTensorField>(mesh, patchi, Zero);
119  setPatchFields<volTensorField>(mesh, patchi, Zero);
120 
121  setPatchFields<surfaceScalarField>(mesh, patchi, Zero);
122  setPatchFields<surfaceVectorField>(mesh, patchi, Zero);
123  setPatchFields<surfaceSphericalTensorField>(mesh, patchi, Zero);
124  setPatchFields<surfaceSymmTensorField>(mesh, patchi, Zero);
125  setPatchFields<surfaceTensorField>(mesh, patchi, Zero);
126 
127  if (mesh.foundObject<pointMesh>(pointMesh::typeName))
128  {
129  pointMesh& pm = const_cast<pointMesh&>(pointMesh::New(mesh));
130  setPatchFields<pointScalarField>(pm, patchi, Zero);
131  setPatchFields<pointVectorField>(pm, patchi, Zero);
132  setPatchFields<pointSphericalTensorField>(pm, patchi, Zero);
133  setPatchFields<pointSymmTensorField>(pm, patchi, Zero);
134  setPatchFields<pointTensorField>(pm, patchi, Zero);
135  }
136 }
137 
138 
140 (
141  fvMesh& mesh,
142  const labelList& oldToNew,
143  const label nNewPatches,
144  const bool validBoundary
145 )
146 {
147  // Note: oldToNew might have entries beyond nNewPatches so
148  // cannot use invert
149  //const labelList newToOld(invert(nNewPatches, oldToNew));
150  labelList newToOld(nNewPatches, -1);
151  forAll(oldToNew, i)
152  {
153  label newi = oldToNew[i];
154 
155  if (newi >= 0 && newi < nNewPatches)
156  {
157  newToOld[newi] = i;
158  }
159  }
160  mesh.reorderPatches(newToOld, validBoundary);
161 }
162 
163 
164 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
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:140
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
const word & name() const
Return name.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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:855
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: fvMesh.C:987
label findPatchID(const word &patchName) const
Find patch index given a name.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
static void zeroPatchFields(fvMesh &mesh, const label patchi)
Change patchField to zero on registered fields.
Definition: fvMeshTools.C:113
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
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:33
static const zero Zero
Definition: zero.H:97
Foam::polyBoundaryMesh.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66