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 -------------------------------------------------------------------------------
10 License
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 
35 void 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",
47  timeName,
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 
73 void Foam::domainDecompositionDryRun::writeVTK
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 // ************************************************************************* //
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
volFields.H
Foam::vtk::fileWriter::output
const fileName & output() const noexcept
Definition: foamVtkFileWriterI.H:85
Foam::domainDecompositionDryRun::mesh
const fvMesh & mesh() const noexcept
Definition: domainDecompositionDryRun.H:104
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::TimePaths::relativePath
fileName relativePath(const fileName &input, const bool caseTag=false) const
Definition: TimePathsI.H:80
domainDecompositionDryRun.H
Foam::Info
messageStream Info
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:36
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::vtk::internalMeshWriter::writeCellData
void writeCellData(const word &fieldName, const UList< Type > &field)
Definition: foamVtkInternalMeshWriterTemplates.C:65
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::vtk::internalMeshWriter::beginCellData
virtual bool beginCellData(label nFields=0)
Definition: foamVtkInternalMeshWriter.C:589
foamVtkInternalMeshWriter.H
Foam::fvMesh::time
const Time & time() const
Definition: fvMesh.H:276
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:81
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:184
writer
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
zeroGradientFvPatchFields.H
Foam::vtk::internalMeshWriter::writeGeometry
virtual bool writeGeometry()
Definition: foamVtkInternalMeshWriter.C:552
Foam::dimless
const dimensionSet dimless
Definition: dimensionSets.C:182