53 #include "Teuchos_XMLParameterListHelpers.hpp"
56 #include "create_inline_mesh.h"
71 int main(
int narg,
char *arg[]) {
73 Tpetra::ScopeGuard tscope(&narg, &arg);
74 Teuchos::RCP<const Teuchos::Comm<int> > CommT = Tpetra::getDefaultComm();
76 int me = CommT->getRank();
77 int numProcs = CommT->getSize();
84 std::string xmlMeshInFileName(
"Poisson.xml");
87 Teuchos::CommandLineProcessor cmdp (
false,
false);
88 cmdp.setOption(
"xmlfile", &xmlMeshInFileName,
89 "XML file with PamGen specifications");
90 cmdp.parse(narg, arg);
93 Teuchos::ParameterList inputMeshList;
95 if(xmlMeshInFileName.length()) {
97 std::cout <<
"\nReading parameter list from the XML file \""
98 <<xmlMeshInFileName<<
"\" ...\n\n";
100 Teuchos::updateParametersFromXmlFile(xmlMeshInFileName,
101 Teuchos::inoutArg(inputMeshList));
103 inputMeshList.print(std::cout,2,
true,
true);
108 std::cout <<
"Cannot read input file: " << xmlMeshInFileName <<
"\n";
113 std::string meshInput = Teuchos::getParameter<std::string>(inputMeshList,
127 if (me == 0) std::cout <<
"Generating mesh ... \n\n";
130 long long maxInt = std::numeric_limits<long long>::max();
131 Create_Pamgen_Mesh(meshInput.c_str(), dim, me, numProcs, maxInt);
134 if (me == 0) std::cout <<
"Creating mesh adapter ... \n\n";
139 inputAdapter_t ia(*CommT,
"region");
140 inputAdapter_t ia2(*CommT,
"vertex");
141 inputAdapter_t::gno_t
const *adjacencyIds=NULL;
142 inputAdapter_t::offset_t
const *offsets=NULL;
143 Teuchos::ArrayRCP<inputAdapter_t::offset_t> moffsets;
144 Teuchos::ArrayRCP<inputAdapter_t::gno_t> madjacencyIds;
147 if (me == 0) std::cout <<
"REGION-BASED TEST" << std::endl;
151 RCP<const base_adapter_t> baseInputAdapter;
153 std::bitset<Zoltan2::NUM_MODEL_FLAGS> modelFlags;
155 if (ia.availAdjs(primaryEType, adjEType)) {
158 if (me == 0) std::cout <<
" Creating GraphModel" << std::endl;
160 baseInputAdapter = (rcp(dynamic_cast<const base_adapter_t *>(&ia),
false));
166 std::cout <<
" Calling get2ndAdjsViewFromAdjs" << std::endl;
169 secondAdjEType, moffsets, madjacencyIds);
171 if (me == 0) std::cout <<
" Checking results" << std::endl;
172 if (ia.avail2ndAdjs(primaryEType, secondAdjEType)) {
173 ia.get2ndAdjsView(primaryEType, secondAdjEType, offsets, adjacencyIds);
176 std::cout <<
"2nd adjacencies not available; cannot check results"
181 for (
size_t telct = 0; telct < ia.getLocalNumOf(primaryEType); telct++) {
182 if (offsets[telct+1]-offsets[telct]!=moffsets[telct+1]-moffsets[telct]) {
183 std::cout <<
"Number of adjacencies do not match" << std::endl;
187 for (inputAdapter_t::offset_t j=moffsets[telct]; j<moffsets[telct+1]; j++) {
188 ssize_t in_list = -1;
190 for (inputAdapter_t::offset_t k=offsets[telct]; k<offsets[telct+1]; k++) {
191 if (adjacencyIds[k] == madjacencyIds[j]) {
198 std::cout <<
"Adjacency missing" << std::endl;
205 std::cout <<
"Adjacencies not available" << std::endl;
209 if (me == 0) std::cout <<
"VERTEX-BASED TEST" << std::endl;
210 primaryEType = ia2.getPrimaryEntityType();
211 adjEType = ia2.getAdjacencyEntityType();
212 secondAdjEType = ia2.getSecondAdjacencyEntityType();
214 if (ia2.availAdjs(primaryEType, adjEType)) {
215 if (ia2.avail2ndAdjs(primaryEType, secondAdjEType)) {
216 ia2.get2ndAdjsView(primaryEType, secondAdjEType, offsets, adjacencyIds);
219 std::cout <<
"2nd adjacencies not available" << std::endl;
225 if (me == 0) std::cout <<
" Creating GraphModel" << std::endl;
227 baseInputAdapter = (rcp(dynamic_cast<const base_adapter_t *>(&ia2),
false));
233 std::cout <<
" Calling get2ndAdjsViewFromAdjs" << std::endl;
236 secondAdjEType, moffsets, madjacencyIds);
238 if (me == 0) std::cout <<
" Checking results" << std::endl;
240 for (
size_t tnoct = 0; tnoct < ia2.getLocalNumOf(primaryEType); tnoct++) {
241 if (offsets[tnoct+1]-offsets[tnoct]!=moffsets[tnoct+1]-moffsets[tnoct]) {
242 std::cout <<
"Number of adjacencies do not match" << std::endl;
246 for (inputAdapter_t::offset_t j=moffsets[tnoct]; j<moffsets[tnoct+1]; j++) {
247 ssize_t in_list = -1;
249 for (inputAdapter_t::offset_t k=offsets[tnoct]; k<offsets[tnoct+1]; k++) {
250 if (adjacencyIds[k] == madjacencyIds[j]) {
257 std::cout <<
"Adjacency missing" << std::endl;
264 std::cout <<
"Adjacencies not available" << std::endl;
269 if (me == 0) std::cout <<
"Deleting the mesh ... \n\n";
271 Delete_Pamgen_Mesh();
274 std::cout <<
"PASS" << std::endl;
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
void get2ndAdjsViewFromAdjs(const Teuchos::RCP< const MeshAdapter< User > > &ia, const RCP< const Comm< int > > comm, Zoltan2::MeshEntityType sourcetarget, Zoltan2::MeshEntityType through, Teuchos::ArrayRCP< typename MeshAdapter< User >::offset_t > &offsets, Teuchos::ArrayRCP< typename MeshAdapter< User >::gno_t > &adjacencyIds)
int main(int narg, char *arg[])
A simple class that can be the User template argument for an InputAdapter.
Defines the PamgenMeshAdapter class.
Defines helper functions for use in the models.
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
GraphModel defines the interface required for graph models.
MeshEntityType
Enumerate entity types for meshes: Regions, Faces, Edges, or Vertices.
Defines the GraphModel interface.
const RCP< const Comm< int > > getComm()
Return the communicator used by the model.
This class represents a mesh.