scotch decomposition method

This commit is contained in:
mattijs
2009-04-23 14:44:23 +01:00
parent 6c1fe9b4b6
commit 2c48a69e50
6 changed files with 586 additions and 5 deletions

View File

@ -2,7 +2,10 @@ decompositionMethod/decompositionMethod.C
geomDecomp/geomDecomp.C geomDecomp/geomDecomp.C
simpleGeomDecomp/simpleGeomDecomp.C simpleGeomDecomp/simpleGeomDecomp.C
hierarchGeomDecomp/hierarchGeomDecomp.C hierarchGeomDecomp/hierarchGeomDecomp.C
metisDecomp/metisDecomp.C
manualDecomp/manualDecomp.C manualDecomp/manualDecomp.C
scotchDecomp/scotchDecomp.C
metisDecomp/metisDecomp.C
LIB = $(FOAM_LIBBIN)/libdecompositionMethods LIB = $(FOAM_LIBBIN)/libdecompositionMethods

View File

@ -1,6 +1,8 @@
EXE_INC = \ EXE_INC = \
-I$(WM_THIRD_PARTY_DIR)/scotch_5.1/src/libscotch/lnInclude \
-I$(WM_THIRD_PARTY_DIR)/metis-5.0pre2/include -I$(WM_THIRD_PARTY_DIR)/metis-5.0pre2/include
LIB_LIBS = \ LIB_LIBS = \
-lscotch \
-lmetis \ -lmetis \
-lGKlib -lGKlib

View File

