Zoltan2
Loading...
Searching...
No Matches
HyperGraphModel.cpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Zoltan2: A package of combinatorial algorithms for scientific computing
4//
5// Copyright 2012 NTESS and the Zoltan2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10//
11// Testing of HyperGraphModel built from APF mesh adapters.
12//
13
20#include <Zoltan2_Standards.hpp>
24
25// Teuchos includes
26#include "Teuchos_RCP.hpp"
27#include "Teuchos_XMLParameterListHelpers.hpp"
28
29// SCOREC includes
30#ifdef HAVE_ZOLTAN2_PARMA
31#include <parma.h>
32#include <apf.h>
33#include <apfMesh.h>
34#include <apfMDS.h>
35#include <apfMesh2.h>
36#include <PCU.h>
37#include <gmi_mesh.h>
38#endif
39
40using Teuchos::ParameterList;
41using Teuchos::RCP;
42
43int main(int narg, char *arg[]) {
44
45 Tpetra::ScopeGuard tscope(&narg, &arg);
46 Teuchos::RCP<const Teuchos::Comm<int> > CommT = Tpetra::getDefaultComm();
47
48#ifdef HAVE_ZOLTAN2_PARMA
49 //Setup for SCOREC
50 PCU_Comm_Init();
51
52 // Generate mesh with MDS
53 gmi_register_mesh();
54 apf::Mesh2* m = apf::loadMdsMesh("../partition/pumiTri14/plate.dmg","../partition/pumiTri14/2/");
55
56 typedef Zoltan2::APFMeshAdapter<apf::Mesh2*> inputAdapter_t;
57 typedef Zoltan2::MeshAdapter<apf::Mesh2*> baseMeshAdapter_t;
58 Teuchos::ParameterList params("test params");
59 params.set("timer_output_stream" , "std::cout");
60 params.set("debug_level", "verbose_detailed_status");
61 params.set("hypergraph_model_type","ghosting");
62
63 RCP<Zoltan2::Environment> env;
64 try{
65 env = rcp(new Zoltan2::Environment(params, Tpetra::getDefaultComm()));
66 }
68
69 RCP<const Zoltan2::Environment> envConst = Teuchos::rcp_const_cast<const Zoltan2::Environment>(env);
70
71
72 inputAdapter_t* ia = new inputAdapter_t(*CommT, m,"vertex","edge",false);
73 inputAdapter_t::scalar_t* arr = new inputAdapter_t::scalar_t[ia->getLocalNumOf(ia->getPrimaryEntityType())];
74 for (size_t i=0;i<ia->getLocalNumOf(ia->getPrimaryEntityType());i++) {
75 arr[i]=PCU_Comm_Self();
76 }
77 const inputAdapter_t::scalar_t* weights=arr;
78 ia->setWeights(ia->getPrimaryEntityType(),weights,1);
79
80 const baseMeshAdapter_t *base_ia = dynamic_cast<const baseMeshAdapter_t*>(ia);
81 Zoltan2::modelFlag_t graphFlags_;
82 RCP<const baseMeshAdapter_t> baseInputAdapter_(base_ia,false);
83
84 Zoltan2::HyperGraphModel<inputAdapter_t> model(baseInputAdapter_,envConst,CommT,
85 graphFlags_,Zoltan2::HYPEREDGE_CENTRIC);
86 ia->destroy();
87 delete ia;
88
89 //Delete APF Mesh;
90 m->destroyNative();
91 apf::destroyMesh(m);
92 //End communications
93 PCU_Comm_Free();
94
95#endif
96
97 std::cout<<"PASS\n";
98 return 0;
99}
Defines the APFMeshAdapter class.
Defines the Environment class.
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Defines the HyperGraphModel interface.
Gathering definitions used in software development.
int main()
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
HyperGraphModel defines the interface required for hyper graph models.
MeshAdapter defines the interface for mesh input.
std::bitset< NUM_MODEL_FLAGS > modelFlag_t
static ArrayRCP< ArrayRCP< zscalar_t > > weights