Loading...
Searching...
No Matches
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,2022 OpenFOAM Foundation
9 Copyright (C) 2016-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "fvMeshLduAddressing.H"
36#include "mapPolyMesh.H"
37#include "MapFvFields.H"
38#include "fvMeshMapper.H"
39#include "mapClouds.H"
40#include "MeshObject.H"
41#include "fvMatrix.H"
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
45namespace Foam
48}
49
50
51// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52
54{
56 <
57 fvMesh,
60 >(*this);
61
63 <
64 lduMesh,
67 >(*this);
68
69 VPtr_.reset(nullptr);
70 SfPtr_.reset(nullptr);
71 magSfPtr_.reset(nullptr);
72 CPtr_.reset(nullptr);
73 CfPtr_.reset(nullptr);
74}
75
76
78{
79 const bool haveV = (VPtr_ != nullptr);
80 const bool haveSf = (SfPtr_ != nullptr);
81 const bool haveMagSf = (magSfPtr_ != nullptr);
82 const bool haveCP = (CPtr_ != nullptr);
83 const bool haveCf = (CfPtr_ != nullptr);
84
85 clearGeomNotOldVol();
86
87 // Now recreate the fields
88 if (haveV)
89 {
90 (void)V();
91 }
92
93 if (haveSf)
94 {
95 (void)Sf();
96 }
97
98 if (haveMagSf)
99 {
100 (void)magSf();
101 }
102
103 if (haveCP)
104 {
105 (void)C();
106 }
107
108 if (haveCf)
109 {
110 (void)Cf();
111 }
112}
113
114
116{
117 clearGeomNotOldVol();
118
119 V0Ptr_.reset(nullptr);
120 V00Ptr_.reset(nullptr);
121
122 // Mesh motion flux cannot be deleted here because the old-time flux
123 // needs to be saved.
124}
125
126
127void Foam::fvMesh::clearAddressing(const bool isMeshUpdate)
128{
129 DebugInFunction << "isMeshUpdate: " << isMeshUpdate << endl;
130
131 if (isMeshUpdate)
132 {
133 // Part of a mesh update. Keep meshObjects that have an updateMesh
134 // callback
136 <
137 fvMesh,
140 >
141 (
142 *this
143 );
145 <
146 lduMesh,
149 >
150 (
151 *this
152 );
153 }
154 else
155 {
158 }
159
160 lduPtr_.reset(nullptr);
161}
162
163
165{
166 if (curTimeIndex_ < time().timeIndex())
167 {
169 << " Storing old time volumes since from time " << curTimeIndex_
170 << " and time now " << time().timeIndex()
171 << " V:" << V.size() << endl;
172
173 if (V00Ptr_ && V0Ptr_)
174 {
175 // Copy V0 into V00 storage
176 *V00Ptr_ = *V0Ptr_;
177 }
178
179 if (!V0Ptr_)
180 {
181 V0Ptr_ = std::make_unique<DimensionedField<scalar, volMesh>>
182 (
184 (
185 "V0",
186 this->time().timeName(),
187 *this,
191 ),
192 *this,
194 );
195 // Note: V0 now sized with current mesh, not with (potentially
196 // different size) V.
197 V0Ptr_->scalarField::resize_nocopy(V.size());
198 }
199
200 // Copy V into V0 storage
201 V0Ptr_->scalarField::operator=(V);
202
203 curTimeIndex_ = time().timeIndex();
204
205 if (debug)
206 {
208 << " Stored old time volumes V0:" << V0Ptr_->size()
209 << endl;
210
211 if (V00Ptr_)
212 {
214 << " Stored oldold time volumes V00:" << V00Ptr_->size()
216 }
217 }
218 }
219}
220
221
222void Foam::fvMesh::clearOutLocal(const bool isMeshUpdate)
223{
224 clearGeom();
226
227 clearAddressing(isMeshUpdate);
228
229 // Clear mesh motion flux
230 phiPtr_.reset(nullptr);
231}
232
233
234void Foam::fvMesh::clearOut(const bool isMeshUpdate)
236 clearOutLocal(isMeshUpdate);
237
238 polyMesh::clearOut(isMeshUpdate);
239}
240
241
243{
244 // Clear mesh motion flux
245 phiPtr_.reset(nullptr);
246}
247
248
249// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
250
251Foam::fvMesh::fvMesh(const IOobject& io, const bool doInit)
252:
253 polyMesh(io, doInit),
254 fvSchemes(static_cast<const objectRegistry&>(*this)),
256 fvSolution(static_cast<const objectRegistry&>(*this)),
257 boundary_(*this, boundaryMesh()),
258 curTimeIndex_(time().timeIndex())
259{
260 DebugInFunction << "Constructing fvMesh from IOobject" << endl;
261
262 if (doInit)
263 {
264 fvMesh::init(false); // do not initialise lower levels
265 }
266}
267
268
269bool Foam::fvMesh::init(const bool doInit)
270{
271 if (doInit)
272 {
273 // Construct basic geometry calculation engine. Note: do before
274 // doing anything with primitiveMesh::cellCentres etc.
275 (void)geometry();
276
277 // Initialise my data
278 polyMesh::init(doInit);
279 }
280
281 // Read some optional fields
282 // ~~~~~~~~~~~~~~~~~~~~~~~~~
283 // For redistributePar + fileHandler can have master processor
284 // find the file but not the subprocs.
285
286 IOobject rio
287 (
288 "none",
289 this->time().timeName(),
290 *this,
294 );
295
296
297 // Read old cell volumes (if present) and set the storage of V00
298 rio.resetHeader("V0");
299 if (returnReduceOr(rio.typeHeaderOk<regIOobject>(false)))
300 {
302 << "Detected V0: " << rio.objectRelPath() << nl;
303
304 V0Ptr_ = std::make_unique<DimensionedField<scalar, volMesh>>
305 (
306 rio,
307 *this,
308 dimensionedScalar(dimVol, Foam::zero{})
309 );
310
311 // Set the moving flag early in case demand-driven geometry
312 // construction checks for it
313 moving(true);
314
315 (void)V00();
316 }
317
318
319 // Read mesh fluxes (if present) and set the mesh to be moving
320 rio.resetHeader("meshPhi");
321 if (returnReduceOr(rio.typeHeaderOk<regIOobject>(false)))
322 {
324 << "Detected meshPhi: " << rio.objectRelPath() << nl;
325
326 // Clear mesh motion flux
327 phiPtr_.reset(nullptr);
328
329 phiPtr_ = std::make_unique<surfaceScalarField>
330 (
331 rio,
332 *this,
333 dimensionedScalar(dimVol/dimTime, Foam::zero{})
334 );
335
336 // Set the moving flag early in case demand-driven geometry
337 // construction checks for it
338 moving(true);
339
340 // The mesh is now considered moving so the old-time cell volumes
341 // will be required for the time derivatives so if they haven't been
342 // read initialise to the current cell volumes
343 if (!V0Ptr_)
344 {
345 V0Ptr_ = std::make_unique<DimensionedField<scalar, volMesh>>
346 (
348 (
349 "V0",
350 this->time().timeName(),
351 *this,
355 ),
356 V()
357 );
358 }
360
361 // Assume something changed
362 return true;
363}
364
365
367(
368 const IOobject& io,
370 faceList&& faces,
371 labelList&& allOwner,
372 labelList&& allNeighbour,
373 const bool syncPar
374)
375:
376 polyMesh
377 (
378 io,
379 std::move(points),
380 std::move(faces),
381 std::move(allOwner),
382 std::move(allNeighbour),
383 syncPar
384 ),
385 fvSchemes(static_cast<const objectRegistry&>(*this)),
386 surfaceInterpolation(*this),
387 fvSolution(static_cast<const objectRegistry&>(*this)),
388 boundary_(*this),
390{
391 DebugInFunction << "Constructing fvMesh from components" << endl;
392}
393
394
396(
397 const IOobject& io,
399 faceList&& faces,
400 cellList&& cells,
401 const bool syncPar
402)
403:
405 (
406 io,
407 std::move(points),
408 std::move(faces),
409 std::move(cells),
410 syncPar
411 ),
412 fvSchemes(static_cast<const objectRegistry&>(*this)),
414 fvSolution(static_cast<const objectRegistry&>(*this)),
415 boundary_(*this),
417{
418 DebugInFunction << "Constructing fvMesh from components" << endl;
419}
420
422Foam::fvMesh::fvMesh(const IOobject& io, const Foam::zero, const bool syncPar)
423:
424 fvMesh(io, pointField(), faceList(), labelList(), labelList(), syncPar)
425{}
426
427
429(
430 const IOobject& io,
431 const fvMesh& baseMesh,
432 const Foam::zero,
433 const bool syncPar
434)
435:
436 fvMesh
437 (
438 io,
439 baseMesh,
440 pointField(),
441 faceList(),
442 labelList(), // owner
443 labelList(), // neighbour
444 syncPar
445 )
446{}
447
448
450(
451 const IOobject& io,
452 const fvMesh& baseMesh,
454 faceList&& faces,
455 labelList&& allOwner,
456 labelList&& allNeighbour,
457 const bool syncPar
458)
459:
461 (
462 io,
463 std::move(points),
464 std::move(faces),
465 std::move(allOwner),
466 std::move(allNeighbour),
467 syncPar
468 ),
470 (
471 static_cast<const objectRegistry&>(*this),
472 io.readOpt(),
473 static_cast<const dictionary*>(baseMesh.hasSchemes())
474 ),
477 (
478 static_cast<const objectRegistry&>(*this),
479 io.readOpt(),
480 static_cast<const dictionary*>(baseMesh.hasSolution())
481 ),
482 boundary_(*this),
483 curTimeIndex_(time().timeIndex())
484{
485 DebugInFunction << "Constructing fvMesh as copy and primitives" << endl;
486
487 // Reset mesh data
488 data().reset(baseMesh.data());
489}
490
491
493(
494 const IOobject& io,
495 const fvMesh& baseMesh,
497 faceList&& faces,
498 cellList&& cells,
499 const bool syncPar
500)
501:
503 (
504 io,
505 std::move(points),
506 std::move(faces),
507 std::move(cells),
508 syncPar
509 ),
511 (
512 static_cast<const objectRegistry&>(*this),
513 io.readOpt(),
514 static_cast<const dictionary*>(baseMesh.hasSchemes())
515 ),
518 (
519 static_cast<const objectRegistry&>(*this),
520 io.readOpt(),
521 static_cast<const dictionary*>(baseMesh.hasSolution())
522 ),
523 boundary_(*this),
524 curTimeIndex_(time().timeIndex())
525{
526 DebugInFunction << "Constructing fvMesh as copy and primitives" << endl;
527
528 // Reset mesh data
529 data().reset(baseMesh.data());
530}
531
532
533// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
534
537 clearOut();
538}
539
540
541// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
544{
545 return static_cast<const fvSchemes*>(this);
546}
547
550{
551 return static_cast<const fvSolution*>(this);
552}
553
556{
557 return static_cast<const fvSchemes&>(*this);
558}
559
562{
563 return static_cast<fvSchemes&>(*this);
564}
565
568{
569 return static_cast<const fvSolution&>(*this);
570}
571
582 const dictionary& dict
583) const
584{
585 // Redirect to fvMatrix solver
587}
588
589
591(
593 const dictionary& dict
594) const
595{
596 // Redirect to fvMatrix solver
598}
599
600
602(
604 const dictionary& dict
605) const
606{
607 // Redirect to fvMatrix solver
609}
610
611
613(
615 const dictionary& dict
616) const
617{
618 // Redirect to fvMatrix solver
620}
621
622
624(
626 const dictionary& dict
627) const
628{
629 // Redirect to fvMatrix solver
631}
632
633
635(
636 polyPatchList& plist,
637 const bool validBoundary
638)
639{
640 if (boundary().size())
641 {
643 << " boundary already exists"
644 << abort(FatalError);
646
647 addPatches(plist, validBoundary);
648 boundary_.addPatches(boundaryMesh());
649}
650
651
653(
654 const List<polyPatch*>& p,
655 const bool validBoundary
656)
657{
658 // Acquire ownership of the pointers
659 polyPatchList plist(const_cast<List<polyPatch*>&>(p));
660
661 addFvPatches(plist, validBoundary);
662}
663
664
666{
667 DebugInFunction << "Removing boundary patches." << endl;
668
669 // Remove fvBoundaryMesh data first.
670 boundary_.clear();
672
673 clearOut();
674}
675
676
678{
679 DebugInFunction << "Updating fvMesh";
680
682
683 if (state == polyMesh::TOPO_PATCH_CHANGE)
684 {
685 DebugInfo << "Boundary and topological update" << endl;
686
687 boundary_.readUpdate(boundaryMesh());
688
689 clearOut();
690
691 }
692 else if (state == polyMesh::TOPO_CHANGE)
693 {
694 DebugInfo << "Topological update" << endl;
695
696 // fvMesh::clearOut() but without the polyMesh::clearOut
697 clearOutLocal();
698 }
699 else if (state == polyMesh::POINTS_MOVED)
700 {
701 DebugInfo << "Point motion update" << endl;
702
703 clearGeom();
704 }
705 else
706 {
707 DebugInfo << "No update" << endl;
708 }
709
710 return state;
711}
712
713
715{
716 if (!lduPtr_)
717 {
719 << "Calculating fvMeshLduAddressing from nFaces:"
720 << nFaces() << endl;
721
722 lduPtr_ = std::make_unique<fvMeshLduAddressing>(*this);
723
724 return *lduPtr_;
725 }
726
727 return *lduPtr_;
728}
729
732{
733 return boundary().interfaces();
734}
735
736
737void Foam::fvMesh::mapFields(const mapPolyMesh& meshMap)
738{
740 << " nOldCells:" << meshMap.nOldCells()
741 << " nCells:" << nCells()
742 << " nOldFaces:" << meshMap.nOldFaces()
743 << " nFaces:" << nFaces()
744 << endl;
745
746 // We require geometric properties valid for the old mesh
747 if
748 (
749 meshMap.cellMap().size() != nCells()
750 || meshMap.faceMap().size() != nFaces()
751 )
752 {
754 << "mapPolyMesh does not correspond to the old mesh."
755 << " nCells:" << nCells()
756 << " cellMap:" << meshMap.cellMap().size()
757 << " nOldCells:" << meshMap.nOldCells()
758 << " nFaces:" << nFaces()
759 << " faceMap:" << meshMap.faceMap().size()
760 << " nOldFaces:" << meshMap.nOldFaces()
761 << exit(FatalError);
762 }
763
764 // Create a mapper
765 const fvMeshMapper mapper(*this, meshMap);
766
767 // Map all the volFields in the objectRegistry
769 (mapper);
771 (mapper);
773 (mapper);
775 (mapper);
777 (mapper);
778
779 // Map all the surfaceFields in the objectRegistry
781 (mapper);
783 (mapper);
785 <
786 sphericalTensor, fvsPatchField, fvMeshMapper, surfaceMesh
787 >
788 (mapper);
790 (mapper);
792 (mapper);
793
794 // Map all the dimensionedFields in the objectRegistry
800
801 // Map all the clouds in the objectRegistry
802 mapClouds(*this, meshMap);
803
804
805 const labelList& cellMap = meshMap.cellMap();
806
807 // Map the old volume. Just map to new cell labels.
808 if (V0Ptr_)
809 {
810 scalarField& V0 = *V0Ptr_;
811
812 scalarField savedV0(V0);
813 V0.resize_nocopy(nCells());
814
815 forAll(V0, i)
816 {
817 if (cellMap[i] > -1)
818 {
819 V0[i] = savedV0[cellMap[i]];
820 }
821 else
822 {
823 V0[i] = 0.0;
824 }
825 }
826
827 // Inject volume of merged cells
828 label nMerged = 0;
829 forAll(meshMap.reverseCellMap(), oldCelli)
830 {
831 label index = meshMap.reverseCellMap()[oldCelli];
832
833 if (index < -1)
834 {
835 label celli = -index-2;
836
837 V0[celli] += savedV0[oldCelli];
838
839 nMerged++;
840 }
841 }
842
844 << "Mapping old time volume V0. Merged "
845 << nMerged << " out of " << nCells() << " cells" << endl;
846 }
847
848
849 // Map the old-old volume. Just map to new cell labels.
850 if (V00Ptr_)
851 {
852 scalarField& V00 = *V00Ptr_;
853
854 scalarField savedV00(V00);
855 V00.resize_nocopy(nCells());
856
857 forAll(V00, i)
858 {
859 if (cellMap[i] > -1)
860 {
861 V00[i] = savedV00[cellMap[i]];
862 }
863 else
864 {
865 V00[i] = 0.0;
866 }
867 }
868
869 // Inject volume of merged cells
870 label nMerged = 0;
871 forAll(meshMap.reverseCellMap(), oldCelli)
872 {
873 label index = meshMap.reverseCellMap()[oldCelli];
874
875 if (index < -1)
876 {
877 label celli = -index-2;
878
879 V00[celli] += savedV00[oldCelli];
880 nMerged++;
881 }
882 }
883
885 << "Mapping old time volume V00. Merged "
886 << nMerged << " out of " << nCells() << " cells" << endl;
887 }
888}
889
890
892{
894
895 // Grab old time volumes if the time has been incremented
896 // This will update V0, V00
897 if (curTimeIndex_ < time().timeIndex())
898 {
899 storeOldVol(V());
900 }
901
902
903 // Move the polyMesh and initialise the mesh motion fluxes field
904 // Note: mesh flux updated by the fvGeometryScheme
905
906 if (!phiPtr_)
907 {
908 DebugInFunction<< "Creating initial meshPhi field" << endl;
909
910 // Create mesh motion flux
911 phiPtr_ = std::make_unique<surfaceScalarField>
912 (
914 (
915 "meshPhi",
916 this->time().timeName(),
917 *this,
921 ),
922 *this,
923 dimensionedScalar(dimVolume/dimTime, Foam::zero{})
924 );
925 }
926 else
927 {
928 // Grab old time mesh motion fluxes if the time has been incremented
929 if (phiPtr_->timeIndex() != time().timeIndex())
930 {
931 DebugInFunction<< "Accessing old-time meshPhi field" << endl;
932 phiPtr_->oldTime();
933 }
934 }
935
937
938 // Update or delete the local geometric properties as early as possible so
939 // they can be used if necessary. These get recreated here instead of
940 // demand driven since they might do parallel transfers which can conflict
941 // with when they're actually being used.
942 // Note that between above "polyMesh::movePoints(p)" and here nothing
943 // should use the local geometric properties.
944 updateGeomNotOldVol();
945
946 // Update other local data
947 boundary_.movePoints();
948
949 // Clear weights, deltaCoeffs, nonOrthoDeltaCoeffs, nonOrthCorrectionVectors
958{
961 // Let surfaceInterpolation handle geometry calculation. Note: this does
962 // lower levels updateGeom
964}
965
966
968{
970
971 // Update polyMesh. This needs to keep volume existent!
973
974 // Our slice of the addressing is no longer valid
975 lduPtr_.reset(nullptr);
976
977 if (VPtr_)
978 {
979 // Grab old time volumes if the time has been incremented
980 // This will update V0, V00
981 storeOldVol(mpm.oldCellVolumes());
982
983 // Few checks
984 if (VPtr_ && (VPtr_->size() != mpm.nOldCells()))
985 {
987 << "V:" << VPtr_->size()
988 << " not equal to the number of old cells "
989 << mpm.nOldCells()
990 << exit(FatalError);
991 }
992 if (V0Ptr_ && (V0Ptr_->size() != mpm.nOldCells()))
993 {
995 << "V0:" << V0Ptr_->size()
996 << " not equal to the number of old cells "
997 << mpm.nOldCells()
998 << exit(FatalError);
999 }
1000 if (V00Ptr_ && (V00Ptr_->size() != mpm.nOldCells()))
1001 {
1003 << "V0:" << V00Ptr_->size()
1004 << " not equal to the number of old cells "
1005 << mpm.nOldCells()
1006 << exit(FatalError);
1007 }
1008 }
1009
1010
1011 // Clear mesh motion flux (note: could instead save & map like volumes)
1012 if (phiPtr_)
1013 {
1014 // Mesh moving and topology change. Recreate meshPhi
1015 phiPtr_.reset(nullptr);
1016
1017 // Create mesh motion flux
1018 phiPtr_ = std::make_unique<surfaceScalarField>
1019 (
1020 IOobject
1021 (
1022 "meshPhi",
1023 this->time().timeName(),
1024 *this,
1028 ),
1029 *this,
1030 dimensionedScalar(dimVolume/dimTime, Foam::zero{})
1031 );
1032 }
1033
1034 // Clear the sliced fields
1035 clearGeomNotOldVol();
1036
1037 // Map all fields
1038 mapFields(mpm);
1039
1040 // Clear the current volume and other geometry factors
1042
1043 // Clear any non-updateable addressing
1045
1048}
1049
1050
1052(
1053 IOstreamOption streamOpt,
1054 const bool writeOnProc
1055) const
1056{
1057 bool ok = true;
1058 if (phiPtr_)
1059 {
1060 ok = phiPtr_->write(writeOnProc);
1061 // NOTE: The old old time mesh phi might be necessary for certain
1062 // solver smooth restart using second order time schemes.
1063 //ok = phiPtr_->oldTime().write();
1064 }
1065 if (V0Ptr_ && V0Ptr_->writeOpt() == IOobject::AUTO_WRITE)
1066 {
1067 // For second order restarts we need to write V0
1068 ok = V0Ptr_->write(writeOnProc);
1069 }
1070
1071 return ok && polyMesh::writeObject(streamOpt, writeOnProc);
1072}
1073
1074
1075bool Foam::fvMesh::write(const bool writeOnProc) const
1076{
1077 return polyMesh::write(writeOnProc);
1078}
1079
1080
1081template<>
1082typename Foam::pTraits<Foam::sphericalTensor>::labelType
1083Foam::fvMesh::validComponents<Foam::sphericalTensor>() const
1088
1089// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1091bool Foam::fvMesh::operator!=(const fvMesh& rhs) const
1092{
1093 return &rhs != this;
1094}
1095
1096
1097bool Foam::fvMesh::operator==(const fvMesh& rhs) const
1098{
1099 return &rhs == this;
1100}
1101
1102
1103// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
label size() const noexcept
Definition HashTable.H:358
@ NO_REGISTER
Do not request registration (bool: false).
readOption readOpt() const noexcept
Get the read option.
@ NO_READ
Nothing to be read.
@ LAZY_READ
Reading is optional [identical to READ_IF_PRESENT].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info. A void type suppresses trait and t...
void resetHeader(const word &newName=word::null)
Clear various bits (headerClassName, note, sizeof...) that would be obtained when reading from a file...
Definition IOobject.C:644
fileName objectRelPath() const
The object path relative to the case.
Definition IOobject.C:581
A simple container for options an IOstream can normally have.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Class holds all the necessary information for mapping fields associated with fvMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
std::unique_ptr< slicedSurfaceVectorField > SfPtr_
Face area vectors.
Definition fvMesh.H:127
virtual bool movePoints()
Avoid masking surfaceInterpolation method.
pTraits< Type >::labelType validComponents() const
Return a labelType of valid component indicators.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write the underlying polyMesh and other data.
Definition fvMesh.C:1045
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
std::unique_ptr< surfaceScalarField > magSfPtr_
Mag face area vectors.
Definition fvMesh.H:132
virtual SolverPerformance< scalar > solve(fvMatrix< scalar > &, const dictionary &) const
Solve returning the solution statistics given convergence tolerance. Use the given solver controls.
Definition fvMesh.C:573
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition fvMesh.C:724
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
void addFvPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition fvMesh.C:628
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition fvMesh.C:707
std::unique_ptr< fvMeshLduAddressing > lduPtr_
Definition fvMesh.H:98
std::unique_ptr< DimensionedField< scalar, volMesh > > V00Ptr_
Cell volumes old-old time level.
Definition fvMesh.H:122
void clearMeshPhi()
Clear cell face motion fluxes.
Definition fvMesh.C:235
const fvSchemes & schemes() const
Read-access to the fvSchemes controls.
Definition fvMesh.C:548
std::unique_ptr< DimensionedField< scalar, volMesh > > V0Ptr_
Cell volumes old time level.
Definition fvMesh.H:117
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition fvMesh.C:262
fvBoundaryMesh boundary_
Boundary mesh.
Definition fvMesh.H:93
const fvSolution * hasSolution() const
Non-null if fvSolution exists (can test as bool).
Definition fvMesh.C:542
void clearGeom()
Clear local geometry.
Definition fvMesh.C:108
void clearOutLocal(const bool isMeshUpdate=false)
Clear local-only storage (geometry, addressing etc).
Definition fvMesh.C:215
fvMesh(const fvMesh &)=delete
No copy construct.
virtual ~fvMesh()
Destructor.
Definition fvMesh.C:528
void updateGeomNotOldVol()
Clear geometry like clearGeomNotOldVol but recreate any.
Definition fvMesh.C:70
void clearGeomNotOldVol()
Clear geometry but not the old-time cell volumes.
Definition fvMesh.C:46
const fvSolution & solution() const
Read-access to the fvSolution controls.
Definition fvMesh.C:560
void storeOldVol(const scalarField &)
Preserve old volume(s).
Definition fvMesh.C:157
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
void clearAddressing(const bool isMeshUpdate=false)
Clear local addressing.
Definition fvMesh.C:120
bool operator!=(const fvMesh &rhs) const
Compares addresses.
Definition fvMesh.C:1084
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition fvMesh.C:960
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition fvMesh.C:670
const surfaceVectorField & Sf() const
Return cell face area vectors.
std::unique_ptr< slicedVolVectorField > CPtr_
Cell centres.
Definition fvMesh.H:137
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
bool operator==(const fvMesh &rhs) const
Compares addresses.
Definition fvMesh.C:1090
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition fvMesh.C:1068
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to these fvPatches.
Definition fvMesh.C:658
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition fvMesh.C:730
std::unique_ptr< slicedSurfaceVectorField > CfPtr_
Face centres.
Definition fvMesh.H:142
label curTimeIndex_
Current time index for cell volumes.
Definition fvMesh.H:107
const fvSchemes * hasSchemes() const
Non-null if fvSchemes exists (can test as bool).
Definition fvMesh.C:536
void clearOut(const bool isMeshUpdate=false)
Clear all geometry and addressing.
Definition fvMesh.C:227
std::unique_ptr< surfaceScalarField > phiPtr_
Face motion fluxes.
Definition fvMesh.H:147
virtual void updateGeom()
Update all geometric data. This gets redirected up from primitiveMesh level.
Definition fvMesh.C:950
std::unique_ptr< SlicedDimensionedField< scalar, volMesh > > VPtr_
Cell volumes.
Definition fvMesh.H:112
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition fvSchemes.H:54
fvSchemes(const fvSchemes &)=delete
No copy construct.
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition fvSolution.H:54
fvSolution(const fvSolution &)=delete
No copy construct.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition lduMesh.H:54
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reverseCellMap() const noexcept
Reverse cell map.
label nOldCells() const noexcept
Number of old cells.
const labelList & faceMap() const noexcept
Old face map.
label nOldFaces() const noexcept
Number of old faces.
const labelList & cellMap() const noexcept
Old cell map.
const scalarField & oldCellVolumes() const
static void clearUpto(objectRegistry &obr)
Clear all meshObject derived from FromType up to (but not including) ToType.
static void clear(objectRegistry &obr)
Clear/remove all meshObject of MeshObjectType via objectRegistry::checkOut().
static void updateMesh(objectRegistry &obr, const mapPolyMesh &mpm)
Update topology using the given map.
static void movePoints(objectRegistry &obr)
Update for mesh motion.
void reset(const meshState &ms)
Reset the dictionary.
Definition meshState.C:95
Registry of regIOobjects.
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
void addPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches.
Definition polyMesh.C:956
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc=true) const
Write items held in the objectRegistry. Normally includes mesh components (points,...
Definition polyMesh.C:1577
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
bool moving() const noexcept
Is mesh moving.
Definition polyMesh.H:732
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition polyMesh.C:341
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition polyMesh.H:92
virtual const meshState & data() const noexcept
Const reference to the mesh and solver state data.
Definition polyMesh.H:559
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition polyMeshIO.C:68
virtual void movePoints(const pointField &)
Move points.
Definition polyMesh.C:1153
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
void removeBoundary()
Remove boundary patches.
label nCells() const noexcept
Number of mesh cells.
void clearAddressing()
Clear topological data.
label nFaces() const noexcept
Number of mesh faces.
const cellList & cells() const
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition regIOobject.H:71
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Cell to surface interpolation scheme. Included in fvMesh.
virtual const fvGeometryScheme & geometry() const
Return reference to geometry calculation scheme.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh for topology changes.
surfaceInterpolation(const fvMesh &)
Construct given an fvMesh.
void clearOut()
Clear all geometry and addressing.
virtual void updateGeom()
Update all geometric data.
Mesh data needed to do the Finite Volume discretisation.
Definition surfaceMesh.H:47
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
faceListList boundary
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
const pointField & points
const cellShapeList & cells
Generic Geometric field mapper. For "real" mapping, add template specialisations for mapping of inter...
word timeName
Definition getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
const expr V(m.psi().mesh().V())
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition polyPatch.H:61
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
void mapClouds(const objectRegistry &db, const mapPolyMesh &mapper)
Generic Geometric field mapper.
Definition mapClouds.H:48
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
errorManip< error > abort(error &err)
Definition errorManip.H:139
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
void MapGeometricFields(const MeshMapper &mapper)
Generic Geometric field mapper.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
vectorField pointField
pointField is a vectorField.
void MapDimensionedFields(const MeshMapper &mapper)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label timeIndex
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.