Skip to content
Snippets Groups Projects
Commit 128a04a2 authored by Simon Praetorius's avatar Simon Praetorius
Browse files

Add dedicated include for metis and parmetis headers

parent 9b924939
Branches
Tags
1 merge request!1315Provide dedicated include for metis and parmetis headers
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
#ifndef DUNE_METIS_HH
#define DUNE_METIS_HH
#if HAVE_METIS
#if HAVE_SCOTCH_METIS
extern "C" {
#include <scotch.h>
}
#endif
extern "C" {
#include <metis.h>
}
#if HAVE_SCOTCH_METIS && !defined(SCOTCH_METIS_RETURN) && !defined(METIS_OK)
// NOTE: scotchmetis does not define a return type for METIS functions
#define METIS_OK 1
#endif
namespace Dune::Metis {
#if defined(REALTYPEWIDTH) || defined(SCOTCH_METIS_DATATYPES)
using real_t = ::real_t;
#else
using real_t = double;
#endif
#if defined(IDXTYPEWIDTH) || defined(SCOTCH_METIS_DATATYPES)
using idx_t = ::idx_t;
#elif HAVE_SCOTCH_METIS
using idx_t = SCOTCH_Num;
#else
using idx_t = int;
#endif
} // end namespace Dune::Metis
#endif // HAVE_METIS
#endif // DUNE_METIS_HH
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
#ifndef DUNE_PARMETIS_HH
#define DUNE_PARMETIS_HH
#if HAVE_PARMETIS
#include <dune/common/parallel/mpi.hh>
#if HAVE_PTSCOTCH_PARMETIS
extern "C" {
#include <ptscotch.h>
}
#endif
extern "C" {
#include <parmetis.h>
}
namespace Dune::ParMetis {
#if defined(REALTYPEWIDTH)
using real_t = ::real_t;
#else
using real_t = float;
#endif
#if defined(IDXTYPEWIDTH)
using idx_t = ::idx_t;
#elif HAVE_PTSCOTCH_PARMETIS
using idx_t = SCOTCH_Num;
#else
using idx_t = int;
#endif
} // end namespace Dune::ParMetis
#endif // HAVE_PARMETIS
#endif // DUNE_PARMETIS_HH
......@@ -12,6 +12,14 @@ dune_add_test(SOURCES indexsettest.cc
LINK_LIBRARIES dunecommon
LABELS quick)
dune_add_test(SOURCES parmetistest.cc
MPI_RANKS 3
TIMEOUT 300
CMAKE_GUARD "(HAVE_MPI && ParMETIS_FOUND)"
LABELS quick)
add_dune_mpi_flags(parmetistest)
add_dune_parmetis_flags(parmetistest)
dune_add_test(SOURCES remoteindicestest.cc
LINK_LIBRARIES dunecommon
MPI_RANKS 1 2 4
......
......@@ -9,33 +9,13 @@
#error "ParMETIS is required for this test."
#endif
#include <mpi.h>
#if HAVE_PTSCOTCH_PARMETIS
extern "C" {
#include <ptscotch.h>
}
#endif
extern "C" {
#include <parmetis.h>
}
#include <dune/common/parallel/mpi.hh>
#include <dune/common/parallel/parmetis.hh>
int main(int argc, char **argv)
{
#if defined(REALTYPEWIDTH)
using real_t = ::real_t;
#else
using real_t = float;
#endif
#if defined(IDXTYPEWIDTH)
using idx_t = ::idx_t;
#elif HAVE_PTSCOTCH_PARMETIS
using idx_t = SCOTCH_Num;
#else
using idx_t = int;
#endif
using real_t = Dune::ParMetis::real_t;
using idx_t = Dune::ParMetis::idx_t;
MPI_Init(&argc, &argv);
......
......@@ -275,13 +275,6 @@ dune_add_test(SOURCES parametertreetest.cc
LINK_LIBRARIES dunecommon
LABELS quick)
dune_add_test(SOURCES parmetistest.cc
MPI_RANKS 3
TIMEOUT 300
CMAKE_GUARD ParMETIS_FOUND
LABELS quick)
add_dune_parmetis_flags(parmetistest)
dune_add_test(SOURCES pathtest.cc
LINK_LIBRARIES dunecommon
LABELS quick)
......
......@@ -8,36 +8,12 @@
#error "METIS is required for this test"
#endif
#if HAVE_SCOTCH_METIS
extern "C" {
#include <scotch.h>
}
#endif
extern "C" {
#include <metis.h>
}
#if HAVE_SCOTCH_METIS && !defined(SCOTCH_METIS_RETURN)
// NOTE: scotchmetis does not define a return type for METIS functions
#define METIS_OK 1
#endif
#include <dune/common/metis.hh>
int main()
{
#if defined(REALTYPEWIDTH) || defined(SCOTCH_METIS_DATATYPES)
using real_t = ::real_t;
#else
using real_t = double;
#endif
#if defined(IDXTYPEWIDTH) || defined(SCOTCH_METIS_DATATYPES)
using idx_t = ::idx_t;
#elif HAVE_SCOTCH_METIS
using idx_t = SCOTCH_Num;
#else
using idx_t = int;
#endif
using real_t = Dune::Metis::real_t;
using idx_t = Dune::Metis::idx_t;
idx_t nVertices = 6; // number of vertices
idx_t nCon = 1; // number of constraints
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment