Skip to content

Commit

Permalink
MueLu: fixing CoordinateTransferFactory_kokkos for structured aggrega…
Browse files Browse the repository at this point in the history
…tion

This adds logic in the CoordinateTransferFactory_kokkos to allow user to perform structured aggregation and hand-off the aggregates to tentative P factory.
Also added xmls and tribits test to check that it all works fine.
  • Loading branch information
lucbv committed Feb 20, 2019
1 parent 761131d commit f793a80
Show file tree
Hide file tree
Showing 5 changed files with 440 additions and 10 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -70,6 +70,9 @@ namespace MueLu {
validParamList->set< RCP<const FactoryBase> >("CoarseMap", Teuchos::null, "Generating factory of the coarse map");
validParamList->set< int > ("write start", -1, "First level at which coordinates should be written to file");
validParamList->set< int > ("write end", -1, "Last level at which coordinates should be written to file");
validParamList->set< bool > ("structured aggregation", false, "Flag specifying that the geometric data is transferred for StructuredAggregationFactory");
validParamList->set<RCP<const FactoryBase> > ("lCoarseNodesPerDim", Teuchos::null, "Factory providing the local number of nodes per spatial dimensions of the mesh");
validParamList->set<RCP<const FactoryBase> > ("numDimensions", Teuchos::null, "Factory providing the number of spatial dimensions of the mesh");

return validParamList;
}
Expand All @@ -78,13 +81,19 @@ namespace MueLu {
void CoordinatesTransferFactory_kokkos<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType>>::DeclareInput(Level& fineLevel, Level& coarseLevel) const {
static bool isAvailableCoords = false;

if (coarseLevel.GetRequestMode() == Level::REQUEST)
isAvailableCoords = coarseLevel.IsAvailable("Coordinates", this);
const ParameterList& pL = GetParameterList();
if(pL.get<bool>("structured aggregation") == true) {
Input(fineLevel, "lCoarseNodesPerDim");
Input(fineLevel, "numDimensions");
} else {
if (coarseLevel.GetRequestMode() == Level::REQUEST)
isAvailableCoords = coarseLevel.IsAvailable("Coordinates", this);

if (isAvailableCoords == false) {
Input(fineLevel, "Coordinates");
Input(fineLevel, "Aggregates");
Input(fineLevel, "CoarseMap");
if (isAvailableCoords == false) {
Input(fineLevel, "Coordinates");
Input(fineLevel, "Aggregates");
Input(fineLevel, "CoarseMap");
}
}
}

Expand All @@ -101,6 +110,17 @@ namespace MueLu {
return;
}


const ParameterList& pL = GetParameterList();
if(pL.get<bool>("structured aggregation") == true) {
Array<LO> lCoarseNodesPerDir = Get<Array<LO> >(fineLevel, "lCoarseNodesPerDim");
Set<Array<LO> >(coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
int numDimensions = Get<int>(fineLevel, "numDimensions");
Set<int>(coarseLevel, "numDimensions", numDimensions);

return;
}

auto aggregates = Get<RCP<Aggregates_kokkos> >(fineLevel, "Aggregates");
auto fineCoords = Get<RCP<doubleMultiVector> >(fineLevel, "Coordinates");
auto coarseMap = Get<RCP<const Map> > (fineLevel, "CoarseMap");
Expand Down Expand Up @@ -184,7 +204,6 @@ namespace MueLu {

Set<RCP<doubleMultiVector> >(coarseLevel, "Coordinates", coarseCoords);

const ParameterList& pL = GetParameterList();
int writeStart = pL.get<int>("write start"), writeEnd = pL.get<int>("write end");
if (writeStart == 0 && fineLevel.GetLevelID() == 0 && writeStart <= writeEnd) {
std::string fileName = "coordinates_before_rebalance_level_" + toString(fineLevel.GetLevelID()) + ".m";
Expand Down
22 changes: 19 additions & 3 deletions packages/muelu/test/structured/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -40,8 +40,24 @@ IF (${PACKAGE_NAME}_ENABLE_Belos AND ${PACKAGE_NAME}_ENABLE_Amesos2)
IF (${PACKAGE_NAME}_ENABLE_Kokkos_Refactor)
TRIBITS_ADD_TEST(
Structured
NAME "Structured_Laplace2D_kokkos"
ARGS "--linAlgebra=Tpetra --xml=structured_1dof_kokkos.xml --matrixType=Laplace2D --nx=25 --ny=25 --nz=25"
NAME "Structured_Interp_Laplace2D_kokkos"
ARGS "--linAlgebra=Tpetra --xml=structured_interp_kokkos.xml --matrixType=Laplace2D --nx=25 --ny=25"
COMM serial mpi
NUM_MPI_PROCS 4
)

TRIBITS_ADD_TEST(
Structured
NAME "Structured_Interp_SA_Laplace2D_kokkos"
ARGS "--linAlgebra=Tpetra --xml=structured_interp_sa_kokkos.xml --matrixType=Laplace2D --nx=25 --ny=25"
COMM serial mpi
NUM_MPI_PROCS 4
)

TRIBITS_ADD_TEST(
Structured
NAME "Structured_Tentative_Laplace2D_kokkos"
ARGS "--linAlgebra=Tpetra --xml=structured_tentative_kokkos.xml --matrixType=Laplace2D --nx=25 --ny=25"
COMM serial mpi
NUM_MPI_PROCS 4
)
Expand Down Expand Up @@ -69,7 +85,7 @@ IF (${PACKAGE_NAME}_ENABLE_Belos AND ${PACKAGE_NAME}_ENABLE_Amesos2)


TRIBITS_COPY_FILES_TO_BINARY_DIR(Structured_cp
SOURCE_FILES structured_1dof.xml structured_2dof.xml structured_3dof.xml structured_1dof_kokkos.xml
SOURCE_FILES structured_1dof.xml structured_2dof.xml structured_3dof.xml structured_interp_kokkos.xml structured_interp_sa_kokkos.xml structured_tentative_kokkos.xml
)


Expand Down
130 changes: 130 additions & 0 deletions packages/muelu/test/structured/structured_interp_kokkos.xml
Original file line number Diff line number Diff line change
@@ -0,0 +1,130 @@
<ParameterList name="MueLu">

<!-- Configuration of the Xpetra operator (fine level) -->
<ParameterList name="Matrix">
<Parameter name="PDE equations" type="int" value="1"/> <!-- Number of PDE equations at each grid node.-->
</ParameterList>

<!-- Factory collection -->
<ParameterList name="Factories">

<ParameterList name="myCoalesceDropFact">
<Parameter name="factory" type="string" value="CoalesceDropFactory_kokkos"/>
<Parameter name="lightweight wrap" type="bool" value="true"/>
<Parameter name="aggregation: drop tol" type="double" value="0.00"/>
</ParameterList>

<ParameterList name="myAggregationFact">
<Parameter name="factory" type="string" value="StructuredAggregationFactory_kokkos"/>
<Parameter name="aggregation: coupling" type="string" value="uncoupled"/>
<Parameter name="aggregation: output type" type="string" value="CrsGraph"/>
<Parameter name="aggregation: coarsening order" type="int" value="0"/>
<Parameter name="aggregation: coarsening rate" type="string" value="{2}"/>
<Parameter name="Graph" type="string" value="myCoalesceDropFact"/>
</ParameterList>


<ParameterList name="myCoarseMapFact">
<Parameter name="factory" type="string" value="CoarseMapFactory_kokkos"/>
<Parameter name="Aggregates" type="string" value="myAggregationFact"/>
</ParameterList>

<!-- Note that ParameterLists must be defined prior to being used -->
<ParameterList name="myProlongatorFact">
<Parameter name="factory" type="string" value="GeometricInterpolationPFactory"/>
<Parameter name="interp: build coarse coordinates" type="bool" value="true"/>
<Parameter name="interp: interpolation order" type="int" value="0"/>
<Parameter name="prolongatorGraph" type="string" value="myAggregationFact"/>
<Parameter name="coarseCoordinatesFineMap" type="string" value="myAggregationFact"/>
<Parameter name="coarseCoordinatesMap" type="string" value="myAggregationFact"/>
</ParameterList>

<ParameterList name="myCoordTransferFact">
<Parameter name="factory" type="string" value="CoordinatesTransferFactory_kokkos"/>
<Parameter name="structured aggregation" type="bool" value="true"/>
<Parameter name="numDimensions" type="string" value="myAggregationFact"/>
<Parameter name="lCoarseNodesPerDim" type="string" value="myAggregationFact"/>
</ParameterList>

<ParameterList name="myNullspaceFact">
<Parameter name="factory" type="string" value="NullspaceFactory_kokkos"/>
<Parameter name="Nullspace" type="string" value="myProlongatorFact"/>
<Parameter name="CoarseMap" type="string" value="myCoarseMapFact"/>
</ParameterList>

<ParameterList name="myRestrictorFact">
<Parameter name="factory" type="string" value="TransPFactory"/>
</ParameterList>

<!-- <ParameterList name="myAggExport"> -->
<!-- <Parameter name="factory" type="string" value="AggregationExportFactory"/> -->
<!-- <Parameter name="Aggregates" type="string" value="myAggregationFact"/> -->
<!-- <Parameter name="aggregation: output filename" type="string" value="structured_aggs"/> -->
<!-- <Parameter name="aggregation: output file: agg style" type="string" value="Jacks"/> -->
<!-- <Parameter name="aggregation: output file: agg style" type="string" value="Convex Hulls"/> -->
<!-- </ParameterList> -->

<ParameterList name="myRAPFact">
<Parameter name="factory" type="string" value="RAPFactory"/>
<Parameter name="P" type="string" value="myProlongatorFact"/>
<Parameter name="R" type="string" value="myRestrictorFact"/>
<ParameterList name="TransferFactories">
<Parameter name="CoordinateTransfer" type="string" value="myCoordTransferFact"/>
<!-- <Parameter name="AggregationExportFactory" type="string" value="myAggExport"/> -->
</ParameterList>
</ParameterList>

<ParameterList name="myILU">
<Parameter name="factory" type="string" value="TrilinosSmoother"/>
<Parameter name="type" type="string" value="RILUK"/>
<ParameterList name="ParameterList">
<Parameter name="schwarz: overlap level" type="int" value="1"/>
<Parameter name="schwarz: combine mode" type="string" value="Zero"/>
<Parameter name="schwarz: use reordering" type="bool" value="false"/>
<Parameter name="fact: iluk level-of-fill" type="int" value="0"/>
<Parameter name="fact: absolute threshold" type="double" value="0."/>
<Parameter name="fact: relative threshold" type="double" value="1."/>
<Parameter name="fact: relax value" type="double" value="0."/>
</ParameterList>
</ParameterList>

<ParameterList name="myMTSGS">
<Parameter name="factory" type="string" value="TrilinosSmoother"/>
<Parameter name="type" type="string" value="RELAXATION"/>
<ParameterList name="ParameterList">
<Parameter name="relaxation: type" type="string" value="MT Symmetric Gauss-Seidel"/>
<Parameter name="relaxation: symmetric matrix structure" type="bool" value="true"/>
<Parameter name="relaxation: sweeps" type="int" value="2"/>
<Parameter name="relaxation: use l1" type="bool" value="true"/>
</ParameterList>
</ParameterList>

</ParameterList>


<!-- Definition of the multigrid preconditioner -->
<ParameterList name="Hierarchy">

<Parameter name="max levels" type="int" value="3"/> <!-- Max number of levels -->
<Parameter name="cycle type" type="string" value="V"/>
<Parameter name="coarse: max size" type="int" value="20"/> <!-- Min number of rows on coarsest level -->
<Parameter name="verbosity" type="string" value="High"/>

<ParameterList name="All">
<Parameter name="PreSmoother" type="string" value="myMTSGS"/>
<Parameter name="PostSmoother" type="string" value="myMTSGS"/>
<Parameter name="Graph" type="string" value="myCoalesceDropFact"/>
<Parameter name="Nullspace" type="string" value="myNullspaceFact"/>
<Parameter name="Aggregates" type="string" value="myAggregationFact"/>
<Parameter name="lCoarseNodesPerDim" type="string" value="myAggregationFact"/>
<Parameter name="P" type="string" value="myProlongatorFact"/>
<Parameter name="R" type="string" value="myRestrictorFact"/>
<Parameter name="A" type="string" value="myRAPFact"/>
<Parameter name="CoarseSolver" type="string" value="DirectSolver"/>
<Parameter name="Coordinates" type="string" value="myProlongatorFact"/>
<Parameter name="lNodesPerDim" type="string" value="myCoordTransferFact"/>
<Parameter name="numDimensions" type="string" value="myCoordTransferFact"/>
</ParameterList>
</ParameterList>

</ParameterList>
Loading

0 comments on commit f793a80

Please sign in to comment.