Zoltan2
GraphModel2ndAdjsFromAdjs.cpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Zoltan2: A package of combinatorial algorithms for scientific computing
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Karen Devine (kddevin@sandia.gov)
39 // Erik Boman (egboman@sandia.gov)
40 // Siva Rajamanickam (srajama@sandia.gov)
41 //
42 // ***********************************************************************
43 //
44 // @HEADER
45 //
46 // Basic testing of Zoltan2::PamgenMeshAdapter
47 
48 #include <Zoltan2_GraphModel.hpp>
49 #include <Zoltan2_ModelHelpers.hpp>
51 
52 // Teuchos includes
53 #include "Teuchos_XMLParameterListHelpers.hpp"
54 
55 // Pamgen includes
56 #include "create_inline_mesh.h"
57 
58 using Teuchos::RCP;
59 
60 /*********************************************************/
61 /* Typedefs */
62 /*********************************************************/
64 
65 
66 
67 /*****************************************************************************/
68 /******************************** MAIN ***************************************/
69 /*****************************************************************************/
70 
71 int main(int narg, char *arg[]) {
72 
73  Tpetra::ScopeGuard tscope(&narg, &arg);
74  Teuchos::RCP<const Teuchos::Comm<int> > CommT = Tpetra::getDefaultComm();
75 
76  int me = CommT->getRank();
77  int numProcs = CommT->getSize();
78 
79  /***************************************************************************/
80  /*************************** GET XML INPUTS ********************************/
81  /***************************************************************************/
82 
83  // default values for command-line arguments
84  std::string xmlMeshInFileName("Poisson.xml");
85 
86  // Read run-time options.
87  Teuchos::CommandLineProcessor cmdp (false, false);
88  cmdp.setOption("xmlfile", &xmlMeshInFileName,
89  "XML file with PamGen specifications");
90  cmdp.parse(narg, arg);
91 
92  // Read xml file into parameter list
93  Teuchos::ParameterList inputMeshList;
94 
95  if(xmlMeshInFileName.length()) {
96  if (me == 0) {
97  std::cout << "\nReading parameter list from the XML file \""
98  <<xmlMeshInFileName<<"\" ...\n\n";
99  }
100  Teuchos::updateParametersFromXmlFile(xmlMeshInFileName,
101  Teuchos::inoutArg(inputMeshList));
102  if (me == 0) {
103  inputMeshList.print(std::cout,2,true,true);
104  std::cout << "\n";
105  }
106  }
107  else {
108  std::cout << "Cannot read input file: " << xmlMeshInFileName << "\n";
109  return 5;
110  }
111 
112  // Get pamgen mesh definition
113  std::string meshInput = Teuchos::getParameter<std::string>(inputMeshList,
114  "meshInput");
115 
116  /***************************************************************************/
117  /********************** GET CELL TOPOLOGY **********************************/
118  /***************************************************************************/
119 
120  // Get dimensions
121  int dim = 3;
122 
123  /***************************************************************************/
124  /***************************** GENERATE MESH *******************************/
125  /***************************************************************************/
126 
127  if (me == 0) std::cout << "Generating mesh ... \n\n";
128 
129  // Generate mesh with Pamgen
130  long long maxInt = std::numeric_limits<long long>::max();
131  Create_Pamgen_Mesh(meshInput.c_str(), dim, me, numProcs, maxInt);
132 
133  // Creating mesh adapter
134  if (me == 0) std::cout << "Creating mesh adapter ... \n\n";
135 
136  typedef Zoltan2::PamgenMeshAdapter<basic_user_t> inputAdapter_t;
138 
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;
145  ia.print(me);
146 
147  if (me == 0) std::cout << "REGION-BASED TEST" << std::endl;
148  Zoltan2::MeshEntityType primaryEType = ia.getPrimaryEntityType();
149  Zoltan2::MeshEntityType adjEType = ia.getAdjacencyEntityType();
150  Zoltan2::MeshEntityType secondAdjEType = ia.getSecondAdjacencyEntityType();
151  RCP<const base_adapter_t> baseInputAdapter;
152  RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(CommT));
153  std::bitset<Zoltan2::NUM_MODEL_FLAGS> modelFlags;
154 
155  if (ia.availAdjs(primaryEType, adjEType)) {
156  // Create a GraphModel based on this input data.
157 
158  if (me == 0) std::cout << " Creating GraphModel" << std::endl;
159 
160  baseInputAdapter = (rcp(dynamic_cast<const base_adapter_t *>(&ia), false));
161 
162  Zoltan2::GraphModel<base_adapter_t> graphModel(baseInputAdapter, env,
163  CommT, modelFlags);
164 
165  if (me == 0)
166  std::cout << " Calling get2ndAdjsViewFromAdjs" << std::endl;
167  Zoltan2::get2ndAdjsViewFromAdjs(baseInputAdapter, graphModel.getComm(),
168  primaryEType,
169  secondAdjEType, moffsets, madjacencyIds);
170 
171  if (me == 0) std::cout << " Checking results" << std::endl;
172  if (ia.avail2ndAdjs(primaryEType, secondAdjEType)) {
173  ia.get2ndAdjsView(primaryEType, secondAdjEType, offsets, adjacencyIds);
174  }
175  else{
176  std::cout << "2nd adjacencies not available; cannot check results"
177  << std::endl;
178  return 2;
179  }
180 
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;
184  return 3;
185  }
186 
187  for (inputAdapter_t::offset_t j=moffsets[telct]; j<moffsets[telct+1]; j++) {
188  ssize_t in_list = -1;
189 
190  for (inputAdapter_t::offset_t k=offsets[telct]; k<offsets[telct+1]; k++) {
191  if (adjacencyIds[k] == madjacencyIds[j]) {
192  in_list = k;
193  break;
194  }
195  }
196 
197  if (in_list < 0) {
198  std::cout << "Adjacency missing" << std::endl;
199  return 4;
200  }
201  }
202  }
203  }
204  else{
205  std::cout << "Adjacencies not available" << std::endl;
206  return 1;
207  }
208 
209  if (me == 0) std::cout << "VERTEX-BASED TEST" << std::endl;
210  primaryEType = ia2.getPrimaryEntityType();
211  adjEType = ia2.getAdjacencyEntityType();
212  secondAdjEType = ia2.getSecondAdjacencyEntityType();
213 
214  if (ia2.availAdjs(primaryEType, adjEType)) {
215  if (ia2.avail2ndAdjs(primaryEType, secondAdjEType)) {
216  ia2.get2ndAdjsView(primaryEType, secondAdjEType, offsets, adjacencyIds);
217  }
218  else{
219  std::cout << "2nd adjacencies not available" << std::endl;
220  return 2;
221  }
222 
223  // Create a GraphModel based on this input data.
224 
225  if (me == 0) std::cout << " Creating GraphModel" << std::endl;
226 
227  baseInputAdapter = (rcp(dynamic_cast<const base_adapter_t *>(&ia2),false));
228 
229  Zoltan2::GraphModel<base_adapter_t> graphModel2(baseInputAdapter, env,
230  CommT, modelFlags);
231 
232  if (me == 0)
233  std::cout << " Calling get2ndAdjsViewFromAdjs" << std::endl;
234  Zoltan2::get2ndAdjsViewFromAdjs(baseInputAdapter, graphModel2.getComm(),
235  primaryEType,
236  secondAdjEType, moffsets, madjacencyIds);
237 
238  if (me == 0) std::cout << " Checking results" << std::endl;
239 
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;
243  return 3;
244  }
245 
246  for (inputAdapter_t::offset_t j=moffsets[tnoct]; j<moffsets[tnoct+1]; j++) {
247  ssize_t in_list = -1;
248 
249  for (inputAdapter_t::offset_t k=offsets[tnoct]; k<offsets[tnoct+1]; k++) {
250  if (adjacencyIds[k] == madjacencyIds[j]) {
251  in_list = k;
252  break;
253  }
254  }
255 
256  if (in_list < 0) {
257  std::cout << "Adjacency missing" << std::endl;
258  return 4;
259  }
260  }
261  }
262  }
263  else{
264  std::cout << "Adjacencies not available" << std::endl;
265  return 1;
266  }
267 
268  // delete mesh
269  if (me == 0) std::cout << "Deleting the mesh ... \n\n";
270 
271  Delete_Pamgen_Mesh();
272 
273  if (me == 0)
274  std::cout << "PASS" << std::endl;
275 
276  return 0;
277 }
278 /*****************************************************************************/
279 /********************************* END MAIN **********************************/
280 /*****************************************************************************/
basic_user_t
Zoltan2::BasicUserTypes< double > basic_user_t
Definition: GraphModel2ndAdjsFromAdjs.cpp:63
Zoltan2_ModelHelpers.hpp
Defines helper functions for use in the models.
Zoltan2_TestingFramework::base_adapter_t
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
Definition: Zoltan2_Typedefs.hpp:168
Zoltan2::MeshEntityType
MeshEntityType
Enumerate entity types for meshes: Regions, Faces, Edges, or Vertices.
Definition: Zoltan2_MeshAdapter.hpp:63
main
int main(int narg, char *arg[])
Definition: GraphModel2ndAdjsFromAdjs.cpp:71
Zoltan2::PamgenMeshAdapter
This class represents a mesh.
Definition: Zoltan2_PamgenMeshAdapter.hpp:91
Zoltan2::get2ndAdjsViewFromAdjs
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)
Definition: Zoltan2_ModelHelpers.hpp:194
Zoltan2::Environment
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
Definition: Zoltan2_Environment.hpp:83
Zoltan2::GraphModel::getComm
const RCP< const Comm< int > > getComm()
Return the communicator used by the model.
Definition: Zoltan2_GraphModel.hpp:137
Zoltan2::BasicUserTypes
A simple class that can be the User template argument for an InputAdapter.
Definition: Zoltan2_InputTraits.hpp:139
Zoltan2_PamgenMeshAdapter.hpp
Defines the PamgenMeshAdapter class.
Zoltan2_GraphModel.hpp
Defines the GraphModel interface.
Zoltan2::GraphModel
GraphModel defines the interface required for graph models.
Definition: Zoltan2_GraphModel.hpp:80