xref: /petsc/src/dm/impls/moab/dmmbutil.cxx (revision ce1fd00907e16585510fb4fe8f5a2f0f252d7469)
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 
49e5136372SVijay Mahadevan static void set_structured_coordinates(PetscInt i, PetscInt j, PetscInt k, PetscReal hxyz, PetscInt vcount, std::vector<double*>& vcoords)
50e5136372SVijay Mahadevan {
51e5136372SVijay Mahadevan   vcoords[0][vcount] = i*hxyz;
52e5136372SVijay Mahadevan   vcoords[1][vcount] = j*hxyz;
53e5136372SVijay Mahadevan   vcoords[2][vcount] = k*hxyz;
54e5136372SVijay Mahadevan }
55e5136372SVijay Mahadevan 
56e5136372SVijay Mahadevan static void set_element_connectivity(PetscInt dim, moab::EntityType etype, PetscInt offset, PetscInt nele, PetscInt i, PetscInt j, PetscInt k, PetscInt vfirst, moab::EntityHandle *connectivity)
57e5136372SVijay Mahadevan {
58e5136372SVijay Mahadevan   std::vector<int>    subent_conn(pow(2,dim));  /* only linear edge, quad, hex supported now */
59e5136372SVijay Mahadevan 
60e5136372SVijay Mahadevan   moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data());
61e5136372SVijay Mahadevan 
62e5136372SVijay Mahadevan   switch(dim) {
63e5136372SVijay Mahadevan     case 1:
64e5136372SVijay Mahadevan       connectivity[offset+subent_conn[0]] = vfirst+i;
65e5136372SVijay Mahadevan       connectivity[offset+subent_conn[1]] = vfirst+(i+1);
66e5136372SVijay Mahadevan       break;
67e5136372SVijay Mahadevan     case 2:
68e5136372SVijay Mahadevan       connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1);
69e5136372SVijay Mahadevan       connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1);
70e5136372SVijay Mahadevan       connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1);
71e5136372SVijay Mahadevan       connectivity[offset+subent_conn[3]] = vfirst+i+(j+1)*(nele+1);
72e5136372SVijay Mahadevan       break;
73e5136372SVijay Mahadevan     case 3:
74e5136372SVijay Mahadevan     default:
75e5136372SVijay Mahadevan       connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
76e5136372SVijay Mahadevan       connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k);
77e5136372SVijay Mahadevan       connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
78e5136372SVijay Mahadevan       connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k);
79e5136372SVijay Mahadevan       connectivity[offset+subent_conn[4]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1));
80e5136372SVijay Mahadevan       connectivity[offset+subent_conn[5]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
81e5136372SVijay Mahadevan       connectivity[offset+subent_conn[6]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1));
82e5136372SVijay Mahadevan       connectivity[offset+subent_conn[7]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
83e5136372SVijay Mahadevan       break;
84e5136372SVijay Mahadevan   }
85e5136372SVijay Mahadevan }
8651d15aeeSVijay Mahadevan 
8751d15aeeSVijay Mahadevan #undef __FUNCT__
8851d15aeeSVijay Mahadevan #define __FUNCT__ "DMMoabCreateBoxMesh"
893a4aeca1SVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt nele, PetscInt user_nghost, DM *dm)
9051d15aeeSVijay Mahadevan {
9151d15aeeSVijay Mahadevan   PetscErrorCode  ierr;
9251d15aeeSVijay Mahadevan   moab::ErrorCode merr;
933a4aeca1SVijay Mahadevan   PetscInt        i,j,k,n,nprocs;
9451d15aeeSVijay Mahadevan   DM_Moab        *dmmoab;
9551d15aeeSVijay Mahadevan   moab::Interface *mbiface;
9651d15aeeSVijay Mahadevan   moab::ParallelComm *pcomm;
97a4717116SVijay Mahadevan   moab::ReadUtilIface* readMeshIface;
98a4717116SVijay Mahadevan 
9951d15aeeSVijay Mahadevan   moab::Tag  id_tag=PETSC_NULL;
100a4717116SVijay Mahadevan   moab::Range         ownedvtx,ownedelms;
1013a4aeca1SVijay Mahadevan   moab::EntityHandle  vfirst,efirst,regionset,faceset,edgeset,vtxset;
102a4717116SVijay Mahadevan   std::vector<double*> vcoords;
103a4717116SVijay Mahadevan   moab::EntityHandle  *connectivity = 0;
10451d15aeeSVijay Mahadevan   moab::EntityType etype;
105c8d5365dSVijay Mahadevan   PetscInt    ise[6];
10651d15aeeSVijay Mahadevan   PetscReal   xse[6];
1073a4aeca1SVijay Mahadevan   /* TODO: Fix nghost > 0 - now relying on exchange_ghost_cells */
1083a4aeca1SVijay Mahadevan   const PetscInt nghost=0;
1093a4aeca1SVijay Mahadevan 
1103a4aeca1SVijay Mahadevan   moab::Tag geom_tag;
1113a4aeca1SVijay Mahadevan 
1123a4aeca1SVijay Mahadevan   moab::Range adj,dim3,dim2;
1133a4aeca1SVijay Mahadevan   bool build_adjacencies=false;
11451d15aeeSVijay Mahadevan 
115a4717116SVijay Mahadevan   const PetscInt npts=nele+1;        /* Number of points in every dimension */
116e5136372SVijay Mahadevan   PetscInt vpere,locnele,locnpts,ghnele,ghnpts;    /* Number of verts/element, vertices, elements owned by this process */
11751d15aeeSVijay Mahadevan 
11851d15aeeSVijay Mahadevan   PetscFunctionBegin;
119e5136372SVijay Mahadevan   if(dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid dimension argument for mesh: dim=[1,3].\n");
120e5136372SVijay Mahadevan 
121c8d5365dSVijay Mahadevan   ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr);
122e5136372SVijay Mahadevan   /* total number of vertices in all dimensions */
123e5136372SVijay Mahadevan   n=pow(npts,dim);
124e5136372SVijay Mahadevan 
125e5136372SVijay Mahadevan   /* do some error checking */
126e5136372SVijay Mahadevan   if(n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2.\n");
127e5136372SVijay 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");
128e5136372SVijay Mahadevan   if(nghost < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of ghost layers cannot be negative.\n");
129e5136372SVijay Mahadevan 
130a4717116SVijay Mahadevan   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
13151d15aeeSVijay Mahadevan   ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr);
13251d15aeeSVijay Mahadevan 
133a4717116SVijay Mahadevan   /* get all the necessary handles from the private DM object */
13451d15aeeSVijay Mahadevan   dmmoab = (DM_Moab*)(*dm)->data;
13551d15aeeSVijay Mahadevan   mbiface = dmmoab->mbiface;
13651d15aeeSVijay Mahadevan   pcomm = dmmoab->pcomm;
13751d15aeeSVijay Mahadevan   id_tag = dmmoab->ltog_tag;
13851d15aeeSVijay Mahadevan   nprocs = pcomm->size();
139c8d5365dSVijay Mahadevan   dmmoab->dim = dim;
14051d15aeeSVijay Mahadevan 
141b5410836SVijay Mahadevan   /* create a file set to associate all entities in current mesh */
142b5410836SVijay Mahadevan   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr);
143b5410836SVijay Mahadevan 
144a4717116SVijay Mahadevan   /* No errors yet; proceed with building the mesh */
145a4717116SVijay Mahadevan   merr = mbiface->query_interface(readMeshIface);MBERRNM(merr);
14651d15aeeSVijay Mahadevan 
14751d15aeeSVijay Mahadevan   ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr);
148a4717116SVijay Mahadevan 
149a4717116SVijay Mahadevan   /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */
15051d15aeeSVijay Mahadevan   ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr);
15151d15aeeSVijay Mahadevan 
152a4717116SVijay Mahadevan   /* set some variables based on dimension */
15351d15aeeSVijay Mahadevan   switch(dim) {
15451d15aeeSVijay Mahadevan    case 1:
15551d15aeeSVijay Mahadevan     vpere = 2;
156a4717116SVijay Mahadevan     locnele = (ise[1]-ise[0]);
157a4717116SVijay Mahadevan     locnpts = (ise[1]-ise[0]+1);
158c8d5365dSVijay Mahadevan     ghnele = (nghost > 0 ? (ise[0] > nghost ? 1 : 0) + (ise[1] < nele - nghost ? 1 : 0) : 0 );
159c8d5365dSVijay Mahadevan     ghnpts = (nghost > 0 ? (ise[0] > 0 ? 1 : 0) + (ise[1] < nele ? 1 : 0) : 0);
16051d15aeeSVijay Mahadevan     etype = moab::MBEDGE;
16151d15aeeSVijay Mahadevan     break;
16251d15aeeSVijay Mahadevan    case 2:
16351d15aeeSVijay Mahadevan     vpere = 4;
164a4717116SVijay Mahadevan     locnele = (ise[1]-ise[0])*(ise[3]-ise[2]);
165a4717116SVijay Mahadevan     locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1);
166c8d5365dSVijay Mahadevan     ghnele = (nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0);
167c8d5365dSVijay Mahadevan     ghnpts = (nghost > 0 ? (ise[2] > 0 ? npts : 0) + (ise[3] < nele ? npts : 0) : 0);
16851d15aeeSVijay Mahadevan     etype = moab::MBQUAD;
16951d15aeeSVijay Mahadevan     break;
17051d15aeeSVijay Mahadevan    case 3:
17151d15aeeSVijay Mahadevan     vpere = 8;
172a4717116SVijay Mahadevan     locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]);
173a4717116SVijay Mahadevan     locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1);
174c8d5365dSVijay Mahadevan     ghnele = (nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0);
175c8d5365dSVijay Mahadevan     ghnpts = (nghost > 0 ? (ise[4] > 0 ? npts*npts : 0) + (ise[5] < nele ? npts*npts : 0) : 0);
17651d15aeeSVijay Mahadevan     etype = moab::MBHEX;
17751d15aeeSVijay Mahadevan     break;
17851d15aeeSVijay Mahadevan   }
17951d15aeeSVijay Mahadevan 
18051d15aeeSVijay Mahadevan   /* we have a domain of size [1,1,1] - now compute local co-ordinate box */
18151d15aeeSVijay Mahadevan   ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr);
182c8d5365dSVijay Mahadevan   for(i=0; i<6; ++i) {
18351d15aeeSVijay Mahadevan     xse[i]=(PetscReal)ise[i]/nele;
18451d15aeeSVijay Mahadevan   }
18551d15aeeSVijay Mahadevan 
186a4717116SVijay Mahadevan   /* Create vertexes and set the coodinate of each vertex */
187c8d5365dSVijay Mahadevan   merr = readMeshIface->get_node_coords(3,locnpts+ghnpts,0,vfirst,vcoords,n);MBERRNM(merr);
18851d15aeeSVijay Mahadevan 
189a4717116SVijay Mahadevan   /* Compute the co-ordinates of vertices and global IDs */
190c8d5365dSVijay Mahadevan   std::vector<int>    vgid(locnpts+ghnpts);
19151d15aeeSVijay Mahadevan   int vcount=0;
19251d15aeeSVijay Mahadevan   double hxyz=1.0/nele;
193e5136372SVijay Mahadevan 
194c8d5365dSVijay Mahadevan   /* create all the owned vertices */
195c8d5365dSVijay Mahadevan   for (k = ise[4]; k <= ise[5]; k++) {
196c8d5365dSVijay Mahadevan     for (j = ise[2]; j <= ise[3]; j++) {
197c8d5365dSVijay Mahadevan       for (i = ise[0]; i <= ise[1]; i++, vcount++) {
198c8d5365dSVijay Mahadevan         set_structured_coordinates(i,j,k,hxyz,vcount,vcoords);
199c8d5365dSVijay Mahadevan         vgid[vcount] = (k*npts+j)*npts+i+1;
200c8d5365dSVijay Mahadevan       }
201c8d5365dSVijay Mahadevan     }
202c8d5365dSVijay Mahadevan   }
203c8d5365dSVijay Mahadevan 
204e5136372SVijay Mahadevan   /* create ghosted vertices requested by user - below the current plane */
205e5136372SVijay Mahadevan   if (ise[2*dim-2] > 0) {
206e5136372SVijay Mahadevan     for (k = (dim==3?ise[4]-nghost:ise[4]); k <= (dim==3?ise[4]-1:ise[5]); k++) {
207e5136372SVijay Mahadevan       for (j = (dim==2?ise[2]-nghost:ise[2]); j <= (dim==2?ise[2]-1:ise[3]); j++) {
208e5136372SVijay Mahadevan         for (i = (dim>1?ise[0]:ise[0]-nghost); i <= (dim>1?ise[1]:ise[0]-1); i++, vcount++) {
209e5136372SVijay Mahadevan           set_structured_coordinates(i,j,k,hxyz,vcount,vcoords);
210e5136372SVijay Mahadevan           vgid[vcount] = (k*npts+j)*npts+i+1;
211e5136372SVijay Mahadevan         }
212e5136372SVijay Mahadevan       }
213e5136372SVijay Mahadevan     }
214e5136372SVijay Mahadevan   }
215e5136372SVijay Mahadevan 
216e5136372SVijay Mahadevan   /* create ghosted vertices requested by user - above the current plane */
217e5136372SVijay Mahadevan   if (ise[2*dim-1] < nele) {
218e5136372SVijay Mahadevan     for (k = (dim==3?ise[5]+1:ise[4]); k <= (dim==3?ise[5]+nghost:ise[5]); k++) {
219e5136372SVijay Mahadevan       for (j = (dim==2?ise[3]+1:ise[2]); j <= (dim==2?ise[3]+nghost:ise[3]); j++) {
220e5136372SVijay Mahadevan         for (i = (dim>1?ise[0]:ise[1]+1); i <= (dim>1?ise[1]:ise[1]+nghost); i++, vcount++) {
221e5136372SVijay Mahadevan           set_structured_coordinates(i,j,k,hxyz,vcount,vcoords);
222e5136372SVijay Mahadevan           vgid[vcount] = (k*npts+j)*npts+i+1;
223e5136372SVijay Mahadevan         }
22451d15aeeSVijay Mahadevan       }
22551d15aeeSVijay Mahadevan     }
22651d15aeeSVijay Mahadevan   }
22751d15aeeSVijay Mahadevan 
228c8d5365dSVijay Mahadevan   if (locnpts+ghnpts != vcount) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,vcount);
229c8d5365dSVijay Mahadevan 
23051d15aeeSVijay Mahadevan   merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr);
23151d15aeeSVijay Mahadevan 
232a4717116SVijay Mahadevan   /* The global ID tag is applied to each owned
233a4717116SVijay Mahadevan      vertex. It acts as an global identifier which MOAB uses to
234a4717116SVijay Mahadevan      assemble the individual pieces of the mesh:
235a4717116SVijay Mahadevan      Set the global ID indices */
23651d15aeeSVijay Mahadevan   merr = mbiface->tag_set_data(id_tag,ownedvtx,vgid.data());MBERRNM(merr);
23751d15aeeSVijay Mahadevan 
238a4717116SVijay Mahadevan   /* Create elements between mesh points using the ReadUtilInterface
239a4717116SVijay Mahadevan      get the reference to element connectivities for all local elements from the ReadUtilInterface */
240e5136372SVijay Mahadevan   merr = readMeshIface->get_element_connect (locnele+ghnele,vpere,etype,1,efirst,connectivity);MBERRNM(merr);
24151d15aeeSVijay Mahadevan 
242a4717116SVijay Mahadevan   /* offset appropriately so that only local ID and not global ID numbers are set for connectivity array */
243c8d5365dSVijay Mahadevan //  PetscPrintf(PETSC_COMM_SELF, "[%D] first local handle %D\n", rank, vfirst);
244e5136372SVijay Mahadevan   vfirst-=vgid[0]-1;
24551d15aeeSVijay Mahadevan 
246a4717116SVijay Mahadevan    /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k):
247a4717116SVijay Mahadevan          and then set the connectivity for each element appropriately */
24851d15aeeSVijay Mahadevan   int ecount=0;
24951d15aeeSVijay Mahadevan 
250e5136372SVijay Mahadevan   /* create ghosted elements requested by user - below the current plane */
251e5136372SVijay Mahadevan   if (ise[2*dim-2] >= nghost) {
252e5136372SVijay Mahadevan     for (k = (dim==3?ise[4]-nghost:ise[4]); k < (dim==3?ise[4]:std::max(ise[5],1)); k++) {
253e5136372SVijay Mahadevan       for (j = (dim==2?ise[2]-nghost:ise[2]); j < (dim==2?ise[2]:std::max(ise[3],1)); j++) {
254e5136372SVijay Mahadevan         for (i = (dim>1?ise[0]:ise[0]-nghost); i < (dim>1?std::max(ise[1],1):ise[0]); i++, ecount++) {
255c8d5365dSVijay Mahadevan           set_element_connectivity(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity);
25651d15aeeSVijay Mahadevan         }
25751d15aeeSVijay Mahadevan       }
25851d15aeeSVijay Mahadevan     }
25951d15aeeSVijay Mahadevan   }
260e5136372SVijay Mahadevan 
261e5136372SVijay Mahadevan   /* create owned elements requested by user */
262e5136372SVijay Mahadevan   for (k = ise[4]; k < std::max(ise[5],1); k++) {
263e5136372SVijay Mahadevan     for (j = ise[2]; j < std::max(ise[3],1); j++) {
264e5136372SVijay Mahadevan       for (i = ise[0]; i < std::max(ise[1],1); i++,ecount++) {
265c8d5365dSVijay Mahadevan         set_element_connectivity(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity);
266e5136372SVijay Mahadevan       }
267e5136372SVijay Mahadevan     }
268e5136372SVijay Mahadevan   }
269e5136372SVijay Mahadevan 
270e5136372SVijay Mahadevan   /* create ghosted elements requested by user - above the current plane */
271e5136372SVijay Mahadevan   if (ise[2*dim-1] <= nele-nghost) {
272e5136372SVijay Mahadevan     for (k = (dim==3?ise[5]:ise[4]); k < (dim==3?ise[5]+nghost:std::max(ise[5],1)); k++) {
273e5136372SVijay Mahadevan       for (j = (dim==2?ise[3]:ise[2]); j < (dim==2?ise[3]+nghost:std::max(ise[3],1)); j++) {
274e5136372SVijay Mahadevan         for (i = (dim>1?ise[0]:ise[1]); i < (dim>1?std::max(ise[1],1):ise[1]+nghost); i++, ecount++) {
275c8d5365dSVijay Mahadevan           set_element_connectivity(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity);
276e5136372SVijay Mahadevan         }
277e5136372SVijay Mahadevan       }
278e5136372SVijay Mahadevan     }
279e5136372SVijay Mahadevan   }
280e5136372SVijay Mahadevan 
281e5136372SVijay Mahadevan   merr = readMeshIface->update_adjacencies(efirst,locnele+ghnele,vpere,connectivity);MBERRNM(merr);
28251d15aeeSVijay Mahadevan 
283a4717116SVijay 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.
284a4717116SVijay Mahadevan         first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */
28551d15aeeSVijay Mahadevan   merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr);
286a4717116SVijay Mahadevan 
287e5136372SVijay Mahadevan   if (locnele+ghnele != (int) ownedelms.size())
288e5136372SVijay Mahadevan     SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele+ghnele,ownedelms.size());
289e5136372SVijay Mahadevan   else if(locnpts+ghnpts != (int) ownedvtx.size())
290e5136372SVijay Mahadevan     SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,ownedvtx.size());
29151d15aeeSVijay Mahadevan   else
292a4717116SVijay Mahadevan     PetscInfo2(*dm, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size());
29351d15aeeSVijay Mahadevan 
2943a4aeca1SVijay Mahadevan   /* lets create some sets */
2953a4aeca1SVijay 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);
2963a4aeca1SVijay Mahadevan 
2973a4aeca1SVijay Mahadevan   merr = mbiface->create_meshset(moab::MESHSET_SET, regionset);MBERRNM(merr);
2983a4aeca1SVijay Mahadevan   merr = mbiface->add_entities(regionset, ownedelms);MBERRNM(merr);
2993a4aeca1SVijay Mahadevan   merr = mbiface->tag_set_data(geom_tag, &regionset, 1, &dmmoab->dim);MBERRNM(merr);
300b5410836SVijay Mahadevan   merr = mbiface->add_parent_child(dmmoab->fileset,regionset);MBERRNM(merr);
3013a4aeca1SVijay Mahadevan   merr = mbiface->unite_meshset(dmmoab->fileset, regionset);MBERRNM(merr);
3023a4aeca1SVijay Mahadevan 
3033a4aeca1SVijay Mahadevan   merr = mbiface->create_meshset(moab::MESHSET_SET, vtxset);MBERRNM(merr);
3043a4aeca1SVijay Mahadevan   merr = mbiface->add_entities(vtxset, ownedvtx);MBERRNM(merr);
3053a4aeca1SVijay Mahadevan   merr = mbiface->add_parent_child(dmmoab->fileset,vtxset);MBERRNM(merr);
3063a4aeca1SVijay Mahadevan   merr = mbiface->unite_meshset(dmmoab->fileset, vtxset);MBERRNM(merr);
3073a4aeca1SVijay Mahadevan 
3083a4aeca1SVijay Mahadevan   if (build_adjacencies) {
3093a4aeca1SVijay Mahadevan     // generate all lower dimensional adjacencies
3103a4aeca1SVijay Mahadevan     merr = mbiface->get_adjacencies( ownedelms, dim-1, true, adj, moab::Interface::UNION );MBERRNM(merr);
3113a4aeca1SVijay Mahadevan     merr = dmmoab->pcomm->get_part_entities(dim2, dim-1);MBERRNM(merr);
3123a4aeca1SVijay Mahadevan     adj.merge(dim2);
3133a4aeca1SVijay Mahadevan 
3143a4aeca1SVijay Mahadevan     /* create face sets */
3153a4aeca1SVijay Mahadevan     merr = mbiface->create_meshset(moab::MESHSET_SET, faceset);MBERRNM(merr);
3163a4aeca1SVijay Mahadevan     merr = mbiface->add_entities(faceset, adj);MBERRNM(merr);
3173a4aeca1SVijay Mahadevan     merr = mbiface->add_parent_child(dmmoab->fileset,faceset);MBERRNM(merr);
3183a4aeca1SVijay Mahadevan     i=dim-1;
3193a4aeca1SVijay Mahadevan     merr = mbiface->tag_set_data(geom_tag, &faceset, 1, &i);MBERRNM(merr);
3203a4aeca1SVijay Mahadevan     merr = mbiface->unite_meshset(dmmoab->fileset, faceset);MBERRNM(merr);
3213a4aeca1SVijay Mahadevan     PetscInfo2(dm, "Found %d %d-Dim quantities.\n", adj.size(), dim-1);
3223a4aeca1SVijay Mahadevan 
3233a4aeca1SVijay Mahadevan     if (dim > 2) {
3243a4aeca1SVijay Mahadevan       dim2.clear();
3253a4aeca1SVijay Mahadevan       /* create edge sets, if appropriate i.e., if dim=3 */
3263a4aeca1SVijay Mahadevan       merr = mbiface->create_meshset(moab::MESHSET_SET, edgeset);MBERRNM(merr);
3273a4aeca1SVijay Mahadevan       merr = mbiface->get_adjacencies(adj, dim-1, true, dim2, moab::Interface::UNION );MBERRNM(merr);
3283a4aeca1SVijay Mahadevan       merr = mbiface->add_entities(edgeset, dim2);MBERRNM(merr);
3293a4aeca1SVijay Mahadevan       merr = mbiface->add_parent_child(dmmoab->fileset,edgeset);MBERRNM(merr);
3303a4aeca1SVijay Mahadevan       i=dim-2;
3313a4aeca1SVijay Mahadevan       merr = mbiface->tag_set_data(geom_tag, &edgeset, 1, &i);MBERRNM(merr);
3323a4aeca1SVijay Mahadevan       merr = mbiface->unite_meshset(dmmoab->fileset, edgeset);MBERRNM(merr);
3333a4aeca1SVijay Mahadevan       PetscInfo2(dm, "Found %d %d-Dim quantities.\n", adj.size(), dim-2);
3343a4aeca1SVijay Mahadevan     }
3353a4aeca1SVijay Mahadevan   }
3363a4aeca1SVijay Mahadevan 
337a4717116SVijay Mahadevan   /* check the handles */
338a4717116SVijay Mahadevan   merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr);
33951d15aeeSVijay Mahadevan 
340c8d5365dSVijay Mahadevan #if 0
341e5136372SVijay Mahadevan   std::stringstream sstr;
342e5136372SVijay Mahadevan   sstr << "test_" << pcomm->rank() << ".vtk";
343e5136372SVijay Mahadevan   mbiface->write_mesh(sstr.str().c_str());
344e5136372SVijay Mahadevan #endif
345e5136372SVijay Mahadevan 
346a4717116SVijay Mahadevan   /* resolve the shared entities by exchanging information to adjacent processors */
347a4717116SVijay Mahadevan   merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr);
348*ce1fd009SVijay Mahadevan   merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,NULL,&id_tag);MBERRV(mbiface,merr);
34951d15aeeSVijay Mahadevan 
350a4717116SVijay Mahadevan   /* Reassign global IDs on all entities. */
351e5136372SVijay Mahadevan   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr);
35251d15aeeSVijay Mahadevan 
3533a4aeca1SVijay Mahadevan   merr = pcomm->exchange_ghost_cells(dim,0,user_nghost,dim,true,false,&dmmoab->fileset);MBERRV(mbiface,merr);
354c8d5365dSVijay Mahadevan 
355a4717116SVijay Mahadevan   /* Everything is set up, now just do a tag exchange to update tags
356a4717116SVijay Mahadevan      on all of the ghost vertexes */
357c8d5365dSVijay Mahadevan   merr = mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr);
358c8d5365dSVijay Mahadevan   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, ownedelms);MBERRNM(merr);
359c8d5365dSVijay Mahadevan 
360a4717116SVijay Mahadevan   merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr);
361a4717116SVijay Mahadevan   merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr);
36251d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
36351d15aeeSVijay Mahadevan }
36451d15aeeSVijay Mahadevan 
36551d15aeeSVijay Mahadevan 
366a4717116SVijay Mahadevan #undef __FUNCT__
367a4717116SVijay Mahadevan #define __FUNCT__ "DMMoab_GetReadOptions_Private"
3687023aa44SVijay Mahadevan PetscErrorCode DMMoab_GetReadOptions_Private(PetscBool by_rank, PetscInt numproc, PetscInt dim, MoabReadMode mode, PetscInt dbglevel, const char* extra_opts, const char** read_opts)
36951d15aeeSVijay Mahadevan {
370a4717116SVijay Mahadevan   std::ostringstream str;
37151d15aeeSVijay Mahadevan 
372a4717116SVijay Mahadevan   PetscFunctionBegin;
373e23c60ebSVijay Mahadevan   /* do parallel read unless using only one processor */
374a4717116SVijay Mahadevan   if (numproc > 1) {
3757023aa44SVijay Mahadevan     str << "PARALLEL=" << mode << ";PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;";
3767023aa44SVijay Mahadevan     str << "PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=" << dim << ".0.1;";
377a4717116SVijay Mahadevan     if (by_rank)
378a4717116SVijay Mahadevan       str << "PARTITION_BY_RANK;";
379a4717116SVijay Mahadevan   }
38051d15aeeSVijay Mahadevan 
381c8d5365dSVijay Mahadevan   if (dbglevel) {
382c8d5365dSVijay Mahadevan     str << "CPUTIME;DEBUG_IO=" << dbglevel << ";";
383c8d5365dSVijay Mahadevan     if (numproc>1)
384c8d5365dSVijay Mahadevan       str << "DEBUG_PIO=" << dbglevel << ";";
385c8d5365dSVijay Mahadevan   }
38651d15aeeSVijay Mahadevan 
387c8d5365dSVijay Mahadevan   if (extra_opts)
388c8d5365dSVijay Mahadevan     str << extra_opts;
38951d15aeeSVijay Mahadevan 
3907023aa44SVijay Mahadevan   *read_opts = str.str().c_str();
39151d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
39251d15aeeSVijay Mahadevan }
39351d15aeeSVijay Mahadevan 
39451d15aeeSVijay Mahadevan 
39551d15aeeSVijay Mahadevan #undef __FUNCT__
396a4717116SVijay Mahadevan #define __FUNCT__ "DMMoabLoadFromFile"
397a4717116SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm)
39851d15aeeSVijay Mahadevan {
399a4717116SVijay Mahadevan   moab::ErrorCode merr;
4007023aa44SVijay Mahadevan   PetscInt        nprocs,dbglevel=0;
4017023aa44SVijay Mahadevan   PetscBool       part_by_rank=PETSC_FALSE;
402a4717116SVijay Mahadevan   DM_Moab        *dmmoab;
403a4717116SVijay Mahadevan   moab::Interface *mbiface;
404a4717116SVijay Mahadevan   moab::ParallelComm *pcomm;
405a4717116SVijay Mahadevan   moab::Range verts,elems;
4067023aa44SVijay Mahadevan   const char *readopts;
407a4717116SVijay Mahadevan   PetscErrorCode ierr;
40851d15aeeSVijay Mahadevan 
40951d15aeeSVijay Mahadevan   PetscFunctionBegin;
410a4717116SVijay Mahadevan   PetscValidPointer(dm,5);
41151d15aeeSVijay Mahadevan 
412a4717116SVijay Mahadevan   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
413a4717116SVijay Mahadevan   ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr);
41451d15aeeSVijay Mahadevan 
415a4717116SVijay Mahadevan   /* get all the necessary handles from the private DM object */
416a4717116SVijay Mahadevan   dmmoab = (DM_Moab*)(*dm)->data;
417a4717116SVijay Mahadevan   mbiface = dmmoab->mbiface;
418a4717116SVijay Mahadevan   pcomm = dmmoab->pcomm;
419a4717116SVijay Mahadevan   nprocs = pcomm->size();
420c8d5365dSVijay Mahadevan   dmmoab->dim = dim;
421a4717116SVijay Mahadevan 
422b5410836SVijay Mahadevan   /* create a file set to associate all entities in current mesh */
423b5410836SVijay Mahadevan   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr);
424b5410836SVijay Mahadevan 
4257023aa44SVijay Mahadevan   /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */
4267023aa44SVijay Mahadevan   ierr  = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab");
4277023aa44SVijay Mahadevan   ierr  = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr);
4287023aa44SVijay Mahadevan   ierr  = PetscOptionsBool("-dmmb_part_by_rank", "Use partition by rank when reading MOAB meshes from file", "dmmbutil.cxx", part_by_rank, &part_by_rank, NULL);CHKERRQ(ierr);
4297023aa44SVijay Mahadevan   ierr  = PetscOptionsEnd();
4307023aa44SVijay Mahadevan 
431a4717116SVijay Mahadevan   /* add mesh loading options specific to the DM */
4327023aa44SVijay Mahadevan   ierr = DMMoab_GetReadOptions_Private(part_by_rank, nprocs, dim, MOAB_PARROPTS_READ_PART, dbglevel, usrreadopts, &readopts);CHKERRQ(ierr);
4337023aa44SVijay Mahadevan 
434e5136372SVijay Mahadevan   PetscInfo2(*dm, "Reading file %s with options: %s\n",filename,readopts);
435a4717116SVijay Mahadevan 
436a4717116SVijay Mahadevan   /* Load the mesh from a file. */
4377023aa44SVijay Mahadevan   merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr);
438a4717116SVijay Mahadevan 
4396e40195eSVijay Mahadevan   /* Reassign global IDs on all entities. */
440e5136372SVijay Mahadevan   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr);
441e5136372SVijay Mahadevan 
442e5136372SVijay Mahadevan   /* load the local vertices */
443e5136372SVijay Mahadevan   merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr);
444e5136372SVijay Mahadevan   /* load the local elements */
445e5136372SVijay Mahadevan   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr);
446e5136372SVijay Mahadevan 
447e5136372SVijay Mahadevan   /* Everything is set up, now just do a tag exchange to update tags
448e5136372SVijay Mahadevan      on all of the ghost vertexes */
449e5136372SVijay Mahadevan   merr = pcomm->exchange_tags(dmmoab->ltog_tag,verts);MBERRV(mbiface,merr);
450e5136372SVijay Mahadevan   merr = pcomm->exchange_tags(dmmoab->ltog_tag,elems);MBERRV(mbiface,merr);
451e5136372SVijay Mahadevan 
452e5136372SVijay Mahadevan   merr = pcomm->exchange_ghost_cells(dim,0,1,0,true,true,&dmmoab->fileset);MBERRV(mbiface,merr);
4536e40195eSVijay Mahadevan 
454a4717116SVijay Mahadevan   merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr);
455a4717116SVijay Mahadevan 
456e5136372SVijay Mahadevan   PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size());
457e5136372SVijay Mahadevan 
458e5136372SVijay Mahadevan #if 1
459e5136372SVijay Mahadevan   std::stringstream sstr;
460e5136372SVijay Mahadevan   sstr << "test_" << pcomm->rank() << ".vtk";
461e5136372SVijay Mahadevan   mbiface->write_mesh(sstr.str().c_str());
462e5136372SVijay Mahadevan #endif
46351d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
46451d15aeeSVijay Mahadevan }
46551d15aeeSVijay Mahadevan 
466