Skip to content

Commit

Permalink
Merge Pull Request #4057 from cgcgcg/Trilinos/tpetraTAFC
Browse files Browse the repository at this point in the history
Automatically Merged using Trilinos Pull Request AutoTester
PR Title: Tpetra: Fixes to TAFC
PR Author: cgcgcg
  • Loading branch information
trilinos-autotester authored Dec 16, 2018
2 parents 833114f + f6ad883 commit f94ba63
Show file tree
Hide file tree
Showing 5 changed files with 41 additions and 12 deletions.
15 changes: 10 additions & 5 deletions packages/muelu/src/Misc/MueLu_RAPFactory_def.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -134,9 +134,9 @@ namespace MueLu {
// Reuse pattern if available (multiple solve)
RCP<ParameterList> APparams;
if(pL.isSublist("matrixmatrix: kernel params"))
APparams=rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
APparams = rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
else
APparams= rcp(new ParameterList);
APparams = rcp(new ParameterList);

// By default, we don't need global constants for A*P
APparams->set("compute global constants: temporaries",APparams->get("compute global constants: temporaries",false));
Expand All @@ -160,8 +160,10 @@ namespace MueLu {

// Reuse coarse matrix memory if available (multiple solve)
RCP<ParameterList> RAPparams;
if(pL.isSublist("matrixmatrix: kernel params")) RAPparams=rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
else RAPparams= rcp(new ParameterList);
if(pL.isSublist("matrixmatrix: kernel params"))
RAPparams = rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
else
RAPparams = rcp(new ParameterList);



Expand Down Expand Up @@ -220,7 +222,10 @@ namespace MueLu {
Set(coarseLevel, "RAP reuse data", RAPparams);
} else {
RCP<ParameterList> RAPparams;
RAPparams= rcp(new ParameterList);
if(pL.isSublist("matrixmatrix: kernel params"))
RAPparams = rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
else
RAPparams = rcp(new ParameterList);

// We *always* need global constants for the RAP, but not for the temps
RAPparams->set("compute global constants: temporaries",RAPparams->get("compute global constants: temporaries",false));
Expand Down
15 changes: 13 additions & 2 deletions packages/tpetra/core/ext/TpetraExt_MatrixMatrix_def.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -374,6 +374,8 @@ void Jacobi(Scalar omega,
int mm_optimization_core_count=0;
auto slist = params->sublist("matrixmatrix: kernel params",false);
mm_optimization_core_count = slist.get("MM_TAFC_OptimizationCoreCount",::Tpetra::Details::Behavior::TAFC_OptimizationCoreCount ());
int foo = params->get("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
if(foo<mm_optimization_core_count) mm_optimization_core_count=foo;
bool isMM = slist.get("isMatrixMatrix_TransferAndFillComplete",false);
auto & ip1slist = importParams1->sublist("matrixmatrix: kernel params",false);
ip1slist.set("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
Expand All @@ -400,6 +402,8 @@ void Jacobi(Scalar omega,
int mm_optimization_core_count = slist.get("MM_TAFC_OptimizationCoreCount",::Tpetra::Details::Behavior::TAFC_OptimizationCoreCount () );
bool isMM = slist.get("isMatrixMatrix_TransferAndFillComplete",false);
auto & ip2slist = importParams2->sublist("matrixmatrix: kernel params",false);
int foo = params->get("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
if(foo<mm_optimization_core_count) mm_optimization_core_count=foo;
ip2slist.set("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
ip2slist.set("isMatrixMatrix_TransferAndFillComplete",isMM);
}
Expand Down Expand Up @@ -1504,6 +1508,8 @@ void mult_AT_B_newmatrix(
bool isMM = slist.get("isMatrixMatrix_TransferAndFillComplete",false);
int mm_optimization_core_count=3000; // ~3000 for serrano
mm_optimization_core_count = slist.get("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
int foo = params->get("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
if(foo<mm_optimization_core_count) mm_optimization_core_count=foo;
auto & sip1 = importParams1->sublist("matrixmatrix: kernel params",false);
sip1.set("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
sip1.set("isMatrixMatrix_TransferAndFillComplete",isMM);
Expand All @@ -1518,6 +1524,8 @@ void mult_AT_B_newmatrix(
bool isMM = slist.get("isMatrixMatrix_TransferAndFillComplete",false);
int mm_optimization_core_count=3000; // ~3000 for serrano
mm_optimization_core_count = slist.get("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
int foo = params->get("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
if(foo<mm_optimization_core_count) mm_optimization_core_count=foo;
auto & sip2 = importParams2->sublist("matrixmatrix: kernel params",false);
sip2.set("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
sip2.set("isMatrixMatrix_TransferAndFillComplete",isMM);
Expand Down Expand Up @@ -1563,8 +1571,11 @@ void mult_AT_B_newmatrix(
if(!params.is_null()) {
Teuchos::ParameterList& params_sublist = params->sublist("matrixmatrix: kernel params",false);
Teuchos::ParameterList& labelList_subList = labelList.sublist("matrixmatrix: kernel params",false);
int foundcount = params_sublist.get("MM_TAFC_OptimizationCoreCount",3000);
labelList_subList.set("MM_TAFC_OptimizationCoreCount",foundcount,"Core Count above which the optimized neighbor discovery is used");
int mm_optimization_core_count = 3000;
mm_optimization_core_count = params_sublist.get("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
int foo = params->get("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
if(foo<mm_optimization_core_count) mm_optimization_core_count=foo;
labelList_subList.set("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count,"Core Count above which the optimized neighbor discovery is used");

labelList_subList.set("isMatrixMatrix_TransferAndFillComplete",true,
"This parameter should be set to true only for MatrixMatrix operations: the optimization in Epetra that was ported to Tpetra does _not_ take into account the possibility that for any given source PID, a particular GID may not exist on the target PID: i.e. a transfer operation. A fix for this general case is in development.");
Expand Down
14 changes: 13 additions & 1 deletion packages/tpetra/core/ext/TpetraExt_TripleMatrixMultiply_def.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -298,7 +298,19 @@ namespace Tpetra {
Teuchos::ParameterList labelList;
labelList.set("Timer Label", label);
RCP<crs_matrix_type> Acprime = rcpFromRef(Ac);
if(!params.is_null()) labelList.set("compute global constants",params->get("compute global constants",true));
if(!params.is_null()) {
Teuchos::ParameterList& params_sublist = params->sublist("matrixmatrix: kernel params",false);
Teuchos::ParameterList& labelList_subList = labelList.sublist("matrixmatrix: kernel params",false);
int mm_optimization_core_count = 3000;
mm_optimization_core_count = params_sublist.get("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
int foo = params->get("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count);
if(foo<mm_optimization_core_count) mm_optimization_core_count=foo;
labelList_subList.set("MM_TAFC_OptimizationCoreCount",mm_optimization_core_count,"Core Count above which the optimized neighbor discovery is used");

labelList_subList.set("isMatrixMatrix_TransferAndFillComplete",true,
"This parameter should be set to true only for MatrixMatrix operations: the optimization in Epetra that was ported to Tpetra does _not_ take into account the possibility that for any given source PID, a particular GID may not exist on the target PID: i.e. a transfer operation. A fix for this general case is in development.");
labelList.set("compute global constants",params->get("compute global constants",true));
}
export_type exporter = export_type(*Pprime->getGraph()->getImporter());
Actemp->exportAndFillComplete(Acprime,
exporter,
Expand Down
4 changes: 2 additions & 2 deletions packages/tpetra/core/src/Tpetra_CrsMatrix_def.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -90,8 +90,8 @@ namespace { // (anonymous)
const ::Teuchos::EReductionType op,
const Teuchos::Comm<int>& comm)
{
Kokkos::View<const int*, Kokkos::HostSpace> localView (&localValue);
Kokkos::View<int*, Kokkos::HostSpace> globalView (&globalValue);
Kokkos::View<const int*, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > localView (&localValue,1);
Kokkos::View<int*, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > globalView (&globalValue, 1);
return ::Tpetra::Details::iallreduce (localView, globalView, op, comm);
}

Expand Down
5 changes: 3 additions & 2 deletions packages/tpetra/core/src/Tpetra_Import_Util2.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -54,6 +54,7 @@
#include "Tpetra_Util.hpp"
#include "Tpetra_Distributor.hpp"
#include "Tpetra_Details_reallocDualViewIfNeeded.hpp"
#include "Tpetra_Details_MpiTypeTraits.hpp"
#include "Tpetra_Vector.hpp"
#include "Kokkos_DualView.hpp"
#include <Teuchos_Array.hpp>
Expand Down Expand Up @@ -383,7 +384,7 @@ reverseNeighborDiscovery(const CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, No
offset+=ReverseRecvSizes[i];
MPI_Irecv(rec_bptr,
recv_data_size,
MPI_LONG_LONG_INT,
::Tpetra::Details::MpiTypeTraits<GO>::getType(rec_bptr[0]),
ProcsTo[i],
recvData_MPI_Tag,
rawComm,
Expand All @@ -397,7 +398,7 @@ reverseNeighborDiscovery(const CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, No
int sendData_MPI_Tag = mpi_tag_base_*2+MyPID;
MPI_Isend(send_bptr,
send_data_size,
MPI_LONG_LONG_INT,
::Tpetra::Details::MpiTypeTraits<GO>::getType(send_bptr[0]),
ProcsFrom[ii],
sendData_MPI_Tag,
rawComm,
Expand Down

0 comments on commit f94ba63

Please sign in to comment.