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-2024 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 
33 (
34  fvMesh& mesh,
35  const polyPatch& patch
36 )
37 {
38  polyBoundaryMesh& polyPatches =
39  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
40 
41  // Check if it is already there
42  label patchi = polyPatches.findIndex(patch.name());
43  if (patchi != -1)
44  {
45  return patchi;
46  }
47 
48  // Append at the end ...
49  label insertPatchi = polyPatches.size();
50 
51  // ... unless there are processor patches, in which case append after the
52  // last global patch
53  if (!isA<processorPolyPatch>(patch))
54  {
55  forAll(polyPatches, patchi)
56  {
57  const polyPatch& pp = polyPatches[patchi];
58 
59  if (isA<processorPolyPatch>(pp))
60  {
61  insertPatchi = patchi;
62  break;
63  }
64  }
65  }
66 
67  mesh.addPatch(insertPatchi, patch);
68 
69  return insertPatchi;
70 }
71 
72 
74 {
76 }
77 
78 
79 void Foam::fvMeshTools::setPatchFields
80 (
81  fvMesh& mesh,
82  const label patchi,
83  const dictionary& patchFieldDict
84 )
85 {
86  #define SetPatchFieldsType(Type, FieldType, Mesh) \
87  setPatchFields<FieldType<Type>>(Mesh, patchi, patchFieldDict);
90  if (mesh.foundObject<pointMesh>(pointMesh::typeName))
91  {
93  (
95  PointField,
97  );
98  }
99  #undef SetPatchFieldsType
100 }
101 
102 
104 (
105  fvMesh& mesh,
106  const labelList& oldToNew,
107  const label nNewPatches,
108  const bool validBoundary
109 )
110 {
111  // Note: oldToNew might have entries beyond nNewPatches so
112  // cannot use `invert(nNewPatches, oldToNew)`
113 
114  labelList newToOld(nNewPatches, -1);
115  forAll(oldToNew, i)
116  {
117  label newi = oldToNew[i];
118  if (newi >= 0 && newi < nNewPatches)
119  {
120  newToOld[newi] = i;
121  }
122  }
123 
124  mesh.reorderPatches(newToOld, validBoundary);
125 }
126 
127 
128 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Generic GeometricField class.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches. If validBoundary call is.
Definition: fvMeshTools.C:104
static label addPatch(fvMesh &mesh, const polyPatch &patch)
Add patch. Inserts patch before all processor patches. Returns the.
Definition: fvMeshTools.C:33
static void addedPatches(fvMesh &mesh)
Complete adding patches.
Definition: fvMeshTools.C:73
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: fvMesh.C:1702
virtual void addPatch(const label insertPatchi, const polyPatch &patch)
Add/insert single patch.
Definition: fvMesh.C:1645
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:53
Foam::polyBoundaryMesh.
label findIndex(const word &patchName) const
Find patch index given a name.
void addedPatches()
Complete addition of single patches.
Definition: polyMesh.C:1319
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define SetPatchFieldsType(Type, FieldType, Mesh)
label patchi
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
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)