openfoam/applications/utilities/parallelProcessing/decomposePar/domainDecompositionDryRunWrite.C
Mark Olesen 61aaacd088 ENH: adjust renumbering methods, extend renumberMesh options
- renumberMesh now has -dry-run, -write-maps, -no-fields,
  -renumber-method, -renumber-coeffs options.

  * Use -dry-run with -write-maps to visualize the before/after
    effects of renumbering (creates a VTK file).

  * -no-fields to renumber the mesh only.
    This is useful and faster when the input fields are uniform
    and the -overwrite option is specified.

  * -renumber-method allows a quick means of specifying a different
    default renumber method (instead of Cuthill-McKee).

    The -renumber-coeffs option allows passing of dictionary content
    for the method.

    Examples,

       // Different ways to specify reverse Cuthill-McKee

       *  -renumber-method RCM
       *  -renumber-coeffs 'reverse true;'
       *  -renumber-method CuthillMcKee
       *  -renumber-coeffs 'reverse true;'
       *  -renumber-coeffs 'method CuthillMcKee; reverse true;'

       // Other (without dictionary coefficients)
       *  renumberMesh -renumber-method random

       // Other (with dictionary coefficients)
       renumberMesh \
           -renumber-method spring \
           -renumber-coeffs 'maxCo 0.1; maxIter 1000; freezeFraction 0.99;'

       // Other (with additional libraries)
       renumberMesh -renumber-method zoltan -lib zoltanRenumber

COMP: build zoltan renumbering to MPI-specific location

- zoltan and Sloan renumbering are now longer automatically linked to
  the renumberMesh utility but must be separately loaded by a
  command-line option or through a dictionary "libs" entry.

ENH: add output cellID for decomposePar -dry-run -cellDist
2024-03-06 17:58:47 +01:00

96 lines
2.8 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2021-2023 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "domainDecompositionDryRun.H"
#include "foamVtkInternalMeshWriter.H"
#include "volFields.H"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::domainDecompositionDryRun::writeVolField
(
const word& timeName,
const labelUList& procIds
) const
{
// Write decomposition as volScalarField for visualization
volScalarField cellDist
(
IOobject
(
"cellDist",
timeName,
this->mesh(),
IOobject::NO_READ,
IOobject::NO_WRITE,
IOobject::NO_REGISTER
),
this->mesh(),
dimensionedScalar(word::null, dimless, -1),
fvPatchFieldBase::zeroGradientType()
);
forAll(procIds, celli)
{
cellDist[celli] = procIds[celli];
}
cellDist.correctBoundaryConditions();
cellDist.write();
Info<< nl << "Wrote decomposition to "
<< cellDist.objectRelPath()
<< " (volScalarField) for visualization."
<< endl;
}
void Foam::domainDecompositionDryRun::writeVTK
(
const fileName& file,
const labelUList& cellToProc
) const
{
const vtk::vtuCells cells(this->mesh());
// not parallel
vtk::internalMeshWriter writer(this->mesh(), cells, file, false);
writer.writeGeometry();
writer.beginCellData();
writer.writeCellData("procID", cellToProc);
writer.writeCellIDs();
Info<< "Wrote decomposition to "
<< this->mesh().time().relativePath(writer.output())
<< " for visualization."
<< endl;
}
// ************************************************************************* //