isoSurfaceTemplates.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 | Copyright (C) 2011-2016 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 "isoSurface.H"
27 #include "polyMesh.H"
28 #include "syncTools.H"
29 #include "surfaceFields.H"
30 #include "OFstream.H"
31 #include "meshTools.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 template<class Type>
37 <
38  Type,
42 >>
43 Foam::isoSurface::adaptPatchFields
44 (
46 ) const
47 {
48  typedef SlicedGeometricField
49  <
50  Type,
53  volMesh
54  > FieldType;
55 
56  tmp<FieldType> tsliceFld
57  (
58  new FieldType
59  (
60  IOobject
61  (
62  fld.name(),
63  fld.instance(),
64  fld.db(),
65  IOobject::NO_READ,
66  IOobject::NO_WRITE,
67  false
68  ),
69  fld, // internal field
70  true // preserveCouples
71  )
72  );
73  FieldType& sliceFld = tsliceFld.ref();
74 
75  const fvMesh& mesh = fld.mesh();
76 
77  const polyBoundaryMesh& patches = mesh.boundaryMesh();
78 
79  typename FieldType::Boundary& sliceFldBf =
80  sliceFld.boundaryFieldRef();
81 
82  forAll(patches, patchi)
83  {
84  const polyPatch& pp = patches[patchi];
85 
86  if
87  (
88  isA<emptyPolyPatch>(pp)
89  && pp.size() != sliceFldBf[patchi].size()
90  )
91  {
92  // Clear old value. Cannot resize it since is a slice.
93  sliceFldBf.set(patchi, NULL);
94 
95  // Set new value we can change
96  sliceFldBf.set
97  (
98  patchi,
100  (
101  mesh.boundary()[patchi],
102  sliceFld
103  )
104  );
105 
106  // Note: cannot use patchInternalField since uses emptyFvPatch::size
107  // Do our own internalField instead.
108  const labelUList& faceCells =
109  mesh.boundary()[patchi].patch().faceCells();
110 
111  Field<Type>& pfld = sliceFldBf[patchi];
112  pfld.setSize(faceCells.size());
113  forAll(faceCells, i)
114  {
115  pfld[i] = sliceFld[faceCells[i]];
116  }
117  }
118  else if (isA<cyclicPolyPatch>(pp))
119  {
120  // Already has interpolate as value
121  }
122  else if (isA<processorPolyPatch>(pp))
123  {
124  fvPatchField<Type>& pfld = const_cast<fvPatchField<Type>&>
125  (
126  sliceFldBf[patchi]
127  );
128 
129  const scalarField& w = mesh.weights().boundaryField()[patchi];
130 
131  tmp<Field<Type>> f =
132  w*pfld.patchInternalField()
133  + (1.0-w)*pfld.patchNeighbourField();
134 
135  PackedBoolList isCollocated
136  (
137  collocatedFaces(refCast<const processorPolyPatch>(pp))
138  );
139 
140  forAll(isCollocated, i)
141  {
142  if (!isCollocated[i])
143  {
144  pfld[i] = f()[i];
145  }
146  }
147  }
148  }
149  return tsliceFld;
150 }
151 
152 
153 template<class Type>
154 Type Foam::isoSurface::generatePoint
155 (
156  const scalar s0,
157  const Type& p0,
158  const bool hasSnap0,
159  const Type& snapP0,
160 
161  const scalar s1,
162  const Type& p1,
163  const bool hasSnap1,
164  const Type& snapP1
165 ) const
166 {
167  scalar d = s1-s0;
168 
169  if (mag(d) > VSMALL)
170  {
171  scalar s = (iso_-s0)/d;
172 
173  if (hasSnap1 && s >= 0.5 && s <= 1)
174  {
175  return snapP1;
176  }
177  else if (hasSnap0 && s >= 0.0 && s <= 0.5)
178  {
179  return snapP0;
180  }
181  else
182  {
183  return s*p1 + (1.0-s)*p0;
184  }
185  }
186  else
187  {
188  scalar s = 0.4999;
189 
190  return s*p1 + (1.0-s)*p0;
191  }
192 }
193 
194 
195 template<class Type>
196 void Foam::isoSurface::generateTriPoints
197 (
198  const scalar s0,
199  const Type& p0,
200  const bool hasSnap0,
201  const Type& snapP0,
202 
203  const scalar s1,
204  const Type& p1,
205  const bool hasSnap1,
206  const Type& snapP1,
207 
208  const scalar s2,
209  const Type& p2,
210  const bool hasSnap2,
211  const Type& snapP2,
212 
213  const scalar s3,
214  const Type& p3,
215  const bool hasSnap3,
216  const Type& snapP3,
217 
218  DynamicList<Type>& points
219 ) const
220 {
221  int triIndex = 0;
222  if (s0 < iso_)
223  {
224  triIndex |= 1;
225  }
226  if (s1 < iso_)
227  {
228  triIndex |= 2;
229  }
230  if (s2 < iso_)
231  {
232  triIndex |= 4;
233  }
234  if (s3 < iso_)
235  {
236  triIndex |= 8;
237  }
238 
239  /* Form the vertices of the triangles for each case */
240  switch (triIndex)
241  {
242  case 0x00:
243  case 0x0F:
244  break;
245 
246  case 0x01:
247  case 0x0E:
248  points.append
249  (
250  generatePoint(s0,p0,hasSnap0,snapP0,s1,p1,hasSnap1,snapP1)
251  );
252  points.append
253  (
254  generatePoint(s0,p0,hasSnap0,snapP0,s2,p2,hasSnap2,snapP2)
255  );
256  points.append
257  (
258  generatePoint(s0,p0,hasSnap0,snapP0,s3,p3,hasSnap3,snapP3)
259  );
260  if (triIndex == 0x0E)
261  {
262  // Flip normals
263  label sz = points.size();
264  Swap(points[sz-2], points[sz-1]);
265  }
266  break;
267 
268  case 0x02:
269  case 0x0D:
270  points.append
271  (
272  generatePoint(s1,p1,hasSnap1,snapP1,s0,p0,hasSnap0,snapP0)
273  );
274  points.append
275  (
276  generatePoint(s1,p1,hasSnap1,snapP1,s3,p3,hasSnap3,snapP3)
277  );
278  points.append
279  (
280  generatePoint(s1,p1,hasSnap1,snapP1,s2,p2,hasSnap2,snapP2)
281  );
282  if (triIndex == 0x0D)
283  {
284  // Flip normals
285  label sz = points.size();
286  Swap(points[sz-2], points[sz-1]);
287  }
288  break;
289 
290  case 0x03:
291  case 0x0C:
292  {
293  Type p0p2 =
294  generatePoint(s0,p0,hasSnap0,snapP0,s2,p2,hasSnap2,snapP2);
295  Type p1p3 =
296  generatePoint(s1,p1,hasSnap1,snapP1,s3,p3,hasSnap3,snapP3);
297 
298  points.append
299  (
300  generatePoint(s0,p0,hasSnap0,snapP0,s3,p3,hasSnap3,snapP3)
301  );
302  points.append(p1p3);
303  points.append(p0p2);
304 
305  points.append(p1p3);
306  points.append
307  (
308  generatePoint(s1,p1,hasSnap1,snapP1,s2,p2,hasSnap2,snapP2)
309  );
310  points.append(p0p2);
311  }
312  if (triIndex == 0x0C)
313  {
314  // Flip normals
315  label sz = points.size();
316  Swap(points[sz-5], points[sz-4]);
317  Swap(points[sz-2], points[sz-1]);
318  }
319  break;
320 
321  case 0x04:
322  case 0x0B:
323  {
324  points.append
325  (
326  generatePoint(s2,p2,hasSnap2,snapP2,s0,p0,hasSnap0,snapP0)
327  );
328  points.append
329  (
330  generatePoint(s2,p2,hasSnap2,snapP2,s1,p1,hasSnap1,snapP1)
331  );
332  points.append
333  (
334  generatePoint(s2,p2,hasSnap2,snapP2,s3,p3,hasSnap3,snapP3)
335  );
336  }
337  if (triIndex == 0x0B)
338  {
339  // Flip normals
340  label sz = points.size();
341  Swap(points[sz-2], points[sz-1]);
342  }
343  break;
344 
345  case 0x05:
346  case 0x0A:
347  {
348  Type p0p1 =
349  generatePoint(s0,p0,hasSnap0,snapP0,s1,p1,hasSnap1,snapP1);
350  Type p2p3 =
351  generatePoint(s2,p2,hasSnap2,snapP2,s3,p3,hasSnap3,snapP3);
352 
353  points.append(p0p1);
354  points.append(p2p3);
355  points.append
356  (
357  generatePoint(s0,p0,hasSnap0,snapP0,s3,p3,hasSnap3,snapP3)
358  );
359 
360  points.append(p0p1);
361  points.append
362  (
363  generatePoint(s1,p1,hasSnap1,snapP1,s2,p2,hasSnap2,snapP2)
364  );
365  points.append(p2p3);
366  }
367  if (triIndex == 0x0A)
368  {
369  // Flip normals
370  label sz = points.size();
371  Swap(points[sz-5], points[sz-4]);
372  Swap(points[sz-2], points[sz-1]);
373  }
374  break;
375 
376  case 0x06:
377  case 0x09:
378  {
379  Type p0p1 =
380  generatePoint(s0,p0,hasSnap0,snapP0,s1,p1,hasSnap1,snapP1);
381  Type p2p3 =
382  generatePoint(s2,p2,hasSnap2,snapP2,s3,p3,hasSnap3,snapP3);
383 
384  points.append(p0p1);
385  points.append
386  (
387  generatePoint(s1,p1,hasSnap1,snapP1,s3,p3,hasSnap3,snapP3)
388  );
389  points.append(p2p3);
390 
391  points.append(p0p1);
392  points.append(p2p3);
393  points.append
394  (
395  generatePoint(s0,p0,hasSnap0,snapP0,s2,p2,hasSnap2,snapP2)
396  );
397  }
398  if (triIndex == 0x09)
399  {
400  // Flip normals
401  label sz = points.size();
402  Swap(points[sz-5], points[sz-4]);
403  Swap(points[sz-2], points[sz-1]);
404  }
405  break;
406 
407  case 0x08:
408  case 0x07:
409  points.append
410  (
411  generatePoint(s3,p3,hasSnap3,snapP3,s0,p0,hasSnap0,snapP0)
412  );
413  points.append
414  (
415  generatePoint(s3,p3,hasSnap3,snapP3,s2,p2,hasSnap2,snapP2)
416  );
417  points.append
418  (
419  generatePoint(s3,p3,hasSnap3,snapP3,s1,p1,hasSnap1,snapP1)
420  );
421  if (triIndex == 0x07)
422  {
423  // Flip normals
424  label sz = points.size();
425  Swap(points[sz-2], points[sz-1]);
426  }
427  break;
428  }
429 }
430 
431 
432 template<class Type>
433 Foam::label Foam::isoSurface::generateFaceTriPoints
434 (
435  const volScalarField& cVals,
436  const scalarField& pVals,
437 
439  const Field<Type>& pCoords,
440 
441  const DynamicList<Type>& snappedPoints,
442  const labelList& snappedCc,
443  const labelList& snappedPoint,
444  const label facei,
445 
446  const scalar neiVal,
447  const Type& neiPt,
448  const bool hasNeiSnap,
449  const Type& neiSnapPt,
450 
451  DynamicList<Type>& triPoints,
452  DynamicList<label>& triMeshCells
453 ) const
454 {
455  label own = mesh_.faceOwner()[facei];
456 
457  label oldNPoints = triPoints.size();
458 
459  const face& f = mesh_.faces()[facei];
460 
461  forAll(f, fp)
462  {
463  label pointi = f[fp];
464  label nextPointi = f[f.fcIndex(fp)];
465 
466  generateTriPoints
467  (
468  pVals[pointi],
469  pCoords[pointi],
470  snappedPoint[pointi] != -1,
471  (
472  snappedPoint[pointi] != -1
473  ? snappedPoints[snappedPoint[pointi]]
474  : Type(Zero)
475  ),
476 
477  pVals[nextPointi],
478  pCoords[nextPointi],
479  snappedPoint[nextPointi] != -1,
480  (
481  snappedPoint[nextPointi] != -1
482  ? snappedPoints[snappedPoint[nextPointi]]
483  : Type(Zero)
484  ),
485 
486  cVals[own],
487  cCoords[own],
488  snappedCc[own] != -1,
489  (
490  snappedCc[own] != -1
491  ? snappedPoints[snappedCc[own]]
492  : Type(Zero)
493  ),
494 
495  neiVal,
496  neiPt,
497  hasNeiSnap,
498  neiSnapPt,
499 
500  triPoints
501  );
502  }
503 
504  // Every three triPoints is a triangle
505  label nTris = (triPoints.size()-oldNPoints)/3;
506  for (label i = 0; i < nTris; i++)
507  {
508  triMeshCells.append(own);
509  }
510 
511  return nTris;
512 }
513 
514 
515 template<class Type>
516 void Foam::isoSurface::generateTriPoints
517 (
518  const volScalarField& cVals,
519  const scalarField& pVals,
520 
522  const Field<Type>& pCoords,
523 
524  const DynamicList<Type>& snappedPoints,
525  const labelList& snappedCc,
526  const labelList& snappedPoint,
527 
528  DynamicList<Type>& triPoints,
529  DynamicList<label>& triMeshCells
530 ) const
531 {
532  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
533  const labelList& own = mesh_.faceOwner();
534  const labelList& nei = mesh_.faceNeighbour();
535 
536  if
537  (
538  (cVals.size() != mesh_.nCells())
539  || (pVals.size() != mesh_.nPoints())
540  || (cCoords.size() != mesh_.nCells())
541  || (pCoords.size() != mesh_.nPoints())
542  || (snappedCc.size() != mesh_.nCells())
543  || (snappedPoint.size() != mesh_.nPoints())
544  )
545  {
547  << "Incorrect size." << endl
548  << "mesh: nCells:" << mesh_.nCells()
549  << " points:" << mesh_.nPoints() << endl
550  << "cVals:" << cVals.size() << endl
551  << "cCoords:" << cCoords.size() << endl
552  << "snappedCc:" << snappedCc.size() << endl
553  << "pVals:" << pVals.size() << endl
554  << "pCoords:" << pCoords.size() << endl
555  << "snappedPoint:" << snappedPoint.size() << endl
556  << abort(FatalError);
557  }
558 
559 
560  // Generate triangle points
561 
562  triPoints.clear();
563  triMeshCells.clear();
564 
565  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
566  {
567  if (faceCutType_[facei] != NOTCUT)
568  {
569  generateFaceTriPoints
570  (
571  cVals,
572  pVals,
573 
574  cCoords,
575  pCoords,
576 
577  snappedPoints,
578  snappedCc,
579  snappedPoint,
580  facei,
581 
582  cVals[nei[facei]],
583  cCoords[nei[facei]],
584  snappedCc[nei[facei]] != -1,
585  (
586  snappedCc[nei[facei]] != -1
587  ? snappedPoints[snappedCc[nei[facei]]]
588  : Type(Zero)
589  ),
590 
591  triPoints,
592  triMeshCells
593  );
594  }
595  }
596 
597 
598  // Determine neighbouring snap status
599  boolList neiSnapped(mesh_.nFaces()-mesh_.nInternalFaces(), false);
600  List<Type> neiSnappedPoint(neiSnapped.size(), Type(Zero));
601  forAll(patches, patchi)
602  {
603  const polyPatch& pp = patches[patchi];
604 
605  if (pp.coupled())
606  {
607  label facei = pp.start();
608  forAll(pp, i)
609  {
610  label bFacei = facei-mesh_.nInternalFaces();
611  label snappedIndex = snappedCc[own[facei]];
612 
613  if (snappedIndex != -1)
614  {
615  neiSnapped[bFacei] = true;
616  neiSnappedPoint[bFacei] = snappedPoints[snappedIndex];
617  }
618  facei++;
619  }
620  }
621  }
622  syncTools::swapBoundaryFaceList(mesh_, neiSnapped);
623  syncTools::swapBoundaryFaceList(mesh_, neiSnappedPoint);
624 
625 
626 
627  forAll(patches, patchi)
628  {
629  const polyPatch& pp = patches[patchi];
630 
631  if (isA<processorPolyPatch>(pp))
632  {
633  const processorPolyPatch& cpp =
634  refCast<const processorPolyPatch>(pp);
635 
636  PackedBoolList isCollocated(collocatedFaces(cpp));
637 
638  forAll(isCollocated, i)
639  {
640  label facei = pp.start()+i;
641 
642  if (faceCutType_[facei] != NOTCUT)
643  {
644  if (isCollocated[i])
645  {
646  generateFaceTriPoints
647  (
648  cVals,
649  pVals,
650 
651  cCoords,
652  pCoords,
653 
654  snappedPoints,
655  snappedCc,
656  snappedPoint,
657  facei,
658 
659  cVals.boundaryField()[patchi][i],
660  cCoords.boundaryField()[patchi][i],
661  neiSnapped[facei-mesh_.nInternalFaces()],
662  neiSnappedPoint[facei-mesh_.nInternalFaces()],
663 
664  triPoints,
665  triMeshCells
666  );
667  }
668  else
669  {
670  generateFaceTriPoints
671  (
672  cVals,
673  pVals,
674 
675  cCoords,
676  pCoords,
677 
678  snappedPoints,
679  snappedCc,
680  snappedPoint,
681  facei,
682 
683  cVals.boundaryField()[patchi][i],
684  cCoords.boundaryField()[patchi][i],
685  false,
686  Type(Zero),
687 
688  triPoints,
689  triMeshCells
690  );
691  }
692  }
693  }
694  }
695  else
696  {
697  label facei = pp.start();
698 
699  forAll(pp, i)
700  {
701  if (faceCutType_[facei] != NOTCUT)
702  {
703  generateFaceTriPoints
704  (
705  cVals,
706  pVals,
707 
708  cCoords,
709  pCoords,
710 
711  snappedPoints,
712  snappedCc,
713  snappedPoint,
714  facei,
715 
716  cVals.boundaryField()[patchi][i],
717  cCoords.boundaryField()[patchi][i],
718  false, // fc not snapped
719  Type(Zero),
720 
721  triPoints,
722  triMeshCells
723  );
724  }
725  facei++;
726  }
727  }
728  }
729 
730  triPoints.shrink();
731  triMeshCells.shrink();
732 }
733 
734 
735 template<class Type>
737 Foam::isoSurface::interpolate
738 (
739  const label nPoints,
740  const labelList& triPointMergeMap,
741  const DynamicList<Type>& unmergedValues
742 )
743 {
744  // One value per point
745  tmp<Field<Type> > tvalues(new Field<Type>(nPoints, Type(Zero)));
746  Field<Type>& values = tvalues.ref();
747  labelList nValues(values.size(), 0);
748 
749  forAll(unmergedValues, i)
750  {
751  label mergedPointi = triPointMergeMap[i];
752 
753  if (mergedPointi >= 0)
754  {
755  values[mergedPointi] += unmergedValues[i];
756  nValues[mergedPointi]++;
757  }
758  }
759 
760  forAll(values, i)
761  {
762  if (nValues[i] > 0)
763  {
764  values[i] /= scalar(nValues[i]);
765  }
766  }
767 
768  return tvalues;
769 }
770 
771 
772 template<class Type>
774 Foam::isoSurface::interpolate
775 (
777  const Field<Type>& pCoords
778 ) const
779 {
780  // Recalculate boundary values
782  <
783  Type,
784  fvPatchField,
786  volMesh
787  >> c2(adaptPatchFields(cCoords));
788 
789 
790  DynamicList<Type> triPoints(3*nCutCells_);
791  DynamicList<label> triMeshCells(nCutCells_);
792 
793  // Dummy snap data
794  DynamicList<Type> snappedPoints;
795  labelList snappedCc(mesh_.nCells(), -1);
796  labelList snappedPoint(mesh_.nPoints(), -1);
797 
798  generateTriPoints
799  (
800  cValsPtr_(),
801  pVals_,
802 
803  c2(),
804  pCoords,
805 
806  snappedPoints,
807  snappedCc,
808  snappedPoint,
809 
810  triPoints,
811  triMeshCells
812  );
813 
814  return interpolate
815  (
816  points().size(),
817  triPointMergeMap_,
818  triPoints
819  );
820 }
821 
822 
823 // ************************************************************************* //
Foam::surfaceFields.
Specialization of fvPatchField which creates the underlying fvPatchField as a slice of the given comp...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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 objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:221
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:225
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const Field< point > & points() const
Return reference to global points.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Generic GeometricField class.
patches[0]
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
label nCells() const
Neighbour processor patch.
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:431
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:53
dynamicFvMesh & mesh
void Swap(T &a, T &b)
Definition: Swap.H:43
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Pre-declare SubField and related Field type.
Definition: Field.H:57
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:310
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:314
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
static const zero Zero
Definition: zero.H:91
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
Foam::polyBoundaryMesh.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:240
labelList f(nPoints)
const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
label size() const
Return the number of elements in the UList.
label nFaces() const
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
void setSize(const label)
Reset size of List.
A bit-packed bool list.
label patchi
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
const Mesh & mesh() const
Return mesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label nPoints() const
dimensioned< scalar > mag(const dimensioned< Type > &)
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:430
A class for managing temporary objects.
Definition: PtrList.H:54
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:224
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const fileName & instance() const
Definition: IOobject.H:337
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545
label nInternalFaces() const
const word & name() const
Return name.
Definition: IOobject.H:260