Skip to content

add support for Zoltan with PT-Scotch #478

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Draft
wants to merge 2 commits into
base: develop
Choose a base branch
from
Draft
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
47 changes: 42 additions & 5 deletions cmake/FindZoltan.cmake
Original file line number Diff line number Diff line change
Expand Up @@ -10,14 +10,51 @@ if(ZOLTAN_PREFIX)
message(STATUS "ZOLTAN_PREFIX ${ZOLTAN_PREFIX}")
endif()

find_path(ZOLTAN_INCLUDE_DIR zoltan.h PATHS "${ZOLTAN_PREFIX}/include")

find_library(ZOLTAN_LIBRARY zoltan PATHS "${ZOLTAN_PREFIX}/lib")
if(ZOLTAN_PREFIX)
find_path(ZOLTAN_INCLUDE_DIR zoltan.h
PATHS "${ZOLTAN_PREFIX}/include"
NO_DEFAULT_PATH
)
find_library(ZOLTAN_LIBRARY zoltan
PATHS "${ZOLTAN_PREFIX}/lib"
NO_DEFAULT_PATH
)
else()
find_path(ZOLTAN_INCLUDE_DIR zoltan.h)
find_library(ZOLTAN_LIBRARY zoltan)
endif()
if (ZOLTAN_INCLUDE_DIR)
message(DEBUG "ZOLTAN_INCLUDE_DIR ${ZOLTAN_INCLUDE_DIR}")
else()
message(FATAL_ERROR "Zoltan not found.")
endif()

set(ZOLTAN_LIBRARIES ${ZOLTAN_LIBRARY} )
set(ZOLTAN_INCLUDE_DIRS ${ZOLTAN_INCLUDE_DIR} )

find_package(Parmetis MODULE REQUIRED)
find_file(ZOLTAN_CONFIG_FILE Zoltan_config.h
"${ZOLTAN_INCLUDE_DIRS}"
NO_DEFAULT_PATH
)
message(DEBUG "ZOLTAN_CONFIG_FILE ${ZOLTAN_CONFIG_FILE}")
file(READ "${ZOLTAN_CONFIG_FILE}" ZOLTAN_CONFIG_TEXT)
string(FIND
"${ZOLTAN_CONFIG_TEXT}" "#define HAVE_PARMETIS"
ZOLTAN_FIND_PARMETIS
)
string(FIND
"${ZOLTAN_CONFIG_TEXT}" "#define HAVE_SCOTCH"
ZOLTAN_FIND_PTSCOTCH
)

if(NOT ZOLTAN_FIND_PARMETIS EQUAL -1)
find_package(Parmetis MODULE REQUIRED)
set(PUMI_HAS_PARMETIS TRUE)
endif()
if(NOT ZOLTAN_FIND_PTSCOTCH EQUAL -1)
find_package(SCOTCH CONFIG REQUIRED)
set(PUMI_HAS_PTSCOTCH TRUE)
endif()

include(FindPackageHandleStandardArgs)
# handle the QUIETLY and REQUIRED arguments and set ZOLTAN_FOUND to TRUE
Expand All @@ -28,4 +65,4 @@ find_package_handle_standard_args(
ZOLTAN_LIBRARY ZOLTAN_INCLUDE_DIR
)

mark_as_advanced(ZOLTAN_INCLUDE_DIR ZOLTAN_LIBRARY )
mark_as_advanced(ZOLTAN_INCLUDE_DIR ZOLTAN_LIBRARY ZOLTAN_CONFIG_FILE)
22 changes: 13 additions & 9 deletions zoltan/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -6,7 +6,6 @@ endif()
# Package options
option(ENABLE_ZOLTAN "Enable Zoltan interface [ON|OFF]" OFF)
xsdk_add_tpl(ZOLTAN)
xsdk_add_tpl(PARMETIS)
message(STATUS "ENABLE_ZOLTAN: " ${ENABLE_ZOLTAN})

