Loading...
Searching...
No Matches
isoAdvection.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) 2016-2017 DHI
9 Modified code Copyright (C) 2016-2022 OpenCFD Ltd.
10 Modified code Copyright (C) 2019-2020 DLR
11 Modified code Copyright (C) 2018, 2021 Johan Roenby
12-------------------------------------------------------------------------------
13License
14 This file is part of OpenFOAM.
15
16 OpenFOAM is free software: you can redistribute it and/or modify it
17 under the terms of the GNU General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24 for more details.
25
26 You should have received a copy of the GNU General Public License
27 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
30
31#include "isoAdvection.H"
32#include "volFields.H"
35#include "fvcSurfaceIntegrate.H"
36#include "fvcGrad.H"
37#include "upwind.H"
38#include "cellSet.H"
39#include "meshTools.H"
40#include "OBJstream.H"
41#include "syncTools.H"
42#include "profiling.H"
45
46// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47
48namespace Foam
49{
50 defineTypeNameAndDebug(isoAdvection, 0);
51}
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
55Foam::isoAdvection::isoAdvection
56(
59 const volVectorField& U
60)
61:
62 mesh_(alpha1.mesh()),
63 dict_(mesh_.solverDict(alpha1.name())),
64 alpha1_(alpha1),
65 alpha1In_(alpha1.primitiveFieldRef()),
66 phi_(phi),
67 U_(U),
68 dVf_
69 (
71 (
72 "dVf_",
73 mesh_.time().timeName(),
74 mesh_,
75 IOobject::NO_READ,
76 IOobject::NO_WRITE
77 ),
78 mesh_,
80 ),
81 alphaPhi_
82 (
84 (
85 "alphaPhi_",
86 mesh_.time().timeName(),
87 mesh_,
88 IOobject::NO_READ,
89 IOobject::NO_WRITE
90 ),
91 mesh_,
93 ),
94 advectionTime_(0),
95 timeIndex_(-1),
96
97 // Tolerances and solution controls
98 nAlphaBounds_(dict_.getOrDefault<label>("nAlphaBounds", 3)),
99 isoFaceTol_(dict_.getOrDefault<scalar>("isoFaceTol", 1e-8)),
100 surfCellTol_(dict_.getOrDefault<scalar>("surfCellTol", 1e-8)),
101 writeIsoFacesToFile_(dict_.getOrDefault("writeIsoFaces", false)),
102
103 // Cell cutting data
104 surfCells_(label(0.2*mesh_.nCells())),
105 surf_(reconstructionSchemes::New(alpha1_, phi_, U_, dict_)),
106 advectFace_(alpha1.mesh(), alpha1),
107 bsFaces_(label(0.2*mesh_.nBoundaryFaces())),
108 bsx0_(bsFaces_.size()),
109 bsn0_(bsFaces_.size()),
110 bsUn0_(bsFaces_.size()),
111
112 // Porosity
113 porosityEnabled_(dict_.getOrDefault<bool>("porosityEnabled", false)),
114 porosityPtr_(nullptr),
115
116 // Parallel run data
117 procPatchLabels_(mesh_.boundary().size()),
118 surfaceCellFacesOnProcPatches_(0)
119{
121
122 // Prepare lists used in parallel runs
123 if (Pstream::parRun())
124 {
125 // Force calculation of required demand driven data (else parallel
126 // communication may crash)
127 mesh_.cellCentres();
128 mesh_.cellVolumes();
129 mesh_.faceCentres();
130 mesh_.faceAreas();
131 mesh_.magSf();
132 mesh_.boundaryMesh().patchID();
133 mesh_.cellPoints();
134 mesh_.cellCells();
135 mesh_.cells();
136
137 // Get boundary mesh and resize the list for parallel comms
138 setProcessorPatches();
139 }
140
141 // Reading porosity properties from constant directory
142 IOdictionary porosityProperties
143 (
144 IOobject
145 (
146 "porosityProperties",
147 mesh_.time().constant(),
148 mesh_,
151 )
152 );
153
154 porosityEnabled_ =
155 porosityProperties.getOrDefault<bool>("porosityEnabled", false);
156
157 if (porosityEnabled_)
158 {
159 porosityPtr_ = mesh_.getObjectPtr<volScalarField>("porosity");
160
161 if (porosityPtr_)
162 {
163 auto limits = gMinMax(porosityPtr_->primitiveField());
164
165 if (limits.min() <= 0 || limits.max() > 1 + SMALL)
166 {
168 << "Porosity field has values <= 0 or > 1"
169 << exit(FatalError);
170 }
171 }
172 else
173 {
175 << "Porosity enabled in constant/porosityProperties "
176 << "but no porosity field is found in object registry."
177 << exit(FatalError);
178 }
179 }
180}
181
182
183// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
184
185void Foam::isoAdvection::setProcessorPatches()
186{
187 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
188 surfaceCellFacesOnProcPatches_.clear();
189 surfaceCellFacesOnProcPatches_.resize(patches.size());
190
191 // Append all processor patch labels to the list
192 procPatchLabels_.clear();
193 forAll(patches, patchi)
194 {
195 if
196 (
198 && !patches[patchi].empty()
199 )
200 {
201 procPatchLabels_.append(patchi);
202 }
203 }
204}
205
206
207void Foam::isoAdvection::extendMarkedCells
208(
209 bitSet& markedCell
210) const
211{
212 // Mark faces using any marked cell
213 bitSet markedFace(mesh_.nFaces());
214
215 for (const label celli : markedCell)
216 {
217 markedFace.set(mesh_.cells()[celli]); // set multiple faces
218 }
219
220 syncTools::syncFaceList(mesh_, markedFace, orEqOp<unsigned int>());
221
222 // Update cells using any markedFace
223 for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
224 {
225 if (markedFace.test(facei))
226 {
227 markedCell.set(mesh_.faceOwner()[facei]);
228 markedCell.set(mesh_.faceNeighbour()[facei]);
229 }
230 }
231 for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); ++facei)
232 {
233 if (markedFace.test(facei))
234 {
235 markedCell.set(mesh_.faceOwner()[facei]);
236 }
237 }
238}
239
240
241void Foam::isoAdvection::timeIntegratedFlux()
242{
243 addProfilingInFunction(geometricVoF);
244 // Get time step
245 const scalar dt = mesh_.time().deltaTValue();
246
247 // Create object for interpolating velocity to isoface centres
249
250 // For each downwind face of each surface cell we "isoadvect" to find dVf
251 label nSurfaceCells = 0;
252
253 // Clear out the data for re-use and reset list containing information
254 // whether cells could possibly need bounding
255 clearIsoFaceData();
256
257 // Get necessary references
258 const scalarField& phiIn = phi_.primitiveField();
259 const scalarField& magSfIn = mesh_.magSf().primitiveField();
260 scalarField& dVfIn = dVf_.primitiveFieldRef();
261
262 // Get necessary mesh data
263 const cellList& cellFaces = mesh_.cells();
264 const labelList& own = mesh_.faceOwner();
265
266
267 // Storage for isoFace points. Only used if writeIsoFacesToFile_
268 DynamicList<List<point>> isoFacePts;
269 const labelUList& interfaceLabels = surf_->interfaceLabels();
270
271 // Calculating isoface normal velocity
272 scalarField Un0(interfaceLabels.size());
273 forAll(Un0, i)
274 {
275 const label celli = interfaceLabels[i];
276 const point x0(surf_->centre()[celli]);
277 const vector n0(normalised(-surf_->normal()[celli]));
278 Un0[i] = UInterp.interpolate(x0, celli) & n0;
279 }
280
281 // Taking acount of porosity if enabled
282 if (porosityEnabled_)
283 {
284 forAll(Un0, i)
285 {
286 const label celli = interfaceLabels[i];
287 Un0[i] /= porosityPtr_->primitiveField()[celli];
288 }
289 }
290
291 // Loop through cells
292 forAll(interfaceLabels, i)
293 {
294 const label celli = interfaceLabels[i];
295 if (mag(surf_->normal()[celli]) != 0)
296 {
297
298 // This is a surface cell, increment counter, append and mark cell
299 nSurfaceCells++;
300 surfCells_.append(celli);
301 const point x0(surf_->centre()[celli]);
302 const vector n0(normalised(-surf_->normal()[celli]));
303
305 << "\n------------ Cell " << celli << " with alpha1 = "
306 << alpha1In_[celli] << " and 1-alpha1 = "
307 << 1.0 - alpha1In_[celli] << " ------------"
308 << endl;
309
310 // Estimate time integrated flux through each downwind face
311 // Note: looping over all cell faces - in reduced-D, some of
312 // these faces will be on empty patches
313 const cell& celliFaces = cellFaces[celli];
314 for (const label facei : celliFaces)
315 {
316 if (mesh_.isInternalFace(facei))
317 {
318 bool isDownwindFace = false;
319
320 if (celli == own[facei])
321 {
322 if (phiIn[facei] >= 0)
323 {
324 isDownwindFace = true;
325 }
326 }
327 else
328 {
329 if (phiIn[facei] < 0)
330 {
331 isDownwindFace = true;
332 }
333 }
334
335 if (isDownwindFace)
336 {
337 dVfIn[facei] = advectFace_.timeIntegratedFaceFlux
338 (
339 facei,
340 x0,
341 n0,
342 Un0[i],
343 dt,
344 phiIn[facei],
345 magSfIn[facei]
346 );
347 }
348
349 }
350 else
351 {
352 bsFaces_.append(facei);
353 bsx0_.append(x0);
354 bsn0_.append(n0);
355 bsUn0_.append(Un0[i]);
356
357 // Note: we must not check if the face is on the
358 // processor patch here.
359 }
360 }
361 }
362 }
363
364 // Get references to boundary fields
365 const polyBoundaryMesh& boundaryMesh = mesh_.boundaryMesh();
366 const surfaceScalarField::Boundary& phib = phi_.boundaryField();
367 const surfaceScalarField::Boundary& magSfb = mesh_.magSf().boundaryField();
368 surfaceScalarField::Boundary& dVfb = dVf_.boundaryFieldRef();
369
370 // Loop through boundary surface faces
371 forAll(bsFaces_, i)
372 {
373 // Get boundary face index (in the global list)
374 const label facei = bsFaces_[i];
375 const label patchi = boundaryMesh.patchID(facei);
376
377 if (!phib[patchi].empty())
378 {
379 const label patchFacei = boundaryMesh[patchi].whichFace(facei);
380
381 const scalar phiP = phib[patchi][patchFacei];
382
383 if (phiP >= 0)
384 {
385 const scalar magSf = magSfb[patchi][patchFacei];
386
387 dVfb[patchi][patchFacei] = advectFace_.timeIntegratedFaceFlux
388 (
389 facei,
390 bsx0_[i],
391 bsn0_[i],
392 bsUn0_[i],
393 dt,
394 phiP,
395 magSf
396 );
397
398 // Handling upwind cyclic boundary patches
399 const polyPatch& pp = boundaryMesh[patchi];
401 if (cpp)
402 {
403 const label neiPatchID(cpp->neighbPolyPatchID());
404 dVfb[neiPatchID][patchFacei] = -dVfb[patchi][patchFacei];
405 }
406
407 // Check if the face is on processor patch and append it to
408 // the list if necessary
409 checkIfOnProcPatch(facei);
410 }
411 }
412 }
413
414 // Synchronize processor patches
415 syncProcPatches(dVf_, phi_);
416
417 writeIsoFaces(isoFacePts);
418
419 DebugInfo << "Number of isoAdvector surface cells = "
420 << returnReduce(nSurfaceCells, sumOp<label>()) << endl;
421}
422
423
424void Foam::isoAdvection::setDownwindFaces
425(
426 const label celli,
427 DynamicLabelList& downwindFaces
428) const
429{
431
432 // Get necessary mesh data and cell information
433 const labelList& own = mesh_.faceOwner();
434 const cellList& cells = mesh_.cells();
435 const cell& c = cells[celli];
436
437 downwindFaces.clear();
438
439 // Check all faces of the cell
440 for (const label facei: c)
441 {
442 // Get face and corresponding flux
443 const scalar phi = faceValue(phi_, facei);
444
445 if (own[facei] == celli)
446 {
447 if (phi >= 0)
448 {
449 downwindFaces.append(facei);
450 }
451 }
452 else if (phi < 0)
453 {
454 downwindFaces.append(facei);
455 }
456 }
457
458 downwindFaces.shrink();
459}
460
461
462Foam::scalar Foam::isoAdvection::netFlux
463(
464 const surfaceScalarField& dVf,
465 const label celli
466) const
467{
468 scalar dV = 0;
469
470 // Get face indices
471 const cell& c = mesh_.cells()[celli];
472
473 // Get mesh data
474 const labelList& own = mesh_.faceOwner();
475
476 for (const label facei : c)
477 {
478 const scalar dVff = faceValue(dVf, facei);
479
480 if (own[facei] == celli)
481 {
482 dV += dVff;
483 }
484 else
485 {
486 dV -= dVff;
487 }
488 }
489
490 return dV;
491}
492
493
494Foam::DynamicList<Foam::label> Foam::isoAdvection::syncProcPatches
495(
497 const surfaceScalarField& phi,
498 bool returnSyncedFaces
499)
500{
501 DynamicLabelList syncedFaces(0);
502 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
503
504 if (Pstream::parRun())
505 {
506 DynamicList<label> neighProcs;
507 PstreamBuffers pBufs;
508
509 // Send
510 for (const label patchi : procPatchLabels_)
511 {
512 const processorPolyPatch& procPatch =
514 const label nbrProci = procPatch.neighbProcNo();
515
516 // Neighbour connectivity
517 neighProcs.push_uniq(nbrProci);
518
519 const scalarField& pFlux = dVf.boundaryField()[patchi];
520 const List<label>& surfCellFacesOnProcPatch =
521 surfaceCellFacesOnProcPatches_[patchi];
522
523 const UIndirectList<scalar> dVfPatch
524 (
525 pFlux,
526 surfCellFacesOnProcPatch
527 );
528
529 UOPstream toNbr(nbrProci, pBufs);
530 toNbr << surfCellFacesOnProcPatch << dVfPatch;
531 }
532
533 // Limited to involved neighbour procs
534 pBufs.finishedNeighbourSends(neighProcs);
535
536
537 // Receive and combine
538 for (const label patchi : procPatchLabels_)
539 {
540 const processorPolyPatch& procPatch =
542 const label nbrProci = procPatch.neighbProcNo();
543
544 List<label> faceIDs;
545 List<scalar> nbrdVfs;
546
547 {
548 UIPstream fromNbr(nbrProci, pBufs);
549 fromNbr >> faceIDs >> nbrdVfs;
550 }
551
552 if (returnSyncedFaces)
553 {
554 List<label> syncedFaceI(faceIDs);
555 for (label& faceI : syncedFaceI)
556 {
557 faceI += procPatch.start();
558 }
559 syncedFaces.append(syncedFaceI);
560 }
561
562 if (debug)
563 {
564 Pout<< "Received at time = " << mesh_.time().value()
565 << ": surfCellFacesOnProcPatch = " << faceIDs << nl
566 << "Received at time = " << mesh_.time().value()
567 << ": dVfPatch = " << nbrdVfs << endl;
568 }
569
570 // Combine fluxes
571 scalarField& localFlux = dVf.boundaryFieldRef()[patchi];
572
573 forAll(faceIDs, i)
574 {
575 const label facei = faceIDs[i];
576 localFlux[facei] = - nbrdVfs[i];
577 if (debug && mag(localFlux[facei] + nbrdVfs[i]) > ROOTVSMALL)
578 {
579 Pout<< "localFlux[facei] = " << localFlux[facei]
580 << " and nbrdVfs[i] = " << nbrdVfs[i]
581 << " for facei = " << facei << endl;
582 }
583 }
584 }
585
586 if (debug)
587 {
588 // Write out results for checking
589 forAll(procPatchLabels_, patchLabeli)
590 {
591 const label patchi = procPatchLabels_[patchLabeli];
592 const scalarField& localFlux = dVf.boundaryField()[patchi];
593 Pout<< "time = " << mesh_.time().value() << ": localFlux = "
594 << localFlux << endl;
595 }
596 }
597
598 // Reinitialising list used for minimal parallel communication
599 forAll(surfaceCellFacesOnProcPatches_, patchi)
600 {
601 surfaceCellFacesOnProcPatches_[patchi].clear();
602 }
603 }
604
605 return syncedFaces;
606}
607
608
609void Foam::isoAdvection::checkIfOnProcPatch(const label facei)
610{
611 if (!mesh_.isInternalFace(facei))
612 {
613 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
614 const label patchi = pbm.patchID(facei);
615
616 if (isA<processorPolyPatch>(pbm[patchi]) && !pbm[patchi].empty())
617 {
618 const label patchFacei = pbm[patchi].whichFace(facei);
619 surfaceCellFacesOnProcPatches_[patchi].append(patchFacei);
620 }
621 }
622}
623
624
625void Foam::isoAdvection::applyBruteForceBounding()
626{
627 addProfilingInFunction(geometricVoF);
628 bool alpha1Changed = false;
629
630 const scalar snapAlphaTol = dict_.getOrDefault<scalar>("snapTol", 0);
631 if (snapAlphaTol > 0)
632 {
633 alpha1_ =
634 alpha1_
635 *pos0(alpha1_ - snapAlphaTol)
636 *neg0(alpha1_ - (1.0 - snapAlphaTol))
637 + pos0(alpha1_ - (1.0 - snapAlphaTol));
638
639 alpha1Changed = true;
640 }
641
642 if (dict_.getOrDefault("clip", true))
643 {
644 alpha1_.clamp_range(zero_one{});
645 alpha1Changed = true;
646 }
647
648 if (alpha1Changed)
649 {
650 alpha1_.correctBoundaryConditions();
651 }
652}
653
654
656{
657 if (!mesh_.time().writeTime()) return;
658
659 if (dict_.getOrDefault("writeSurfCells", false))
660 {
661 cellSet cSet
662 (
663 IOobject
664 (
665 "surfCells",
666 mesh_.time().timeName(),
667 mesh_,
669 )
670 );
671
672 cSet.insert(surfCells_);
673
674 cSet.write();
675 }
676}
677
678
680(
681 const UList<List<point>>& faces
682) const
683{
684 if (!writeIsoFacesToFile_ || !mesh_.time().writeTime()) return;
685
686 // Writing isofaces to obj file for inspection, e.g. in paraview
687 const fileName outputFile
688 (
689 mesh_.time().globalPath()
690 / "isoFaces"
691 / word::printf("isoFaces_%012d.obj", mesh_.time().timeIndex())
692 );
693
694 if (Pstream::parRun())
695 {
696 // Collect points from all the processors
697 List<List<List<point>>> allProcFaces(Pstream::nProcs());
698 allProcFaces[Pstream::myProcNo()] = faces;
699 Pstream::gatherList(allProcFaces);
700
701 if (Pstream::master())
702 {
703 mkDir(outputFile.path());
704 OBJstream os(outputFile);
705 Info<< nl << "isoAdvection: writing iso faces to file: "
706 << os.name() << nl << endl;
707
708 for (const auto& procFacePts : allProcFaces)
709 {
710 for (const auto& facePts : procFacePts)
711 {
712 os.writeFace(facePts, false);
713 }
714 }
715 }
716 }
717 else
718 {
719 mkDir(outputFile.path());
720 OBJstream os(outputFile);
721 Info<< nl << "isoAdvection: writing iso faces to file: "
722 << os.name() << nl << endl;
723
724 for (const auto& facePts : faces)
725 {
726 os.writeFace(facePts, false);
727 }
728 }
729}
730
731
732// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
const volScalarField & alpha1
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
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
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
Ostream & writeFace(const UList< point > &points, const bool lines=true)
Write face loop points with lines/filled-polygon.
Definition OBJstream.C:279
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void finishedNeighbourSends(const labelUList &neighProcs, const bool wait=true)
Mark the send phase as being finished, with communication being limited to a known subset of send/rec...
void append(autoPtr< T > &ptr)
Move append an element to the end of the list.
Definition PtrList.H:364
const word & constant() const noexcept
Return constant name.
Definition TimePathsI.H:131
Input inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UIPstream.H:313
A List with indirect addressing. Like IndirectList but does not store addressing.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
Output inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UOPstream.H:408
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
@ gatherList
gatherList [manual algorithm]
Definition UPstream.H:194
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A collection of cell labels.
Definition cellSet.H:50
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
static int debug
Definition cutFace.H:143
Cyclic plane patch.
A class for handling file names.
Definition fileName.H:75
static std::string path(const std::string &str)
Return directory path name (part before last /).
Definition fileNameI.H:169
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
An implementation of the isoAdvector geometric Volume-of-Fluid method advancing the provided volume f...
void writeIsoFaces(const UList< List< point > > &isoFacePts) const
Write isoface points to .obj file.
void writeSurfaceCells() const
Return cellSet of surface cells.
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const labelList & patchID() const
Per boundary face label the patch index.
void clear()
Clear the patch list and all demand-driven data.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const vectorField & faceCentres() const
const scalarField & cellVolumes() const
const vectorField & cellCentres() const
const labelListList & cellCells() const
const labelListList & cellPoints() const
const vectorField & faceAreas() const
const cellList & cells() const
Neighbour processor patch.
Original code supplied by Henning Scheufler, DLR (2019).
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const bool parRun=UPstream::parRun())
Synchronize values on all mesh faces.
Definition syncTools.H:465
static word printf(const char *fmt, const PrimitiveType &val)
Use a printf-style formatter for a primitive.
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition pTraits.H:51
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
auto limits
Definition setRDeltaT.H:186
const polyBoundaryMesh & patches
faceListList boundary
dynamicFvMesh & mesh
IOdictionary porosityProperties(IOobject("porosityProperties", runTime.constant(), runTime, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
Calculate the gradient of the given field.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
const cellShapeList & cells
word timeName
Definition getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
phib
Definition pEqn.H:190
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedScalar pos0(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition POSIX.C:616
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
vector point
Point is a vector.
Definition point.H:37
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar neg0(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define addProfilingInFunction(Name)
Define profiling trigger with specified name and description corresponding to the compiler-defined fu...
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
conserve primitiveFieldRef()+