xref: /petsc/src/dm/impls/moab/dmmbutil.cxx (revision c528d872bb44edeaa2f3a1e7c43d15e14485f957)
1af0996ceSBarry Smith #include <petsc/private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
2af0996ceSBarry Smith #include <petsc/private/vecimpl.h>
351d15aeeSVijay Mahadevan 
451d15aeeSVijay Mahadevan #include <petscdmmoab.h>
551d15aeeSVijay Mahadevan #include <MBTagConventions.hpp>
651d15aeeSVijay Mahadevan #include <moab/ReadUtilIface.hpp>
751d15aeeSVijay Mahadevan #include <moab/ScdInterface.hpp>
851d15aeeSVijay Mahadevan #include <moab/CN.hpp>
951d15aeeSVijay Mahadevan 
1051d15aeeSVijay Mahadevan 
11aa859218SVijay Mahadevan static PetscErrorCode DMMoabComputeDomainBounds_Private(moab::ParallelComm* pcomm, PetscInt dim, PetscInt neleglob, PetscInt *ise)
1251d15aeeSVijay Mahadevan {
1351d15aeeSVijay Mahadevan   PetscInt size,rank;
1451d15aeeSVijay Mahadevan   PetscInt fraction,remainder;
1551d15aeeSVijay Mahadevan   PetscInt starts[3],sizes[3];
1651d15aeeSVijay Mahadevan 
1751d15aeeSVijay Mahadevan   PetscFunctionBegin;
1851d15aeeSVijay Mahadevan   if(dim<1 && dim>3) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The problem dimension is invalid: %D",dim);
1951d15aeeSVijay Mahadevan 
2051d15aeeSVijay Mahadevan   size=pcomm->size();
2151d15aeeSVijay Mahadevan   rank=pcomm->rank();
2251d15aeeSVijay Mahadevan   fraction=neleglob/size;    /* partition only by the largest dimension */
2351d15aeeSVijay Mahadevan   remainder=neleglob%size;   /* remainder after partition which gets evenly distributed by round-robin */
2451d15aeeSVijay Mahadevan 
2551d15aeeSVijay Mahadevan   if(fraction==0) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The leading dimension size should be greater than number of processors: %D < %D",neleglob,size);
2651d15aeeSVijay Mahadevan 
2751d15aeeSVijay Mahadevan   starts[0]=starts[1]=starts[2]=0;       /* default dimensional element = 1 */
2851d15aeeSVijay Mahadevan   sizes[0]=sizes[1]=sizes[2]=neleglob;   /* default dimensional element = 1 */
2951d15aeeSVijay Mahadevan 
30b8ecf6d3SVijay Mahadevan   if(rank < remainder) { /* This process gets "fraction+1" elements */
3151d15aeeSVijay Mahadevan     sizes[dim-1] = (fraction + 1);
3251d15aeeSVijay Mahadevan     starts[dim-1] = rank * (fraction+1);
33b8ecf6d3SVijay Mahadevan   } else { /* This process gets "fraction" elements */
3451d15aeeSVijay Mahadevan     sizes[dim-1] = fraction;
3551d15aeeSVijay Mahadevan     starts[dim-1] = (remainder*(fraction+1) + fraction*(rank-remainder));
3651d15aeeSVijay Mahadevan   }
3751d15aeeSVijay Mahadevan 
3851d15aeeSVijay Mahadevan   for(int i=dim-1; i>=0; --i) {
3951d15aeeSVijay Mahadevan     ise[2*i]=starts[i];ise[2*i+1]=starts[i]+sizes[i];
4051d15aeeSVijay Mahadevan   }
4151d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
4251d15aeeSVijay Mahadevan }
4351d15aeeSVijay Mahadevan 
44aa859218SVijay Mahadevan static void DMMoab_SetStructuredCoords_Private(PetscInt i, PetscInt j, PetscInt k, PetscReal hx, PetscReal hy, PetscReal hz, PetscInt vcount, std::vector<double*>& vcoords)
45e5136372SVijay Mahadevan {
4666f88a78SVijay Mahadevan   vcoords[0][vcount] = i*hx;
4766f88a78SVijay Mahadevan   vcoords[1][vcount] = j*hy;
4866f88a78SVijay Mahadevan   vcoords[2][vcount] = k*hz;
49e5136372SVijay Mahadevan }
50e5136372SVijay Mahadevan 
51c3dd00c7SVijay Mahadevan static void DMMoab_SetTensorElementConnectivity_Private(PetscInt dim, moab::EntityType etype, PetscInt offset, PetscInt nele, PetscInt i, PetscInt j, PetscInt k, PetscInt vfirst, moab::EntityHandle *connectivity)
52e5136372SVijay Mahadevan {
53e5136372SVijay Mahadevan   std::vector<int>    subent_conn(pow(2,dim));  /* only linear edge, quad, hex supported now */
54e5136372SVijay Mahadevan 
55e5136372SVijay Mahadevan   moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data());
56e5136372SVijay Mahadevan 
57e5136372SVijay Mahadevan   switch(dim) {
58e5136372SVijay Mahadevan     case 1:
59e5136372SVijay Mahadevan       connectivity[offset+subent_conn[0]] = vfirst+i;
60e5136372SVijay Mahadevan       connectivity[offset+subent_conn[1]] = vfirst+(i+1);
61e5136372SVijay Mahadevan       break;
62e5136372SVijay Mahadevan     case 2:
63e5136372SVijay Mahadevan       connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1);
64e5136372SVijay Mahadevan       connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1);
65e5136372SVijay Mahadevan       connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1);
66e5136372SVijay Mahadevan       connectivity[offset+subent_conn[3]] = vfirst+i+(j+1)*(nele+1);
67e5136372SVijay Mahadevan       break;
68e5136372SVijay Mahadevan     case 3:
69e5136372SVijay Mahadevan     default:
70e5136372SVijay Mahadevan       connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
71e5136372SVijay Mahadevan       connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k);
72e5136372SVijay Mahadevan       connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
73e5136372SVijay Mahadevan       connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k);
74e5136372SVijay Mahadevan       connectivity[offset+subent_conn[4]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1));
75e5136372SVijay Mahadevan       connectivity[offset+subent_conn[5]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
76e5136372SVijay Mahadevan       connectivity[offset+subent_conn[6]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1));
77e5136372SVijay Mahadevan       connectivity[offset+subent_conn[7]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
78e5136372SVijay Mahadevan       break;
79e5136372SVijay Mahadevan   }
80e5136372SVijay Mahadevan }
8151d15aeeSVijay Mahadevan 
82c3dd00c7SVijay Mahadevan static void DMMoab_SetSimplexElementConnectivity_Private(PetscInt dim, PetscInt subelem, moab::EntityType etype, PetscInt offset, PetscInt nele, PetscInt i, PetscInt j, PetscInt k, PetscInt vfirst, moab::EntityHandle *connectivity)
83c3dd00c7SVijay Mahadevan {
84c3dd00c7SVijay Mahadevan   std::vector<int>    subent_conn(pow(2,dim));  /* only linear edge, quad, hex supported now */
85c3dd00c7SVijay Mahadevan 
86c3dd00c7SVijay Mahadevan   moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data());
87c3dd00c7SVijay Mahadevan 
88c3dd00c7SVijay Mahadevan   switch(dim) {
89c3dd00c7SVijay Mahadevan     case 1:
90c3dd00c7SVijay Mahadevan       connectivity[offset+subent_conn[0]] = vfirst+i;
91c3dd00c7SVijay Mahadevan       connectivity[offset+subent_conn[1]] = vfirst+(i+1);
92c3dd00c7SVijay Mahadevan       break;
93c3dd00c7SVijay Mahadevan     case 2:
94c3dd00c7SVijay Mahadevan       if (subelem) { /* 1 2 3 of a QUAD */
95c3dd00c7SVijay Mahadevan         connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1);
96c3dd00c7SVijay Mahadevan         connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1);
97c3dd00c7SVijay Mahadevan         connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1);
98c3dd00c7SVijay Mahadevan       }
99c3dd00c7SVijay Mahadevan       else {        /* 3 4 1 of a QUAD */
100c3dd00c7SVijay Mahadevan         connectivity[offset+subent_conn[0]] = vfirst+(i+1)+(j+1)*(nele+1);
101c3dd00c7SVijay Mahadevan         connectivity[offset+subent_conn[1]] = vfirst+i+(j+1)*(nele+1);
102c3dd00c7SVijay Mahadevan         connectivity[offset+subent_conn[2]] = vfirst+i+j*(nele+1);
103c3dd00c7SVijay Mahadevan       }
104c3dd00c7SVijay Mahadevan       break;
105c3dd00c7SVijay Mahadevan     case 3:
106c3dd00c7SVijay Mahadevan     default:
107c3dd00c7SVijay Mahadevan       switch(subelem) {
108c3dd00c7SVijay Mahadevan         case 0: /* 0 1 2 5 of a HEX */
109c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
110c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k);
111c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
112c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[3]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
113c3dd00c7SVijay Mahadevan           break;
114c3dd00c7SVijay Mahadevan         case 1: /* 0 2 7 5 of a HEX */
115c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
116c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
117c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[2]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
118c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[3]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
119c3dd00c7SVijay Mahadevan           break;
120c3dd00c7SVijay Mahadevan         case 2: /* 0 2 3 7 of a HEX */
121c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
122c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
123c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[2]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k);
124c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
125c3dd00c7SVijay Mahadevan           break;
126c3dd00c7SVijay Mahadevan         case 3: /* 0 5 7 4 of a HEX */
127c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
128c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
129c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[2]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
130c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1));
131c3dd00c7SVijay Mahadevan           break;
132c3dd00c7SVijay Mahadevan         case 4: /* 2 7 5 6 of a HEX */
133c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[0]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
134c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[1]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
135c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
136c3dd00c7SVijay Mahadevan           connectivity[offset+subent_conn[3]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1));
137c3dd00c7SVijay Mahadevan           break;
138c3dd00c7SVijay Mahadevan       }
139c3dd00c7SVijay Mahadevan       break;
140c3dd00c7SVijay Mahadevan   }
141c3dd00c7SVijay Mahadevan }
142c3dd00c7SVijay Mahadevan 
143c3dd00c7SVijay Mahadevan static void DMMoab_SetElementConnectivity_Private(PetscBool useSimplex, PetscInt dim, moab::EntityType etype, PetscInt *ecount, PetscInt vpere, PetscInt nele, PetscInt i, PetscInt j, PetscInt k, PetscInt vfirst, moab::EntityHandle *connectivity)
144c3dd00c7SVijay Mahadevan {
145c3dd00c7SVijay Mahadevan   PetscInt m,subelem;
146c3dd00c7SVijay Mahadevan   if (useSimplex) {
147a8e33bb7SVijay Mahadevan     subelem=(dim==1 ? 1 : (dim==2 ? 2 : 5));
148c3dd00c7SVijay Mahadevan     for (m=0; m<subelem; m++)
149c3dd00c7SVijay Mahadevan       DMMoab_SetSimplexElementConnectivity_Private(dim, m, etype, (*ecount+m)*vpere, nele, i, j, k, vfirst, connectivity);
150a8e33bb7SVijay Mahadevan 
151c3dd00c7SVijay Mahadevan   }
152c3dd00c7SVijay Mahadevan   else {
153c3dd00c7SVijay Mahadevan     subelem=1;
154c3dd00c7SVijay Mahadevan     DMMoab_SetTensorElementConnectivity_Private(dim, etype, (*ecount)*vpere, nele, i, j, k, vfirst, connectivity);
155c3dd00c7SVijay Mahadevan   }
156c3dd00c7SVijay Mahadevan   *ecount+=subelem;
157c3dd00c7SVijay Mahadevan }
158c3dd00c7SVijay Mahadevan 
159c3dd00c7SVijay Mahadevan 
160aa859218SVijay Mahadevan /*@
161aa859218SVijay Mahadevan   DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with user specified bounds.
162aa859218SVijay Mahadevan 
163aa859218SVijay Mahadevan   Collective on MPI_Comm
164aa859218SVijay Mahadevan 
165aa859218SVijay Mahadevan   Input Parameters:
166aa859218SVijay Mahadevan + comm - The communicator for the DM object
167aa859218SVijay Mahadevan . dim - The spatial dimension
168b8ecf6d3SVijay Mahadevan . bounds - The bounds of the box specified with [x-left, x-right, y-bottom, y-top, z-bottom, z-top] depending on the spatial dimension
169b8ecf6d3SVijay Mahadevan . nele - The number of discrete elements in each direction
170b8ecf6d3SVijay Mahadevan . user_nghost - The number of ghosted layers needed in the partitioned mesh
171aa859218SVijay Mahadevan 
172aa859218SVijay Mahadevan   Output Parameter:
173aa859218SVijay Mahadevan . dm  - The DM object
174aa859218SVijay Mahadevan 
175aa859218SVijay Mahadevan   Level: beginner
176aa859218SVijay Mahadevan 
177aa859218SVijay Mahadevan .keywords: DM, create
178aa859218SVijay Mahadevan .seealso: DMSetType(), DMCreate(), DMMoabLoadFromFile()
179aa859218SVijay Mahadevan @*/
180c3dd00c7SVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool useSimplex, const PetscReal* bounds, PetscInt nele, PetscInt user_nghost, DM *dm)
18151d15aeeSVijay Mahadevan {
18251d15aeeSVijay Mahadevan   PetscErrorCode  ierr;
18351d15aeeSVijay Mahadevan   moab::ErrorCode merr;
1843a4aeca1SVijay Mahadevan   PetscInt        i,j,k,n,nprocs;
18551d15aeeSVijay Mahadevan   DM_Moab        *dmmoab;
18651d15aeeSVijay Mahadevan   moab::Interface *mbiface;
18751d15aeeSVijay Mahadevan   moab::ParallelComm *pcomm;
188a4717116SVijay Mahadevan   moab::ReadUtilIface* readMeshIface;
189a4717116SVijay Mahadevan 
190*c528d872SBarry Smith   moab::Tag  id_tag=NULL;
191a4717116SVijay Mahadevan   moab::Range         ownedvtx,ownedelms;
1923a4aeca1SVijay Mahadevan   moab::EntityHandle  vfirst,efirst,regionset,faceset,edgeset,vtxset;
193a4717116SVijay Mahadevan   std::vector<double*> vcoords;
194a4717116SVijay Mahadevan   moab::EntityHandle  *connectivity = 0;
1955e168b13SVijay Mahadevan   moab::EntityType etype=moab::MBHEX;
196c8d5365dSVijay Mahadevan   PetscInt    ise[6];
19766f88a78SVijay Mahadevan   PetscReal   xse[6],defbounds[6];
1983a4aeca1SVijay Mahadevan   /* TODO: Fix nghost > 0 - now relying on exchange_ghost_cells */
1993a4aeca1SVijay Mahadevan   const PetscInt nghost=0;
2003a4aeca1SVijay Mahadevan 
2013a4aeca1SVijay Mahadevan   moab::Tag geom_tag;
2023a4aeca1SVijay Mahadevan 
2033a4aeca1SVijay Mahadevan   moab::Range adj,dim3,dim2;
2043a4aeca1SVijay Mahadevan   bool build_adjacencies=false;
20551d15aeeSVijay Mahadevan 
206a4717116SVijay Mahadevan   const PetscInt npts=nele+1;        /* Number of points in every dimension */
2075e168b13SVijay Mahadevan   PetscInt vpere=0,locnele=0,locnpts=0,ghnele,ghnpts;    /* Number of verts/element, vertices, elements owned by this process */
20851d15aeeSVijay Mahadevan 
20951d15aeeSVijay Mahadevan   PetscFunctionBegin;
210e5136372SVijay Mahadevan   if(dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid dimension argument for mesh: dim=[1,3].\n");
211e5136372SVijay Mahadevan 
212c8d5365dSVijay Mahadevan   ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr);
213e5136372SVijay Mahadevan   /* total number of vertices in all dimensions */
214e5136372SVijay Mahadevan   n=pow(npts,dim);
215e5136372SVijay Mahadevan 
216e5136372SVijay Mahadevan   /* do some error checking */
217e5136372SVijay Mahadevan   if(n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2.\n");
218e5136372SVijay Mahadevan   if(nprocs > n) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of processors must be less than or equal to number of elements.\n");
219e5136372SVijay Mahadevan   if(nghost < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of ghost layers cannot be negative.\n");
220e5136372SVijay Mahadevan 
221a4717116SVijay Mahadevan   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
222*c528d872SBarry Smith   ierr = DMMoabCreateMoab(comm, NULL, NULL, NULL, NULL, dm);CHKERRQ(ierr);
22351d15aeeSVijay Mahadevan 
224a4717116SVijay Mahadevan   /* get all the necessary handles from the private DM object */
22551d15aeeSVijay Mahadevan   dmmoab = (DM_Moab*)(*dm)->data;
22651d15aeeSVijay Mahadevan   mbiface = dmmoab->mbiface;
22751d15aeeSVijay Mahadevan   pcomm = dmmoab->pcomm;
22851d15aeeSVijay Mahadevan   id_tag = dmmoab->ltog_tag;
22951d15aeeSVijay Mahadevan   nprocs = pcomm->size();
230c8d5365dSVijay Mahadevan   dmmoab->dim = dim;
23151d15aeeSVijay Mahadevan 
232b5410836SVijay Mahadevan   /* create a file set to associate all entities in current mesh */
233b5410836SVijay Mahadevan   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr);
234b5410836SVijay Mahadevan 
235a4717116SVijay Mahadevan   /* No errors yet; proceed with building the mesh */
236a4717116SVijay Mahadevan   merr = mbiface->query_interface(readMeshIface);MBERRNM(merr);
23751d15aeeSVijay Mahadevan 
23851d15aeeSVijay Mahadevan   ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr);
239a4717116SVijay Mahadevan 
240a4717116SVijay Mahadevan   /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */
24151d15aeeSVijay Mahadevan   ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr);
24251d15aeeSVijay Mahadevan 
243a4717116SVijay Mahadevan   /* set some variables based on dimension */
24451d15aeeSVijay Mahadevan   switch(dim) {
24551d15aeeSVijay Mahadevan    case 1:
24651d15aeeSVijay Mahadevan     vpere = 2;
247a4717116SVijay Mahadevan     locnele = (ise[1]-ise[0]);
248a4717116SVijay Mahadevan     locnpts = (ise[1]-ise[0]+1);
249c8d5365dSVijay Mahadevan     ghnele = (nghost > 0 ? (ise[0] > nghost ? 1 : 0) + (ise[1] < nele - nghost ? 1 : 0) : 0 );
250c8d5365dSVijay Mahadevan     ghnpts = (nghost > 0 ? (ise[0] > 0 ? 1 : 0) + (ise[1] < nele ? 1 : 0) : 0);
25151d15aeeSVijay Mahadevan     etype = moab::MBEDGE;
25251d15aeeSVijay Mahadevan     break;
25351d15aeeSVijay Mahadevan    case 2:
254c3dd00c7SVijay Mahadevan     locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1);
255c3dd00c7SVijay Mahadevan     ghnpts = (nghost > 0 ? (ise[2] > 0 ? npts : 0) + (ise[3] < nele ? npts : 0) : 0);
256c3dd00c7SVijay Mahadevan     if (useSimplex) {
257c3dd00c7SVijay Mahadevan       vpere = 3;
258c3dd00c7SVijay Mahadevan       locnele = 2*(ise[1]-ise[0])*(ise[3]-ise[2]);
259c3dd00c7SVijay Mahadevan       ghnele = 2*(nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0);
260c3dd00c7SVijay Mahadevan       etype = moab::MBTRI;
261c3dd00c7SVijay Mahadevan     }
262c3dd00c7SVijay Mahadevan     else {
26351d15aeeSVijay Mahadevan       vpere = 4;
264a4717116SVijay Mahadevan       locnele = (ise[1]-ise[0])*(ise[3]-ise[2]);
265c8d5365dSVijay Mahadevan       ghnele = (nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0);
26651d15aeeSVijay Mahadevan       etype = moab::MBQUAD;
267c3dd00c7SVijay Mahadevan     }
26851d15aeeSVijay Mahadevan     break;
26951d15aeeSVijay Mahadevan    case 3:
2705e168b13SVijay Mahadevan    default:
271c3dd00c7SVijay Mahadevan     locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1);
272c3dd00c7SVijay Mahadevan     ghnpts = (nghost > 0 ? (ise[4] > 0 ? npts*npts : 0) + (ise[5] < nele ? npts*npts : 0) : 0);
273c3dd00c7SVijay Mahadevan     if (useSimplex) {
274c3dd00c7SVijay Mahadevan       vpere = 4;
275c3dd00c7SVijay Mahadevan       locnele = 5*(ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]);
276c3dd00c7SVijay Mahadevan       ghnele = 5*(nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0);
277c3dd00c7SVijay Mahadevan       etype = moab::MBTET;
278c3dd00c7SVijay Mahadevan     }
279c3dd00c7SVijay Mahadevan     else {
28051d15aeeSVijay Mahadevan       vpere = 8;
281a4717116SVijay Mahadevan       locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]);
282c8d5365dSVijay Mahadevan       ghnele = (nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0);
28351d15aeeSVijay Mahadevan       etype = moab::MBHEX;
284c3dd00c7SVijay Mahadevan     }
28551d15aeeSVijay Mahadevan     break;
28651d15aeeSVijay Mahadevan   }
28751d15aeeSVijay Mahadevan 
28851d15aeeSVijay Mahadevan   /* we have a domain of size [1,1,1] - now compute local co-ordinate box */
28951d15aeeSVijay Mahadevan   ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr);
290c8d5365dSVijay Mahadevan   for(i=0; i<6; ++i) {
29151d15aeeSVijay Mahadevan     xse[i]=(PetscReal)ise[i]/nele;
29251d15aeeSVijay Mahadevan   }
29351d15aeeSVijay Mahadevan 
294a4717116SVijay Mahadevan   /* Create vertexes and set the coodinate of each vertex */
295c8d5365dSVijay Mahadevan   merr = readMeshIface->get_node_coords(3,locnpts+ghnpts,0,vfirst,vcoords,n);MBERRNM(merr);
29651d15aeeSVijay Mahadevan 
297a4717116SVijay Mahadevan   /* Compute the co-ordinates of vertices and global IDs */
298c8d5365dSVijay Mahadevan   std::vector<int>    vgid(locnpts+ghnpts);
29951d15aeeSVijay Mahadevan   int vcount=0;
30066f88a78SVijay Mahadevan 
30166f88a78SVijay Mahadevan   if (!bounds) { /* default box mesh is defined on a unit-cube */
30266f88a78SVijay Mahadevan     defbounds[0]=0.0; defbounds[1]=1.0;
30366f88a78SVijay Mahadevan     defbounds[2]=0.0; defbounds[3]=1.0;
30466f88a78SVijay Mahadevan     defbounds[4]=0.0; defbounds[5]=1.0;
30566f88a78SVijay Mahadevan     bounds=defbounds;
30666f88a78SVijay Mahadevan   }
30766f88a78SVijay Mahadevan   else {
30866f88a78SVijay Mahadevan     /* validate the bounds data */
30966f88a78SVijay Mahadevan     if(bounds[0] >= bounds[1]) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"X-dim: Left boundary cannot be greater than right. [%G >= %G]\n",bounds[0],bounds[1]);
31066f88a78SVijay Mahadevan     if(dim > 1 && (bounds[2] >= bounds[3])) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Y-dim: Left boundary cannot be greater than right. [%G >= %G]\n",bounds[2],bounds[3]);
31166f88a78SVijay Mahadevan     if(dim > 2 && (bounds[4] >= bounds[5])) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Z-dim: Left boundary cannot be greater than right. [%G >= %G]\n",bounds[4],bounds[5]);
31266f88a78SVijay Mahadevan   }
31366f88a78SVijay Mahadevan 
31466f88a78SVijay Mahadevan   const double hx=(bounds[1]-bounds[0])/nele;
31566f88a78SVijay Mahadevan   const double hy=(dim > 1 ? (bounds[3]-bounds[2])/nele : 0.0);
31666f88a78SVijay Mahadevan   const double hz=(dim > 2 ? (bounds[5]-bounds[4])/nele : 0.0);
317e5136372SVijay Mahadevan 
318c8d5365dSVijay Mahadevan   /* create all the owned vertices */
319c8d5365dSVijay Mahadevan   for (k = ise[4]; k <= ise[5]; k++) {
320c8d5365dSVijay Mahadevan     for (j = ise[2]; j <= ise[3]; j++) {
321c8d5365dSVijay Mahadevan       for (i = ise[0]; i <= ise[1]; i++, vcount++) {
322aa859218SVijay Mahadevan         DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords);
323c8d5365dSVijay Mahadevan         vgid[vcount] = (k*npts+j)*npts+i+1;
324c8d5365dSVijay Mahadevan       }
325c8d5365dSVijay Mahadevan     }
326c8d5365dSVijay Mahadevan   }
327c8d5365dSVijay Mahadevan 
328e5136372SVijay Mahadevan   /* create ghosted vertices requested by user - below the current plane */
329e5136372SVijay Mahadevan   if (ise[2*dim-2] > 0) {
330e5136372SVijay Mahadevan     for (k = (dim==3?ise[4]-nghost:ise[4]); k <= (dim==3?ise[4]-1:ise[5]); k++) {
331e5136372SVijay Mahadevan       for (j = (dim==2?ise[2]-nghost:ise[2]); j <= (dim==2?ise[2]-1:ise[3]); j++) {
332e5136372SVijay Mahadevan         for (i = (dim>1?ise[0]:ise[0]-nghost); i <= (dim>1?ise[1]:ise[0]-1); i++, vcount++) {
333aa859218SVijay Mahadevan           DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords);
334e5136372SVijay Mahadevan           vgid[vcount] = (k*npts+j)*npts+i+1;
335e5136372SVijay Mahadevan         }
336e5136372SVijay Mahadevan       }
337e5136372SVijay Mahadevan     }
338e5136372SVijay Mahadevan   }
339e5136372SVijay Mahadevan 
340e5136372SVijay Mahadevan   /* create ghosted vertices requested by user - above the current plane */
341e5136372SVijay Mahadevan   if (ise[2*dim-1] < nele) {
342e5136372SVijay Mahadevan     for (k = (dim==3?ise[5]+1:ise[4]); k <= (dim==3?ise[5]+nghost:ise[5]); k++) {
343e5136372SVijay Mahadevan       for (j = (dim==2?ise[3]+1:ise[2]); j <= (dim==2?ise[3]+nghost:ise[3]); j++) {
344e5136372SVijay Mahadevan         for (i = (dim>1?ise[0]:ise[1]+1); i <= (dim>1?ise[1]:ise[1]+nghost); i++, vcount++) {
345aa859218SVijay Mahadevan           DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords);
346e5136372SVijay Mahadevan           vgid[vcount] = (k*npts+j)*npts+i+1;
347e5136372SVijay Mahadevan         }
34851d15aeeSVijay Mahadevan       }
34951d15aeeSVijay Mahadevan     }
35051d15aeeSVijay Mahadevan   }
35151d15aeeSVijay Mahadevan 
352c8d5365dSVijay Mahadevan   if (locnpts+ghnpts != vcount) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,vcount);
353c8d5365dSVijay Mahadevan 
35451d15aeeSVijay Mahadevan   merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr);
35551d15aeeSVijay Mahadevan 
356a4717116SVijay Mahadevan   /* The global ID tag is applied to each owned
357a4717116SVijay Mahadevan      vertex. It acts as an global identifier which MOAB uses to
358a4717116SVijay Mahadevan      assemble the individual pieces of the mesh:
359a4717116SVijay Mahadevan      Set the global ID indices */
36051d15aeeSVijay Mahadevan   merr = mbiface->tag_set_data(id_tag,ownedvtx,vgid.data());MBERRNM(merr);
36151d15aeeSVijay Mahadevan 
362a4717116SVijay Mahadevan   /* Create elements between mesh points using the ReadUtilInterface
363a4717116SVijay Mahadevan      get the reference to element connectivities for all local elements from the ReadUtilInterface */
364e5136372SVijay Mahadevan   merr = readMeshIface->get_element_connect (locnele+ghnele,vpere,etype,1,efirst,connectivity);MBERRNM(merr);
36551d15aeeSVijay Mahadevan 
366a4717116SVijay Mahadevan   /* offset appropriately so that only local ID and not global ID numbers are set for connectivity array */
367e5136372SVijay Mahadevan   vfirst-=vgid[0]-1;
36851d15aeeSVijay Mahadevan 
369a4717116SVijay Mahadevan    /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k):
370a4717116SVijay Mahadevan          and then set the connectivity for each element appropriately */
37151d15aeeSVijay Mahadevan   int ecount=0;
37251d15aeeSVijay Mahadevan 
373e5136372SVijay Mahadevan   /* create ghosted elements requested by user - below the current plane */
374e5136372SVijay Mahadevan   if (ise[2*dim-2] >= nghost) {
375e5136372SVijay Mahadevan     for (k = (dim==3?ise[4]-nghost:ise[4]); k < (dim==3?ise[4]:std::max(ise[5],1)); k++) {
376e5136372SVijay Mahadevan       for (j = (dim==2?ise[2]-nghost:ise[2]); j < (dim==2?ise[2]:std::max(ise[3],1)); j++) {
377c3dd00c7SVijay Mahadevan         for (i = (dim>1?ise[0]:ise[0]-nghost); i < (dim>1?std::max(ise[1],1):ise[0]); i++) {
378c3dd00c7SVijay Mahadevan           DMMoab_SetElementConnectivity_Private(useSimplex, dim, etype, &ecount, vpere, nele, i, j, k, vfirst, connectivity);
37951d15aeeSVijay Mahadevan         }
38051d15aeeSVijay Mahadevan       }
38151d15aeeSVijay Mahadevan     }
38251d15aeeSVijay Mahadevan   }
383e5136372SVijay Mahadevan 
384e5136372SVijay Mahadevan   /* create owned elements requested by user */
385e5136372SVijay Mahadevan   for (k = ise[4]; k < std::max(ise[5],1); k++) {
386e5136372SVijay Mahadevan     for (j = ise[2]; j < std::max(ise[3],1); j++) {
387c3dd00c7SVijay Mahadevan       for (i = ise[0]; i < std::max(ise[1],1); i++) {
388c3dd00c7SVijay Mahadevan         DMMoab_SetElementConnectivity_Private(useSimplex, dim, etype, &ecount, vpere, nele, i, j, k, vfirst, connectivity);
389e5136372SVijay Mahadevan       }
390e5136372SVijay Mahadevan     }
391e5136372SVijay Mahadevan   }
392e5136372SVijay Mahadevan 
393e5136372SVijay Mahadevan   /* create ghosted elements requested by user - above the current plane */
394e5136372SVijay Mahadevan   if (ise[2*dim-1] <= nele-nghost) {
395e5136372SVijay Mahadevan     for (k = (dim==3?ise[5]:ise[4]); k < (dim==3?ise[5]+nghost:std::max(ise[5],1)); k++) {
396e5136372SVijay Mahadevan       for (j = (dim==2?ise[3]:ise[2]); j < (dim==2?ise[3]+nghost:std::max(ise[3],1)); j++) {
397c3dd00c7SVijay Mahadevan         for (i = (dim>1?ise[0]:ise[1]); i < (dim>1?std::max(ise[1],1):ise[1]+nghost); i++) {
398c3dd00c7SVijay Mahadevan           DMMoab_SetElementConnectivity_Private(useSimplex, dim, etype, &ecount, vpere, nele, i, j, k, vfirst, connectivity);
399e5136372SVijay Mahadevan         }
400e5136372SVijay Mahadevan       }
401e5136372SVijay Mahadevan     }
402e5136372SVijay Mahadevan   }
403e5136372SVijay Mahadevan 
404e5136372SVijay Mahadevan   merr = readMeshIface->update_adjacencies(efirst,locnele+ghnele,vpere,connectivity);MBERRNM(merr);
40551d15aeeSVijay Mahadevan 
406a4717116SVijay Mahadevan   /* 2. Get the vertices and hexes from moab and check their numbers against I*J*K and (I-1)*(J-1)*(K-1), resp.
407a4717116SVijay Mahadevan         first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */
40851d15aeeSVijay Mahadevan   merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr);
409a4717116SVijay Mahadevan 
4106c4ed002SBarry Smith   if (locnele+ghnele != (int) ownedelms.size()) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele+ghnele,ownedelms.size());
4116c4ed002SBarry Smith   else if(locnpts+ghnpts != (int) ownedvtx.size()) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,ownedvtx.size());
4126c4ed002SBarry Smith   else {
4136c4ed002SBarry Smith     ierr = PetscInfo2(NULL, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size());CHKERRQ(ierr);
4146c4ed002SBarry Smith   }
41551d15aeeSVijay Mahadevan 
4163a4aeca1SVijay Mahadevan   /* lets create some sets */
4173a4aeca1SVijay Mahadevan   merr = mbiface->tag_get_handle(GEOM_DIMENSION_TAG_NAME, 1, moab::MB_TYPE_INTEGER, geom_tag, moab::MB_TAG_SPARSE|moab::MB_TAG_CREAT);MBERRNM(merr);
4183a4aeca1SVijay Mahadevan 
4193a4aeca1SVijay Mahadevan   merr = mbiface->create_meshset(moab::MESHSET_SET, regionset);MBERRNM(merr);
4203a4aeca1SVijay Mahadevan   merr = mbiface->add_entities(regionset, ownedelms);MBERRNM(merr);
4213a4aeca1SVijay Mahadevan   merr = mbiface->tag_set_data(geom_tag, &regionset, 1, &dmmoab->dim);MBERRNM(merr);
422b5410836SVijay Mahadevan   merr = mbiface->add_parent_child(dmmoab->fileset,regionset);MBERRNM(merr);
4233a4aeca1SVijay Mahadevan   merr = mbiface->unite_meshset(dmmoab->fileset, regionset);MBERRNM(merr);
4243a4aeca1SVijay Mahadevan 
4253a4aeca1SVijay Mahadevan   merr = mbiface->create_meshset(moab::MESHSET_SET, vtxset);MBERRNM(merr);
4263a4aeca1SVijay Mahadevan   merr = mbiface->add_entities(vtxset, ownedvtx);MBERRNM(merr);
4273a4aeca1SVijay Mahadevan   merr = mbiface->add_parent_child(dmmoab->fileset,vtxset);MBERRNM(merr);
4283a4aeca1SVijay Mahadevan   merr = mbiface->unite_meshset(dmmoab->fileset, vtxset);MBERRNM(merr);
4293a4aeca1SVijay Mahadevan 
4303a4aeca1SVijay Mahadevan   if (build_adjacencies) {
4313a4aeca1SVijay Mahadevan     // generate all lower dimensional adjacencies
4323a4aeca1SVijay Mahadevan     merr = mbiface->get_adjacencies( ownedelms, dim-1, true, adj, moab::Interface::UNION );MBERRNM(merr);
4333a4aeca1SVijay Mahadevan     merr = dmmoab->pcomm->get_part_entities(dim2, dim-1);MBERRNM(merr);
4343a4aeca1SVijay Mahadevan     adj.merge(dim2);
4353a4aeca1SVijay Mahadevan 
4363a4aeca1SVijay Mahadevan     /* create face sets */
4373a4aeca1SVijay Mahadevan     merr = mbiface->create_meshset(moab::MESHSET_SET, faceset);MBERRNM(merr);
4383a4aeca1SVijay Mahadevan     merr = mbiface->add_entities(faceset, adj);MBERRNM(merr);
4393a4aeca1SVijay Mahadevan     merr = mbiface->add_parent_child(dmmoab->fileset,faceset);MBERRNM(merr);
4403a4aeca1SVijay Mahadevan     i=dim-1;
4413a4aeca1SVijay Mahadevan     merr = mbiface->tag_set_data(geom_tag, &faceset, 1, &i);MBERRNM(merr);
4423a4aeca1SVijay Mahadevan     merr = mbiface->unite_meshset(dmmoab->fileset, faceset);MBERRNM(merr);
443b8ecf6d3SVijay Mahadevan     PetscInfo2(NULL, "Found %d %d-Dim quantities.\n", adj.size(), dim-1);
4443a4aeca1SVijay Mahadevan 
4453a4aeca1SVijay Mahadevan     if (dim > 2) {
4463a4aeca1SVijay Mahadevan       dim2.clear();
4473a4aeca1SVijay Mahadevan       /* create edge sets, if appropriate i.e., if dim=3 */
4483a4aeca1SVijay Mahadevan       merr = mbiface->create_meshset(moab::MESHSET_SET, edgeset);MBERRNM(merr);
4493a4aeca1SVijay Mahadevan       merr = mbiface->get_adjacencies(adj, dim-1, true, dim2, moab::Interface::UNION );MBERRNM(merr);
4503a4aeca1SVijay Mahadevan       merr = mbiface->add_entities(edgeset, dim2);MBERRNM(merr);
4513a4aeca1SVijay Mahadevan       merr = mbiface->add_parent_child(dmmoab->fileset,edgeset);MBERRNM(merr);
4523a4aeca1SVijay Mahadevan       i=dim-2;
4533a4aeca1SVijay Mahadevan       merr = mbiface->tag_set_data(geom_tag, &edgeset, 1, &i);MBERRNM(merr);
4543a4aeca1SVijay Mahadevan       merr = mbiface->unite_meshset(dmmoab->fileset, edgeset);MBERRNM(merr);
455b8ecf6d3SVijay Mahadevan       PetscInfo2(NULL, "Found %d %d-Dim quantities.\n", adj.size(), dim-2);
4563a4aeca1SVijay Mahadevan     }
4573a4aeca1SVijay Mahadevan   }
4583a4aeca1SVijay Mahadevan 
459a4717116SVijay Mahadevan   /* check the handles */
460a4717116SVijay Mahadevan   merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr);
46151d15aeeSVijay Mahadevan 
462a4717116SVijay Mahadevan   /* resolve the shared entities by exchanging information to adjacent processors */
463a4717116SVijay Mahadevan   merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr);
464ce1fd009SVijay Mahadevan   merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,NULL,&id_tag);MBERRV(mbiface,merr);
46551d15aeeSVijay Mahadevan 
4663a4aeca1SVijay Mahadevan   merr = pcomm->exchange_ghost_cells(dim,0,user_nghost,dim,true,false,&dmmoab->fileset);MBERRV(mbiface,merr);
467c8d5365dSVijay Mahadevan 
468e427d9c9SVijay Mahadevan   /* Reassign global IDs on all entities. */
469e427d9c9SVijay Mahadevan   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,false,true,false);MBERRNM(merr);
470e427d9c9SVijay Mahadevan 
471a4717116SVijay Mahadevan   /* Everything is set up, now just do a tag exchange to update tags
472a4717116SVijay Mahadevan      on all of the ghost vertexes */
473c8d5365dSVijay Mahadevan   merr = mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr);
474c8d5365dSVijay Mahadevan   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, ownedelms);MBERRNM(merr);
475c8d5365dSVijay Mahadevan 
476a4717116SVijay Mahadevan   merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr);
477a4717116SVijay Mahadevan   merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr);
47851d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
47951d15aeeSVijay Mahadevan }
48051d15aeeSVijay Mahadevan 
48151d15aeeSVijay Mahadevan 
4822e4e7c01SVijay Mahadevan PetscErrorCode DMMoab_GetReadOptions_Private(PetscBool by_rank, PetscInt numproc, PetscInt dim, MoabReadMode mode, PetscInt dbglevel, const char* dm_opts, const char* extra_opts, const char** read_opts)
48351d15aeeSVijay Mahadevan {
4846acfe860SVijay Mahadevan   char           *ropts;
48561eb6e27SVijay Mahadevan   char           ropts_par[PETSC_MAX_PATH_LEN];
48661eb6e27SVijay Mahadevan   char           ropts_dbg[PETSC_MAX_PATH_LEN];
487f90c3b0eSVijay Mahadevan   PetscErrorCode ierr;
48851d15aeeSVijay Mahadevan 
489a4717116SVijay Mahadevan   PetscFunctionBegin;
49095dccacaSBarry Smith   ierr = PetscMalloc1(PETSC_MAX_PATH_LEN,&ropts);CHKERRQ(ierr);
49161eb6e27SVijay Mahadevan   ierr = PetscMemzero(&ropts_par,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
49261eb6e27SVijay Mahadevan   ierr = PetscMemzero(&ropts_dbg,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
49361eb6e27SVijay Mahadevan 
494e23c60ebSVijay Mahadevan   /* do parallel read unless using only one processor */
495a4717116SVijay Mahadevan   if (numproc > 1) {
4966acfe860SVijay Mahadevan     ierr = PetscSNPrintf(ropts_par, PETSC_MAX_PATH_LEN, "PARALLEL=%s;PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=%d.0.1%s;",MoabReadModes[mode],dim,(by_rank ? ";PARTITION_BY_RANK":""));CHKERRQ(ierr);
4972e4e7c01SVijay Mahadevan   }
4982e4e7c01SVijay Mahadevan 
499c8d5365dSVijay Mahadevan   if (dbglevel) {
500f90c3b0eSVijay Mahadevan     if (numproc>1) {
5016acfe860SVijay Mahadevan       ierr = PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "%sCPUTIME;DEBUG_IO=%d;DEBUG_PIO=%d;",dbglevel,dbglevel);CHKERRQ(ierr);
50261eb6e27SVijay Mahadevan     }
50361eb6e27SVijay Mahadevan     else {
5046acfe860SVijay Mahadevan       ierr = PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "%sCPUTIME;DEBUG_IO=%d;",dbglevel);CHKERRQ(ierr);
505f90c3b0eSVijay Mahadevan     }
506c8d5365dSVijay Mahadevan   }
50751d15aeeSVijay Mahadevan 
5086acfe860SVijay Mahadevan   ierr = PetscSNPrintf(ropts, PETSC_MAX_PATH_LEN, "%s%s%s%s",ropts_par,ropts_dbg,(extra_opts?extra_opts:""),(dm_opts?dm_opts:""));CHKERRQ(ierr);
509f90c3b0eSVijay Mahadevan   *read_opts = ropts;
51051d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
51151d15aeeSVijay Mahadevan }
51251d15aeeSVijay Mahadevan 
51351d15aeeSVijay Mahadevan 
514aa859218SVijay Mahadevan /*@
515aa859218SVijay Mahadevan   DMMoabLoadFromFile - Creates a DM object by loading the mesh from a user specified file.
516aa859218SVijay Mahadevan 
517aa859218SVijay Mahadevan   Collective on MPI_Comm
518aa859218SVijay Mahadevan 
519aa859218SVijay Mahadevan   Input Parameters:
520aa859218SVijay Mahadevan + comm - The communicator for the DM object
521aa859218SVijay Mahadevan . dim - The spatial dimension
522b8ecf6d3SVijay Mahadevan . filename - The name of the mesh file to be loaded
523b8ecf6d3SVijay Mahadevan . usrreadopts - The options string to read a MOAB mesh.
524aa859218SVijay Mahadevan   Reference (Parallel Mesh Initialization: http://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo)
525aa859218SVijay Mahadevan 
526aa859218SVijay Mahadevan   Output Parameter:
527aa859218SVijay Mahadevan . dm  - The DM object
528aa859218SVijay Mahadevan 
529aa859218SVijay Mahadevan   Level: beginner
530aa859218SVijay Mahadevan 
531aa859218SVijay Mahadevan .keywords: DM, create
532b8ecf6d3SVijay Mahadevan 
533aa859218SVijay Mahadevan .seealso: DMSetType(), DMCreate(), DMMoabCreateBoxMesh()
534aa859218SVijay Mahadevan @*/
535a4717116SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm)
53651d15aeeSVijay Mahadevan {
537a4717116SVijay Mahadevan   moab::ErrorCode merr;
5382e4e7c01SVijay Mahadevan   PetscInt        nprocs;
539a4717116SVijay Mahadevan   DM_Moab        *dmmoab;
540a4717116SVijay Mahadevan   moab::Interface *mbiface;
541a4717116SVijay Mahadevan   moab::ParallelComm *pcomm;
542a4717116SVijay Mahadevan   moab::Range verts,elems;
5437023aa44SVijay Mahadevan   const char *readopts;
544a4717116SVijay Mahadevan   PetscErrorCode ierr;
54551d15aeeSVijay Mahadevan 
54651d15aeeSVijay Mahadevan   PetscFunctionBegin;
547a4717116SVijay Mahadevan   PetscValidPointer(dm,5);
54851d15aeeSVijay Mahadevan 
549a4717116SVijay Mahadevan   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
550*c528d872SBarry Smith   ierr = DMMoabCreateMoab(comm, NULL, NULL, NULL, NULL, dm);CHKERRQ(ierr);
55151d15aeeSVijay Mahadevan 
552a4717116SVijay Mahadevan   /* get all the necessary handles from the private DM object */
553a4717116SVijay Mahadevan   dmmoab = (DM_Moab*)(*dm)->data;
554a4717116SVijay Mahadevan   mbiface = dmmoab->mbiface;
555a4717116SVijay Mahadevan   pcomm = dmmoab->pcomm;
556a4717116SVijay Mahadevan   nprocs = pcomm->size();
557aa859218SVijay Mahadevan   /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */
558c8d5365dSVijay Mahadevan   dmmoab->dim = dim;
559a4717116SVijay Mahadevan 
560b5410836SVijay Mahadevan   /* create a file set to associate all entities in current mesh */
561b5410836SVijay Mahadevan   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr);
562b5410836SVijay Mahadevan 
563a4717116SVijay Mahadevan   /* add mesh loading options specific to the DM */
5642e4e7c01SVijay Mahadevan   ierr = DMMoab_GetReadOptions_Private(dmmoab->partition_by_rank, nprocs, dim, dmmoab->read_mode,
5652e4e7c01SVijay Mahadevan                                         dmmoab->rw_dbglevel, dmmoab->extra_read_options, usrreadopts, &readopts);CHKERRQ(ierr);
5667023aa44SVijay Mahadevan 
567e5136372SVijay Mahadevan   PetscInfo2(*dm, "Reading file %s with options: %s\n",filename,readopts);
568a4717116SVijay Mahadevan 
569a4717116SVijay Mahadevan   /* Load the mesh from a file. */
5707023aa44SVijay Mahadevan   merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr);
571a4717116SVijay Mahadevan 
5726e40195eSVijay Mahadevan   /* Reassign global IDs on all entities. */
573e5136372SVijay Mahadevan   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr);
574e5136372SVijay Mahadevan 
575e5136372SVijay Mahadevan   /* load the local vertices */
576e5136372SVijay Mahadevan   merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr);
577e5136372SVijay Mahadevan   /* load the local elements */
578e5136372SVijay Mahadevan   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr);
579e5136372SVijay Mahadevan 
580e5136372SVijay Mahadevan   /* Everything is set up, now just do a tag exchange to update tags
581e5136372SVijay Mahadevan      on all of the ghost vertexes */
582e5136372SVijay Mahadevan   merr = pcomm->exchange_tags(dmmoab->ltog_tag,verts);MBERRV(mbiface,merr);
583e5136372SVijay Mahadevan   merr = pcomm->exchange_tags(dmmoab->ltog_tag,elems);MBERRV(mbiface,merr);
584e5136372SVijay Mahadevan 
585e5136372SVijay Mahadevan   merr = pcomm->exchange_ghost_cells(dim,0,1,0,true,true,&dmmoab->fileset);MBERRV(mbiface,merr);
5866e40195eSVijay Mahadevan 
587a4717116SVijay Mahadevan   merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr);
588a4717116SVijay Mahadevan 
589e5136372SVijay Mahadevan   PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size());
5906acfe860SVijay Mahadevan   ierr = PetscFree(readopts);CHKERRQ(ierr);
59151d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
59251d15aeeSVijay Mahadevan }
59351d15aeeSVijay Mahadevan 
594