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