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 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"
32
33// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34
35void Foam::domainDecompositionDryRun::writeVolField
36(
37 const word& timeName,
38 const labelUList& procIds
39) const
40{
41 // Write decomposition as volScalarField for visualization
42 volScalarField cellDist
43 (
44 IOobject
45 (
46 "cellDist",
48 this->mesh(),
51 false
52 ),
53 this->mesh(),
54 dimensionedScalar("cellDist", dimless, -1),
55 zeroGradientFvPatchScalarField::typeName
56 );
57
58 forAll(procIds, celli)
59 {
60 cellDist[celli] = procIds[celli];
61 }
62
63 cellDist.correctBoundaryConditions();
64 cellDist.write();
65
66 Info<< nl << "Wrote decomposition to "
67 << cellDist.objectRelPath()
68 << " (volScalarField) for visualization."
69 << endl;
70}
71
72
74(
75 const fileName& file,
76 const labelUList& cellToProc
77) const
78{
79 const vtk::vtuCells cells(this->mesh());
80
81 // not parallel
82 vtk::internalMeshWriter writer(this->mesh(), cells, file, false);
83
86 writer.writeCellData("procID", cellToProc);
87
88 Info<< "Wrote decomposition to "
89 << this->mesh().time().relativePath(writer.output())
90 << " for visualization."
91 << endl;
92}
93
94
95// ************************************************************************* //
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
fileName relativePath(const fileName &input, const bool caseTag=false) const
Definition: TimePathsI.H:87
const fvMesh & mesh() const noexcept
The mesh.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
void writeVTK() const
Write VTK freeSurface mesh.
const fileName & output() const noexcept
The current output file name.
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.
virtual bool writeGeometry()
Write mesh topology.
dynamicFvMesh & mesh
const cellShapeList & cells
word timeName
Definition: getTimeIndex.H:3
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333