Zoltan2
Loading...
Searching...
No Matches
PamgenMeshInput.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// Basic testing of Zoltan2::PamgenMeshAdapter
12
15
16// Teuchos includes
17#include "Teuchos_XMLParameterListHelpers.hpp"
18
19// Pamgen includes
20#include "create_inline_mesh.h"
21
22/*********************************************************/
23/* Typedefs */
24/*********************************************************/
26
27/*****************************************************************************/
28/******************************** MAIN ***************************************/
29/*****************************************************************************/
30
31int main(int narg, char *arg[]) {
32
33 Tpetra::ScopeGuard tscope(&narg, &arg);
34 Teuchos::RCP<const Teuchos::Comm<int> > CommT = Tpetra::getDefaultComm();
35
36 int me = CommT->getRank();
37 int numProcs = CommT->getSize();
38
39 /***************************************************************************/
40 /*************************** GET XML INPUTS ********************************/
41 /***************************************************************************/
42
43 // default values for command-line arguments
44 std::string xmlMeshInFileName("Poisson.xml");
45
46 // Read run-time options.
47 Teuchos::CommandLineProcessor cmdp (false, false);
48 cmdp.setOption("xmlfile", &xmlMeshInFileName,
49 "XML file with PamGen specifications");
50 cmdp.parse(narg, arg);
51
52 // Read xml file into parameter list
53 Teuchos::ParameterList inputMeshList;
54
55 if(xmlMeshInFileName.length()) {
56 if (me == 0) {
57 std::cout << "\nReading parameter list from the XML file \""
58 <<xmlMeshInFileName<<"\" ...\n\n";
59 }
60 Teuchos::updateParametersFromXmlFile(xmlMeshInFileName,
61 Teuchos::inoutArg(inputMeshList));
62 if (me == 0) {
63 inputMeshList.print(std::cout,2,true,true);
64 std::cout << "\n";
65 }
66 }
67 else {
68 std::cout << "Cannot read input file: " << xmlMeshInFileName << "\n";
69 return 5;
70 }
71
72 // Get pamgen mesh definition
73 std::string meshInput = Teuchos::getParameter<std::string>(inputMeshList,
74 "meshInput");
75
76 /***************************************************************************/
77 /********************** GET CELL TOPOLOGY **********************************/
78 /***************************************************************************/
79
80 // Get dimensions
81 int dim = 3;
82
83 /***************************************************************************/
84 /***************************** GENERATE MESH *******************************/
85 /***************************************************************************/
86
87 if (me == 0) std::cout << "Generating mesh ... \n\n";
88
89 // Generate mesh with Pamgen
90 long long maxInt = 9223372036854775807LL;
91 Create_Pamgen_Mesh(meshInput.c_str(), dim, me, numProcs, maxInt);
92
93 // Creating mesh adapter
94 if (me == 0) std::cout << "Creating mesh adapter ... \n\n";
95
96 typedef Zoltan2::PamgenMeshAdapter<basic_user_t> inputAdapter_t;
97
98 inputAdapter_t ia(*CommT, "region");
99 inputAdapter_t ia2(*CommT, "vertex");
100 inputAdapter_t::gno_t const *adjacencyIds=NULL;
101 inputAdapter_t::offset_t const *offsets=NULL;
102 ia.print(me);
103
104 // Exercise the componentMetrics on the input; make sure the adapter works
105 {
107 std::cout << me << " Region-based: Number of components on processor = "
108 << cc.getNumComponents() << std::endl;
109 std::cout << me << " Region-based: Max component size on processor = "
110 << cc.getMaxComponentSize() << std::endl;
111 std::cout << me << " Region-based: Min component size on processor = "
112 << cc.getMinComponentSize() << std::endl;
113 std::cout << me << " Region-based: Avg component size on processor = "
114 << cc.getAvgComponentSize() << std::endl;
115 }
116
117 Zoltan2::MeshEntityType primaryEType = ia.getPrimaryEntityType();
118 Zoltan2::MeshEntityType adjEType = ia.getAdjacencyEntityType();
119
120 int dimension, num_nodes, num_elem;
121 int error = 0;
122 char title[100];
123 int exoid = 0;
124 int num_elem_blk, num_node_sets, num_side_sets;
125 error += im_ex_get_init(exoid, title, &dimension, &num_nodes, &num_elem,
126 &num_elem_blk, &num_node_sets, &num_side_sets);
127
128 int *element_num_map = new int [num_elem];
129 error += im_ex_get_elem_num_map(exoid, element_num_map);
130
131 int *node_num_map = new int [num_nodes];
132 error += im_ex_get_node_num_map(exoid, node_num_map);
133
134 int *elem_blk_ids = new int [num_elem_blk];
135 error += im_ex_get_elem_blk_ids(exoid, elem_blk_ids);
136
137 int *num_nodes_per_elem = new int [num_elem_blk];
138 int *num_attr = new int [num_elem_blk];
139 int *num_elem_this_blk = new int [num_elem_blk];
140 char **elem_type = new char * [num_elem_blk];
141 int **connect = new int * [num_elem_blk];
142
143 for(int i = 0; i < num_elem_blk; i++){
144 elem_type[i] = new char [MAX_STR_LENGTH + 1];
145 error += im_ex_get_elem_block(exoid, elem_blk_ids[i], elem_type[i],
146 (int*)&(num_elem_this_blk[i]),
147 (int*)&(num_nodes_per_elem[i]),
148 (int*)&(num_attr[i]));
149 delete[] elem_type[i];
150 }
151
152 delete[] elem_type;
153 elem_type = NULL;
154 delete[] num_attr;
155 num_attr = NULL;
156
157 for(int b = 0; b < num_elem_blk; b++) {
158 connect[b] = new int [num_nodes_per_elem[b]*num_elem_this_blk[b]];
159 error += im_ex_get_elem_conn(exoid, elem_blk_ids[b], connect[b]);
160 }
161
162 delete[] elem_blk_ids;
163 elem_blk_ids = NULL;
164
165 int telct = 0;
166
167 if (ia.availAdjs(primaryEType, adjEType)) {
168 ia.getAdjsView(primaryEType, adjEType, offsets, adjacencyIds);
169
170 if ((int)ia.getLocalNumOf(primaryEType) != num_elem) {
171 std::cout << "Number of elements do not match\n";
172 return 2;
173 }
174
175 for (int b = 0; b < num_elem_blk; b++) {
176 for (int i = 0; i < num_elem_this_blk[b]; i++) {
177 if (offsets[telct + 1] - offsets[telct] != num_nodes_per_elem[b]) {
178 std::cout << "Number of adjacencies do not match" << std::endl;
179 return 3;
180 }
181
182 for (int j = 0; j < num_nodes_per_elem[b]; j++) {
183 ssize_t in_list = -1;
184
185 for(inputAdapter_t::offset_t k=offsets[telct];k<offsets[telct+1];k++){
186 if(adjacencyIds[k] ==
187 node_num_map[connect[b][i*num_nodes_per_elem[b]+j]-1]) {
188 in_list = k;
189 break;
190 }
191 }
192
193 if (in_list < 0) {
194 std::cout << "Adjacency missing" << std::endl;
195 return 4;
196 }
197 }
198
199 ++telct;
200 }
201 }
202
203 if (telct != num_elem) {
204 std::cout << "Number of elements do not match\n";
205 return 2;
206 }
207 }
208 else{
209 std::cout << "Adjacencies not available" << std::endl;
210 return 1;
211 }
212
213 ia2.print(me);
214 {
216 std::cout << me << " Vertex-based: Number of components on processor = "
217 << cc.getNumComponents() << std::endl;
218 std::cout << me << " Vertex-based: Max component size on processor = "
219 << cc.getMaxComponentSize() << std::endl;
220 std::cout << me << " Vertex-based: Min component size on processor = "
221 << cc.getMinComponentSize() << std::endl;
222 std::cout << me << " Vertex-based: Avg component size on processor = "
223 << cc.getAvgComponentSize() << std::endl;
224 }
225
226 primaryEType = ia2.getPrimaryEntityType();
227 adjEType = ia2.getAdjacencyEntityType();
228
229 if (ia2.availAdjs(primaryEType, adjEType)) {
230 ia2.getAdjsView(primaryEType, adjEType, offsets, adjacencyIds);
231
232 if ((int)ia2.getLocalNumOf(primaryEType) != num_nodes) {
233 std::cout << "Number of nodes do not match\n";
234 return 2;
235 }
236
237 telct = 0;
238 int *num_adj = new int[num_nodes];
239
240 for (int i = 0; i < num_nodes; i++) {
241 num_adj[i] = 0;
242 }
243
244 for (int b = 0; b < num_elem_blk; b++) {
245 for (int i = 0; i < num_elem_this_blk[b]; i++) {
246 for (int j = 0; j < num_nodes_per_elem[b]; j++) {
247 ssize_t in_list = -1;
248 ++num_adj[connect[b][i * num_nodes_per_elem[b] + j] - 1];
249
250 for(inputAdapter_t::lno_t k =
251 offsets[connect[b][i * num_nodes_per_elem[b] + j] - 1];
252 k < offsets[connect[b][i * num_nodes_per_elem[b] + j]]; k++) {
253 if(adjacencyIds[k] == element_num_map[telct]) {
254 in_list = k;
255 break;
256 }
257 }
258
259 if (in_list < 0) {
260 std::cout << "Adjacency missing" << std::endl;
261 return 4;
262 }
263 }
264
265 ++telct;
266 }
267 }
268
269 for (int i = 0; i < num_nodes; i++) {
270 if (offsets[i + 1] - offsets[i] != num_adj[i]) {
271 std::cout << "Number of adjacencies do not match" << std::endl;
272 return 3;
273 }
274 }
275
276 delete[] num_adj;
277 num_adj = NULL;
278 }
279 else{
280 std::cout << "Adjacencies not available" << std::endl;
281 return 1;
282 }
283
284 // delete mesh
285 if (me == 0) std::cout << "Deleting the mesh ... \n\n";
286
287 Delete_Pamgen_Mesh();
288 // clean up - reduce the result codes
289
290 // make sure another process doesn't mangle the PASS output or the test will read as a fail when it should pass
291 std::cout << std::flush;
292 CommT->barrier();
293 if (me == 0)
294 std::cout << "PASS" << std::endl;
295
296 return 0;
297}
298/*****************************************************************************/
299/********************************* END MAIN **********************************/
300/*****************************************************************************/
Zoltan2::BasicUserTypes basic_user_t
Defines the PamgenMeshAdapter class.
Identify and compute the number of connected components in a processor's input Note that this routine...
int main()
A simple class that can be the User template argument for an InputAdapter.
This class represents a mesh.
MeshEntityType
Enumerate entity types for meshes: Regions, Faces, Edges, or Vertices.