if(SCOREC_NO_MPI AND ENABLE_ZOLTAN)
Expand Down Expand Up @@ -58,15 +57,20 @@ target_link_libraries(apf_zoltan PUBLIC pcu apf)

# Do extra work if zoltan is enabled
if(ENABLE_ZOLTAN)
target_include_directories(apf_zoltan PRIVATE
${ZOLTAN_INCLUDE_DIRS}
${PARMETIS_INCLUDE_DIRS}
)
target_link_libraries(apf_zoltan PUBLIC
${ZOLTAN_LIBRARIES}
${PARMETIS_LIBRARIES}
)
target_include_directories(apf_zoltan PRIVATE ${ZOLTAN_INCLUDE_DIRS})
target_link_libraries(apf_zoltan PUBLIC ${ZOLTAN_LIBRARIES})
target_compile_definitions(apf_zoltan PUBLIC PUMI_HAS_ZOLTAN)
if(PUMI_HAS_PTSCOTCH)
target_link_libraries(apf_zoltan PUBLIC
SCOTCH::ptscotch SCOTCH::ptscotcherr
)
target_compile_definitions(apf_zoltan PUBLIC PUMI_HAS_PTSCOTCH)
endif()
if(PUMI_HAS_PARMETIS)
target_include_directories(apf_zoltan PRIVATE ${PARMETIS_INCLUDE_DIRS})
target_link_libraries(apf_zoltan PUBLIC ${PARMETIS_LIBRARIES})
target_compile_definitions(apf_zoltan PUBLIC PUMI_HAS_PARMETIS)
endif()
endif()

scorec_export_library(apf_zoltan)
Expand Down
1 change: 1 addition & 0 deletions zoltan/apfZoltan.h
Original file line number Diff line number Diff line change
Expand Up @@ -42,6 +42,7 @@ enum ZoltanMethod {
HYPERGRAPH,
/** \brief Use ParMetis */
PARMETIS,
PTSCOTCH, /**< Use PT-Scotch */
/** \brief General graph partitionig */
GRAPH
};
Expand Down
40 changes: 38 additions & 2 deletions zoltan/apfZoltanCallbacks.cc
Original file line number Diff line number Diff line change
Expand Up @@ -9,7 +9,9 @@
#include "apfZoltanMesh.h"
#include "apfZoltan.h"
#include "apfShape.h"
#ifdef PUMI_HAS_PARMETIS
#include <metis.h>
#endif
#include <pcu_util.h>
#include <lionPrint.h>
#include <cstdlib>
Expand All @@ -29,9 +31,33 @@ static int setZoltanLbMethod(struct Zoltan_Struct* ztn, ZoltanMesh* zb)
case HYPERGRAPH:
lbMethod = "HYPERGRAPH"; break;
case PARMETIS: //fall into GRAPH settings
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

we should remove the 'fall into GRAPH settings' comment

lbMethod = "GRAPH";
#ifdef PUMI_HAS_PARMETIS
Zoltan_Set_Param(ztn, "GRAPH_PACKAGE", "PARMETIS");
#else
lion_oprint(1, "WARNING: ParMETIS ZoltanMethod requested but ParMETIS"
" was not enabled at build time.\n");
Comment on lines +38 to +39
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think we should make this an error (i.e., stop execution) or, indicate to the user that PHG will be used instead (IIUC):
https://sandialabs.github.io/Zoltan/ug_html/ug_alg_graph.html

Also, lion_oprint will write on all ranks; we should probably only print from rank 0 here.

#endif
break;
case PTSCOTCH:
lbMethod = "GRAPH";
#ifdef PUMI_HAS_PTSCOTCH
Zoltan_Set_Param(ztn, "GRAPH_PACKAGE", "Scotch");
#else
lion_oprint(1, "WARNING: PT-Scotch ZoltanMethod requested but PT-Scotch"
" was not enabled at build time.\n");
Comment on lines +47 to +48
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think we should make this an error (i.e., stop execution) or, indicate to the user that PHG will be used instead (IIUC):
https://sandialabs.github.io/Zoltan/ug_html/ug_alg_graph.html

