wallBoundedStreamLine.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-2014 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 "Pstream.H"
27 #include "functionObjectList.H"
28 #include "wallBoundedStreamLine.H"
29 #include "fvMesh.H"
31 #include "ReadFields.H"
32 #include "meshSearch.H"
33 #include "sampledSet.H"
34 #include "globalIndex.H"
35 #include "mapDistribute.H"
36 #include "interpolationCellPoint.H"
37 #include "PatchTools.H"
38 #include "meshSearchMeshObject.H"
39 #include "faceSet.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 defineTypeNameAndDebug(wallBoundedStreamLine, 0);
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
52 Foam::wallBoundedStreamLine::wallPatch() const
53 {
54  const fvMesh& mesh = dynamic_cast<const fvMesh&>(obr_);
55 
56  const polyBoundaryMesh& patches = mesh.boundaryMesh();
57 
58  label nFaces = 0;
59 
60  forAll(patches, patchI)
61  {
62  //if (!polyPatch::constraintType(patches[patchI].type()))
63  if (isA<wallPolyPatch>(patches[patchI]))
64  {
65  nFaces += patches[patchI].size();
66  }
67  }
68 
69  labelList addressing(nFaces);
70 
71  nFaces = 0;
72 
73  forAll(patches, patchI)
74  {
75  //if (!polyPatch::constraintType(patches[patchI].type()))
76  if (isA<wallPolyPatch>(patches[patchI]))
77  {
78  const polyPatch& pp = patches[patchI];
79 
80  forAll(pp, i)
81  {
82  addressing[nFaces++] = pp.start()+i;
83  }
84  }
85  }
86 
87  return autoPtr<indirectPrimitivePatch>
88  (
90  (
91  IndirectList<face>
92  (
93  mesh.faces(),
94  addressing
95  ),
96  mesh.points()
97  )
98  );
99 }
100 
101 
102 Foam::tetIndices Foam::wallBoundedStreamLine::findNearestTet
103 (
104  const PackedBoolList& isWallPatch,
105  const point& seedPt,
106  const label cellI
107 ) const
108 {
109  const fvMesh& mesh = dynamic_cast<const fvMesh&>(obr_);
110 
111  const cell& cFaces = mesh.cells()[cellI];
112 
113  label minFaceI = -1;
114  label minTetPtI = -1;
115  scalar minDistSqr = sqr(GREAT);
116 
117  forAll(cFaces, cFaceI)
118  {
119  label faceI = cFaces[cFaceI];
120 
121  if (isWallPatch[faceI])
122  {
123  const face& f = mesh.faces()[faceI];
124  const label fp0 = mesh.tetBasePtIs()[faceI];
125  const point& basePoint = mesh.points()[f[fp0]];
126 
127  label fp = f.fcIndex(fp0);
128  for (label i = 2; i < f.size(); i++)
129  {
130  const point& thisPoint = mesh.points()[f[fp]];
131  label nextFp = f.fcIndex(fp);
132  const point& nextPoint = mesh.points()[f[nextFp]];
133 
134  const triPointRef tri(basePoint, thisPoint, nextPoint);
135 
136  scalar d2 = magSqr(tri.centre() - seedPt);
137  if (d2 < minDistSqr)
138  {
139  minDistSqr = d2;
140  minFaceI = faceI;
141  minTetPtI = i-1;
142  }
143  fp = nextFp;
144  }
145  }
146  }
147 
148  // Put particle in tet
149  return tetIndices
150  (
151  cellI,
152  minFaceI,
153  minTetPtI,
154  mesh
155  );
156 }
157 
158 
159 void Foam::wallBoundedStreamLine::track()
160 {
161  const Time& runTime = obr_.time();
162  const fvMesh& mesh = dynamic_cast<const fvMesh&>(obr_);
163 
164 
165  // Determine the 'wall' patches
166  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
167  // These are the faces that need to be followed
168 
169  autoPtr<indirectPrimitivePatch> boundaryPatch(wallPatch());
170  PackedBoolList isWallPatch(mesh.nFaces());
171  forAll(boundaryPatch().addressing(), i)
172  {
173  isWallPatch[boundaryPatch().addressing()[i]] = 1;
174  }
175 
176 
177 
178  // Find nearest wall particle for the seedPoints
179  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
180 
181  IDLList<wallBoundedStreamLineParticle> initialParticles;
182  wallBoundedStreamLineParticleCloud particles
183  (
184  mesh,
185  cloudName_,
186  initialParticles
187  );
188 
189  {
190  // Get the seed points
191  // ~~~~~~~~~~~~~~~~~~~
192 
193  const sampledSet& seedPoints = sampledSetPtr_();
194 
195 
196  forAll(seedPoints, i)
197  {
198  const point& seedPt = seedPoints[i];
199  label cellI = seedPoints.cells()[i];
200 
201  tetIndices ids(findNearestTet(isWallPatch, seedPt, cellI));
202 
203  if (ids.face() != -1 && isWallPatch[ids.face()])
204  {
205  //Pout<< "Seeding particle :" << nl
206  // << " seedPt:" << seedPt << nl
207  // << " face :" << ids.face() << nl
208  // << " at :" << mesh.faceCentres()[ids.face()] << nl
209  // << " cell :" << mesh.cellCentres()[ids.cell()] << nl
210  // << endl;
211 
212  particles.addParticle
213  (
214  new wallBoundedStreamLineParticle
215  (
216  mesh,
217  ids.faceTri(mesh).centre(),
218  ids.cell(),
219  ids.face(), // tetFace
220  ids.tetPt(),
221  -1, // not on a mesh edge
222  -1, // not on a diagonal edge
223  lifeTime_ // lifetime
224  )
225  );
226  }
227  else
228  {
229  Pout<< type() << " : ignoring seed " << seedPt
230  << " since not in wall cell." << endl;
231  }
232  }
233  }
234 
235  label nSeeds = returnReduce(particles.size(), sumOp<label>());
236 
237  Info<< type() << " : seeded " << nSeeds << " particles." << endl;
238 
239 
240 
241  // Read or lookup fields
242  PtrList<volScalarField> vsFlds;
243  PtrList<interpolation<scalar> > vsInterp;
244  PtrList<volVectorField> vvFlds;
245  PtrList<interpolation<vector> > vvInterp;
246 
247  label UIndex = -1;
248 
249  if (loadFromFiles_)
250  {
251  IOobjectList allObjects(mesh, runTime.timeName());
252 
253  IOobjectList objects(2*fields_.size());
254  forAll(fields_, i)
255  {
256  objects.add(*allObjects[fields_[i]]);
257  }
258 
259  ReadFields(mesh, objects, vsFlds);
260  vsInterp.setSize(vsFlds.size());
261  forAll(vsFlds, i)
262  {
263  vsInterp.set
264  (
265  i,
267  (
268  interpolationScheme_,
269  vsFlds[i]
270  )
271  );
272  }
273  ReadFields(mesh, objects, vvFlds);
274  vvInterp.setSize(vvFlds.size());
275  forAll(vvFlds, i)
276  {
277  vvInterp.set
278  (
279  i,
281  (
282  interpolationScheme_,
283  vvFlds[i]
284  )
285  );
286  }
287  }
288  else
289  {
290  label nScalar = 0;
291  label nVector = 0;
292 
293  forAll(fields_, i)
294  {
295  if (mesh.foundObject<volScalarField>(fields_[i]))
296  {
297  nScalar++;
298  }
299  else if (mesh.foundObject<volVectorField>(fields_[i]))
300  {
301  nVector++;
302  }
303  else
304  {
305  FatalErrorIn("wallBoundedStreamLine::execute()")
306  << "Cannot find field " << fields_[i] << endl
307  << "Valid scalar fields are:"
308  << mesh.names(volScalarField::typeName) << endl
309  << "Valid vector fields are:"
310  << mesh.names(volVectorField::typeName)
311  << exit(FatalError);
312  }
313  }
314  vsInterp.setSize(nScalar);
315  nScalar = 0;
316  vvInterp.setSize(nVector);
317  nVector = 0;
318 
319  forAll(fields_, i)
320  {
321  if (mesh.foundObject<volScalarField>(fields_[i]))
322  {
323  const volScalarField& f = mesh.lookupObject<volScalarField>
324  (
325  fields_[i]
326  );
327  vsInterp.set
328  (
329  nScalar++,
331  (
332  interpolationScheme_,
333  f
334  )
335  );
336  }
337  else if (mesh.foundObject<volVectorField>(fields_[i]))
338  {
339  const volVectorField& f = mesh.lookupObject<volVectorField>
340  (
341  fields_[i]
342  );
343 
344  if (f.name() == UName_)
345  {
346  UIndex = nVector;
347  }
348 
349  vvInterp.set
350  (
351  nVector++,
353  (
354  interpolationScheme_,
355  f
356  )
357  );
358  }
359  }
360  }
361 
362  // Store the names
363  scalarNames_.setSize(vsInterp.size());
364  forAll(vsInterp, i)
365  {
366  scalarNames_[i] = vsInterp[i].psi().name();
367  }
368  vectorNames_.setSize(vvInterp.size());
369  forAll(vvInterp, i)
370  {
371  vectorNames_[i] = vvInterp[i].psi().name();
372  }
373 
374  // Check that we know the index of U in the interpolators.
375 
376  if (UIndex == -1)
377  {
378  FatalErrorIn("wallBoundedStreamLine::execute()")
379  << "Cannot find field to move particles with : " << UName_
380  << endl
381  << "This field has to be present in the sampled fields "
382  << fields_
383  << " and in the objectRegistry." << endl
384  << exit(FatalError);
385  }
386 
387  // Sampled data
388  // ~~~~~~~~~~~~
389 
390  // Size to maximum expected sizes.
391  allTracks_.clear();
392  allTracks_.setCapacity(nSeeds);
393  allScalars_.setSize(vsInterp.size());
394  forAll(allScalars_, i)
395  {
396  allScalars_[i].clear();
397  allScalars_[i].setCapacity(nSeeds);
398  }
399  allVectors_.setSize(vvInterp.size());
400  forAll(allVectors_, i)
401  {
402  allVectors_[i].clear();
403  allVectors_[i].setCapacity(nSeeds);
404  }
405 
406 
407  // additional particle info
408  wallBoundedStreamLineParticle::trackingData td
409  (
410  particles,
411  vsInterp,
412  vvInterp,
413  UIndex, // index of U in vvInterp
414  trackForward_, // track in +u direction?
415  trackLength_, // fixed track length
416  isWallPatch, // which faces are to follow
417 
418  allTracks_,
419  allScalars_,
420  allVectors_
421  );
422 
423 
424  // Set very large dt. Note: cannot use GREAT since 1/GREAT is SMALL
425  // which is a trigger value for the tracking...
426  const scalar trackTime = Foam::sqrt(GREAT);
427 
428  // Track
429  particles.move(td, trackTime);
430 }
431 
432 
433 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
434 
435 Foam::wallBoundedStreamLine::wallBoundedStreamLine
436 (
437  const word& name,
438  const objectRegistry& obr,
439  const dictionary& dict,
440  const bool loadFromFiles
441 )
442 :
443  dict_(dict),
444  name_(name),
445  obr_(obr),
446  loadFromFiles_(loadFromFiles),
447  active_(true)
448 {
449  // Only active if a fvMesh is available
450  if (isA<fvMesh>(obr_))
451  {
452  read(dict_);
453  }
454  else
455  {
456  active_ = false;
457  WarningIn
458  (
459  "wallBoundedStreamLine::wallBoundedStreamLine\n"
460  "("
461  "const word&, "
462  "const objectRegistry&, "
463  "const dictionary&, "
464  "const bool "
465  ")"
466  ) << "No fvMesh available, deactivating " << name_
467  << nl << endl;
468  }
469 }
470 
471 
472 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
473 
475 {}
476 
477 
478 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
479 
481 {
482  if (active_)
483  {
484  //dict_ = dict;
485  dict.lookup("fields") >> fields_;
486  if (dict.found("UName"))
487  {
488  dict.lookup("UName") >> UName_;
489  }
490  else
491  {
492  UName_ = "U";
493  if (dict.found("U"))
494  {
496  (
497  "wallBoundedStreamLine::read(const dictionary&)",
498  dict
499  ) << "Using deprecated entry \"U\"."
500  << " Please use \"UName\" instead."
501  << endl;
502  dict.lookup("U") >> UName_;
503  }
504  }
505 
506  if (findIndex(fields_, UName_) == -1)
507  {
509  (
510  "wallBoundedStreamLine::read(const dictionary&)",
511  dict
512  ) << "Velocity field for tracking " << UName_
513  << " should be present in the list of fields " << fields_
514  << exit(FatalIOError);
515  }
516 
517 
518  dict.lookup("trackForward") >> trackForward_;
519  dict.lookup("lifeTime") >> lifeTime_;
520  if (lifeTime_ < 1)
521  {
522  FatalErrorIn(":wallBoundedStreamLine::read(const dictionary&)")
523  << "Illegal value " << lifeTime_ << " for lifeTime"
524  << exit(FatalError);
525  }
526  trackLength_ = VGREAT;
527  if (dict.found("trackLength"))
528  {
529  dict.lookup("trackLength") >> trackLength_;
530 
531  Info<< type() << " : fixed track length specified : "
532  << trackLength_ << nl << endl;
533  }
534 
535 
536  interpolationScheme_ = dict.lookupOrDefault
537  (
538  "interpolationScheme",
540  );
541 
542  //Info<< typeName << " using interpolation " << interpolationScheme_
543  // << endl;
544 
545  cloudName_ = dict.lookupOrDefault<word>
546  (
547  "cloudName",
548  "wallBoundedStreamLine"
549  );
550  dict.lookup("seedSampleSet") >> seedSet_;
551 
552  const fvMesh& mesh = dynamic_cast<const fvMesh&>(obr_);
553 
554  const dictionary& coeffsDict = dict.subDict(seedSet_ + "Coeffs");
555 
556  sampledSetPtr_ = sampledSet::New
557  (
558  seedSet_,
559  mesh,
561  coeffsDict
562  );
563  coeffsDict.lookup("axis") >> sampledSetAxis_;
564 
565  scalarFormatterPtr_ = writer<scalar>::New(dict.lookup("setFormat"));
566  vectorFormatterPtr_ = writer<vector>::New(dict.lookup("setFormat"));
567 
568 
569  // Make sure that the mesh is trackable
570  if (debug)
571  {
572  // 1. positive volume decomposition tets
573  faceSet faces(mesh, "lowQualityTetFaces", mesh.nFaces()/100+1);
574  if
575  (
577  (
578  mesh,
580  true,
581  &faces
582  )
583  )
584  {
585  label nFaces = returnReduce(faces.size(), sumOp<label>());
586 
587  WarningIn("wallBoundedStreamLine::read(const dictionary&)")
588  << "Found " << nFaces
589  <<" faces with low quality or negative volume "
590  << "decomposition tets. Writing to faceSet " << faces.name()
591  << endl;
592  }
593 
594  // 2. all edges on a cell having two faces
595  EdgeMap<label> numFacesPerEdge;
596  forAll(mesh.cells(), cellI)
597  {
598  const cell& cFaces = mesh.cells()[cellI];
599 
600  numFacesPerEdge.clear();
601 
602  forAll(cFaces, cFaceI)
603  {
604  label faceI = cFaces[cFaceI];
605  const face& f = mesh.faces()[faceI];
606  forAll(f, fp)
607  {
608  const edge e(f[fp], f.nextLabel(fp));
610  numFacesPerEdge.find(e);
611  if (eFnd != numFacesPerEdge.end())
612  {
613  eFnd()++;
614  }
615  else
616  {
617  numFacesPerEdge.insert(e, 1);
618  }
619  }
620  }
621 
622  forAllConstIter(EdgeMap<label>, numFacesPerEdge, iter)
623  {
624  if (iter() != 2)
625  {
627  (
628  "wallBoundedStreamLine::read(const dictionary&)"
629  ) << "problem cell:" << cellI
630  << abort(FatalError);
631  }
632  }
633  }
634  }
635  }
636 }
637 
638 
640 {}
641 
642 
644 {}
645 
646 
648 {}
649 
650 
652 {
653  if (active_)
654  {
655  const Time& runTime = obr_.time();
656  const fvMesh& mesh = dynamic_cast<const fvMesh&>(obr_);
657 
658 
659  // Do all injection and tracking
660  track();
661 
662 
663  if (Pstream::parRun())
664  {
665  // Append slave tracks to master ones
666  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
667 
668  globalIndex globalTrackIDs(allTracks_.size());
669 
670  // Construct a distribution map to pull all to the master.
671  labelListList sendMap(Pstream::nProcs());
672  labelListList recvMap(Pstream::nProcs());
673 
674  if (Pstream::master())
675  {
676  // Master: receive all. My own first, then consecutive
677  // processors.
678  label trackI = 0;
679 
680  forAll(recvMap, procI)
681  {
682  labelList& fromProc = recvMap[procI];
683  fromProc.setSize(globalTrackIDs.localSize(procI));
684  forAll(fromProc, i)
685  {
686  fromProc[i] = trackI++;
687  }
688  }
689  }
690 
691  labelList& toMaster = sendMap[0];
692  toMaster.setSize(globalTrackIDs.localSize());
693  forAll(toMaster, i)
694  {
695  toMaster[i] = i;
696  }
697 
698  const mapDistribute distMap
699  (
700  globalTrackIDs.size(),
701  sendMap.xfer(),
702  recvMap.xfer()
703  );
704 
705 
706  // Distribute the track positions. Note: use scheduled comms
707  // to prevent buffering.
709  (
711  distMap.schedule(),
712  distMap.constructSize(),
713  distMap.subMap(),
714  distMap.constructMap(),
715  allTracks_
716  );
717 
718  // Distribute the scalars
719  forAll(allScalars_, scalarI)
720  {
722  (
724  distMap.schedule(),
725  distMap.constructSize(),
726  distMap.subMap(),
727  distMap.constructMap(),
728  allScalars_[scalarI]
729  );
730  }
731  // Distribute the vectors
732  forAll(allVectors_, vectorI)
733  {
735  (
737  distMap.schedule(),
738  distMap.constructSize(),
739  distMap.subMap(),
740  distMap.constructMap(),
741  allVectors_[vectorI]
742  );
743  }
744  }
745 
746 
747  label n = 0;
748  forAll(allTracks_, trackI)
749  {
750  n += allTracks_[trackI].size();
751  }
752 
753  Info<< " Tracks:" << allTracks_.size() << nl
754  << " Total samples:" << n << endl;
755 
756 
757  // Massage into form suitable for writers
758  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
759 
760  if (Pstream::master() && allTracks_.size())
761  {
762  // Make output directory
763 
764  fileName vtkPath
765  (
767  ? runTime.path()/".."/"postProcessing"/"sets"/name()
768  : runTime.path()/"postProcessing"/"sets"/name()
769  );
770  if (mesh.name() != fvMesh::defaultRegion)
771  {
772  vtkPath = vtkPath/mesh.name();
773  }
774  vtkPath = vtkPath/mesh.time().timeName();
775 
776  mkDir(vtkPath);
777 
778  // Convert track positions
779 
780  PtrList<coordSet> tracks(allTracks_.size());
781  forAll(allTracks_, trackI)
782  {
783  tracks.set
784  (
785  trackI,
786  new coordSet
787  (
788  "track" + Foam::name(trackI),
789  sampledSetAxis_ //"xyz"
790  )
791  );
792  tracks[trackI].transfer(allTracks_[trackI]);
793  }
794 
795  // Convert scalar values
796 
797  if (allScalars_.size() > 0)
798  {
799  List<List<scalarField> > scalarValues(allScalars_.size());
800 
801  forAll(allScalars_, scalarI)
802  {
803  DynamicList<scalarList>& allTrackVals =
804  allScalars_[scalarI];
805  scalarValues[scalarI].setSize(allTrackVals.size());
806 
807  forAll(allTrackVals, trackI)
808  {
809  scalarList& trackVals = allTrackVals[trackI];
810  scalarValues[scalarI][trackI].transfer(trackVals);
811  }
812  }
813 
814  fileName vtkFile
815  (
816  vtkPath
817  / scalarFormatterPtr_().getFileName
818  (
819  tracks[0],
820  scalarNames_
821  )
822  );
823 
824  Info<< "Writing data to " << vtkFile.path() << endl;
825 
826  scalarFormatterPtr_().write
827  (
828  true, // writeTracks
829  tracks,
830  scalarNames_,
831  scalarValues,
832  OFstream(vtkFile)()
833  );
834  }
835 
836  // Convert vector values
837 
838  if (allVectors_.size() > 0)
839  {
840  List<List<vectorField> > vectorValues(allVectors_.size());
841 
842  forAll(allVectors_, vectorI)
843  {
844  DynamicList<vectorList>& allTrackVals =
845  allVectors_[vectorI];
846  vectorValues[vectorI].setSize(allTrackVals.size());
847 
848  forAll(allTrackVals, trackI)
849  {
850  vectorList& trackVals = allTrackVals[trackI];
851  vectorValues[vectorI][trackI].transfer(trackVals);
852  }
853  }
854 
855  fileName vtkFile
856  (
857  vtkPath
858  / vectorFormatterPtr_().getFileName
859  (
860  tracks[0],
861  vectorNames_
862  )
863  );
864 
865  //Info<< "Writing vector data to " << vtkFile << endl;
866 
867  vectorFormatterPtr_().write
868  (
869  true, // writeTracks
870  tracks,
871  vectorNames_,
872  vectorValues,
873  OFstream(vtkFile)()
874  );
875  }
876  }
877  }
878 }
879 
880 
882 {
883  read(dict_);
884 }
885 
886 
888 {
889  // Moving mesh affects the search tree
890  read(dict_);
891 }
892 
893 
894 //void Foam::wallBoundedStreamLine::readUpdate
895 //(const polyMesh::readUpdateState state)
896 //{
897 // if (state != UNCHANGED)
898 // {
899 // read(dict_);
900 // }
901 //}
902 
903 
904 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
Output to file stream.
Definition: OFstream.H:81
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:380
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.
Definition: writer.C:35
void transfer(List< T > &)
Transfer contents of the argument List into this.
Definition: DynamicListI.H:277
static const char *const typeName
Definition: Field.H:94
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
label nFaces() const
A list of face labels.
Definition: faceSet.H:48
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
vector point
Point is a vector.
Definition: point.H:41
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
fileName path() const
Return path.
Definition: Time.H:281
static void distribute(const Pstream::commsTypes commsType, const List< labelPair > &schedule, const label constructSize, const labelListList &subMap, const labelListList &constructMap, List< T > &, const int tag=UPstream::msgType())
Distribute data. Note:schedule only used for Pstream::scheduled.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Holds list of sampling positions.
Definition: coordSet.H:49
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
static autoPtr< sampledSet > New(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Return a reference to the selected sampledSet.
Definition: sampledSet.C:448
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const meshSearchMeshObject & New(const polyMesh &mesh)
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
Xfer< List< T > > xfer()
Transfer contents to the Xfer container.
Definition: ListI.H:90
Class containing processor-to-processor mapping information.
virtual void execute()
Execute the averaging.
static const scalar minTetQuality
Minimum tetrahedron quality.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
A class for handling words, derived from string.
Definition: word.H:59
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
static bool checkFaceTets(const polyMesh &mesh, scalar tol=minTetQuality, const bool report=false, labelHashSet *setPtr=NULL)
Check face-decomposition tet volume.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const cellList & cells() const
PrimitivePatch< face, IndirectList, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
messageStream Info
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:117
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:287
virtual ~wallBoundedStreamLine()
Destructor.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:638
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Namespace for OpenFOAM.
#define IOWarningIn(functionName, ios)
Report an IO warning using Foam::Warning.
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
virtual void write()
Calculate the field average data and write.
label n
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:386
static const char nl
Definition: Ostream.H:260
const double e
Elementary charge.
Definition: doubleFloat.H:78
void setSize(const label)
Reset size of List.
Definition: List.C:318
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOerror FatalIOError
#define WarningIn(functionName)
Report a warning using Foam::Warning.
void clear()
Clear all entries from table.
Definition: HashTable.C:473
wordList ReadFields(const Mesh &mesh, const IOobjectList &objects, PtrList< GeoField > &fields, const bool syncPar=true)
Helper routine to read fields.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
#define forAll(list, i)
Definition: UList.H:421
virtual Ostream & write(const token &)=0
Write next token to stream.
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
void setSize(const label)
Alter the addressed list size.
Definition: DynamicListI.H:167
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
Helper routine to read fields.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Registry of regIOobjects.
List< label > labelList
A List of labels.
Definition: labelList.H:56
A class for handling file names.
Definition: fileName.H:69
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:398
bool read(const char *, int32_t &)
Definition: int32IO.C:87
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:106
const Time & time() const
Return time.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1060
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:139
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:420
static autoPtr< interpolation< scalar > > New(const word &interpolationType, const GeometricField< scalar, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:73
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
virtual void end()
Execute the averaging at the final time-loop, currently does nothing.
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
virtual void read(const dictionary &)
Read the field average data.