PatchToolsNormals.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) 2011-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 "PatchTools.H"
27 #include "polyMesh.H"
28 #include "indirectPrimitivePatch.H"
29 #include "globalMeshData.H"
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 template<class FaceList, class PointField>
35 (
36  const polyMesh& mesh,
38 )
39 {
40  const globalMeshData& globalData = mesh.globalData();
41  const indirectPrimitivePatch& coupledPatch = globalData.coupledPatch();
42  const Map<label>& coupledPatchMP = coupledPatch.meshPointMap();
43  const mapDistribute& map = globalData.globalPointSlavesMap();
44  const globalIndexAndTransform& transforms =
45  globalData.globalTransforms();
46 
47 
48 
49 
50  // Combine normals. Note: do on all master points. Cannot just use
51  // patch points since the master point does not have to be on the
52  // patch!
53 
54  pointField coupledPointNormals(map.constructSize(), Zero);
55 
56  {
57  // Collect local pointFaces (sized on patch points only)
58  List<List<point>> pointFaceNormals(map.constructSize());
59  forAll(p.meshPoints(), patchPointi)
60  {
61  label meshPointi = p.meshPoints()[patchPointi];
62  Map<label>::const_iterator fnd = coupledPatchMP.find(meshPointi);
63  if (fnd != coupledPatchMP.end())
64  {
65  label coupledPointi = fnd();
66 
67  List<point>& pNormals = pointFaceNormals[coupledPointi];
68  const labelList& pFaces = p.pointFaces()[patchPointi];
69  pNormals.setSize(pFaces.size());
70  forAll(pFaces, i)
71  {
72  pNormals[i] = p.faceNormals()[pFaces[i]];
73  }
74  }
75  }
76 
77 
78  // Pull remote data into local slots
79  map.distribute
80  (
81  transforms,
82  pointFaceNormals,
84  );
85 
86 
87  // Combine all face normals (-local, -remote,untransformed,
88  // -remote,transformed)
89 
90  const labelListList& slaves = globalData.globalPointSlaves();
91  const labelListList& transformedSlaves =
92  globalData.globalPointTransformedSlaves();
93 
94  forAll(slaves, coupledPointi)
95  {
96  const labelList& slaveSlots = slaves[coupledPointi];
97  const labelList& transformedSlaveSlots =
98  transformedSlaves[coupledPointi];
99 
100  point& n = coupledPointNormals[coupledPointi];
101 
102  // Local entries
103  const List<point>& local = pointFaceNormals[coupledPointi];
104 
105  label nFaces =
106  local.size()
107  + slaveSlots.size()
108  + transformedSlaveSlots.size();
109 
110  n = sum(local);
111 
112  // Add any remote face normals
113  forAll(slaveSlots, i)
114  {
115  n += sum(pointFaceNormals[slaveSlots[i]]);
116  }
117  forAll(transformedSlaveSlots, i)
118  {
119  n += sum(pointFaceNormals[transformedSlaveSlots[i]]);
120  }
121 
122  if (nFaces >= 1)
123  {
124  n /= mag(n)+vSmall;
125  }
126 
127  // Put back into slave slots
128  forAll(slaveSlots, i)
129  {
130  coupledPointNormals[slaveSlots[i]] = n;
131  }
132  forAll(transformedSlaveSlots, i)
133  {
134  coupledPointNormals[transformedSlaveSlots[i]] = n;
135  }
136  }
137 
138 
139  // Send back
141  (
142  transforms,
143  coupledPointNormals.size(),
144  coupledPointNormals,
146  );
147  }
148 
149 
150  // 1. Start off with local normals (note:without calculating pointNormals
151  // to avoid them being stored)
152 
153  tmp<pointField> textrudeN(new pointField(p.nPoints(), Zero));
154  pointField& extrudeN = textrudeN.ref();
155  {
156  const faceList& localFaces = p.localFaces();
157  const vectorField& faceNormals = p.faceNormals();
158 
159  forAll(localFaces, facei)
160  {
161  const face& f = localFaces[facei];
162  const vector& n = faceNormals[facei];
163  forAll(f, fp)
164  {
165  extrudeN[f[fp]] += n;
166  }
167  }
168  extrudeN /= mag(extrudeN)+vSmall;
169  }
170 
171 
172  // 2. Override patch normals on coupled points
173  forAll(p.meshPoints(), patchPointi)
174  {
175  label meshPointi = p.meshPoints()[patchPointi];
176  Map<label>::const_iterator fnd = coupledPatchMP.find(meshPointi);
177  if (fnd != coupledPatchMP.end())
178  {
179  label coupledPointi = fnd();
180  extrudeN[patchPointi] = coupledPointNormals[coupledPointi];
181  }
182  }
183 
184  return textrudeN;
185 }
186 
187 
188 template<class FaceList, class PointField>
190 (
191  const polyMesh& mesh,
193  const labelList& patchEdges,
194  const labelList& coupledEdges
195 )
196 {
197  // 1. Start off with local normals
198 
199  tmp<pointField> tedgeNormals(new pointField(p.nEdges(), Zero));
200  pointField& edgeNormals = tedgeNormals.ref();
201  {
202  const labelListList& edgeFaces = p.edgeFaces();
203  const vectorField& faceNormals = p.faceNormals();
204 
205  forAll(edgeFaces, edgeI)
206  {
207  const labelList& eFaces = edgeFaces[edgeI];
208  forAll(eFaces, i)
209  {
210  edgeNormals[edgeI] += faceNormals[eFaces[i]];
211  }
212  }
213  edgeNormals /= mag(edgeNormals)+vSmall;
214  }
215 
216 
217 
218  const globalMeshData& globalData = mesh.globalData();
219  const mapDistribute& map = globalData.globalEdgeSlavesMap();
220 
221 
222  // Convert patch-edge data into cpp-edge data
223  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
224 
225  //- Construct with all data in consistent orientation
226  pointField cppEdgeData(map.constructSize(), Zero);
227 
228  forAll(patchEdges, i)
229  {
230  label patchEdgeI = patchEdges[i];
231  label coupledEdgeI = coupledEdges[i];
232  cppEdgeData[coupledEdgeI] = edgeNormals[patchEdgeI];
233  }
234 
235 
236  // Synchronise
237  // ~~~~~~~~~~~
238 
239  globalData.syncData
240  (
241  cppEdgeData,
242  globalData.globalEdgeSlaves(),
243  globalData.globalEdgeTransformedSlaves(),
244  map,
245  globalData.globalTransforms(),
246  plusEqOp<point>(), // add since normalised later on
248  );
249  cppEdgeData /= mag(cppEdgeData)+vSmall;
250 
251 
252  // Back from cpp-edge to patch-edge data
253  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
254 
255  forAll(patchEdges, i)
256  {
257  label patchEdgeI = patchEdges[i];
258  label coupledEdgeI = coupledEdges[i];
259  edgeNormals[patchEdgeI] = cppEdgeData[coupledEdgeI];
260  }
261 
262  return tedgeNormals;
263 }
264 
265 
266 // ************************************************************************* //
label nPoints() const
Return number of points supporting patch faces.
const labelListList & globalPointSlaves() const
#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 labelListList & globalPointTransformedSlaves() const
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
const mapDistribute & globalEdgeSlavesMap() const
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
const mapDistribute & globalPointSlavesMap() const
const labelListList & globalEdgeSlaves() const
static tmp< pointField > edgeNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &, const labelList &patchEdges, const labelList &coupledEdges)
Return parallel consistent edge normals for patches using mesh points.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
A list of faces which address into the list of points.
const Map< label > & meshPointMap() const
Mesh point map. Given the global point index find its.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Default transformation behaviour.
const globalIndexAndTransform & globalTransforms() const
Global transforms numbering.
const labelListList & edgeFaces() const
Return edge-face addressing.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
static const zero Zero
Definition: zero.H:97
const Field< PointType > & faceNormals() const
Return face normals for patch.
const labelListList & pointFaces() const
Return point-face addressing.
label nEdges() const
Return number of edges in patch.
labelList f(nPoints)
static tmp< pointField > pointNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &)
Return parallel consistent point normals for patches using mesh points.
void reverseDistribute(const label constructSize, List< T > &, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
static void syncData(List< Type > &pointData, const labelListList &slaves, const labelListList &transformedSlaves, const mapDistribute &slavesMap, const globalIndexAndTransform &, const CombineOp &cop, const TransformOp &top)
Helper: synchronise data with transforms.
label constructSize() const
Constructed data size.
void setSize(const label)
Reset size of List.
Definition: List.C:281
Class containing processor-to-processor mapping information.
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: PtrList.H:53
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< small) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
Determination and storage of the possible independent transforms introduced by coupledPolyPatches, as well as all of the possible permutations of these transforms generated by the presence of multiple coupledPolyPatches, i.e. more than one cyclic boundary. Note that any given point can be on maximum 3 transforms only (and these transforms have to be perpendicular)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
const labelListList & globalEdgeTransformedSlaves() const