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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fvMesh.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "slicedVolFields.H"
33 #include "slicedSurfaceFields.H"
34 #include "SubField.H"
35 #include "demandDrivenData.H"
36 #include "fvMeshLduAddressing.H"
37 #include "mapPolyMesh.H"
38 #include "MapFvFields.H"
39 #include "fvMeshMapper.H"
40 #include "mapClouds.H"
41 #include "MeshObject.H"
42 #include "fvMatrix.H"
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48  defineTypeNameAndDebug(fvMesh, 0);
49 }
50 
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
54 void Foam::fvMesh::clearGeomNotOldVol()
55 {
57  <
58  fvMesh,
59  GeometricMeshObject,
60  MoveableMeshObject
61  >(*this);
62 
64  <
65  lduMesh,
66  GeometricMeshObject,
67  MoveableMeshObject
68  >(*this);
69 
70  slicedVolScalarField::Internal* VPtr =
71  static_cast<slicedVolScalarField::Internal*>(VPtr_);
73  VPtr_ = nullptr;
74 
75  deleteDemandDrivenData(SfPtr_);
76  deleteDemandDrivenData(magSfPtr_);
78  deleteDemandDrivenData(CfPtr_);
79 }
80 
81 
82 void Foam::fvMesh::updateGeomNotOldVol()
83 {
84  bool haveV = (VPtr_ != nullptr);
85  bool haveSf = (SfPtr_ != nullptr);
86  bool haveMagSf = (magSfPtr_ != nullptr);
87  bool haveCP = (CPtr_ != nullptr);
88  bool haveCf = (CfPtr_ != nullptr);
89 
90  clearGeomNotOldVol();
91 
92  // Now recreate the fields
93  if (haveV)
94  {
95  (void)V();
96  }
97  if (haveSf)
98  {
99  (void)Sf();
100  }
101  if (haveMagSf)
102  {
103  (void)magSf();
104  }
105  if (haveCP)
106  {
107  (void)C();
108  }
109  if (haveCf)
110  {
111  (void)Cf();
112  }
113 }
114 
115 
116 void Foam::fvMesh::clearGeom()
117 {
118  clearGeomNotOldVol();
119 
120  deleteDemandDrivenData(V0Ptr_);
121  deleteDemandDrivenData(V00Ptr_);
122 
123  // Mesh motion flux cannot be deleted here because the old-time flux
124  // needs to be saved.
125 }
126 
127 
128 void Foam::fvMesh::clearAddressing(const bool isMeshUpdate)
129 {
130  if (debug)
131  {
132  InfoInFunction << "isMeshUpdate: " << isMeshUpdate << endl;
133  }
134 
135  if (isMeshUpdate)
136  {
137  // Part of a mesh update. Keep meshObjects that have an updateMesh
138  // callback
140  <
141  fvMesh,
142  TopologicalMeshObject,
143  UpdateableMeshObject
144  >
145  (
146  *this
147  );
149  <
150  lduMesh,
151  TopologicalMeshObject,
152  UpdateableMeshObject
153  >
154  (
155  *this
156  );
157  }
158  else
159  {
160  meshObject::clear<fvMesh, TopologicalMeshObject>(*this);
161  meshObject::clear<lduMesh, TopologicalMeshObject>(*this);
162  }
163  deleteDemandDrivenData(lduPtr_);
164 }
165 
166 
167 void Foam::fvMesh::storeOldVol(const scalarField& V)
168 {
169  if (curTimeIndex_ < time().timeIndex())
170  {
171  if (debug)
172  {
174  << " Storing old time volumes since from time " << curTimeIndex_
175  << " and time now " << time().timeIndex()
176  << " V:" << V.size()
177  << endl;
178  }
179 
180 
181  if (V00Ptr_ && V0Ptr_)
182  {
183  // Copy V0 into V00 storage
184  *V00Ptr_ = *V0Ptr_;
185  }
186 
187  if (V0Ptr_)
188  {
189  // Copy V into V0 storage
190  V0Ptr_->scalarField::operator=(V);
191  }
192  else
193  {
194  // Allocate V0 storage, fill with V
195  V0Ptr_ = new DimensionedField<scalar, volMesh>
196  (
197  IOobject
198  (
199  "V0",
200  time().timeName(),
201  *this,
204  false
205  ),
206  *this,
207  dimVolume
208  );
209  scalarField& V0 = *V0Ptr_;
210  // Note: V0 now sized with current mesh, not with (potentially
211  // different size) V.
212  V0.setSize(V.size());
213  V0 = V;
214  }
215 
216  curTimeIndex_ = time().timeIndex();
217 
218  if (debug)
219  {
221  << " Stored old time volumes V0:" << V0Ptr_->size()
222  << endl;
223  if (V00Ptr_)
224  {
226  << " Stored oldold time volumes V00:" << V00Ptr_->size()
227  << endl;
228  }
229  }
230  }
231 }
232 
233 
235 {
236  clearGeom();
238 
239  clearAddressing();
240 
241  // Clear mesh motion flux
242  deleteDemandDrivenData(phiPtr_);
243 
245 }
246 
247 
248 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
249 
250 Foam::fvMesh::fvMesh(const IOobject& io)
251 :
252  polyMesh(io),
253  surfaceInterpolation(*this),
254  fvSchemes(static_cast<const objectRegistry&>(*this)),
255  fvSolution(static_cast<const objectRegistry&>(*this)),
256  data(static_cast<const objectRegistry&>(*this)),
257  boundary_(*this, boundaryMesh()),
258  lduPtr_(nullptr),
259  curTimeIndex_(time().timeIndex()),
260  VPtr_(nullptr),
261  V0Ptr_(nullptr),
262  V00Ptr_(nullptr),
263  SfPtr_(nullptr),
264  magSfPtr_(nullptr),
265  CPtr_(nullptr),
266  CfPtr_(nullptr),
267  phiPtr_(nullptr)
268 {
269  if (debug)
270  {
271  InfoInFunction << "Constructing fvMesh from IOobject" << endl;
272  }
273 
274  // Check the existence of the cell volumes and read if present
275  // and set the storage of V00
276  if (fileHandler().isFile(time().timePath()/dbDir()/"V0"))
277  {
279  (
280  IOobject
281  (
282  "V0",
283  time().timeName(),
284  *this,
287  false
288  ),
289  *this
290  );
291 
292  V00();
293  }
294 
295  // Check the existence of the mesh fluxes, read if present and set the
296  // mesh to be moving
297  if (fileHandler().isFile(time().timePath()/dbDir()/"meshPhi"))
298  {
299  phiPtr_ = new surfaceScalarField
300  (
301  IOobject
302  (
303  "meshPhi",
304  time().timeName(),
305  *this,
308  false
309  ),
310  *this
311  );
312 
313  // The mesh is now considered moving so the old-time cell volumes
314  // will be required for the time derivatives so if they haven't been
315  // read initialise to the current cell volumes
316  if (!V0Ptr_)
317  {
319  (
320  IOobject
321  (
322  "V0",
323  time().timeName(),
324  *this,
327  false
328  ),
329  V()
330  );
331  }
332 
333  moving(true);
334  }
335 }
336 
337 
338 Foam::fvMesh::fvMesh
339 (
340  const IOobject& io,
341  pointField&& points,
342  faceList&& faces,
343  labelList&& allOwner,
344  labelList&& allNeighbour,
345  const bool syncPar
346 )
347 :
348  polyMesh
349  (
350  io,
351  std::move(points),
352  std::move(faces),
353  std::move(allOwner),
354  std::move(allNeighbour),
355  syncPar
356  ),
357  surfaceInterpolation(*this),
358  fvSchemes(static_cast<const objectRegistry&>(*this)),
359  fvSolution(static_cast<const objectRegistry&>(*this)),
360  data(static_cast<const objectRegistry&>(*this)),
361  boundary_(*this, boundaryMesh()),
362  lduPtr_(nullptr),
363  curTimeIndex_(time().timeIndex()),
364  VPtr_(nullptr),
365  V0Ptr_(nullptr),
366  V00Ptr_(nullptr),
367  SfPtr_(nullptr),
368  magSfPtr_(nullptr),
369  CPtr_(nullptr),
370  CfPtr_(nullptr),
371  phiPtr_(nullptr)
372 {
373  if (debug)
374  {
375  InfoInFunction << "Constructing fvMesh from components" << endl;
376  }
377 }
378 
379 
380 Foam::fvMesh::fvMesh
381 (
382  const IOobject& io,
383  pointField&& points,
384  faceList&& faces,
385  cellList&& cells,
386  const bool syncPar
387 )
388 :
389  polyMesh
390  (
391  io,
392  std::move(points),
393  std::move(faces),
394  std::move(cells),
395  syncPar
396  ),
397  surfaceInterpolation(*this),
398  fvSchemes(static_cast<const objectRegistry&>(*this)),
399  fvSolution(static_cast<const objectRegistry&>(*this)),
400  data(static_cast<const objectRegistry&>(*this)),
401  boundary_(*this),
402  lduPtr_(nullptr),
403  curTimeIndex_(time().timeIndex()),
404  VPtr_(nullptr),
405  V0Ptr_(nullptr),
406  V00Ptr_(nullptr),
407  SfPtr_(nullptr),
408  magSfPtr_(nullptr),
409  CPtr_(nullptr),
410  CfPtr_(nullptr),
411  phiPtr_(nullptr)
412 {
413  if (debug)
414  {
415  InfoInFunction << "Constructing fvMesh from components" << endl;
416  }
417 }
418 
419 
420 Foam::fvMesh::fvMesh(const IOobject& io, const zero, const bool syncPar)
421 :
422  fvMesh(io, pointField(), faceList(), labelList(), labelList(), syncPar)
423 {}
424 
425 
426 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
427 
429 {
430  clearOut();
431 }
432 
433 
434 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
435 
437 (
438  fvMatrix<scalar>& m,
439  const dictionary& dict
440 ) const
441 {
442  // Redirect to fvMatrix solver
443  return m.solveSegregatedOrCoupled(dict);
444 }
445 
446 
448 (
449  fvMatrix<vector>& m,
450  const dictionary& dict
451 ) const
452 {
453  // Redirect to fvMatrix solver
454  return m.solveSegregatedOrCoupled(dict);
455 }
456 
457 
459 (
461  const dictionary& dict
462 ) const
463 {
464  // Redirect to fvMatrix solver
465  return m.solveSegregatedOrCoupled(dict);
466 }
467 
468 
470 (
472  const dictionary& dict
473 ) const
474 {
475  // Redirect to fvMatrix solver
476  return m.solveSegregatedOrCoupled(dict);
477 }
478 
479 
481 (
482  fvMatrix<tensor>& m,
483  const dictionary& dict
484 ) const
485 {
486  // Redirect to fvMatrix solver
487  return m.solveSegregatedOrCoupled(dict);
488 }
489 
490 
492 (
493  PtrList<polyPatch>& plist,
494  const bool validBoundary
495 )
496 {
497  if (boundary().size())
498  {
500  << " boundary already exists"
501  << abort(FatalError);
502  }
503 
504  addPatches(plist, validBoundary);
505  boundary_.addPatches(boundaryMesh());
506 }
507 
508 
510 (
511  const List<polyPatch*>& p,
512  const bool validBoundary
513 )
514 {
515  // Acquire ownership of the pointers
516  PtrList<polyPatch> plist(const_cast<List<polyPatch*>&>(p));
517 
518  addFvPatches(plist, validBoundary);
519 }
520 
521 
523 {
524  if (debug)
525  {
526  InfoInFunction << "Removing boundary patches." << endl;
527  }
528 
529  // Remove fvBoundaryMesh data first.
530  boundary_.clear();
531  boundary_.setSize(0);
533 
534  clearOut();
535 }
536 
537 
539 {
540  if (debug)
541  {
542  InfoInFunction << "Updating fvMesh. ";
543  }
544 
546 
547  if (state == polyMesh::TOPO_PATCH_CHANGE)
548  {
549  if (debug)
550  {
551  Info<< "Boundary and topological update" << endl;
552  }
553 
554  boundary_.readUpdate(boundaryMesh());
555 
556  clearOut();
557 
558  }
559  else if (state == polyMesh::TOPO_CHANGE)
560  {
561  if (debug)
562  {
563  Info<< "Topological update" << endl;
564  }
565 
566  clearOut();
567  }
568  else if (state == polyMesh::POINTS_MOVED)
569  {
570  if (debug)
571  {
572  Info<< "Point motion update" << endl;
573  }
574 
575  clearGeom();
576  }
577  else
578  {
579  if (debug)
580  {
581  Info<< "No update" << endl;
582  }
583  }
584 
585  return state;
586 }
587 
588 
590 {
591  return boundary_;
592 }
593 
594 
596 {
597  if (!lduPtr_)
598  {
599  if (debug)
600  {
602  << " calculating fvMeshLduAddressing from nFaces:"
603  << nFaces() << endl;
604  }
605 
606  lduPtr_ = new fvMeshLduAddressing(*this);
607  }
608 
609  return *lduPtr_;
610 }
611 
612 
614 {
615  if (debug)
616  {
618  << " nOldCells:" << meshMap.nOldCells()
619  << " nCells:" << nCells()
620  << " nOldFaces:" << meshMap.nOldFaces()
621  << " nFaces:" << nFaces()
622  << endl;
623  }
624 
625 
626  // We require geometric properties valid for the old mesh
627  if
628  (
629  meshMap.cellMap().size() != nCells()
630  || meshMap.faceMap().size() != nFaces()
631  )
632  {
634  << "mapPolyMesh does not correspond to the old mesh."
635  << " nCells:" << nCells()
636  << " cellMap:" << meshMap.cellMap().size()
637  << " nOldCells:" << meshMap.nOldCells()
638  << " nFaces:" << nFaces()
639  << " faceMap:" << meshMap.faceMap().size()
640  << " nOldFaces:" << meshMap.nOldFaces()
641  << exit(FatalError);
642  }
643 
644  // Create a mapper
645  const fvMeshMapper mapper(*this, meshMap);
646 
647  // Map all the volFields in the objectRegistry
648  MapGeometricFields<scalar, fvPatchField, fvMeshMapper, volMesh>
649  (mapper);
650  MapGeometricFields<vector, fvPatchField, fvMeshMapper, volMesh>
651  (mapper);
652  MapGeometricFields<sphericalTensor, fvPatchField, fvMeshMapper, volMesh>
653  (mapper);
654  MapGeometricFields<symmTensor, fvPatchField, fvMeshMapper, volMesh>
655  (mapper);
656  MapGeometricFields<tensor, fvPatchField, fvMeshMapper, volMesh>
657  (mapper);
658 
659  // Map all the surfaceFields in the objectRegistry
660  MapGeometricFields<scalar, fvsPatchField, fvMeshMapper, surfaceMesh>
661  (mapper);
662  MapGeometricFields<vector, fvsPatchField, fvMeshMapper, surfaceMesh>
663  (mapper);
664  MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
665  (mapper);
666  MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
667  (mapper);
668  MapGeometricFields<tensor, fvsPatchField, fvMeshMapper, surfaceMesh>
669  (mapper);
670 
671  // Map all the dimensionedFields in the objectRegistry
672  MapDimensionedFields<scalar, fvMeshMapper, volMesh>(mapper);
673  MapDimensionedFields<vector, fvMeshMapper, volMesh>(mapper);
674  MapDimensionedFields<sphericalTensor, fvMeshMapper, volMesh>(mapper);
675  MapDimensionedFields<symmTensor, fvMeshMapper, volMesh>(mapper);
676  MapDimensionedFields<tensor, fvMeshMapper, volMesh>(mapper);
677 
678  // Map all the clouds in the objectRegistry
679  mapClouds(*this, meshMap);
680 
681 
682  const labelList& cellMap = meshMap.cellMap();
683 
684  // Map the old volume. Just map to new cell labels.
685  if (V0Ptr_)
686  {
687  scalarField& V0 = *V0Ptr_;
688 
689  scalarField savedV0(V0);
690  V0.setSize(nCells());
691 
692  forAll(V0, i)
693  {
694  if (cellMap[i] > -1)
695  {
696  V0[i] = savedV0[cellMap[i]];
697  }
698  else
699  {
700  V0[i] = 0.0;
701  }
702  }
703 
704  // Inject volume of merged cells
705  label nMerged = 0;
706  forAll(meshMap.reverseCellMap(), oldCelli)
707  {
708  label index = meshMap.reverseCellMap()[oldCelli];
709 
710  if (index < -1)
711  {
712  label celli = -index-2;
713 
714  V0[celli] += savedV0[oldCelli];
715 
716  nMerged++;
717  }
718  }
719 
720  if (debug)
721  {
722  Info<< "Mapping old time volume V0. Merged "
723  << nMerged << " out of " << nCells() << " cells" << endl;
724  }
725  }
726 
727 
728  // Map the old-old volume. Just map to new cell labels.
729  if (V00Ptr_)
730  {
731  scalarField& V00 = *V00Ptr_;
732 
733  scalarField savedV00(V00);
734  V00.setSize(nCells());
735 
736  forAll(V00, i)
737  {
738  if (cellMap[i] > -1)
739  {
740  V00[i] = savedV00[cellMap[i]];
741  }
742  else
743  {
744  V00[i] = 0.0;
745  }
746  }
747 
748  // Inject volume of merged cells
749  label nMerged = 0;
750  forAll(meshMap.reverseCellMap(), oldCelli)
751  {
752  label index = meshMap.reverseCellMap()[oldCelli];
753 
754  if (index < -1)
755  {
756  label celli = -index-2;
757 
758  V00[celli] += savedV00[oldCelli];
759  nMerged++;
760  }
761  }
762 
763  if (debug)
764  {
765  Info<< "Mapping old time volume V00. Merged "
766  << nMerged << " out of " << nCells() << " cells" << endl;
767  }
768  }
769 }
770 
771 
773 {
774  // Grab old time volumes if the time has been incremented
775  // This will update V0, V00
776  if (curTimeIndex_ < time().timeIndex())
777  {
778  storeOldVol(V());
779  }
780 
781  if (!phiPtr_)
782  {
783  // Create mesh motion flux
784  phiPtr_ = new surfaceScalarField
785  (
786  IOobject
787  (
788  "meshPhi",
789  this->time().timeName(),
790  *this,
793  false
794  ),
795  *this,
797  );
798  }
799  else
800  {
801  // Grab old time mesh motion fluxes if the time has been incremented
802  if (phiPtr_->timeIndex() != time().timeIndex())
803  {
804  phiPtr_->oldTime();
805  }
806  }
807 
808  surfaceScalarField& phi = *phiPtr_;
809 
810  // Move the polyMesh and set the mesh motion fluxes to the swept-volumes
811 
812  scalar rDeltaT = 1.0/time().deltaTValue();
813 
815  scalarField& sweptVols = tsweptVols.ref();
816 
818  scalarField::subField(sweptVols, nInternalFaces());
819  phi.primitiveFieldRef() *= rDeltaT;
820 
821  const fvPatchList& patches = boundary();
822 
823  surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
824 
825  forAll(patches, patchi)
826  {
827  phibf[patchi] = patches[patchi].patchSlice(sweptVols);
828  phibf[patchi] *= rDeltaT;
829  }
830 
831  // Update or delete the local geometric properties as early as possible so
832  // they can be used if necessary. These get recreated here instead of
833  // demand driven since they might do parallel transfers which can conflict
834  // with when they're actually being used.
835  // Note that between above "polyMesh::movePoints(p)" and here nothing
836  // should use the local geometric properties.
837  updateGeomNotOldVol();
838 
839 
840  // Update other local data
841  boundary_.movePoints();
843 
844  meshObject::movePoints<fvMesh>(*this);
845  meshObject::movePoints<lduMesh>(*this);
846 
847  return tsweptVols;
848 }
849 
850 
852 {
853  // Update polyMesh. This needs to keep volume existent!
855 
856  // Our slice of the addressing is no longer valid
857  deleteDemandDrivenData(lduPtr_);
858 
859  if (VPtr_)
860  {
861  // Grab old time volumes if the time has been incremented
862  // This will update V0, V00
863  storeOldVol(mpm.oldCellVolumes());
864 
865  // Few checks
866  if (VPtr_ && (V().size() != mpm.nOldCells()))
867  {
869  << "V:" << V().size()
870  << " not equal to the number of old cells "
871  << mpm.nOldCells()
872  << exit(FatalError);
873  }
874  if (V0Ptr_ && (V0Ptr_->size() != mpm.nOldCells()))
875  {
877  << "V0:" << V0Ptr_->size()
878  << " not equal to the number of old cells "
879  << mpm.nOldCells()
880  << exit(FatalError);
881  }
882  if (V00Ptr_ && (V00Ptr_->size() != mpm.nOldCells()))
883  {
885  << "V0:" << V00Ptr_->size()
886  << " not equal to the number of old cells "
887  << mpm.nOldCells()
888  << exit(FatalError);
889  }
890  }
891 
892 
893  // Clear mesh motion flux (note: could instead save & map like volumes)
894  deleteDemandDrivenData(phiPtr_);
895 
896  // Clear the sliced fields
897  clearGeomNotOldVol();
898 
899  // Map all fields
900  mapFields(mpm);
901 
902  // Clear the current volume and other geometry factors
904 
905  // Clear any non-updateable addressing
906  clearAddressing(true);
907 
908  meshObject::updateMesh<fvMesh>(*this, mpm);
909  meshObject::updateMesh<lduMesh>(*this, mpm);
910 }
911 
912 
914 (
918  const bool valid
919 ) const
920 {
921  bool ok = true;
922  if (phiPtr_)
923  {
924  ok = phiPtr_->write(valid);
925  // NOTE: The old old time mesh phi might be necessary for certain
926  // solver smooth restart using second order time schemes.
927  //ok = phiPtr_->oldTime().write();
928  }
929  if (V0Ptr_ && V0Ptr_->writeOpt() == IOobject::AUTO_WRITE)
930  {
931  // For second order restarts we need to write V0
932  ok = V0Ptr_->write(valid);
933  }
934 
935  return ok && polyMesh::writeObject(fmt, ver, cmp, valid);
936 }
937 
938 
939 bool Foam::fvMesh::write(const bool valid) const
940 {
941  return polyMesh::write(valid);
942 }
943 
944 
945 template<>
947 Foam::fvMesh::validComponents<Foam::sphericalTensor>() const
948 {
950 }
951 
952 
953 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
954 
955 bool Foam::fvMesh::operator!=(const fvMesh& rhs) const
956 {
957  return &rhs != this;
958 }
959 
960 
961 bool Foam::fvMesh::operator==(const fvMesh& rhs) const
962 {
963  return &rhs == this;
964 }
965 
966 
967 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::lduAddressing
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Definition: lduAddressing.H:114
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::surfaceInterpolation::clearOut
void clearOut()
Clear all geometry and addressing.
Definition: surfaceInterpolation.C:48
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::fvMesh::~fvMesh
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:428
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
SubField.H
Foam::fvMesh::write
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:939
Foam::polyMesh::POINTS_MOVED
Definition: polyMesh.H:94
Foam::primitiveMesh::clearAddressing
void clearAddressing()
Clear topological data.
Definition: primitiveMeshClear.C:143
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::polyMesh::movePoints
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1107
Foam::fvMatrix::solveSegregatedOrCoupled
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:62
Foam::fvMesh::writeObject
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool valid) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:914
Foam::fvSchemes
Selector class for finite volume differencing schemes. fvMesh is derived from fvShemes so that all fi...
Definition: fvSchemes.H:52
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:513
Foam::polyMesh::dbDir
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:798
mapPolyMesh.H
Foam::fvMesh::mapFields
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: fvMesh.C:613
Foam::mapPolyMesh::nOldFaces
label nOldFaces() const
Number of old faces.
Definition: mapPolyMesh.H:380
Foam::fvMesh::removeFvBoundary
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:522
Foam::isFile
bool isFile(const fileName &name, const bool checkGzip=true, const bool followLink=true)
Does the name exist as a FILE in the file system?
Definition: MSwindows.C:658
Foam::fileHandler
const fileOperation & fileHandler()
Get current file handler.
Definition: fileOperation.C:1181
fvMatrix.H
Foam::polyMesh::clearOut
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
Definition: polyMeshClear.C:232
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
surfaceFields.H
Foam::surfaceFields.
Foam::fvMesh::operator!=
bool operator!=(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:955
Foam::fvMesh::V00
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
Definition: fvMeshGeometry.C:247
Foam::fvMesh::operator==
bool operator==(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:961
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
C
volScalarField & C
Definition: readThermalProperties.H:102
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
mapClouds.H
Generic Geometric field mapper. For "real" mapping, add template specialisations for mapping of inter...
Foam::fvMeshLduAddressing
Foam::fvMeshLduAddressing.
Definition: fvMeshLduAddressing.H:51
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:57
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::fvMesh::readUpdate
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:538
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::polyMesh::TOPO_PATCH_CHANGE
Definition: polyMesh.H:96
Foam::fvSchemes::debug
static int debug
Debug switch.
Definition: fvSchemes.H:105
Foam::fvMesh::solve
virtual SolverPerformance< scalar > solve(fvMatrix< scalar > &, const dictionary &) const
Definition: fvMesh.C:437
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
fvMeshMapper.H
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:62
Foam::IOstreamOption::versionNumber
Representation of a major/minor version number.
Definition: IOstreamOption.H:79
Foam::Field< vector >
Foam::regIOobject::write
virtual bool write(const bool valid=true) const
Write using setting from DB.
Definition: regIOobjectWrite.C:170
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fvMesh::lduAddr
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:595
Foam::polyMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
Definition: polyMeshUpdate.C:41
Foam::fvSolution
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:49
Foam::fvMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:851
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::polyMesh::TOPO_CHANGE
Definition: polyMesh.H:95
Foam::objectRegistry::writeObject
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool valid) const
Write the objects.
Definition: objectRegistry.C:475
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::IOstreamOption::streamFormat
streamFormat
Data format (ascii | binary)
Definition: IOstreamOption.H:64
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::fvMesh::addFvPatches
void addFvPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:492
Foam::surfaceInterpolation::movePoints
bool movePoints()
Do what is necessary 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:84
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::mapPolyMesh::oldCellVolumes
const scalarField & oldCellVolumes() const
Definition: mapPolyMesh.H:650
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::polyMesh::removeBoundary
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:38
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:735
Foam::polyMesh::readUpdateState
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:589
Foam::mapPolyMesh::reverseCellMap
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:531
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::fvMeshMapper
Class holds all the necessary information for mapping fields associated with fvMesh.
Definition: fvMeshMapper.H:57
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:752
Foam::List< face >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:52
Foam::mapClouds
void mapClouds(const objectRegistry &db, const mapPolyMesh &mapper)
Generic Geometric field mapper.
Definition: mapClouds.H:51
Foam::mapPolyMesh::nOldCells
label nOldCells() const
Number of old cells.
Definition: mapPolyMesh.H:386
Foam::meshObject::clearUpto
static void clearUpto(objectRegistry &obr)
Definition: MeshObject.C:217
points
const pointField & points
Definition: gmvOutputHeader.H:1
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::mapPolyMesh::faceMap
const labelList & faceMap() const
Old face map.
Definition: mapPolyMesh.H:409
Foam::surfaceInterpolation
Cell to surface interpolation scheme. Included in fvMesh.
Definition: surfaceInterpolation.H:54
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
slicedSurfaceFields.H
MeshObject.H
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
Foam::IOstreamOption::compressionType
compressionType
Compression treatment (UNCOMPRESSED | COMPRESSED)
Definition: IOstreamOption.H:71
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:61
Foam::Field< scalar >::subField
SubField< scalar > subField
Declare type of subField.
Definition: Field.H:90
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:51
fvMeshLduAddressing.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:54
boundary
faceListList boundary
Definition: createBlockMesh.H:4
Foam::mapPolyMesh::cellMap
const labelList & cellMap() const
Old cell map.
Definition: mapPolyMesh.H:434
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::fvMesh::clearOut
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:234
Foam::polyMesh::readUpdate
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:77
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:190
Foam::zero
A class representing the concept of 0 (zero), which can be used to avoid manipulating objects that ar...
Definition: zero.H:61
Foam::IOobject::MUST_READ
Definition: IOobject.H:120