@ -0,0 +1,425 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
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 2 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, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
S Strat=b{job=t,map=t,poli=S,sep=(m{asc=b{bnd=d{pass=40,dif=1,rem=1}f{move=80,pass=-1,bal=0.005},org=f{move=80,pass=-1,bal=0.005},width=3},low=h{pass=10}f{move=80,pass=-1,bal=0.0005},type=h,vert=80,rat=0.8}|m{asc=b{bnd=d{pass=40,dif=1,rem=1}f{move=80,pass=-1,bal=0.005},org=f{move=80,pass=-1,bal=0.005},width=3},low=h{pass=10}f{move=80,pass=-1,bal=0.0005},type=h,vert=80,rat=0.8})}
\*---------------------------------------------------------------------------*/
#include "scotchDecomp.H"
#include "addToRunTimeSelectionTable.H"
#include "floatScalar.H"
#include "IFstream.H"
#include "Time.H"
#include "cyclicPolyPatch.H"
extern "C"
{
#define OMPI_SKIP_MPICXX
#include "module.h"
#include "common.h"
#include "scotch.h"
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(scotchDecomp, 0);
addToRunTimeSelectionTable
(
decompositionMethod,
scotchDecomp,
dictionaryMesh
);
}
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::scotchDecomp::check(const int retVal, const char* str)
{
if (retVal)
{
FatalErrorIn("scotchDecomp::decompose(..)")
<< "Called to scotch routine " << str << " failed."
<< exit(FatalError);
}
}
// Call Metis with options from dictionary.
Foam::label Foam::scotchDecomp::decompose
(
const List<int>& adjncy,
const List<int>& xadj,
List<int>& finalDecomp
)
{
// Strategy
// ~~~~~~~~
// Default.
SCOTCH_Strat stradat;
check(SCOTCH_stratInit(&stradat), "SCOTCH_stratInit");
//SCOTCH_stratGraphMap(&stradat, &argv[i][2]);
//fprintf(stdout, "S\tStrat=");
//SCOTCH_stratSave(&stradat, stdout);
//fprintf(stdout, "\n");
// Graph
// ~~~~~
SCOTCH_Graph grafdat;
check(SCOTCH_graphInit(&grafdat), "SCOTCH_graphInit");
check
(
SCOTCH_graphBuild
(
&grafdat,
0, // baseval, c-style numbering
xadj.size()-1, // vertnbr, nCells
xadj.begin(), // verttab, start index per cell into adjncy
&xadj[1], // vendtab, end index ,,
NULL, // velotab, vertex weights
NULL, // vlbltab
adjncy.size(), // edgenbr, number of arcs
adjncy.begin(), // edgetab
NULL // edlotab, edge weights
),
"SCOTCH_graphBuild"
);
check(SCOTCH_graphCheck(&grafdat), "SCOTCH_graphCheck");
//// Architecture
//// ~~~~~~~~~~~~
//// (fully connected network topology since using switch)
//
//SCOTCH_Arch archdat;
//check(SCOTCH_archInit(&archdat), "SCOTCH_archInit");
//check
//(
// // SCOTCH_archCmpltw for weighted.
// SCOTCH_archCmplt(&archdat, nProcessors_),
// "SCOTCH_archCmplt"
//);
//SCOTCH_Mapping mapdat;
//SCOTCH_graphMapInit(&grafdat, &mapdat, &archdat, NULL);
//SCOTCH_graphMapCompute(&grafdat, &mapdat, &stradat); /* Perform mapping */
finalDecomp.setSize(xadj.size()-1);
check
(
SCOTCH_graphPart
(
&grafdat,
nProcessors_, // partnbr
&stradat, // const SCOTCH_Strat *
finalDecomp.begin() // parttab
),
"SCOTCH_graphPart"
);
// Release storage for graph
SCOTCH_graphExit(&grafdat);
SCOTCH_stratExit(&stradat);
//// Release storage for network topology
//SCOTCH_archExit(&archdat);
return 0;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::scotchDecomp::scotchDecomp
(
const dictionary& decompositionDict,
const polyMesh& mesh
)
:
decompositionMethod(decompositionDict),
mesh_(mesh)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::labelList Foam::scotchDecomp::decompose(const pointField& points)
{
if (points.size() != mesh_.nCells())
{
FatalErrorIn("scotchDecomp::decompose(const pointField&)")
<< "Can use this decomposition method only for the whole mesh"
<< endl
<< "and supply one coordinate (cellCentre) for every cell." << endl
<< "The number of coordinates " << points.size() << endl
<< "The number of cells in the mesh " << mesh_.nCells()
<< exit(FatalError);
}
// Make Metis CSR (Compressed Storage Format) storage
// adjncy : contains neighbours (= edges in graph)
// xadj(celli) : start of information in adjncy for celli
List<int> xadj(mesh_.nCells()+1);
// Initialise the number of internal faces of the cells to twice the
// number of internal faces
label nInternalFaces = 2*mesh_.nInternalFaces();
// Check the boundary for coupled patches and add to the number of
// internal faces
const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
forAll(pbm, patchi)
{
if (isA<cyclicPolyPatch>(pbm[patchi]))
{
nInternalFaces += pbm[patchi].size();
}
}
// Create the adjncy array the size of the total number of internal and
// coupled faces
List<int> adjncy(nInternalFaces);
// Fill in xadj
// ~~~~~~~~~~~~
label freeAdj = 0;
for (label cellI = 0; cellI < mesh_.nCells(); cellI++)
{
xadj[cellI] = freeAdj;
const labelList& cFaces = mesh_.cells()[cellI];
forAll(cFaces, i)
{
label faceI = cFaces[i];
if
(
mesh_.isInternalFace(faceI)
|| isA<cyclicPolyPatch>(pbm[pbm.whichPatch(faceI)])
)
{
freeAdj++;
}
}
}
xadj[mesh_.nCells()] = freeAdj;
// Fill in adjncy
// ~~~~~~~~~~~~~~
labelList nFacesPerCell(mesh_.nCells(), 0);
// Internal faces
for (label faceI = 0; faceI < mesh_.nInternalFaces(); faceI++)
{
label own = mesh_.faceOwner()[faceI];
label nei = mesh_.faceNeighbour()[faceI];
adjncy[xadj[own] + nFacesPerCell[own]++] = nei;
adjncy[xadj[nei] + nFacesPerCell[nei]++] = own;
}
// Coupled faces. Only cyclics done.
forAll(pbm, patchi)
{
if (isA<cyclicPolyPatch>(pbm[patchi]))
{
const unallocLabelList& faceCells = pbm[patchi].faceCells();
label sizeby2 = faceCells.size()/2;
for (label facei=0; facei<sizeby2; facei++)
{
label own = faceCells[facei];
label nei = faceCells[facei + sizeby2];
adjncy[xadj[own] + nFacesPerCell[own]++] = nei;
adjncy[xadj[nei] + nFacesPerCell[nei]++] = own;
}
}
}
// Decompose using default weights
List<int> finalDecomp;
decompose(adjncy, xadj, finalDecomp);
// Copy back to labelList
labelList decomp(finalDecomp.size());
forAll(decomp, i)
{
decomp[i] = finalDecomp[i];
}
return decomp;
}
// From cell-cell connections to Metis format (like CompactListList)
void Foam::scotchDecomp::calcMetisCSR
(
const labelListList& cellCells,
List<int>& adjncy,
List<int>& xadj
)
{
// Count number of internal faces
label nConnections = 0;
forAll(cellCells, coarseI)
{
nConnections += cellCells[coarseI].size();
}
// Create the adjncy array as twice the size of the total number of
// internal faces
adjncy.setSize(nConnections);
xadj.setSize(cellCells.size()+1);
// Fill in xadj
// ~~~~~~~~~~~~
label freeAdj = 0;
forAll(cellCells, coarseI)
{
xadj[coarseI] = freeAdj;
const labelList& cCells = cellCells[coarseI];
forAll(cCells, i)
{
adjncy[freeAdj++] = cCells[i];
}
}
xadj[cellCells.size()] = freeAdj;
}
Foam::labelList Foam::scotchDecomp::decompose
(
const labelList& agglom,
const pointField& agglomPoints
)
{
if (agglom.size() != mesh_.nCells())
{
FatalErrorIn
(
"parMetisDecomp::decompose(const labelList&, const pointField&)"
) << "Size of cell-to-coarse map " << agglom.size()
<< " differs from number of cells in mesh " << mesh_.nCells()
<< exit(FatalError);
}
// Make Metis CSR (Compressed Storage Format) storage
// adjncy : contains neighbours (= edges in graph)
// xadj(celli) : start of information in adjncy for celli
List<int> adjncy;
List<int> xadj;
{
// Get cellCells on coarse mesh.
labelListList cellCells;
calcCellCells
(
mesh_,
agglom,
agglomPoints.size(),
cellCells
);
calcMetisCSR(cellCells, adjncy, xadj);
}
// Decompose using default weights
List<int> finalDecomp;
decompose(adjncy, xadj, finalDecomp);
// Rework back into decomposition for original mesh_
labelList fineDistribution(agglom.size());
forAll(fineDistribution, i)
{
fineDistribution[i] = finalDecomp[agglom[i]];
}
return fineDistribution;
}
Foam::labelList Foam::scotchDecomp::decompose
(
const labelListList& globalCellCells,
const pointField& cellCentres
)
{
if (cellCentres.size() != globalCellCells.size())
{
FatalErrorIn
(
"scotchDecomp::decompose(const pointField&, const labelListList&)"
) << "Inconsistent number of cells (" << globalCellCells.size()
<< ") and number of cell centres (" << cellCentres.size()
<< ")." << exit(FatalError);
}
// Make Metis CSR (Compressed Storage Format) storage
// adjncy : contains neighbours (= edges in graph)
// xadj(celli) : start of information in adjncy for celli
List<int> adjncy;
List<int> xadj;
calcMetisCSR(globalCellCells, adjncy, xadj);
// Decompose using default weights
List<int> finalDecomp;
decompose(adjncy, xadj, finalDecomp);
// Copy back to labelList
labelList decomp(finalDecomp.size());
forAll(decomp, i)
{
decomp[i] = finalDecomp[i];
}
return decomp;
}
// ************************************************************************* //

View File

@ -0,0 +1,150 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
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 2 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, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Class
Foam::scotchDecomp
Description
Scotch domain decomposition
SourceFiles
scotchDecomp.C
\*---------------------------------------------------------------------------*/
#ifndef scotchDecomp_H
#define scotchDecomp_H
#include "decompositionMethod.H"
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class scotchDecomp Declaration
\*---------------------------------------------------------------------------*/
class scotchDecomp
:
public decompositionMethod
{
// Private data
const polyMesh& mesh_;
// Private Member Functions
//- Check and print error message
static void check(const int, const char*);
label decompose
(
const List<int>& adjncy,
const List<int>& xadj,
List<int>& finalDecomp
);
//- Disallow default bitwise copy construct and assignment
void operator=(const scotchDecomp&);
scotchDecomp(const scotchDecomp&);
public:
//- Runtime type information
TypeName("scotch");
// Constructors
//- Construct given the decomposition dictionary and mesh
scotchDecomp
(
const dictionary& decompositionDict,
const polyMesh& mesh
);
// Destructor
virtual ~scotchDecomp()
{}
// Member Functions
virtual bool parallelAware() const
{
// Metis does not know about proc boundaries
return false;
}
//- Return for every coordinate the wanted processor number. Use the
// mesh connectivity (if needed)
virtual labelList decompose(const pointField&);
//- Return for every coordinate the wanted processor number. Gets
// passed agglomeration map (from fine to coarse cells) and coarse cell
// location. Can be overridden by decomposers that provide this
// functionality natively.
virtual labelList decompose
(
const labelList& agglom,
const pointField&
);
//- Return for every coordinate the wanted processor number. Explicitly
// provided mesh connectivity.
// The connectivity is equal to mesh.cellCells() except for
// - in parallel the cell numbers are global cell numbers (starting
// from 0 at processor0 and then incrementing all through the
// processors)
// - the connections are across coupled patches
virtual labelList decompose
(
const labelListList& globalCellCells,
const pointField& cc
);
//- Helper to convert cellcells into compact row storage
static void calcMetisCSR
(
const labelListList& globalCellCells,
List<int>& adjncy,
List<int>& xadj
);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -1,8 +1,8 @@
.SUFFIXES: .y .Y .SUFFIXES: .y .Y
ytoo = bison -v $(YYPREFIX) -d $$SOURCE ; mv $$SOURCE_DIR/*.tab.c $*.c ; mv $$SOURCE_DIR/*.tab.h $*.h ; $(cc) $(cFLAGS) -c $*.c -o $@ ytoo = bison -v -d -y $$SOURCE ; mv y.tab.c $*.c ; mv y.tab.h $*.h ; $(cc) $(cFLAGS) -c $*.c -o $@
Ytoo = bison -v $(YYPREFIX) -d $$SOURCE ; mv $$SOURCE_DIR/*.tab.c $*.C ; mv $$SOURCE_DIR/*.tab.h $*.H ; $(CC) $(c++FLAGS) -c $*.C -o $@ Ytoo = bison -v -d -y $$SOURCE ; mv y.tab.c $*.C ; mv y.tab.h $*.H ; $(CC) $(c++FLAGS) -c $*.C -o $@
.y.dep: .y.dep:
$(MAKE_DEP) $(MAKE_DEP)

View File

@ -5,6 +5,7 @@ include $(GENERAL_RULES)/sourceToDep
include $(GENERAL_RULES)/java include $(GENERAL_RULES)/java
include $(GENERAL_RULES)/flex include $(GENERAL_RULES)/flex
include $(GENERAL_RULES)/flex++ include $(GENERAL_RULES)/flex++
include $(GENERAL_RULES)/byacc #include $(GENERAL_RULES)/byacc
include $(GENERAL_RULES)/btyacc++ #include $(GENERAL_RULES)/btyacc++
include $(GENERAL_RULES)/bison
include $(GENERAL_RULES)/moc include $(GENERAL_RULES)/moc