{
MPI_Init(&argc, &argv);
int nprocs, rank;
MPI_Comm comm=MPI_COMM_WORLD;
MPI_Comm_size(comm, &nprocs);
MPI_Comm_rank(comm, &rank);
char * filen = "p8ex1.h5m";
if (argc>1)
filen = argv[1];
ErrCode rc = iMOAB_Initialize(argc, argv);
CHECKRC(rc, "failed to initialize MOAB");
int num_global_vertices=0, num_global_elements=0, num_dimension=0, num_parts=0;
rc = iMOAB_ReadHeaderInfo ( filen, &num_global_vertices, &num_global_elements, &num_dimension,
&num_parts, (int)strlen(filen) );
CHECKRC(rc, "failed to read header info");
if (0==rank)
{
printf("file %s has %d vertices, %d elements, %d parts in partition\n", filen,
num_global_vertices, num_global_elements, num_parts);
}
int appID;
iMOAB_AppID pid=&appID;
rc = iMOAB_RegisterApplication( "PROTEUS", &comm, pid);
CHECKRC(rc, "failed to register application");
char *read_opts="PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;PARALLEL_RESOLVE_SHARED_ENTS";
int num_ghost_layers=1;
rc = iMOAB_LoadMesh( pid, filen, read_opts, &num_ghost_layers, strlen(filen), strlen(read_opts) );
CHECKRC(rc, "failed to load mesh");
int nverts[3], nelem[3], nblocks[3], nsbc[3], ndbc[3];
rc = iMOAB_GetMeshInfo( pid, nverts, nelem, nblocks, nsbc, ndbc);
CHECKRC(rc, "failed to get mesh info");
iMOAB_GlobalID * vGlobalID = (iMOAB_GlobalID*)malloc(nverts[2]*sizeof(iMOAB_GlobalID)) ;
rc = iMOAB_GetVertexID(pid, &nverts[2], vGlobalID);
CHECKRC(rc, "failed to get vertex id info");
int * vranks = (int*)malloc(nverts[2]*sizeof(int));
rc = iMOAB_GetVertexOwnership(pid, &nverts[2], vranks );
CHECKRC(rc, "failed to get vertex ranks");
double * coords = (double*) malloc(3*nverts[2]*sizeof(double));
int size_coords= 3*nverts[2];
rc = iMOAB_GetVisibleVerticesCoordinates( pid, &size_coords, coords);
CHECKRC(rc, "failed to get coordinates");
iMOAB_GlobalID * gbIDs = (iMOAB_GlobalID*) malloc(nblocks[2]*sizeof(iMOAB_GlobalID));
rc = iMOAB_GetBlockID(pid, &nblocks[2], gbIDs);
CHECKRC(rc, "failed to get block info");
int tagIndex[2];
int entTypes[2] = {0, 1};
int tagTypes[2] = { DENSE_INTEGER, DENSE_DOUBLE } ;
int num_components = 1;
rc = iMOAB_DefineTagStorage(pid, "INTFIELD", &tagTypes[0], &num_components, &tagIndex[0], strlen("INTFIELD") );
CHECKRC(rc, "failed to get tag INTFIELD ");
rc = iMOAB_DefineTagStorage(pid, "DFIELD", &tagTypes[1], &num_components, &tagIndex[1], strlen("DFIELD") );
CHECKRC(rc, "failed to get tag DFIELD ");
int num_tags_to_sync=1;
rc = iMOAB_SynchronizeTags(pid, &num_tags_to_sync, &tagIndex[0], &tagTypes[0] );
CHECKRC(rc, "failed to sync tag INTFIELD ");
for (int irank=0; irank<nprocs; irank++)
{
if (irank==rank)
{
printf("on rank %d, there are \n"
" %3d visible vertices of which %3d local %3d ghost \n"
" %3d visible elements of which %3d owned %3d ghost \n"
" %3d visible blocks\n"
" %3d visible neumann BCs\n"
" %3d visible dirichlet BCs\n", rank,
nverts[2], nverts[0], nverts[1],
nelem[2], nelem[0], nelem[1],
nblocks[2], nsbc[2], ndbc[2]);
int numToPrint = nverts[2];
printf("on rank %d vertex info:\n", rank);
for (int i=0; i<numToPrint; i++)
printf(" vertex local id: %3d, rank ID:%d global ID: %3d coords: %g, %g, %g\n", i, vranks[i], vGlobalID[i],
coords[3*i], coords[3*i+1], coords[3*i+2]);
iMOAB_GlobalID * element_global_IDs = (iMOAB_GlobalID*)malloc(nelem[2]*sizeof(iMOAB_GlobalID));
iMOAB_GlobalID * block_IDs = (iMOAB_GlobalID*)malloc(nelem[2]*sizeof(iMOAB_GlobalID));
int * ranks = (int*)malloc(nelem[2]*sizeof(int));
rc = iMOAB_GetVisibleElementsInfo(pid, &nelem[2], element_global_IDs, ranks, block_IDs);
CHECKRC(rc, "failed to get all elem info");
for (int i=0; i<nelem[2]; i++)
printf(" element local id: %3d, global ID: %3d rank:%d block ID: %2d \n", i, element_global_IDs[i],
ranks[i], block_IDs[i]);
free(element_global_IDs);
free(ranks);
free(block_IDs);
int conn[27];
int nv = 27;
int eindex = 0;
rc = iMOAB_GetElementConnectivity(pid, &eindex, &nv, conn);
CHECKRC(rc, "failed to get first element connectivity");
printf(" conn for first element: \n");
for (int i=0; i<nv; i++)
printf(" %3d", conn[i]);
printf("\n");
int local_index = 0;
int num_adjacent_elements = 10;
iMOAB_LocalID adjacent_element_IDs[10];
rc = iMOAB_GetNeighborElements(pid, &local_index, &num_adjacent_elements, adjacent_element_IDs);
CHECKRC(rc, "failed to get first element neighbors");
printf(" neighbors for first element:\n");
for (int i=0; i<num_adjacent_elements; i++)
{
printf(" %4d", adjacent_element_IDs[i]);
}
printf("\n");
for (int i=0; i<nblocks[2]; i++)
{
printf(" block index: %3d, block ID: %3d \n", i, gbIDs[i] );
int vertices_per_element, num_elements_in_block;
rc = iMOAB_GetBlockInfo(pid, &gbIDs[i] , &vertices_per_element, &num_elements_in_block);
CHECKRC(rc, "failed to elem block info");
printf(" has %4d elements with %d vertices per element\n", num_elements_in_block, vertices_per_element);
int size_conn= num_elements_in_block*vertices_per_element;
iMOAB_LocalID * element_connectivity = (iMOAB_LocalID*) malloc (sizeof(iMOAB_LocalID)*size_conn);
rc = iMOAB_GetBlockElementConnectivities(pid, &gbIDs[i], &size_conn, element_connectivity);
CHECKRC(rc, "failed to get block elem connectivity");
int * element_ownership = (int*) malloc (sizeof(int)*num_elements_in_block);
rc = iMOAB_GetElementOwnership(pid, &gbIDs[i], &num_elements_in_block, element_ownership);
CHECKRC(rc, "failed to get block elem ownership");
iMOAB_GlobalID* global_element_ID = (iMOAB_GlobalID*)malloc(sizeof(iMOAB_GlobalID)*num_elements_in_block);
iMOAB_LocalID* local_element_ID =(iMOAB_LocalID*)malloc(sizeof(iMOAB_LocalID)*num_elements_in_block);
rc = iMOAB_GetElementID(pid, &gbIDs[i], &num_elements_in_block, global_element_ID, local_element_ID);
CHECKRC(rc, "failed to get block elem IDs");
for (int j=0; j< num_elements_in_block; j++)
{
printf(" elem %3d owned by %d gid: %4d lid: %4d -- ", j, element_ownership[j], global_element_ID[j], local_element_ID[j]);
for (int k=0; k<vertices_per_element; k++)
printf( " %5d", element_connectivity[j*vertices_per_element+k]);
printf("\n");
}
free(global_element_ID);
free(local_element_ID);
free (element_connectivity);
free (element_ownership);
}
int * int_tag_vals = (int *) malloc (sizeof(int) * nverts[2]);
rc = iMOAB_GetIntTagStorage(pid, "INTFIELD", &nverts[2], &entTypes[0],
int_tag_vals, strlen("INTFIELD"));
CHECKRC(rc, "failed to get INTFIELD tag");
printf("INTFIELD tag values:\n");
for (int i=0; i<nverts[2]; i++)
{
printf(" %4d", int_tag_vals[i]);
if (i%20==19)
printf("\n");
}
printf("\n");
free(int_tag_vals);
double * double_tag_vals = (double *) malloc (sizeof(double) * nelem[2]);
rc = iMOAB_GetDoubleTagStorage(pid, "DFIELD", &nelem[2], &entTypes[1],
double_tag_vals, strlen("DFIELD"));
CHECKRC(rc, "failed to get DFIELD tag");
printf("DFIELD tag values: (not exchanged) \n");
for (int i=0; i<nelem[2]; i++)
{
printf(" %f", double_tag_vals[i]);
if (i%8==7)
printf("\n");
}
printf("\n");
free(double_tag_vals);
iMOAB_LocalID * surfBC_ID = (iMOAB_LocalID*) malloc (sizeof(iMOAB_LocalID)*nsbc[2]);
int * ref_surf = (int *) malloc (sizeof(int)*nsbc[2]);
int * bc_value = (int *) malloc (sizeof(int)*nsbc[2]);
rc = iMOAB_GetPointerToSurfaceBC(pid, &nsbc[2], surfBC_ID, ref_surf, bc_value);
CHECKRC(rc, "failed to get surf boundary conditions");
printf(" Surface boundary conditions:\n");
for (int i=0; i<nsbc[2]; i++)
{
printf(" elem_localID %4d side:%d BC:%2d\n",surfBC_ID[i] ,ref_surf[i], bc_value[i] );
}
free(surfBC_ID);
free(ref_surf);
free(bc_value);
iMOAB_LocalID * vertBC_ID = (iMOAB_LocalID*) malloc (sizeof(iMOAB_LocalID)*ndbc[2]);
int * vertBC_value = (int *) malloc (sizeof(int)*ndbc[2]);
rc = iMOAB_GetPointerToVertexBC(pid, &ndbc[2], vertBC_ID, vertBC_value);
CHECKRC(rc, "failed to get vertex boundary conditions");
printf(" Vertex boundary conditions:\n");
for (int i=0; i<ndbc[2]; i++)
{
printf(" vertex %4d BC:%2d\n",vertBC_ID[i], vertBC_value[i] );
}
free(vertBC_ID);
free(vertBC_value);
}
MPI_Barrier(comm);
}
free(coords);
free (vGlobalID);
free (vranks);
char outputFile[] = "fnew.h5m";
char writeOptions[] ="PARALLEL=WRITE_PART";
rc = iMOAB_WriteMesh(pid, outputFile, writeOptions,
strlen(outputFile), strlen(writeOptions) );
rc = iMOAB_DeregisterApplication(pid);
CHECKRC(rc, "failed to de-register application");
rc = iMOAB_Finalize();
CHECKRC(rc, "failed to finalize MOAB");
MPI_Finalize();
return 0;
}