topoCellLooper.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "topoCellLooper.H"
29 #include "cellFeatures.H"
30 #include "polyMesh.H"
31 #include "unitConversion.H"
32 #include "DynamicList.H"
33 #include "ListOps.H"
34 #include "meshTools.H"
35 #include "hexMatcher.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(topoCellLooper, 0);
43  addToRunTimeSelectionTable(cellLooper, topoCellLooper, word);
44 }
45 
46 // Angle for polys to be considered splitHexes.
47 const Foam::scalar Foam::topoCellLooper::featureCos = Foam::cos(degToRad(10.0));
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 // In-memory truncate a list
53 template<class T>
54 void Foam::topoCellLooper::subsetList
55 (
56  const label startI,
57  const label freeI,
58  DynamicList<T>& lst
59 )
60 {
61  if (startI == 0)
62  {
63  // Truncate (setSize decides itself not to do anything if nothing
64  // changed)
65  if (freeI < 0)
66  {
68  << " lst:" << lst << abort(FatalError);
69  }
70  lst.setCapacity(freeI);
71  }
72  else
73  {
74  // Shift elements down
75  label newI = 0;
76  for (label elemI = startI; elemI < freeI; elemI++)
77  {
78  lst[newI++] = lst[elemI];
79  }
80 
81  if ((freeI - startI) < 0)
82  {
84  << " lst:" << lst << abort(FatalError);
85  }
86 
87  lst.setCapacity(freeI - startI);
88  }
89 }
90 
91 
92 // Returns label of edge nFeaturePts away from startEdge (in the direction of
93 // startVertI) and not counting non-featurePoints.
94 //
95 // When stepping to this face it can happen in 3 different ways:
96 //
97 // --|------
98 // |
99 // 1 | 0
100 // |A
101 // |
102 // |
103 // --|------
104 //
105 // A: jumping from face0 to face1 across edge A.
106 // startEdge != -1
107 // startVert == -1
108 //
109 // --|------
110 // |
111 // 1 | 0
112 // +B
113 // |
114 // |
115 // --|------
116 //
117 // B: jumping from face0 to face1 across (non-feature) point B
118 // startEdge == -1
119 // startVert != -1
120 //
121 // --|------
122 // 0 | 1
123 // |C
124 // --+
125 // |
126 // |
127 // --|------
128 //
129 // C: jumping from face0 to face1 across (feature) point C on edge.
130 // startEdge != -1
131 // startVert != -1
132 //
133 void Foam::topoCellLooper::walkFace
134 (
135  const cellFeatures& features,
136  const label facei,
137  const label startEdgeI,
138  const label startVertI,
139  const label nFeaturePts,
140 
141  label& edgeI,
142  label& vertI
143 ) const
144 {
145  const labelList& fEdges = mesh().faceEdges()[facei];
146 
147  edgeI = startEdgeI;
148 
149  vertI = startVertI;
150 
151  // Number of feature points crossed so far
152  label nVisited = 0;
153 
154  if (vertI == -1)
155  {
156  // Started on edge. Go to one of its endpoints.
157  vertI = mesh().edges()[edgeI].start();
158 
159  if (features.isFeatureVertex(facei, vertI))
160  {
161  nVisited++;
162  }
163  }
164 
165  if ((edgeI == -1) || !meshTools::edgeOnFace(mesh(), facei, edgeI))
166  {
167  // Either edge is not set or not on current face. Just take one of
168  // the edges on this face as starting edge.
169  edgeI = getFirstVertEdge(facei, vertI);
170  }
171 
172  // Now we should have starting edge on face and a vertex on that edge.
173 
174  do
175  {
176 
177  edgeI = meshTools::otherEdge(mesh(), fEdges, edgeI, vertI);
178 
179  if (nVisited == nFeaturePts)
180  {
181  break;
182  }
183 
184  vertI = mesh().edges()[edgeI].otherVertex(vertI);
185 
186  if (features.isFeatureVertex(facei, vertI))
187  {
188  nVisited++;
189  }
190  }
191  while (true);
192 }
193 
194 
195 // Returns list of vertices on 'superEdge' i.e. list of edges connected by
196 // non-feature points. First and last are feature points, ones inbetween are
197 // not.
198 Foam::labelList Foam::topoCellLooper::getSuperEdge
199 (
200  const cellFeatures& features,
201  const label facei,
202  const label startEdgeI,
203  const label startVertI
204 ) const
205 {
206  const labelList& fEdges = mesh().faceEdges()[facei];
207 
208  labelList superVerts(fEdges.size());
209  label superVertI = 0;
210 
211 
212  label edgeI = startEdgeI;
213 
214  label vertI = startVertI;
215 
216  superVerts[superVertI++] = vertI;
217 
218  label prevEdgeI = -1;
219 
220  do
221  {
222  vertI = mesh().edges()[edgeI].otherVertex(vertI);
223 
224  superVerts[superVertI++] = vertI;
225 
226  prevEdgeI = edgeI;
227 
228  edgeI = meshTools::otherEdge(mesh(), fEdges, edgeI, vertI);
229  }
230  while (!features.isFeaturePoint(prevEdgeI, edgeI));
231 
232  superVerts.setSize(superVertI);
233 
234  return superVerts;
235 }
236 
237 
238 // Return non-feature edge from cells' edges that is most perpendicular
239 // to refinement direction.
240 Foam::label Foam::topoCellLooper::getAlignedNonFeatureEdge
241 (
242  const vector& refDir,
243  const label celli,
244  const cellFeatures& features
245 ) const
246 {
247  const labelList& cEdges = mesh().cellEdges()[celli];
248 
249  const point& ctr = mesh().cellCentres()[celli];
250 
251  label cutEdgeI = -1;
252  scalar maxCos = -GREAT;
253 
254  forAll(cEdges, cEdgeI)
255  {
256  label edgeI = cEdges[cEdgeI];
257 
258  if (!features.isFeatureEdge(edgeI))
259  {
260  const edge& e = mesh().edges()[edgeI];
261 
262  // Get plane spanned by e.start, e.end and cell centre.
263  vector e0 = mesh().points()[e.start()] - ctr;
264  vector e1 = mesh().points()[e.end()] - ctr;
265 
266  const vector n = normalised(e0 ^ e1);
267 
268  scalar cosAngle = mag(refDir & n);
269 
270  if (cosAngle > maxCos)
271  {
272  maxCos = cosAngle;
273 
274  cutEdgeI = edgeI;
275  }
276  }
277  }
278 
279  return cutEdgeI;
280 }
281 
282 
283 // Starts from edge and vertex on edge on face (or neighbouring face)
284 // and steps either to existing vertex (vertI != -1) or to edge (vertI == -1)
285 // by walking point-edge and crossing nFeats featurePoints.
286 void Foam::topoCellLooper::walkAcrossFace
287 (
288  const cellFeatures& features,
289  const label facei,
290  const label startEdgeI,
291  const label startVertI,
292  const label nFeats,
293 
294  label& edgeI,
295  label& vertI
296 ) const
297 {
298  label oppositeVertI = -1;
299  label oppositeEdgeI = -1;
300 
301  // Go to oppositeEdge and a vertex on that.
302  walkFace
303  (
304  features,
305  facei,
306  startEdgeI,
307  startVertI,
308  nFeats,
309 
310  oppositeEdgeI,
311  oppositeVertI
312  );
313 
314  // Loop over super edge to find internal points if there are any.
315 
316  labelList superEdge =
317  getSuperEdge
318  (
319  features,
320  facei,
321  oppositeEdgeI,
322  oppositeVertI
323  );
324 
325  label sz = superEdge.size();
326 
327  if (sz == 2)
328  {
329  // No non-feature point inbetween feature points.
330  // Mark edge.
331 
332  vertI = -1;
333  edgeI = oppositeEdgeI;
334  }
335  else if (sz == 3)
336  {
337  vertI = superEdge[1];
338  edgeI = -1;
339  }
340  else
341  {
342  //Should choose acc. to geometry!
343  label index = sz/2;
344 
345  if (debug)
346  {
347  Pout<< " Don't know what to do. Stepped to non-feature point "
348  << "at index " << index << " in superEdge:" << superEdge
349  << endl;
350  }
351 
352  vertI = superEdge[index];
353  edgeI = -1;
354  }
355 }
356 
357 
358 // Walks cell circumference. Updates face, edge, vertex.
359 //
360 // Position on face is given by:
361 //
362 // vertI == -1, facei != -1, edgeI != -1
363 // on edge of face. Cross edge to neighbouring face.
364 //
365 // vertI != -1, edgeI != -1, facei == -1
366 // coming from edge onto vertex vertI. Need to step to one
367 // of the faces not using edgeI.
368 //
369 // vertI != -1, edgeI == -1, facei != -1
370 // coming from vertex on side of face. Step to one of the faces
371 // using vertI but not facei
372 void Foam::topoCellLooper::walkSplitHex
373 (
374  const label celli,
375  const cellFeatures& features,
376  const label fromFacei,
377  const label fromEdgeI,
378  const label fromVertI,
379 
380  DynamicList<label>& loop,
381  DynamicList<scalar>& loopWeights
382 ) const
383 {
384  // Work vars giving position on cell
385  label facei = fromFacei;
386  label edgeI = fromEdgeI;
387  label vertI = fromVertI;
388 
389  do
390  {
391  if (debug)
392  {
393  Pout<< "Entering walk with : cell:" << celli << " face:" << facei;
394  if (facei != -1)
395  {
396  Pout<< " verts:" << mesh().faces()[facei];
397  }
398  Pout<< " edge:" << edgeI;
399  if (edgeI != -1)
400  {
401  Pout<< " verts:" << mesh().edges()[edgeI];
402  }
403  Pout<< " vert:" << vertI << endl;
404  }
405 
406  label startLoop = -1;
407 
408  if
409  (
410  vertI != -1
411  && (startLoop = loop.find(vertToEVert(vertI))) != -1
412  )
413  {
414  // Breaking walk since vertI already cut
415  label firstFree = loop.size();
416 
417  subsetList(startLoop, firstFree, loop);
418  subsetList(startLoop, firstFree, loopWeights);
419 
420  break;
421  }
422  if
423  (
424  edgeI != -1
425  && (startLoop = loop.find(edgeToEVert(edgeI))) != -1
426  )
427  {
428  // Breaking walk since edgeI already cut
429  label firstFree = loop.size();
430 
431  subsetList(startLoop, firstFree, loop);
432  subsetList(startLoop, firstFree, loopWeights);
433 
434  break;
435  }
436 
437 
438  if (vertI == -1)
439  {
440  // On edge
441  if (edgeI == -1)
442  {
444  << abort(FatalError);
445  }
446 
447  loop.append(edgeToEVert(edgeI));
448  loopWeights.append(0.5);
449 
450  // Cross edge to next face
451  facei = meshTools::otherFace(mesh(), celli, facei, edgeI);
452 
453  if (debug)
454  {
455  Pout<< " stepped across edge " << mesh().edges()[edgeI]
456  << " to face " << facei << " verts:"
457  << mesh().faces()[facei] << endl;
458  }
459 
460  label nextEdgeI = -1;
461  label nextVertI = -1;
462 
463  // Walk face along its edges
464  walkAcrossFace
465  (
466  features,
467  facei,
468  edgeI,
469  vertI,
470  2,
471 
472  nextEdgeI,
473  nextVertI
474  );
475 
476  edgeI = nextEdgeI;
477  vertI = nextVertI;
478  }
479  else
480  {
481  // On vertex.
482 
483  loop.append(vertToEVert(vertI));
484  loopWeights.append(-GREAT);
485 
486  if (edgeI == -1)
487  {
488  // Normal vertex on edge of face. Get edges connected to it
489  // which are not on facei.
490  labelList nextEdges = getVertEdgesNonFace
491  (
492  celli,
493  facei,
494  vertI
495  );
496 
497  if (nextEdges.empty())
498  {
499  // Cross to other face (there is only one since no edges)
500  const labelList& pFaces = mesh().pointFaces()[vertI];
501 
502  forAll(pFaces, pFacei)
503  {
504  label thisFacei = pFaces[pFacei];
505 
506  if
507  (
508  (thisFacei != facei)
509  && meshTools::faceOnCell(mesh(), celli, thisFacei)
510  )
511  {
512  facei = thisFacei;
513  break;
514  }
515  }
516 
517  if (debug)
518  {
519  Pout<< " stepped from non-edge vertex " << vertI
520  << " to face " << facei << " verts:"
521  << mesh().faces()[facei]
522  << " since candidate edges:" << nextEdges << endl;
523  }
524 
525  label nextEdgeI = -1;
526  label nextVertI = -1;
527 
528  walkAcrossFace
529  (
530  features,
531  facei,
532  edgeI,
533  vertI,
534  2, // 2 vertices to cross
535 
536  nextEdgeI,
537  nextVertI
538  );
539 
540  edgeI = nextEdgeI;
541  vertI = nextVertI;
542  }
543  else if (nextEdges.size() == 1)
544  {
545  // One edge. Go along this one.
546  edgeI = nextEdges[0];
547 
548  if (debug)
549  {
550  Pout<< " stepped from non-edge vertex " << vertI
551  << " along edge " << edgeI << " verts:"
552  << mesh().edges()[edgeI]
553  << " out of candidate edges:"
554  << nextEdges << endl;
555  }
556 
557  vertI = mesh().edges()[edgeI].otherVertex(vertI);
558 
559  facei = -1;
560  }
561  else
562  {
563  // Multiple edges to choose from. Pick any one.
564  // (ideally should be geometric)
565 
566  label index = nextEdges.size()/2;
567 
568  edgeI = nextEdges[index];
569 
570  if (debug)
571  {
572  Pout<< " stepped from non-edge vertex " << vertI
573  << " along edge " << edgeI << " verts:"
574  << mesh().edges()[edgeI]
575  << " out of candidate edges:" << nextEdges << endl;
576  }
577 
578  vertI = mesh().edges()[edgeI].otherVertex(vertI);
579 
580  facei = -1;
581  }
582  }
583  else
584  {
585  // Get faces connected to startVertI but not startEdgeI
586  labelList nextFaces =
587  getVertFacesNonEdge
588  (
589  celli,
590  edgeI,
591  vertI
592  );
593 
594  if (nextFaces.size() == 1)
595  {
596  // Only one face to cross.
597  facei = nextFaces[0];
598 
599  label nextEdgeI = -1;
600  label nextVertI = -1;
601 
602  walkAcrossFace
603  (
604  features,
605  facei,
606  edgeI,
607  vertI,
608  2, // 2 vertices to cross
609 
610  nextEdgeI,
611  nextVertI
612  );
613 
614  edgeI = nextEdgeI;
615  vertI = nextVertI;
616  }
617  else if (nextFaces.size() == 2)
618  {
619  // Split face. Get edge inbetween.
620  facei = -1;
621 
622  edgeI =
624  (
625  mesh(),
626  nextFaces[0],
627  nextFaces[1]
628  );
629 
630  vertI = mesh().edges()[edgeI].otherVertex(vertI);
631  }
632  else
633  {
635  << "Choosing from more than "
636  << "two candidates:" << nextFaces
637  << " when coming from vertex " << vertI << " on cell "
638  << celli << abort(FatalError);
639  }
640  }
641  }
642 
643  if (debug)
644  {
645  Pout<< "Walked to : face:" << facei;
646  if (facei != -1)
647  {
648  Pout<< " verts:" << mesh().faces()[facei];
649  }
650  Pout<< " edge:" << edgeI;
651  if (edgeI != -1)
652  {
653  Pout<< " verts:" << mesh().edges()[edgeI];
654  }
655  Pout<< " vert:" << vertI << endl;
656  }
657  }
658  while (true);
659 }
660 
661 
662 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
663 
664 Foam::topoCellLooper::topoCellLooper(const polyMesh& mesh)
665 :
667 {}
668 
669 
670 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
671 
673 (
674  const vector& refDir,
675  const label celli,
676  const boolList& vertIsCut,
677  const boolList& edgeIsCut,
678  const scalarField& edgeWeight,
679 
680  labelList& loop,
681  scalarField& loopWeights
682 ) const
683 {
684  if (mesh().cellShapes()[celli].model() == hex_)
685  {
686  // Let parent handle hex case.
687  return
689  (
690  refDir,
691  celli,
692  vertIsCut,
693  edgeIsCut,
694  edgeWeight,
695  loop,
696  loopWeights
697  );
698  }
699  else
700  {
701  cellFeatures superCell(mesh(), featureCos, celli);
702 
703  if (hexMatcher().isA(superCell.faces()))
704  {
705  label edgeI =
706  getAlignedNonFeatureEdge
707  (
708  refDir,
709  celli,
710  superCell
711  );
712 
713  label vertI = -1;
714 
715  label facei = -1;
716 
717  if (edgeI != -1)
718  {
719  // Found non-feature edge. Start walking from vertex on edge.
720  vertI = mesh().edges()[edgeI].start();
721  }
722  else
723  {
724  // No 'matching' non-feature edge found on cell. Get starting
725  // normal i.e. feature edge.
726  edgeI = getMisAlignedEdge(refDir, celli);
727 
728  // Get any face using edge
729  label face0;
730  label face1;
731  meshTools::getEdgeFaces(mesh(), celli, edgeI, face0, face1);
732 
733  facei = face0;
734  }
735 
736  label nEstCuts = 2*mesh().cells()[celli].size();
737 
738  DynamicList<label> localLoop(nEstCuts);
739  DynamicList<scalar> localLoopWeights(nEstCuts);
740 
741  walkSplitHex
742  (
743  celli,
744  superCell,
745  facei,
746  edgeI,
747  vertI,
748 
749  localLoop,
750  localLoopWeights
751  );
752 
753  if (localLoop.size() <=2)
754  {
755  return false;
756  }
757  else
758  {
759  loop.transfer(localLoop);
760  loopWeights.transfer(localLoopWeights);
761 
762  return true;
763  }
764  }
765  else
766  {
767  // Let parent handle poly case.
768  return hexCellLooper::cut
769  (
770  refDir,
771  celli,
772  vertIsCut,
773  edgeIsCut,
774  edgeWeight,
775  loop,
776  loopWeights
777  );
778  }
779  }
780 }
781 
782 
784 (
785  const plane& cutPlane,
786  const label celli,
787  const boolList& vertIsCut,
788  const boolList& edgeIsCut,
789  const scalarField& edgeWeight,
790 
791  labelList& loop,
792  scalarField& loopWeights
793 ) const
794 {
795  // Let parent handle cut with plane.
796  return
798  (
799  cutPlane,
800  celli,
801  vertIsCut,
802  edgeIsCut,
803  edgeWeight,
804 
805  loop,
806  loopWeights
807  );
808 }
809 
810 
811 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
meshTools.H
Foam::meshTools::getSharedEdge
label getSharedEdge(const primitiveMesh &mesh, const label f0, const label f1)
Return edge shared by two faces. Throws error if no edge found.
Definition: meshTools.C:408
Foam::topoCellLooper::featureCos
static const scalar featureCos
Cos of angle for feature recognition (of splitHexes)
Definition: topoCellLooper.H:163
Foam::meshTools::otherFace
label otherFace(const primitiveMesh &mesh, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:555
Foam::DynamicList< label >
hexMatcher.H
Foam::primitiveMesh::pointFaces
const labelListList & pointFaces() const
Definition: primitiveMeshPointFaces.C:34
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::meshTools::walkFace
label walkFace(const primitiveMesh &mesh, const label facei, const label startEdgeI, const label startVertI, const label nEdges)
Returns label of edge nEdges away from startEdge (in the direction.
Definition: meshTools.C:603
Foam::plane
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:89
Foam::meshTools::faceOnCell
bool faceOnCell(const primitiveMesh &mesh, const label celli, const label facei)
Is face used by cell.
Definition: meshTools.C:330
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::hexMatcher
A cellMatcher for hex cells (cellModel::HEX).
Definition: hexMatcher.H:53
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::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
Foam::Field< scalar >
Foam::primitiveMesh::cellEdges
const labelListList & cellEdges() const
Definition: primitiveMeshCellEdges.C:120
Foam::cellFeatures::faces
const faceList & faces() const
Definition: cellFeatures.H:141
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:459
Foam::primitiveMesh::faceEdges
const labelListList & faceEdges() const
Definition: primitiveMeshEdges.C:528
Foam::meshTools::getEdgeFaces
void getEdgeFaces(const primitiveMesh &mesh, const label celli, const label edgeI, label &face0, label &face1)
Get faces on cell using edgeI. Throws error if no two found.
Definition: meshTools.C:479
Foam::topoCellLooper::cut
virtual bool cut(const vector &refDir, const label celli, const boolList &vertIsCut, const boolList &edgeIsCut, const scalarField &edgeWeight, labelList &loop, scalarField &loopWeights) const
Create cut along circumference of celli. Gets current mesh cuts.
Definition: topoCellLooper.C:673
Foam::FatalError
error FatalError
Foam::cellFeatures
Cell analysis class.
Definition: cellFeatures.H:63
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::meshTools::otherEdge
label otherEdge(const primitiveMesh &mesh, const labelList &edgeLabels, const label thisEdgeI, const label thisVertI)
Return label of other edge (out of candidates edgeLabels)
Definition: meshTools.C:521
cellFeatures.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:483
Foam::meshTools::edgeOnFace
bool edgeOnFace(const primitiveMesh &mesh, const label facei, const label edgeI)
Is edge used by face.
Definition: meshTools.C:319
Foam::isA
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Definition: typeInfo.H:197
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:175
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:102
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
topoCellLooper.H
DynamicList.H
cellShapes
cellShapeList cellShapes
Definition: createBlockMesh.H:3
ListOps.H
Various functions to operate on Lists.
Foam::hexCellLooper::cut
virtual bool cut(const vector &refDir, const label celli, const boolList &vertIsCut, const boolList &edgeIsCut, const scalarField &edgeWeight, labelList &loop, scalarField &loopWeights) const
Create cut along circumference of celli. Gets current mesh cuts.
Definition: hexCellLooper.C:162
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::hexCellLooper
Implementation of cellLooper.
Definition: hexCellLooper.H:62
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265
Foam::subsetList
List< T > subsetList(const UList< T > &input, const UnaryPredicate &pred, const bool invert=false)
Copy a subset of the input list when predicate is true.