Also, lion_oprint will write on all ranks; we should probably only print from rank 0 here.

#endif
break;
case GRAPH:
lbMethod = "GRAPH";
Zoltan_Set_Param(ztn, "GRAPH_PACKAGE", "PARMETIS"); // instead of PHG
// Prefer ParMETIS, then PT-Scotch, then Zoltan-native PHG.
#if defined(PUMI_HAS_PARMETIS)
Zoltan_Set_Param(ztn, "GRAPH_PACKAGE", "PARMETIS");
#elif defined(PUMI_HAS_PTSCOTCH)
Zoltan_Set_Param(ztn, "GRAPH_PACKAGE", "Scotch");
#else
// Zoltan_Set_Param(ztn, "HYPERGRAPH_PACKAGE", "PHG"); //FIXME: not required?
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

It looks like we are providing at least some of the hypergraph callbacks

Zoltan_Set_Fn(ztn, ZOLTAN_HG_SIZE_CS_FN_TYPE, (void (*)()) getHgSize, (void*) (zb));
Zoltan_Set_Fn(ztn, ZOLTAN_HG_CS_FN_TYPE, (void (*)()) getHg, (void*) (zb));

so I think this should work as a decent fallback.

Copy link
Collaborator Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think even if the hypergraph functions are not fully implemented we can use PHG and the hypergraph will be constructed using the graph query operations.

https://sandialabs.github.io/Zoltan/ug_html/ug_graph_vs_hg.html

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

OK. Thank you.

So, whatever we do for the above 'warnings' (i.e., GRAPH requested but PHG provided) we should do that here in the #else block.
https://sandialabs.github.io/Zoltan/ug_html/ug_alg_hypergraph.html

#endif
break;
default:
lion_oprint(1,"ERROR %s Invalid LB_METHOD %d\n",__func__, zb->method);
Expand Down Expand Up @@ -68,8 +94,14 @@ static int setZoltanLbApproach(struct Zoltan_Struct* ztn, ZoltanMesh* zb)
return 1;
}
Zoltan_Set_Param(ztn, "LB_APPROACH", ptnAp.c_str());
if ( (3 == zb->method) || (4 == zb->method) )
if (zb->method == PARMETIS
#ifdef PUMI_HAS_PARMETIS // in this case GRAPH implies PARMETIS.
|| zb->method == GRAPH
#endif
) {
Zoltan_Set_Param(ztn, "PARMETIS_METHOD", pMethod.c_str());
}
// No zb->method == PTSCOTCH because Zoltan only supports RBISECT.
return 0;
}

Expand Down Expand Up @@ -304,7 +336,9 @@ ZoltanData::~ZoltanData()
void ZoltanData::run()
{
/* ensure Metis indices are the same size as Zoltan indices */
#ifdef PUMI_HAS_PARMETIS
PCU_ALWAYS_ASSERT(IDXTYPEWIDTH == sizeof(ZOLTAN_ID_TYPE)*8);
#endif
setup();
ptn();
}
Expand All @@ -331,7 +365,9 @@ void ZoltanData::setup()
if ( zb->isLocal && 0 != m->getPCU()->Self() )
snprintf(paramStr, 128, "%d", 0); //if local silence all but rank 0
Zoltan_Set_Param(ztn, "debug_level", paramStr);
#ifdef PUMI_HAS_PARMETIS
Zoltan_Set_Param(ztn, "PARMETIS_OUTPUT_LEVEL", paramStr);
#endif
Zoltan_Set_Param(ztn, "CHECK_GRAPH", "0");
Zoltan_Set_Param(ztn, "CHECK_HYPERGRAPH", "0");

Expand Down