52 #include "Teuchos_XMLParameterListHelpers.hpp"
55 #include "create_inline_mesh.h"
66 int main(
int narg,
char *arg[]) {
68 Tpetra::ScopeGuard tscope(&narg, &arg);
69 Teuchos::RCP<const Teuchos::Comm<int> > CommT = Tpetra::getDefaultComm();
71 int me = CommT->getRank();
72 int numProcs = CommT->getSize();
79 std::string xmlMeshInFileName(
"Poisson.xml");
82 Teuchos::CommandLineProcessor cmdp (
false,
false);
83 cmdp.setOption(
"xmlfile", &xmlMeshInFileName,
84 "XML file with PamGen specifications");
85 cmdp.parse(narg, arg);
88 Teuchos::ParameterList inputMeshList;
90 if(xmlMeshInFileName.length()) {
92 std::cout <<
"\nReading parameter list from the XML file \""
93 <<xmlMeshInFileName<<
"\" ...\n\n";
95 Teuchos::updateParametersFromXmlFile(xmlMeshInFileName,
96 Teuchos::inoutArg(inputMeshList));
98 inputMeshList.print(std::cout,2,
true,
true);
103 std::cout <<
"Cannot read input file: " << xmlMeshInFileName <<
"\n";
108 std::string meshInput = Teuchos::getParameter<std::string>(inputMeshList,
122 if (me == 0) std::cout <<
"Generating mesh ... \n\n";
125 long long maxInt = 9223372036854775807LL;
126 Create_Pamgen_Mesh(meshInput.c_str(), dim, me, numProcs, maxInt);
129 if (me == 0) std::cout <<
"Creating mesh adapter ... \n\n";
133 inputAdapter_t ia(*CommT,
"region");
134 inputAdapter_t ia2(*CommT,
"vertex");
135 inputAdapter_t::gno_t
const *adjacencyIds=NULL;
136 inputAdapter_t::offset_t
const *offsets=NULL;
142 std::cout << me <<
" Region-based: Number of components on processor = "
144 std::cout << me <<
" Region-based: Max component size on processor = "
146 std::cout << me <<
" Region-based: Min component size on processor = "
148 std::cout << me <<
" Region-based: Avg component size on processor = "
155 int dimension, num_nodes, num_elem;
159 int num_elem_blk, num_node_sets, num_side_sets;
160 error += im_ex_get_init(exoid, title, &dimension, &num_nodes, &num_elem,
161 &num_elem_blk, &num_node_sets, &num_side_sets);
163 int *element_num_map =
new int [num_elem];
164 error += im_ex_get_elem_num_map(exoid, element_num_map);
166 int *node_num_map =
new int [num_nodes];
167 error += im_ex_get_node_num_map(exoid, node_num_map);
169 int *elem_blk_ids =
new int [num_elem_blk];
170 error += im_ex_get_elem_blk_ids(exoid, elem_blk_ids);
172 int *num_nodes_per_elem =
new int [num_elem_blk];
173 int *num_attr =
new int [num_elem_blk];
174 int *num_elem_this_blk =
new int [num_elem_blk];
175 char **elem_type =
new char * [num_elem_blk];
176 int **connect =
new int * [num_elem_blk];
178 for(
int i = 0; i < num_elem_blk; i++){
179 elem_type[i] =
new char [MAX_STR_LENGTH + 1];
180 error += im_ex_get_elem_block(exoid, elem_blk_ids[i], elem_type[i],
181 (
int*)&(num_elem_this_blk[i]),
182 (
int*)&(num_nodes_per_elem[i]),
183 (
int*)&(num_attr[i]));
184 delete[] elem_type[i];
192 for(
int b = 0; b < num_elem_blk; b++) {
193 connect[b] =
new int [num_nodes_per_elem[b]*num_elem_this_blk[b]];
194 error += im_ex_get_elem_conn(exoid, elem_blk_ids[b], connect[b]);
197 delete[] elem_blk_ids;
202 if (ia.availAdjs(primaryEType, adjEType)) {
203 ia.getAdjsView(primaryEType, adjEType, offsets, adjacencyIds);
205 if ((
int)ia.getLocalNumOf(primaryEType) != num_elem) {
206 std::cout <<
"Number of elements do not match\n";
210 for (
int b = 0; b < num_elem_blk; b++) {
211 for (
int i = 0; i < num_elem_this_blk[b]; i++) {
212 if (offsets[telct + 1] - offsets[telct] != num_nodes_per_elem[b]) {
213 std::cout <<
"Number of adjacencies do not match" << std::endl;
217 for (
int j = 0; j < num_nodes_per_elem[b]; j++) {
218 ssize_t in_list = -1;
220 for(inputAdapter_t::offset_t k=offsets[telct];k<offsets[telct+1];k++){
221 if(adjacencyIds[k] ==
222 node_num_map[connect[b][i*num_nodes_per_elem[b]+j]-1]) {
229 std::cout <<
"Adjacency missing" << std::endl;
238 if (telct != num_elem) {
239 std::cout <<
"Number of elements do not match\n";
244 std::cout <<
"Adjacencies not available" << std::endl;
251 std::cout << me <<
" Vertex-based: Number of components on processor = "
253 std::cout << me <<
" Vertex-based: Max component size on processor = "
255 std::cout << me <<
" Vertex-based: Min component size on processor = "
257 std::cout << me <<
" Vertex-based: Avg component size on processor = "
261 primaryEType = ia2.getPrimaryEntityType();
262 adjEType = ia2.getAdjacencyEntityType();
264 if (ia2.availAdjs(primaryEType, adjEType)) {
265 ia2.getAdjsView(primaryEType, adjEType, offsets, adjacencyIds);
267 if ((
int)ia2.getLocalNumOf(primaryEType) != num_nodes) {
268 std::cout <<
"Number of nodes do not match\n";
273 int *num_adj =
new int[num_nodes];
275 for (
int i = 0; i < num_nodes; i++) {
279 for (
int b = 0; b < num_elem_blk; b++) {
280 for (
int i = 0; i < num_elem_this_blk[b]; i++) {
281 for (
int j = 0; j < num_nodes_per_elem[b]; j++) {
282 ssize_t in_list = -1;
283 ++num_adj[connect[b][i * num_nodes_per_elem[b] + j] - 1];
285 for(inputAdapter_t::lno_t k =
286 offsets[connect[b][i * num_nodes_per_elem[b] + j] - 1];
287 k < offsets[connect[b][i * num_nodes_per_elem[b] + j]]; k++) {
288 if(adjacencyIds[k] == element_num_map[telct]) {
295 std::cout <<
"Adjacency missing" << std::endl;
304 for (
int i = 0; i < num_nodes; i++) {
305 if (offsets[i + 1] - offsets[i] != num_adj[i]) {
306 std::cout <<
"Number of adjacencies do not match" << std::endl;
315 std::cout <<
"Adjacencies not available" << std::endl;
320 if (me == 0) std::cout <<
"Deleting the mesh ... \n\n";
322 Delete_Pamgen_Mesh();
326 std::cout << std::flush;
329 std::cout <<
"PASS" << std::endl;
int main(int narg, char *arg[])
A simple class that can be the User template argument for an InputAdapter.
double getAvgComponentSize()
Defines the PamgenMeshAdapter class.
size_t getMaxComponentSize()
size_t getNumComponents()
MeshEntityType
Enumerate entity types for meshes: Regions, Faces, Edges, or Vertices.
size_t getMinComponentSize()
Identify and compute the number of connected components in a processor's input Note that this routine...
This class represents a mesh.