fvMesh.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 "fvMesh.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "slicedVolFields.H"
30 #include "slicedSurfaceFields.H"
31 #include "SubField.H"
32 #include "demandDrivenData.H"
33 #include "fvMeshLduAddressing.H"
34 #include "mapPolyMesh.H"
35 #include "MapFvFields.H"
36 #include "fvMeshMapper.H"
37 #include "mapClouds.H"
38 #include "MeshObject.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44  defineTypeNameAndDebug(fvMesh, 0);
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::fvMesh::clearGeomNotOldVol()
51 {
53  <
54  fvMesh,
55  GeometricMeshObject,
56  MoveableMeshObject
57  >(*this);
58 
60  <
61  lduMesh,
62  GeometricMeshObject,
63  MoveableMeshObject
64  >(*this);
65 
67  static_cast<slicedVolScalarField::Internal*>(VPtr_);
69  VPtr_ = NULL;
70 
71  deleteDemandDrivenData(SfPtr_);
72  deleteDemandDrivenData(magSfPtr_);
74  deleteDemandDrivenData(CfPtr_);
75 }
76 
77 
78 void Foam::fvMesh::updateGeomNotOldVol()
79 {
80  bool haveV = (VPtr_ != NULL);
81  bool haveSf = (SfPtr_ != NULL);
82  bool haveMagSf = (magSfPtr_ != NULL);
83  bool haveCP = (CPtr_ != NULL);
84  bool haveCf = (CfPtr_ != NULL);
85 
86  clearGeomNotOldVol();
87 
88  // Now recreate the fields
89  if (haveV)
90  {
91  (void)V();
92  }
93  if (haveSf)
94  {
95  (void)Sf();
96  }
97  if (haveMagSf)
98  {
99  (void)magSf();
100  }
101  if (haveCP)
102  {
103  (void)C();
104  }
105  if (haveCf)
106  {
107  (void)Cf();
108  }
109 }
110 
111 
112 void Foam::fvMesh::clearGeom()
113 {
114  clearGeomNotOldVol();
115 
116  deleteDemandDrivenData(V0Ptr_);
117  deleteDemandDrivenData(V00Ptr_);
118 
119  // Mesh motion flux cannot be deleted here because the old-time flux
120  // needs to be saved.
121 }
122 
123 
124 void Foam::fvMesh::clearAddressing(const bool isMeshUpdate)
125 {
126  if (debug)
127  {
128  InfoInFunction << "isMeshUpdate: " << isMeshUpdate << endl;
129  }
130 
131  if (isMeshUpdate)
132  {
133  // Part of a mesh update. Keep meshObjects that have an updateMesh
134  // callback
136  <
137  fvMesh,
138  TopologicalMeshObject,
139  UpdateableMeshObject
140  >
141  (
142  *this
143  );
145  <
146  lduMesh,
147  TopologicalMeshObject,
148  UpdateableMeshObject
149  >
150  (
151  *this
152  );
153  }
154  else
155  {
156  meshObject::clear<fvMesh, TopologicalMeshObject>(*this);
157  meshObject::clear<lduMesh, TopologicalMeshObject>(*this);
158  }
159  deleteDemandDrivenData(lduPtr_);
160 }
161 
162 
163 void Foam::fvMesh::storeOldVol(const scalarField& V)
164 {
165  if (curTimeIndex_ < time().timeIndex())
166  {
167  if (debug)
168  {
170  << " Storing old time volumes since from time " << curTimeIndex_
171  << " and time now " << time().timeIndex()
172  << " V:" << V.size()
173  << endl;
174  }
175 
176 
177  if (V00Ptr_ && V0Ptr_)
178  {
179  // Copy V0 into V00 storage
180  *V00Ptr_ = *V0Ptr_;
181  }
182 
183  if (V0Ptr_)
184  {
185  // Copy V into V0 storage
186  V0Ptr_->scalarField::operator=(V);
187  }
188  else
189  {
190  // Allocate V0 storage, fill with V
191  V0Ptr_ = new DimensionedField<scalar, volMesh>
192  (
193  IOobject
194  (
195  "V0",
196  time().timeName(),
197  *this,
200  false
201  ),
202  *this,
203  dimVolume
204  );
205  scalarField& V0 = *V0Ptr_;
206  // Note: V0 now sized with current mesh, not with (potentially
207  // different size) V.
208  V0.setSize(V.size());
209  V0 = V;
210  }
211 
212  curTimeIndex_ = time().timeIndex();
213 
214  if (debug)
215  {
217  << " Stored old time volumes V0:" << V0Ptr_->size()
218  << endl;
219  if (V00Ptr_)
220  {
222  << " Stored oldold time volumes V00:" << V00Ptr_->size()
223  << endl;
224  }
225  }
226  }
227 }
228 
229 
231 {
232  clearGeom();
234 
235  clearAddressing();
236 
237  // Clear mesh motion flux
238  deleteDemandDrivenData(phiPtr_);
239 
241 }
242 
243 
244 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
245 
246 Foam::fvMesh::fvMesh(const IOobject& io)
247 :
248  polyMesh(io),
249  surfaceInterpolation(*this),
250  fvSchemes(static_cast<const objectRegistry&>(*this)),
251  fvSolution(static_cast<const objectRegistry&>(*this)),
252  data(static_cast<const objectRegistry&>(*this)),
253  boundary_(*this, boundaryMesh()),
254  lduPtr_(NULL),
255  curTimeIndex_(time().timeIndex()),
256  VPtr_(NULL),
257  V0Ptr_(NULL),
258  V00Ptr_(NULL),
259  SfPtr_(NULL),
260  magSfPtr_(NULL),
261  CPtr_(NULL),
262  CfPtr_(NULL),
263  phiPtr_(NULL)
264 {
265  if (debug)
266  {
267  InfoInFunction << "Constructing fvMesh from IOobject" << endl;
268  }
269 
270  // Check the existance of the cell volumes and read if present
271  // and set the storage of V00
272  if (isFile(time().timePath()/"V0"))
273  {
275  (
276  IOobject
277  (
278  "V0",
279  time().timeName(),
280  *this,
283  false
284  ),
285  *this
286  );
287 
288  V00();
289  }
290 
291  // Check the existance of the mesh fluxes, read if present and set the
292  // mesh to be moving
293  if (isFile(time().timePath()/"meshPhi"))
294  {
295  phiPtr_ = new surfaceScalarField
296  (
297  IOobject
298  (
299  "meshPhi",
300  time().timeName(),
301  *this,
304  false
305  ),
306  *this
307  );
308 
309  // The mesh is now considered moving so the old-time cell volumes
310  // will be required for the time derivatives so if they haven't been
311  // read initialise to the current cell volumes
312  if (!V0Ptr_)
313  {
315  (
316  IOobject
317  (
318  "V0",
319  time().timeName(),
320  *this,
323  false
324  ),
325  V()
326  );
327  }
328 
329  moving(true);
330  }
331 }
332 
333 
334 Foam::fvMesh::fvMesh
335 (
336  const IOobject& io,
337  const Xfer<pointField>& points,
338  const cellShapeList& shapes,
339  const faceListList& boundaryFaces,
340  const wordList& boundaryPatchNames,
341  const PtrList<dictionary>& boundaryDicts,
342  const word& defaultBoundaryPatchName,
343  const word& defaultBoundaryPatchType,
344  const bool syncPar
345 )
346 :
347  polyMesh
348  (
349  io,
350  points,
351  shapes,
352  boundaryFaces,
353  boundaryPatchNames,
354  boundaryDicts,
355  defaultBoundaryPatchName,
356  defaultBoundaryPatchType,
357  syncPar
358  ),
359  surfaceInterpolation(*this),
360  fvSchemes(static_cast<const objectRegistry&>(*this)),
361  fvSolution(static_cast<const objectRegistry&>(*this)),
362  data(static_cast<const objectRegistry&>(*this)),
363  boundary_(*this, boundaryMesh()),
364  lduPtr_(NULL),
365  curTimeIndex_(time().timeIndex()),
366  VPtr_(NULL),
367  V0Ptr_(NULL),
368  V00Ptr_(NULL),
369  SfPtr_(NULL),
370  magSfPtr_(NULL),
371  CPtr_(NULL),
372  CfPtr_(NULL),
373  phiPtr_(NULL)
374 {
375  if (debug)
376  {
377  InfoInFunction << "Constructing fvMesh from cellShapes" << endl;
378  }
379 }
380 
381 
382 Foam::fvMesh::fvMesh
383 (
384  const IOobject& io,
385  const Xfer<pointField>& points,
386  const Xfer<faceList>& faces,
387  const Xfer<labelList>& allOwner,
388  const Xfer<labelList>& allNeighbour,
389  const bool syncPar
390 )
391 :
392  polyMesh(io, points, faces, allOwner, allNeighbour, syncPar),
393  surfaceInterpolation(*this),
394  fvSchemes(static_cast<const objectRegistry&>(*this)),
395  fvSolution(static_cast<const objectRegistry&>(*this)),
396  data(static_cast<const objectRegistry&>(*this)),
397  boundary_(*this, boundaryMesh()),
398  lduPtr_(NULL),
399  curTimeIndex_(time().timeIndex()),
400  VPtr_(NULL),
401  V0Ptr_(NULL),
402  V00Ptr_(NULL),
403  SfPtr_(NULL),
404  magSfPtr_(NULL),
405  CPtr_(NULL),
406  CfPtr_(NULL),
407  phiPtr_(NULL)
408 {
409  if (debug)
410  {
411  InfoInFunction << "Constructing fvMesh from components" << endl;
412  }
413 }
414 
415 
416 Foam::fvMesh::fvMesh
417 (
418  const IOobject& io,
419  const Xfer<pointField>& points,
420  const Xfer<faceList>& faces,
421  const Xfer<cellList>& cells,
422  const bool syncPar
423 )
424 :
425  polyMesh(io, points, faces, cells, syncPar),
426  surfaceInterpolation(*this),
427  fvSchemes(static_cast<const objectRegistry&>(*this)),
428  fvSolution(static_cast<const objectRegistry&>(*this)),
429  data(static_cast<const objectRegistry&>(*this)),
430  boundary_(*this),
431  lduPtr_(NULL),
432  curTimeIndex_(time().timeIndex()),
433  VPtr_(NULL),
434  V0Ptr_(NULL),
435  V00Ptr_(NULL),
436  SfPtr_(NULL),
437  magSfPtr_(NULL),
438  CPtr_(NULL),
439  CfPtr_(NULL),
440  phiPtr_(NULL)
441 {
442  if (debug)
443  {
444  InfoInFunction << "Constructing fvMesh from components" << endl;
445  }
446 }
447 
448 
449 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
450 
452 {
453  clearOut();
454 }
455 
456 
457 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
458 
460 (
461  const List<polyPatch*> & p,
462  const bool validBoundary
463 )
464 {
465  if (boundary().size())
466  {
468  << " boundary already exists"
469  << abort(FatalError);
470  }
471 
472  // first add polyPatches
473  addPatches(p, validBoundary);
474  boundary_.addPatches(boundaryMesh());
475 }
476 
477 
479 {
480  if (debug)
481  {
482  InfoInFunction << "Removing boundary patches." << endl;
483  }
484 
485  // Remove fvBoundaryMesh data first.
486  boundary_.clear();
487  boundary_.setSize(0);
489 
490  clearOut();
491 }
492 
493 
495 {
496  if (debug)
497  {
498  InfoInFunction << "Updating fvMesh. ";
499  }
500 
502 
503  if (state == polyMesh::TOPO_PATCH_CHANGE)
504  {
505  if (debug)
506  {
507  Info<< "Boundary and topological update" << endl;
508  }
509 
510  boundary_.readUpdate(boundaryMesh());
511 
512  clearOut();
513 
514  }
515  else if (state == polyMesh::TOPO_CHANGE)
516  {
517  if (debug)
518  {
519  Info<< "Topological update" << endl;
520  }
521 
522  clearOut();
523  }
524  else if (state == polyMesh::POINTS_MOVED)
525  {
526  if (debug)
527  {
528  Info<< "Point motion update" << endl;
529  }
530 
531  clearGeom();
532  }
533  else
534  {
535  if (debug)
536  {
537  Info<< "No update" << endl;
538  }
539  }
540 
541  return state;
542 }
543 
544 
546 {
547  return boundary_;
548 }
549 
550 
552 {
553  if (!lduPtr_)
554  {
555  lduPtr_ = new fvMeshLduAddressing(*this);
556  }
557 
558  return *lduPtr_;
559 }
560 
561 
563 {
564  if (debug)
565  {
567  << " nOldCells:" << meshMap.nOldCells()
568  << " nCells:" << nCells()
569  << " nOldFaces:" << meshMap.nOldFaces()
570  << " nFaces:" << nFaces()
571  << endl;
572  }
573 
574 
575  // We require geometric properties valid for the old mesh
576  if
577  (
578  meshMap.cellMap().size() != nCells()
579  || meshMap.faceMap().size() != nFaces()
580  )
581  {
583  << "mapPolyMesh does not correspond to the old mesh."
584  << " nCells:" << nCells()
585  << " cellMap:" << meshMap.cellMap().size()
586  << " nOldCells:" << meshMap.nOldCells()
587  << " nFaces:" << nFaces()
588  << " faceMap:" << meshMap.faceMap().size()
589  << " nOldFaces:" << meshMap.nOldFaces()
590  << exit(FatalError);
591  }
592 
593  // Create a mapper
594  const fvMeshMapper mapper(*this, meshMap);
595 
596  // Map all the volFields in the objectRegistry
597  MapGeometricFields<scalar, fvPatchField, fvMeshMapper, volMesh>
598  (mapper);
599  MapGeometricFields<vector, fvPatchField, fvMeshMapper, volMesh>
600  (mapper);
601  MapGeometricFields<sphericalTensor, fvPatchField, fvMeshMapper, volMesh>
602  (mapper);
603  MapGeometricFields<symmTensor, fvPatchField, fvMeshMapper, volMesh>
604  (mapper);
605  MapGeometricFields<tensor, fvPatchField, fvMeshMapper, volMesh>
606  (mapper);
607 
608  // Map all the surfaceFields in the objectRegistry
609  MapGeometricFields<scalar, fvsPatchField, fvMeshMapper, surfaceMesh>
610  (mapper);
611  MapGeometricFields<vector, fvsPatchField, fvMeshMapper, surfaceMesh>
612  (mapper);
613  MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
614  (mapper);
615  MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
616  (mapper);
617  MapGeometricFields<tensor, fvsPatchField, fvMeshMapper, surfaceMesh>
618  (mapper);
619 
620  // Map all the dimensionedFields in the objectRegistry
621  MapDimensionedFields<scalar, fvMeshMapper, volMesh>(mapper);
622  MapDimensionedFields<vector, fvMeshMapper, volMesh>(mapper);
623  MapDimensionedFields<sphericalTensor, fvMeshMapper, volMesh>(mapper);
624  MapDimensionedFields<symmTensor, fvMeshMapper, volMesh>(mapper);
625  MapDimensionedFields<tensor, fvMeshMapper, volMesh>(mapper);
626 
627  // Map all the clouds in the objectRegistry
628  mapClouds(*this, meshMap);
629 
630 
631  const labelList& cellMap = meshMap.cellMap();
632 
633  // Map the old volume. Just map to new cell labels.
634  if (V0Ptr_)
635  {
636  scalarField& V0 = *V0Ptr_;
637 
638  scalarField savedV0(V0);
639  V0.setSize(nCells());
640 
641  forAll(V0, i)
642  {
643  if (cellMap[i] > -1)
644  {
645  V0[i] = savedV0[cellMap[i]];
646  }
647  else
648  {
649  V0[i] = 0.0;
650  }
651  }
652 
653  // Inject volume of merged cells
654  label nMerged = 0;
655  forAll(meshMap.reverseCellMap(), oldCelli)
656  {
657  label index = meshMap.reverseCellMap()[oldCelli];
658 
659  if (index < -1)
660  {
661  label celli = -index-2;
662 
663  V0[celli] += savedV0[oldCelli];
664 
665  nMerged++;
666  }
667  }
668 
669  if (debug)
670  {
671  Info<< "Mapping old time volume V0. Merged "
672  << nMerged << " out of " << nCells() << " cells" << endl;
673  }
674  }
675 
676 
677  // Map the old-old volume. Just map to new cell labels.
678  if (V00Ptr_)
679  {
680  scalarField& V00 = *V00Ptr_;
681 
682  scalarField savedV00(V00);
683  V00.setSize(nCells());
684 
685  forAll(V00, i)
686  {
687  if (cellMap[i] > -1)
688  {
689  V00[i] = savedV00[cellMap[i]];
690  }
691  else
692  {
693  V00[i] = 0.0;
694  }
695  }
696 
697  // Inject volume of merged cells
698  label nMerged = 0;
699  forAll(meshMap.reverseCellMap(), oldCelli)
700  {
701  label index = meshMap.reverseCellMap()[oldCelli];
702 
703  if (index < -1)
704  {
705  label celli = -index-2;
706 
707  V00[celli] += savedV00[oldCelli];
708  nMerged++;
709  }
710  }
711 
712  if (debug)
713  {
714  Info<< "Mapping old time volume V00. Merged "
715  << nMerged << " out of " << nCells() << " cells" << endl;
716  }
717  }
718 }
719 
720 
722 {
723  // Grab old time volumes if the time has been incremented
724  // This will update V0, V00
725  if (curTimeIndex_ < time().timeIndex())
726  {
727  storeOldVol(V());
728  }
729 
730  if (!phiPtr_)
731  {
732  // Create mesh motion flux
733  phiPtr_ = new surfaceScalarField
734  (
735  IOobject
736  (
737  "meshPhi",
738  this->time().timeName(),
739  *this,
742  false
743  ),
744  *this,
746  );
747  }
748  else
749  {
750  // Grab old time mesh motion fluxes if the time has been incremented
751  if (phiPtr_->timeIndex() != time().timeIndex())
752  {
753  phiPtr_->oldTime();
754  }
755  }
756 
757  surfaceScalarField& phi = *phiPtr_;
758 
759  // Move the polyMesh and set the mesh motion fluxes to the swept-volumes
760 
761  scalar rDeltaT = 1.0/time().deltaTValue();
762 
763  tmp<scalarField> tsweptVols = polyMesh::movePoints(p);
764  scalarField& sweptVols = tsweptVols.ref();
765 
766  phi.primitiveFieldRef() =
768  phi.primitiveFieldRef() *= rDeltaT;
769 
770  const fvPatchList& patches = boundary();
771 
772  surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
773 
774  forAll(patches, patchi)
775  {
776  phibf[patchi] = patches[patchi].patchSlice(sweptVols);
777  phibf[patchi] *= rDeltaT;
778  }
779 
780  // Update or delete the local geometric properties as early as possible so
781  // they can be used if necessary. These get recreated here instead of
782  // demand driven since they might do parallel transfers which can conflict
783  // with when they're actually being used.
784  // Note that between above "polyMesh::movePoints(p)" and here nothing
785  // should use the local geometric properties.
786  updateGeomNotOldVol();
787 
788 
789  // Update other local data
790  boundary_.movePoints();
792 
793  meshObject::movePoints<fvMesh>(*this);
794  meshObject::movePoints<lduMesh>(*this);
795 
796  return tsweptVols;
797 }
798 
799 
801 {
802  // Update polyMesh. This needs to keep volume existent!
804 
805  if (VPtr_)
806  {
807  // Grab old time volumes if the time has been incremented
808  // This will update V0, V00
809  storeOldVol(mpm.oldCellVolumes());
810 
811  // Few checks
812  if (VPtr_ && (V().size() != mpm.nOldCells()))
813  {
815  << "V:" << V().size()
816  << " not equal to the number of old cells "
817  << mpm.nOldCells()
818  << exit(FatalError);
819  }
820  if (V0Ptr_ && (V0Ptr_->size() != mpm.nOldCells()))
821  {
823  << "V0:" << V0Ptr_->size()
824  << " not equal to the number of old cells "
825  << mpm.nOldCells()
826  << exit(FatalError);
827  }
828  if (V00Ptr_ && (V00Ptr_->size() != mpm.nOldCells()))
829  {
831  << "V0:" << V00Ptr_->size()
832  << " not equal to the number of old cells "
833  << mpm.nOldCells()
834  << exit(FatalError);
835  }
836  }
837 
838 
839  // Clear mesh motion flux (note: could instead save & map like volumes)
840  deleteDemandDrivenData(phiPtr_);
841 
842  // Clear the sliced fields
843  clearGeomNotOldVol();
844 
845  // Map all fields
846  mapFields(mpm);
847 
848  // Clear the current volume and other geometry factors
850 
851  // Clear any non-updateable addressing
852  clearAddressing(true);
853 
854  meshObject::updateMesh<fvMesh>(*this, mpm);
855  meshObject::updateMesh<lduMesh>(*this, mpm);
856 }
857 
858 
860 (
864 ) const
865 {
866  return polyMesh::writeObject(fmt, ver, cmp);
867 }
868 
869 
871 {
872  bool ok = true;
873  if (phiPtr_)
874  {
875  ok = phiPtr_->write();
876  }
877 
878  return ok && polyMesh::write();
879 }
880 
881 
882 template<>
884 Foam::fvMesh::validComponents<Foam::sphericalTensor>() const
885 {
887 }
888 
889 
890 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
891 
892 bool Foam::fvMesh::operator!=(const fvMesh& bm) const
893 {
894  return &bm != this;
895 }
896 
897 
898 bool Foam::fvMesh::operator==(const fvMesh& bm) const
899 {
900  return &bm == this;
901 }
902 
903 
904 // ************************************************************************* //
Foam::surfaceFields.
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:230
virtual bool writeObjects(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:860
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1048
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
void clearAddressing()
Clear topological data.
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:478
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void mapClouds(const objectRegistry &db, const mapPolyMesh &mapper)
Generic Geometric field mapper.
Definition: mapClouds.H:48
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
static int debug
Debug switch.
Definition: fvSchemes.H:103
const surfaceScalarField & phi() const
Return cell face motion fluxes.
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:114
bool operator!=(const fvMesh &) const
Definition: fvMesh.C:892
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:59
bool isFile(const fileName &, const bool checkGzip=true)
Does the name exist as a FILE in the file system?
Definition: POSIX.C:492
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Traits class for primitives.
Definition: pTraits.H:50
bool movePoints()
Do what is neccessary if the mesh has moved.
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
bool operator==(const fvMesh &) const
Definition: fvMesh.C:898
label size() const
Return number of elements in table.
Cell to surface interpolation scheme. Included in fvMesh.
patches[0]
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
Class holds all the necessary information for mapping fields associated with fvMesh.
Definition: fvMeshMapper.H:55
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:460
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
bool moving() const
Is mesh moving.
Definition: polyMesh.H:493
label nCells() const
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:451
void readUpdate(const polyBoundaryMesh &)
Update boundary based on new polyBoundaryMesh.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:800
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:494
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
surfaceInterpolation(const fvMesh &)
Construct given an fvMesh.
A class for handling words, derived from string.
Definition: word.H:59
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp) const
Write the objects.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static void clearUpto(objectRegistry &)
Clear all meshObject derived from FromType up to (but not including)
Definition: MeshObject.C:401
label timeIndex() const
Return the time index of the field.
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
word timeName
Definition: getTimeIndex.H:3
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
Foam::fvMeshLduAddressing.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
defineTypeNameAndDebug(combustionModel, 0)
Database for solution data, solver performance and other reduced data.
Definition: data.H:52
volScalarField & C
label nOldFaces() const
Number of old faces.
Definition: mapPolyMesh.H:375
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
Generic Geometric field mapper. For "real" mapping, add template specialisations for mapping of inter...
label nOldCells() const
Number of old cells.
Definition: mapPolyMesh.H:381
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:878
label nFaces() const
void setSize(const label)
Reset size of List.
Definition: List.C:295
Template functions to aid in the implementation of demand driven data.
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:47
void clearOut()
Clear all geometry and addressing.
const labelList & cellMap() const
Old cell map.
Definition: mapPolyMesh.H:429
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
Foam::fvBoundaryMesh.
Selector class for finite volume differencing schemes. fvMesh is derived from fvShemes so that all fi...
Definition: fvSchemes.H:50
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label timeIndex
Definition: getTimeIndex.H:4
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const scalarField & oldCellVolumes() const
Definition: mapPolyMesh.H:645
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:526
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Version number type.
Definition: IOstream.H:96
virtual bool write() const
Write using setting from DB.
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:551
void movePoints()
Correct patches after moving points.
void clear()
Clear the PtrList, i.e. set size to zero deleting all the.
Definition: PtrList.C:174
messageStream Info
const labelList & faceMap() const
Old face map.
Definition: mapPolyMesh.H:404
SubField< scalar > subField
Declare type of subField.
Definition: Field.H:89
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:88
A class for managing temporary objects.
Definition: PtrList.H:54
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:65
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Registry of regIOobjects.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
void deleteDemandDrivenData(DataPtr &dataPtr)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545
label nInternalFaces() const
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:870
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: fvMesh.C:562
#define InfoInFunction
Report an information message using Foam::Info.