Loading...
Searching...
No Matches
domainDecompositionDryRunWrite.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) 2021-2023 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "volFields.H"
31
32// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33
34void Foam::domainDecompositionDryRun::writeVolField
35(
36 const word& timeName,
37 const labelUList& procIds
38) const
39{
40 // Write decomposition as volScalarField for visualization
41 volScalarField cellDist
42 (
43 IOobject
44 (
45 "cellDist",
47 this->mesh(),
51 ),
52 this->mesh(),
55 );
56
57 forAll(procIds, celli)
58 {
59 cellDist[celli] = procIds[celli];
60 }
61
62 cellDist.correctBoundaryConditions();
63 cellDist.write();
64
65 Info<< nl << "Wrote decomposition to "
66 << cellDist.objectRelPath()
67 << " (volScalarField) for visualization."
68 << endl;
69}
70
71
72void Foam::domainDecompositionDryRun::writeVTK
73(
74 const fileName& file,
75 const labelUList& cellToProc
76) const
77{
78 const vtk::vtuCells cells(this->mesh());
79
80 // not parallel
81 vtk::internalMeshWriter writer(this->mesh(), cells, file, false);
82
85 writer.writeCellData("procID", cellToProc);
86 writer.writeCellIDs();
87
88 Info<< "Wrote decomposition to "
89 << this->mesh().time().relativePath(writer.output())
90 << " for visualization."
91 << endl;
92}
93
94
95// ************************************************************************* //
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edgesCentres")))
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
fileName relativePath(const fileName &input, const bool caseTag=false) const
Return the input relative to the globalPath by stripping off a leading value of the globalPath.
Definition TimePathsI.H:122
const fvMesh & mesh() const noexcept
The mesh.
A class for handling file names.
Definition fileName.H:75
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
const fileName & output() const noexcept
The current output file name.
Write an OpenFOAM volume (internal) geometry and internal fields as a vtu file or a legacy vtk file.
virtual bool writeGeometry()
Write patch topology.
virtual bool beginCellData(label nFields=0)
Begin CellData output section for specified number of fields.
void writeCellData(const word &fieldName, const UList< Type > &field)
Write primitive field of CellData.
A deep-copy description of an OpenFOAM volume mesh in data structures suitable for VTK UnstructuredGr...
static const word null
An empty word.
Definition word.H:84
dynamicFvMesh & mesh
const cellShapeList & cells
word timeName
Definition getTimeIndex.H:3
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299