xref: /petsc/src/dm/impls/moab/dmmbutil.cxx (revision 6e40195ecc592446ee0fa1655e2767ca3348244e)
151d15aeeSVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I  "petscdm.h"   I*/
251d15aeeSVijay Mahadevan #include <petsc-private/vecimpl.h> /*I  "petscdm.h"   I*/
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 
1151d15aeeSVijay Mahadevan #undef __FUNCT__
1251d15aeeSVijay Mahadevan #define __FUNCT__ "DMMoabComputeDomainBounds_Private"
1351d15aeeSVijay Mahadevan PetscErrorCode DMMoabComputeDomainBounds_Private(moab::ParallelComm* pcomm, PetscInt dim, PetscInt neleglob, PetscInt *ise)
1451d15aeeSVijay Mahadevan {
1551d15aeeSVijay Mahadevan   PetscInt size,rank;
1651d15aeeSVijay Mahadevan   PetscInt fraction,remainder;
1751d15aeeSVijay Mahadevan   PetscInt starts[3],sizes[3];
1851d15aeeSVijay Mahadevan 
1951d15aeeSVijay Mahadevan   PetscFunctionBegin;
2051d15aeeSVijay Mahadevan   if(dim<1 && dim>3) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The problem dimension is invalid: %D",dim);
2151d15aeeSVijay Mahadevan 
2251d15aeeSVijay Mahadevan   size=pcomm->size();
2351d15aeeSVijay Mahadevan   rank=pcomm->rank();
2451d15aeeSVijay Mahadevan   fraction=neleglob/size;    /* partition only by the largest dimension */
2551d15aeeSVijay Mahadevan   remainder=neleglob%size;   /* remainder after partition which gets evenly distributed by round-robin */
2651d15aeeSVijay Mahadevan 
2751d15aeeSVijay 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);
2851d15aeeSVijay Mahadevan 
2951d15aeeSVijay Mahadevan   starts[0]=starts[1]=starts[2]=0;       /* default dimensional element = 1 */
3051d15aeeSVijay Mahadevan   sizes[0]=sizes[1]=sizes[2]=neleglob;   /* default dimensional element = 1 */
3151d15aeeSVijay Mahadevan 
3251d15aeeSVijay Mahadevan   if(rank < remainder) {
3351d15aeeSVijay Mahadevan     /* This process gets "fraction+1" elements */
3451d15aeeSVijay Mahadevan     sizes[dim-1] = (fraction + 1);
3551d15aeeSVijay Mahadevan     starts[dim-1] = rank * (fraction+1);
3651d15aeeSVijay Mahadevan   } else {
3751d15aeeSVijay Mahadevan     /* This process gets "fraction" elements */
3851d15aeeSVijay Mahadevan     sizes[dim-1] = fraction;
3951d15aeeSVijay Mahadevan     starts[dim-1] = (remainder*(fraction+1) + fraction*(rank-remainder));
4051d15aeeSVijay Mahadevan   }
4151d15aeeSVijay Mahadevan 
4251d15aeeSVijay Mahadevan   for(int i=dim-1; i>=0; --i) {
4351d15aeeSVijay Mahadevan     ise[2*i]=starts[i];ise[2*i+1]=starts[i]+sizes[i];
4451d15aeeSVijay Mahadevan   }
4551d15aeeSVijay Mahadevan 
4651d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
4751d15aeeSVijay Mahadevan }
4851d15aeeSVijay Mahadevan 
4951d15aeeSVijay Mahadevan 
5051d15aeeSVijay Mahadevan #undef __FUNCT__
5151d15aeeSVijay Mahadevan #define __FUNCT__ "DMMoabCreateBoxMesh"
5251d15aeeSVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt nele, PetscInt nghost, DM *dm)
5351d15aeeSVijay Mahadevan {
5451d15aeeSVijay Mahadevan   PetscErrorCode ierr;
5551d15aeeSVijay Mahadevan   moab::ErrorCode merr;
56a4717116SVijay Mahadevan   PetscInt            nprocs;
5751d15aeeSVijay Mahadevan   DM_Moab        *dmmoab;
5851d15aeeSVijay Mahadevan   moab::Interface *mbiface;
5951d15aeeSVijay Mahadevan   moab::ParallelComm *pcomm;
60a4717116SVijay Mahadevan   moab::ReadUtilIface* readMeshIface;
61a4717116SVijay Mahadevan 
6251d15aeeSVijay Mahadevan   moab::Tag  id_tag=PETSC_NULL;
63a4717116SVijay Mahadevan   moab::Range         ownedvtx,ownedelms;
64a4717116SVijay Mahadevan   moab::EntityHandle  vfirst,efirst;
65a4717116SVijay Mahadevan   std::vector<double*> vcoords;
66a4717116SVijay Mahadevan   moab::EntityHandle  *connectivity = 0;
67a4717116SVijay Mahadevan   std::vector<int>    subent_conn;
6851d15aeeSVijay Mahadevan   moab::EntityType etype;
6951d15aeeSVijay Mahadevan   PetscInt    ise[6];
7051d15aeeSVijay Mahadevan   PetscReal   xse[6];
7151d15aeeSVijay Mahadevan 
72a4717116SVijay Mahadevan   const PetscInt npts=nele+1;        /* Number of points in every dimension */
73a4717116SVijay Mahadevan   PetscInt vpere,locnele,locnpts;    /* Number of verts/element, vertices, elements owned by this process */
7451d15aeeSVijay Mahadevan 
7551d15aeeSVijay Mahadevan   PetscFunctionBegin;
76a4717116SVijay Mahadevan   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
7751d15aeeSVijay Mahadevan   ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr);
7851d15aeeSVijay Mahadevan 
79a4717116SVijay Mahadevan   /* get all the necessary handles from the private DM object */
8051d15aeeSVijay Mahadevan   dmmoab = (DM_Moab*)(*dm)->data;
8151d15aeeSVijay Mahadevan   mbiface = dmmoab->mbiface;
8251d15aeeSVijay Mahadevan   pcomm = dmmoab->pcomm;
8351d15aeeSVijay Mahadevan   id_tag = dmmoab->ltog_tag;
8451d15aeeSVijay Mahadevan   nprocs = pcomm->size();
8551d15aeeSVijay Mahadevan 
86a4717116SVijay Mahadevan   /* do some error checking */
87a4717116SVijay Mahadevan   if(pow(npts,dim) < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2");
88a4717116SVijay Mahadevan   if(nprocs > pow(nele,dim)) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of processors must be less than or equal to number of elements");
8951d15aeeSVijay Mahadevan 
90a4717116SVijay Mahadevan   /* No errors yet; proceed with building the mesh */
91a4717116SVijay Mahadevan   merr = mbiface->query_interface(readMeshIface);MBERRNM(merr);
9251d15aeeSVijay Mahadevan 
9351d15aeeSVijay Mahadevan   ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr);
94a4717116SVijay Mahadevan 
95a4717116SVijay Mahadevan   /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */
9651d15aeeSVijay Mahadevan   ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr);
9751d15aeeSVijay Mahadevan 
98a4717116SVijay Mahadevan   /* set some variables based on dimension */
9951d15aeeSVijay Mahadevan   switch(dim) {
10051d15aeeSVijay Mahadevan    case 1:
10151d15aeeSVijay Mahadevan     vpere = 2;
102a4717116SVijay Mahadevan     locnele = (ise[1]-ise[0]);
103a4717116SVijay Mahadevan     locnpts = (ise[1]-ise[0]+1);
10451d15aeeSVijay Mahadevan     etype = moab::MBEDGE;
10551d15aeeSVijay Mahadevan     break;
10651d15aeeSVijay Mahadevan    case 2:
10751d15aeeSVijay Mahadevan     vpere = 4;
108a4717116SVijay Mahadevan     locnele = (ise[1]-ise[0])*(ise[3]-ise[2]);
109a4717116SVijay Mahadevan     locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1);
11051d15aeeSVijay Mahadevan     etype = moab::MBQUAD;
11151d15aeeSVijay Mahadevan     break;
11251d15aeeSVijay Mahadevan    case 3:
11351d15aeeSVijay Mahadevan     vpere = 8;
114a4717116SVijay Mahadevan     locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]);
115a4717116SVijay Mahadevan     locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1);
11651d15aeeSVijay Mahadevan     etype = moab::MBHEX;
11751d15aeeSVijay Mahadevan     break;
11851d15aeeSVijay Mahadevan   }
11951d15aeeSVijay Mahadevan 
12051d15aeeSVijay Mahadevan   /* we have a domain of size [1,1,1] - now compute local co-ordinate box */
12151d15aeeSVijay Mahadevan   ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr);
12251d15aeeSVijay Mahadevan   for(int i=0; i<6; ++i) {
12351d15aeeSVijay Mahadevan     xse[i]=(PetscReal)ise[i]/nele;
12451d15aeeSVijay Mahadevan   }
12551d15aeeSVijay Mahadevan 
126a4717116SVijay Mahadevan   /* Create vertexes and set the coodinate of each vertex */
127a4717116SVijay Mahadevan   const int sequence_size = (locnele + vpere) + 1;
128a4717116SVijay Mahadevan   merr = readMeshIface->get_node_coords(3,locnpts,0,vfirst,vcoords,sequence_size);MBERRNM(merr);
12951d15aeeSVijay Mahadevan 
130a4717116SVijay Mahadevan   /* Compute the co-ordinates of vertices and global IDs */
131a4717116SVijay Mahadevan   std::vector<int>    vgid(locnpts);
13251d15aeeSVijay Mahadevan   int vcount=0;
13351d15aeeSVijay Mahadevan   double hxyz=1.0/nele;
13451d15aeeSVijay Mahadevan   for (int k = ise[4]; k <= ise[5]; k++) {
13551d15aeeSVijay Mahadevan     for (int j = ise[2]; j <= ise[3]; j++) {
13651d15aeeSVijay Mahadevan       for (int i = ise[0]; i <= ise[1]; i++, vcount++) {
13751d15aeeSVijay Mahadevan         vcoords[0][vcount] = i*hxyz;
13851d15aeeSVijay Mahadevan         vcoords[1][vcount] = j*hxyz;
13951d15aeeSVijay Mahadevan         vcoords[2][vcount] = k*hxyz;
14051d15aeeSVijay Mahadevan         vgid[vcount] = (k*nele+j)*nele+i;
14151d15aeeSVijay Mahadevan       }
14251d15aeeSVijay Mahadevan     }
14351d15aeeSVijay Mahadevan   }
14451d15aeeSVijay Mahadevan 
14551d15aeeSVijay Mahadevan   merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr);
146a4717116SVijay Mahadevan   merr = mbiface->add_entities (dmmoab->fileset, ownedvtx);MBERRNM(merr);
14751d15aeeSVijay Mahadevan 
148a4717116SVijay Mahadevan   /* The global ID tag is applied to each owned
149a4717116SVijay Mahadevan      vertex. It acts as an global identifier which MOAB uses to
150a4717116SVijay Mahadevan      assemble the individual pieces of the mesh:
151a4717116SVijay Mahadevan      Set the global ID indices */
15251d15aeeSVijay Mahadevan   merr = mbiface->tag_set_data(id_tag,ownedvtx,vgid.data());MBERRNM(merr);
15351d15aeeSVijay Mahadevan 
154a4717116SVijay Mahadevan   /* Create elements between mesh points using the ReadUtilInterface
155a4717116SVijay Mahadevan      get the reference to element connectivities for all local elements from the ReadUtilInterface */
156a4717116SVijay Mahadevan   merr = readMeshIface->get_element_connect (locnele,vpere,etype,1,efirst,connectivity);MBERRNM(merr);
15751d15aeeSVijay Mahadevan 
158a4717116SVijay Mahadevan   /* offset appropriately so that only local ID and not global ID numbers are set for connectivity array */
159a4717116SVijay Mahadevan   vfirst-=vgid[0];
16051d15aeeSVijay Mahadevan 
161a4717116SVijay Mahadevan    /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k):
162a4717116SVijay Mahadevan          and then set the connectivity for each element appropriately */
16351d15aeeSVijay Mahadevan   int ecount=0;
164a4717116SVijay Mahadevan   subent_conn.resize(vpere);
16551d15aeeSVijay Mahadevan   for (int k = ise[4]; k < std::max(ise[5],1); k++) {
16651d15aeeSVijay Mahadevan     for (int j = ise[2]; j < std::max(ise[3],1); j++) {
16751d15aeeSVijay Mahadevan       for (int i = ise[0]; i < std::max(ise[1],1); i++,ecount++) {
16851d15aeeSVijay Mahadevan         const int offset = ecount*vpere;
16951d15aeeSVijay Mahadevan         moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data());
17051d15aeeSVijay Mahadevan 
17151d15aeeSVijay Mahadevan         switch(dim) {
17251d15aeeSVijay Mahadevan           case 1:
17351d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[0]] = vfirst+i;
17451d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[1]] = vfirst+(i+1);
175a4717116SVijay Mahadevan             PetscPrintf(PETSC_COMM_WORLD, "ELEMENT[%D,%D,%D]: CONNECTIVITY = %D, %D\n", i,j,k,connectivity[offset+subent_conn[0]], connectivity[offset+subent_conn[1]]);
17651d15aeeSVijay Mahadevan             break;
17751d15aeeSVijay Mahadevan           case 2:
17851d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1);
17951d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1);
18051d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1);
18151d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[3]] = vfirst+i+(j+1)*(nele+1);
182a4717116SVijay Mahadevan             PetscPrintf(PETSC_COMM_WORLD, "ELEMENT[%D,%D,%D]: CONNECTIVITY = %D, %D, %D, %D\n", i,j,k,connectivity[offset+subent_conn[0]], connectivity[offset+subent_conn[1]], connectivity[offset+subent_conn[2]], connectivity[offset+subent_conn[3]]);
18351d15aeeSVijay Mahadevan             break;
18451d15aeeSVijay Mahadevan           case 3:
18551d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
18651d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k);
18751d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
18851d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k);
18951d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[4]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1));
19051d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[5]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
19151d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[6]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1));
19251d15aeeSVijay Mahadevan             connectivity[offset+subent_conn[7]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
19351d15aeeSVijay Mahadevan             break;
19451d15aeeSVijay Mahadevan         }
19551d15aeeSVijay Mahadevan       }
19651d15aeeSVijay Mahadevan     }
19751d15aeeSVijay Mahadevan   }
198a4717116SVijay Mahadevan   merr = readMeshIface->update_adjacencies(efirst,locnele,vpere,connectivity);MBERRNM(merr);
19951d15aeeSVijay Mahadevan 
200a4717116SVijay 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.
201a4717116SVijay Mahadevan         first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */
20251d15aeeSVijay Mahadevan   merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr);
203a4717116SVijay Mahadevan   merr = mbiface->add_entities (dmmoab->fileset, ownedelms);MBERRNM(merr);
20451d15aeeSVijay Mahadevan 
205a4717116SVijay Mahadevan //  merr = mbiface->unite_meshset(dmmoab->fileset, 0);MBERRV(mbiface,merr);
206a4717116SVijay Mahadevan 
207a4717116SVijay Mahadevan   if (locnele != (int) ownedelms.size())
208a4717116SVijay Mahadevan     SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele,ownedelms.size());
209a4717116SVijay Mahadevan   else if(locnpts != (int) ownedvtx.size())
210a4717116SVijay Mahadevan     SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts,ownedvtx.size());
21151d15aeeSVijay Mahadevan   else
212a4717116SVijay Mahadevan     PetscInfo2(*dm, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size());
21351d15aeeSVijay Mahadevan 
214a4717116SVijay Mahadevan   /* check the handles */
215a4717116SVijay Mahadevan   merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr);
21651d15aeeSVijay Mahadevan 
217a4717116SVijay Mahadevan   /* resolve the shared entities by exchanging information to adjacent processors */
218a4717116SVijay Mahadevan   merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr);
219a4717116SVijay Mahadevan   merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,&id_tag);MBERRV(mbiface,merr);
22051d15aeeSVijay Mahadevan 
221a4717116SVijay Mahadevan   /* Reassign global IDs on all entities. */
222a4717116SVijay Mahadevan   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,false,true);MBERRNM(merr);
223a4717116SVijay Mahadevan   merr = pcomm->exchange_ghost_cells(dim,0,nghost,0,true);MBERRV(mbiface,merr);
22451d15aeeSVijay Mahadevan 
225a4717116SVijay Mahadevan   /* Everything is set up, now just do a tag exchange to update tags
226a4717116SVijay Mahadevan      on all of the ghost vertexes */
227a4717116SVijay Mahadevan   merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr);
228a4717116SVijay Mahadevan   merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr);
22951d15aeeSVijay Mahadevan 
230a4717116SVijay Mahadevan //  ierr = DMMoabSetLocalVertices(*dm, &ownedvtx);CHKERRQ(ierr);
231a4717116SVijay Mahadevan //  ierr = DMMoabSetLocalElements(*dm, &ownedelms);CHKERRQ(ierr);
23251d15aeeSVijay Mahadevan 
23351d15aeeSVijay Mahadevan   merr = mbiface->set_dimension(dim);MBERRNM(merr);
23451d15aeeSVijay Mahadevan 
23551d15aeeSVijay Mahadevan   std::stringstream sstr;
236a4717116SVijay Mahadevan   sstr << "test_" << pcomm->rank() << ".vtk";
23751d15aeeSVijay Mahadevan   mbiface->write_mesh(sstr.str().c_str());
23851d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
23951d15aeeSVijay Mahadevan }
24051d15aeeSVijay Mahadevan 
24151d15aeeSVijay Mahadevan 
242a4717116SVijay Mahadevan #undef __FUNCT__
243a4717116SVijay Mahadevan #define __FUNCT__ "DMMoab_GetReadOptions_Private"
244a4717116SVijay Mahadevan PetscErrorCode DMMoab_GetReadOptions_Private(PetscBool by_rank, PetscInt numproc, PetscInt dim, MoabReadMode mode, PetscInt dbglevel, const char* extra_opts, const char* read_opts)
24551d15aeeSVijay Mahadevan {
246a4717116SVijay Mahadevan   std::ostringstream str;
24751d15aeeSVijay Mahadevan 
248a4717116SVijay Mahadevan   PetscFunctionBegin;
249a4717116SVijay Mahadevan   // do parallel read unless only one processor
250a4717116SVijay Mahadevan   if (numproc > 1) {
251a4717116SVijay Mahadevan     str << "PARALLEL=" << mode << ";PARTITION=PARALLEL_PARTITION;";
252a4717116SVijay Mahadevan     str << "PARTITION_DISTRIBUTE;PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=" << dim << ".0.1;";
253a4717116SVijay Mahadevan     if (by_rank)
254a4717116SVijay Mahadevan       str << "PARTITION_BY_RANK;";
255a4717116SVijay Mahadevan   }
25651d15aeeSVijay Mahadevan 
257a4717116SVijay Mahadevan   if (extra_opts)
258a4717116SVijay Mahadevan     str << extra_opts << ";";
25951d15aeeSVijay Mahadevan 
260a4717116SVijay Mahadevan   if (dbglevel)
261a4717116SVijay Mahadevan     str << "DEBUG_IO=" << dbglevel << ";CPUTIME;";
26251d15aeeSVijay Mahadevan 
263a4717116SVijay Mahadevan   read_opts = str.str().c_str();
26451d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
26551d15aeeSVijay Mahadevan }
26651d15aeeSVijay Mahadevan 
26751d15aeeSVijay Mahadevan 
26851d15aeeSVijay Mahadevan #undef __FUNCT__
269a4717116SVijay Mahadevan #define __FUNCT__ "DMMoabLoadFromFile"
270a4717116SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm)
27151d15aeeSVijay Mahadevan {
272a4717116SVijay Mahadevan   moab::ErrorCode merr;
273a4717116SVijay Mahadevan   PetscInt            nprocs;
274a4717116SVijay Mahadevan   DM_Moab        *dmmoab;
275a4717116SVijay Mahadevan   moab::Interface *mbiface;
276a4717116SVijay Mahadevan   moab::ParallelComm *pcomm;
277a4717116SVijay Mahadevan   moab::Range verts,elems;
278a4717116SVijay Mahadevan   const char* readopts=0;
279a4717116SVijay Mahadevan   PetscErrorCode ierr;
28051d15aeeSVijay Mahadevan 
28151d15aeeSVijay Mahadevan   PetscFunctionBegin;
282a4717116SVijay Mahadevan   PetscValidPointer(dm,5);
28351d15aeeSVijay Mahadevan 
284a4717116SVijay Mahadevan   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
285a4717116SVijay Mahadevan   ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr);
28651d15aeeSVijay Mahadevan 
287a4717116SVijay Mahadevan   /* get all the necessary handles from the private DM object */
288a4717116SVijay Mahadevan   dmmoab = (DM_Moab*)(*dm)->data;
289a4717116SVijay Mahadevan   mbiface = dmmoab->mbiface;
290a4717116SVijay Mahadevan   pcomm = dmmoab->pcomm;
291a4717116SVijay Mahadevan   nprocs = pcomm->size();
292a4717116SVijay Mahadevan 
293a4717116SVijay Mahadevan   /* add mesh loading options specific to the DM */
294a4717116SVijay Mahadevan   ierr = DMMoab_GetReadOptions_Private(PETSC_TRUE, nprocs, dim, MOAB_PARROPTS_READ_PART, 0, usrreadopts, readopts);CHKERRQ(ierr);
295a4717116SVijay Mahadevan 
296a4717116SVijay Mahadevan   /* Load the mesh from a file. */
297a4717116SVijay Mahadevan   merr = mbiface->load_file(filename, &dmmoab->fileset, ""/*readopts*/);MBERRVM(mbiface,"Reading MOAB file failed.", merr);
298a4717116SVijay Mahadevan 
299*6e40195eSVijay Mahadevan   /* Reassign global IDs on all entities. */
300*6e40195eSVijay Mahadevan   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,false,true);MBERRNM(merr);
301*6e40195eSVijay Mahadevan 
302a4717116SVijay Mahadevan   merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr);
303a4717116SVijay Mahadevan 
304a4717116SVijay Mahadevan   /* load the local vertices */
305a4717116SVijay Mahadevan   merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr);
306a4717116SVijay Mahadevan   ierr = DMMoabSetLocalVertices(*dm, &verts);CHKERRQ(ierr);
307a4717116SVijay Mahadevan 
308a4717116SVijay Mahadevan   /* load the local elements */
309a4717116SVijay Mahadevan   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr);
310a4717116SVijay Mahadevan //  ierr = DMMoabSetLocalElements(*dm, &elems);CHKERRQ(ierr);
311a4717116SVijay Mahadevan 
312a4717116SVijay Mahadevan   merr = mbiface->set_dimension(dim);MBERRNM(merr);
313a4717116SVijay Mahadevan 
314a4717116SVijay Mahadevan   PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size());
31551d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
31651d15aeeSVijay Mahadevan }
31751d15aeeSVijay Mahadevan 
318