xref: /petsc/src/dm/impls/moab/dmmbutil.cxx (revision 2e4e7c01f2fc776e49aa0e3bac11897d06784bfb)
1b8ecf6d3SVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
2b8ecf6d3SVijay Mahadevan #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 
1151d15aeeSVijay Mahadevan #undef __FUNCT__
1251d15aeeSVijay Mahadevan #define __FUNCT__ "DMMoabComputeDomainBounds_Private"
13aa859218SVijay Mahadevan static 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 
32b8ecf6d3SVijay Mahadevan   if(rank < remainder) { /* This process gets "fraction+1" elements */
3351d15aeeSVijay Mahadevan     sizes[dim-1] = (fraction + 1);
3451d15aeeSVijay Mahadevan     starts[dim-1] = rank * (fraction+1);
35b8ecf6d3SVijay Mahadevan   } else { /* This process gets "fraction" elements */
3651d15aeeSVijay Mahadevan     sizes[dim-1] = fraction;
3751d15aeeSVijay Mahadevan     starts[dim-1] = (remainder*(fraction+1) + fraction*(rank-remainder));
3851d15aeeSVijay Mahadevan   }
3951d15aeeSVijay Mahadevan 
4051d15aeeSVijay Mahadevan   for(int i=dim-1; i>=0; --i) {
4151d15aeeSVijay Mahadevan     ise[2*i]=starts[i];ise[2*i+1]=starts[i]+sizes[i];
4251d15aeeSVijay Mahadevan   }
4351d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
4451d15aeeSVijay Mahadevan }
4551d15aeeSVijay Mahadevan 
46aa859218SVijay Mahadevan #undef __FUNCT__
47aa859218SVijay Mahadevan #define __FUNCT__ "DMMoab_SetStructuredCoords_Private"
48aa859218SVijay Mahadevan static void DMMoab_SetStructuredCoords_Private(PetscInt i, PetscInt j, PetscInt k, PetscReal hx, PetscReal hy, PetscReal hz, PetscInt vcount, std::vector<double*>& vcoords)
49e5136372SVijay Mahadevan {
5066f88a78SVijay Mahadevan   vcoords[0][vcount] = i*hx;
5166f88a78SVijay Mahadevan   vcoords[1][vcount] = j*hy;
5266f88a78SVijay Mahadevan   vcoords[2][vcount] = k*hz;
53e5136372SVijay Mahadevan }
54e5136372SVijay Mahadevan 
55aa859218SVijay Mahadevan #undef __FUNCT__
56aa859218SVijay Mahadevan #define __FUNCT__ "DMMoab_SetElementConnectivity_Private"
57aa859218SVijay Mahadevan static void DMMoab_SetElementConnectivity_Private(PetscInt dim, moab::EntityType etype, PetscInt offset, PetscInt nele, PetscInt i, PetscInt j, PetscInt k, PetscInt vfirst, moab::EntityHandle *connectivity)
58e5136372SVijay Mahadevan {
59e5136372SVijay Mahadevan   std::vector<int>    subent_conn(pow(2,dim));  /* only linear edge, quad, hex supported now */
60e5136372SVijay Mahadevan 
61e5136372SVijay Mahadevan   moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data());
62e5136372SVijay Mahadevan 
63e5136372SVijay Mahadevan   switch(dim) {
64e5136372SVijay Mahadevan     case 1:
65e5136372SVijay Mahadevan       connectivity[offset+subent_conn[0]] = vfirst+i;
66e5136372SVijay Mahadevan       connectivity[offset+subent_conn[1]] = vfirst+(i+1);
67e5136372SVijay Mahadevan       break;
68e5136372SVijay Mahadevan     case 2:
69e5136372SVijay Mahadevan       connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1);
70e5136372SVijay Mahadevan       connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1);
71e5136372SVijay Mahadevan       connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1);
72e5136372SVijay Mahadevan       connectivity[offset+subent_conn[3]] = vfirst+i+(j+1)*(nele+1);
73e5136372SVijay Mahadevan       break;
74e5136372SVijay Mahadevan     case 3:
75e5136372SVijay Mahadevan     default:
76e5136372SVijay Mahadevan       connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
77e5136372SVijay Mahadevan       connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k);
78e5136372SVijay Mahadevan       connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
79e5136372SVijay Mahadevan       connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k);
80e5136372SVijay Mahadevan       connectivity[offset+subent_conn[4]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1));
81e5136372SVijay Mahadevan       connectivity[offset+subent_conn[5]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
82e5136372SVijay Mahadevan       connectivity[offset+subent_conn[6]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1));
83e5136372SVijay Mahadevan       connectivity[offset+subent_conn[7]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
84e5136372SVijay Mahadevan       break;
85e5136372SVijay Mahadevan   }
86e5136372SVijay Mahadevan }
8751d15aeeSVijay Mahadevan 
8851d15aeeSVijay Mahadevan #undef __FUNCT__
8951d15aeeSVijay Mahadevan #define __FUNCT__ "DMMoabCreateBoxMesh"
90aa859218SVijay Mahadevan /*@
91aa859218SVijay Mahadevan   DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with user specified bounds.
92aa859218SVijay Mahadevan 
93aa859218SVijay Mahadevan   Collective on MPI_Comm
94aa859218SVijay Mahadevan 
95aa859218SVijay Mahadevan   Input Parameters:
96aa859218SVijay Mahadevan + comm - The communicator for the DM object
97aa859218SVijay Mahadevan . dim - The spatial dimension
98b8ecf6d3SVijay 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
99b8ecf6d3SVijay Mahadevan . nele - The number of discrete elements in each direction
100b8ecf6d3SVijay Mahadevan . user_nghost - The number of ghosted layers needed in the partitioned mesh
101aa859218SVijay Mahadevan 
102aa859218SVijay Mahadevan   Output Parameter:
103aa859218SVijay Mahadevan . dm  - The DM object
104aa859218SVijay Mahadevan 
105aa859218SVijay Mahadevan   Level: beginner
106aa859218SVijay Mahadevan 
107aa859218SVijay Mahadevan .keywords: DM, create
108aa859218SVijay Mahadevan .seealso: DMSetType(), DMCreate(), DMMoabLoadFromFile()
109aa859218SVijay Mahadevan @*/
11066f88a78SVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, const PetscReal* bounds, PetscInt nele, PetscInt user_nghost, DM *dm)
11151d15aeeSVijay Mahadevan {
11251d15aeeSVijay Mahadevan   PetscErrorCode  ierr;
11351d15aeeSVijay Mahadevan   moab::ErrorCode merr;
1143a4aeca1SVijay Mahadevan   PetscInt        i,j,k,n,nprocs;
11551d15aeeSVijay Mahadevan   DM_Moab        *dmmoab;
11651d15aeeSVijay Mahadevan   moab::Interface *mbiface;
11751d15aeeSVijay Mahadevan   moab::ParallelComm *pcomm;
118a4717116SVijay Mahadevan   moab::ReadUtilIface* readMeshIface;
119a4717116SVijay Mahadevan 
12051d15aeeSVijay Mahadevan   moab::Tag  id_tag=PETSC_NULL;
121a4717116SVijay Mahadevan   moab::Range         ownedvtx,ownedelms;
1223a4aeca1SVijay Mahadevan   moab::EntityHandle  vfirst,efirst,regionset,faceset,edgeset,vtxset;
123a4717116SVijay Mahadevan   std::vector<double*> vcoords;
124a4717116SVijay Mahadevan   moab::EntityHandle  *connectivity = 0;
12551d15aeeSVijay Mahadevan   moab::EntityType etype;
126c8d5365dSVijay Mahadevan   PetscInt    ise[6];
12766f88a78SVijay Mahadevan   PetscReal   xse[6],defbounds[6];
1283a4aeca1SVijay Mahadevan   /* TODO: Fix nghost > 0 - now relying on exchange_ghost_cells */
1293a4aeca1SVijay Mahadevan   const PetscInt nghost=0;
1303a4aeca1SVijay Mahadevan 
1313a4aeca1SVijay Mahadevan   moab::Tag geom_tag;
1323a4aeca1SVijay Mahadevan 
1333a4aeca1SVijay Mahadevan   moab::Range adj,dim3,dim2;
1343a4aeca1SVijay Mahadevan   bool build_adjacencies=false;
13551d15aeeSVijay Mahadevan 
136a4717116SVijay Mahadevan   const PetscInt npts=nele+1;        /* Number of points in every dimension */
137e5136372SVijay Mahadevan   PetscInt vpere,locnele,locnpts,ghnele,ghnpts;    /* Number of verts/element, vertices, elements owned by this process */
13851d15aeeSVijay Mahadevan 
13951d15aeeSVijay Mahadevan   PetscFunctionBegin;
140e5136372SVijay Mahadevan   if(dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid dimension argument for mesh: dim=[1,3].\n");
141e5136372SVijay Mahadevan 
142c8d5365dSVijay Mahadevan   ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr);
143e5136372SVijay Mahadevan   /* total number of vertices in all dimensions */
144e5136372SVijay Mahadevan   n=pow(npts,dim);
145e5136372SVijay Mahadevan 
146e5136372SVijay Mahadevan   /* do some error checking */
147e5136372SVijay Mahadevan   if(n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2.\n");
148e5136372SVijay 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");
149e5136372SVijay Mahadevan   if(nghost < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of ghost layers cannot be negative.\n");
150e5136372SVijay Mahadevan 
151a4717116SVijay Mahadevan   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
15251d15aeeSVijay Mahadevan   ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr);
15351d15aeeSVijay Mahadevan 
154a4717116SVijay Mahadevan   /* get all the necessary handles from the private DM object */
15551d15aeeSVijay Mahadevan   dmmoab = (DM_Moab*)(*dm)->data;
15651d15aeeSVijay Mahadevan   mbiface = dmmoab->mbiface;
15751d15aeeSVijay Mahadevan   pcomm = dmmoab->pcomm;
15851d15aeeSVijay Mahadevan   id_tag = dmmoab->ltog_tag;
15951d15aeeSVijay Mahadevan   nprocs = pcomm->size();
160c8d5365dSVijay Mahadevan   dmmoab->dim = dim;
16151d15aeeSVijay Mahadevan 
162b5410836SVijay Mahadevan   /* create a file set to associate all entities in current mesh */
163b5410836SVijay Mahadevan   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr);
164b5410836SVijay Mahadevan 
165a4717116SVijay Mahadevan   /* No errors yet; proceed with building the mesh */
166a4717116SVijay Mahadevan   merr = mbiface->query_interface(readMeshIface);MBERRNM(merr);
16751d15aeeSVijay Mahadevan 
16851d15aeeSVijay Mahadevan   ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr);
169a4717116SVijay Mahadevan 
170a4717116SVijay Mahadevan   /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */
17151d15aeeSVijay Mahadevan   ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr);
17251d15aeeSVijay Mahadevan 
173a4717116SVijay Mahadevan   /* set some variables based on dimension */
17451d15aeeSVijay Mahadevan   switch(dim) {
17551d15aeeSVijay Mahadevan    case 1:
17651d15aeeSVijay Mahadevan     vpere = 2;
177a4717116SVijay Mahadevan     locnele = (ise[1]-ise[0]);
178a4717116SVijay Mahadevan     locnpts = (ise[1]-ise[0]+1);
179c8d5365dSVijay Mahadevan     ghnele = (nghost > 0 ? (ise[0] > nghost ? 1 : 0) + (ise[1] < nele - nghost ? 1 : 0) : 0 );
180c8d5365dSVijay Mahadevan     ghnpts = (nghost > 0 ? (ise[0] > 0 ? 1 : 0) + (ise[1] < nele ? 1 : 0) : 0);
18151d15aeeSVijay Mahadevan     etype = moab::MBEDGE;
18251d15aeeSVijay Mahadevan     break;
18351d15aeeSVijay Mahadevan    case 2:
18451d15aeeSVijay Mahadevan     vpere = 4;
185a4717116SVijay Mahadevan     locnele = (ise[1]-ise[0])*(ise[3]-ise[2]);
186a4717116SVijay Mahadevan     locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1);
187c8d5365dSVijay Mahadevan     ghnele = (nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0);
188c8d5365dSVijay Mahadevan     ghnpts = (nghost > 0 ? (ise[2] > 0 ? npts : 0) + (ise[3] < nele ? npts : 0) : 0);
18951d15aeeSVijay Mahadevan     etype = moab::MBQUAD;
19051d15aeeSVijay Mahadevan     break;
19151d15aeeSVijay Mahadevan    case 3:
19251d15aeeSVijay Mahadevan     vpere = 8;
193a4717116SVijay Mahadevan     locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]);
194a4717116SVijay Mahadevan     locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1);
195c8d5365dSVijay Mahadevan     ghnele = (nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0);
196c8d5365dSVijay Mahadevan     ghnpts = (nghost > 0 ? (ise[4] > 0 ? npts*npts : 0) + (ise[5] < nele ? npts*npts : 0) : 0);
19751d15aeeSVijay Mahadevan     etype = moab::MBHEX;
19851d15aeeSVijay Mahadevan     break;
19951d15aeeSVijay Mahadevan   }
20051d15aeeSVijay Mahadevan 
20151d15aeeSVijay Mahadevan   /* we have a domain of size [1,1,1] - now compute local co-ordinate box */
20251d15aeeSVijay Mahadevan   ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr);
203c8d5365dSVijay Mahadevan   for(i=0; i<6; ++i) {
20451d15aeeSVijay Mahadevan     xse[i]=(PetscReal)ise[i]/nele;
20551d15aeeSVijay Mahadevan   }
20651d15aeeSVijay Mahadevan 
207a4717116SVijay Mahadevan   /* Create vertexes and set the coodinate of each vertex */
208c8d5365dSVijay Mahadevan   merr = readMeshIface->get_node_coords(3,locnpts+ghnpts,0,vfirst,vcoords,n);MBERRNM(merr);
20951d15aeeSVijay Mahadevan 
210a4717116SVijay Mahadevan   /* Compute the co-ordinates of vertices and global IDs */
211c8d5365dSVijay Mahadevan   std::vector<int>    vgid(locnpts+ghnpts);
21251d15aeeSVijay Mahadevan   int vcount=0;
21366f88a78SVijay Mahadevan 
21466f88a78SVijay Mahadevan   if (!bounds) { /* default box mesh is defined on a unit-cube */
21566f88a78SVijay Mahadevan     defbounds[0]=0.0; defbounds[1]=1.0;
21666f88a78SVijay Mahadevan     defbounds[2]=0.0; defbounds[3]=1.0;
21766f88a78SVijay Mahadevan     defbounds[4]=0.0; defbounds[5]=1.0;
21866f88a78SVijay Mahadevan     bounds=defbounds;
21966f88a78SVijay Mahadevan   }
22066f88a78SVijay Mahadevan   else {
22166f88a78SVijay Mahadevan     /* validate the bounds data */
22266f88a78SVijay 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]);
22366f88a78SVijay 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]);
22466f88a78SVijay 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]);
22566f88a78SVijay Mahadevan   }
22666f88a78SVijay Mahadevan 
22766f88a78SVijay Mahadevan   const double hx=(bounds[1]-bounds[0])/nele;
22866f88a78SVijay Mahadevan   const double hy=(dim > 1 ? (bounds[3]-bounds[2])/nele : 0.0);
22966f88a78SVijay Mahadevan   const double hz=(dim > 2 ? (bounds[5]-bounds[4])/nele : 0.0);
230e5136372SVijay Mahadevan 
231c8d5365dSVijay Mahadevan   /* create all the owned vertices */
232c8d5365dSVijay Mahadevan   for (k = ise[4]; k <= ise[5]; k++) {
233c8d5365dSVijay Mahadevan     for (j = ise[2]; j <= ise[3]; j++) {
234c8d5365dSVijay Mahadevan       for (i = ise[0]; i <= ise[1]; i++, vcount++) {
235aa859218SVijay Mahadevan         DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords);
236c8d5365dSVijay Mahadevan         vgid[vcount] = (k*npts+j)*npts+i+1;
237c8d5365dSVijay Mahadevan       }
238c8d5365dSVijay Mahadevan     }
239c8d5365dSVijay Mahadevan   }
240c8d5365dSVijay Mahadevan 
241e5136372SVijay Mahadevan   /* create ghosted vertices requested by user - below the current plane */
242e5136372SVijay Mahadevan   if (ise[2*dim-2] > 0) {
243e5136372SVijay Mahadevan     for (k = (dim==3?ise[4]-nghost:ise[4]); k <= (dim==3?ise[4]-1:ise[5]); k++) {
244e5136372SVijay Mahadevan       for (j = (dim==2?ise[2]-nghost:ise[2]); j <= (dim==2?ise[2]-1:ise[3]); j++) {
245e5136372SVijay Mahadevan         for (i = (dim>1?ise[0]:ise[0]-nghost); i <= (dim>1?ise[1]:ise[0]-1); i++, vcount++) {
246aa859218SVijay Mahadevan           DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords);
247e5136372SVijay Mahadevan           vgid[vcount] = (k*npts+j)*npts+i+1;
248e5136372SVijay Mahadevan         }
249e5136372SVijay Mahadevan       }
250e5136372SVijay Mahadevan     }
251e5136372SVijay Mahadevan   }
252e5136372SVijay Mahadevan 
253e5136372SVijay Mahadevan   /* create ghosted vertices requested by user - above the current plane */
254e5136372SVijay Mahadevan   if (ise[2*dim-1] < nele) {
255e5136372SVijay Mahadevan     for (k = (dim==3?ise[5]+1:ise[4]); k <= (dim==3?ise[5]+nghost:ise[5]); k++) {
256e5136372SVijay Mahadevan       for (j = (dim==2?ise[3]+1:ise[2]); j <= (dim==2?ise[3]+nghost:ise[3]); j++) {
257e5136372SVijay Mahadevan         for (i = (dim>1?ise[0]:ise[1]+1); i <= (dim>1?ise[1]:ise[1]+nghost); i++, vcount++) {
258aa859218SVijay Mahadevan           DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords);
259e5136372SVijay Mahadevan           vgid[vcount] = (k*npts+j)*npts+i+1;
260e5136372SVijay Mahadevan         }
26151d15aeeSVijay Mahadevan       }
26251d15aeeSVijay Mahadevan     }
26351d15aeeSVijay Mahadevan   }
26451d15aeeSVijay Mahadevan 
265c8d5365dSVijay Mahadevan   if (locnpts+ghnpts != vcount) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,vcount);
266c8d5365dSVijay Mahadevan 
26751d15aeeSVijay Mahadevan   merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr);
26851d15aeeSVijay Mahadevan 
269a4717116SVijay Mahadevan   /* The global ID tag is applied to each owned
270a4717116SVijay Mahadevan      vertex. It acts as an global identifier which MOAB uses to
271a4717116SVijay Mahadevan      assemble the individual pieces of the mesh:
272a4717116SVijay Mahadevan      Set the global ID indices */
27351d15aeeSVijay Mahadevan   merr = mbiface->tag_set_data(id_tag,ownedvtx,vgid.data());MBERRNM(merr);
27451d15aeeSVijay Mahadevan 
275a4717116SVijay Mahadevan   /* Create elements between mesh points using the ReadUtilInterface
276a4717116SVijay Mahadevan      get the reference to element connectivities for all local elements from the ReadUtilInterface */
277e5136372SVijay Mahadevan   merr = readMeshIface->get_element_connect (locnele+ghnele,vpere,etype,1,efirst,connectivity);MBERRNM(merr);
27851d15aeeSVijay Mahadevan 
279a4717116SVijay Mahadevan   /* offset appropriately so that only local ID and not global ID numbers are set for connectivity array */
280c8d5365dSVijay Mahadevan //  PetscPrintf(PETSC_COMM_SELF, "[%D] first local handle %D\n", rank, vfirst);
281e5136372SVijay Mahadevan   vfirst-=vgid[0]-1;
28251d15aeeSVijay Mahadevan 
283a4717116SVijay Mahadevan    /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k):
284a4717116SVijay Mahadevan          and then set the connectivity for each element appropriately */
28551d15aeeSVijay Mahadevan   int ecount=0;
28651d15aeeSVijay Mahadevan 
287e5136372SVijay Mahadevan   /* create ghosted elements requested by user - below the current plane */
288e5136372SVijay Mahadevan   if (ise[2*dim-2] >= nghost) {
289e5136372SVijay Mahadevan     for (k = (dim==3?ise[4]-nghost:ise[4]); k < (dim==3?ise[4]:std::max(ise[5],1)); k++) {
290e5136372SVijay Mahadevan       for (j = (dim==2?ise[2]-nghost:ise[2]); j < (dim==2?ise[2]:std::max(ise[3],1)); j++) {
291e5136372SVijay Mahadevan         for (i = (dim>1?ise[0]:ise[0]-nghost); i < (dim>1?std::max(ise[1],1):ise[0]); i++, ecount++) {
292aa859218SVijay Mahadevan           DMMoab_SetElementConnectivity_Private(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity);
29351d15aeeSVijay Mahadevan         }
29451d15aeeSVijay Mahadevan       }
29551d15aeeSVijay Mahadevan     }
29651d15aeeSVijay Mahadevan   }
297e5136372SVijay Mahadevan 
298e5136372SVijay Mahadevan   /* create owned elements requested by user */
299e5136372SVijay Mahadevan   for (k = ise[4]; k < std::max(ise[5],1); k++) {
300e5136372SVijay Mahadevan     for (j = ise[2]; j < std::max(ise[3],1); j++) {
301e5136372SVijay Mahadevan       for (i = ise[0]; i < std::max(ise[1],1); i++,ecount++) {
302aa859218SVijay Mahadevan         DMMoab_SetElementConnectivity_Private(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity);
303e5136372SVijay Mahadevan       }
304e5136372SVijay Mahadevan     }
305e5136372SVijay Mahadevan   }
306e5136372SVijay Mahadevan 
307e5136372SVijay Mahadevan   /* create ghosted elements requested by user - above the current plane */
308e5136372SVijay Mahadevan   if (ise[2*dim-1] <= nele-nghost) {
309e5136372SVijay Mahadevan     for (k = (dim==3?ise[5]:ise[4]); k < (dim==3?ise[5]+nghost:std::max(ise[5],1)); k++) {
310e5136372SVijay Mahadevan       for (j = (dim==2?ise[3]:ise[2]); j < (dim==2?ise[3]+nghost:std::max(ise[3],1)); j++) {
311e5136372SVijay Mahadevan         for (i = (dim>1?ise[0]:ise[1]); i < (dim>1?std::max(ise[1],1):ise[1]+nghost); i++, ecount++) {
312aa859218SVijay Mahadevan           DMMoab_SetElementConnectivity_Private(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity);
313e5136372SVijay Mahadevan         }
314e5136372SVijay Mahadevan       }
315e5136372SVijay Mahadevan     }
316e5136372SVijay Mahadevan   }
317e5136372SVijay Mahadevan 
318e5136372SVijay Mahadevan   merr = readMeshIface->update_adjacencies(efirst,locnele+ghnele,vpere,connectivity);MBERRNM(merr);
31951d15aeeSVijay Mahadevan 
320a4717116SVijay 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.
321a4717116SVijay Mahadevan         first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */
32251d15aeeSVijay Mahadevan   merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr);
323a4717116SVijay Mahadevan 
324e5136372SVijay Mahadevan   if (locnele+ghnele != (int) ownedelms.size())
325e5136372SVijay Mahadevan     SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele+ghnele,ownedelms.size());
326e5136372SVijay Mahadevan   else if(locnpts+ghnpts != (int) ownedvtx.size())
327e5136372SVijay Mahadevan     SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,ownedvtx.size());
32851d15aeeSVijay Mahadevan   else
329e427d9c9SVijay Mahadevan     PetscInfo2(NULL, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size());
33051d15aeeSVijay Mahadevan 
3313a4aeca1SVijay Mahadevan   /* lets create some sets */
3323a4aeca1SVijay 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);
3333a4aeca1SVijay Mahadevan 
3343a4aeca1SVijay Mahadevan   merr = mbiface->create_meshset(moab::MESHSET_SET, regionset);MBERRNM(merr);
3353a4aeca1SVijay Mahadevan   merr = mbiface->add_entities(regionset, ownedelms);MBERRNM(merr);
3363a4aeca1SVijay Mahadevan   merr = mbiface->tag_set_data(geom_tag, &regionset, 1, &dmmoab->dim);MBERRNM(merr);
337b5410836SVijay Mahadevan   merr = mbiface->add_parent_child(dmmoab->fileset,regionset);MBERRNM(merr);
3383a4aeca1SVijay Mahadevan   merr = mbiface->unite_meshset(dmmoab->fileset, regionset);MBERRNM(merr);
3393a4aeca1SVijay Mahadevan 
3403a4aeca1SVijay Mahadevan   merr = mbiface->create_meshset(moab::MESHSET_SET, vtxset);MBERRNM(merr);
3413a4aeca1SVijay Mahadevan   merr = mbiface->add_entities(vtxset, ownedvtx);MBERRNM(merr);
3423a4aeca1SVijay Mahadevan   merr = mbiface->add_parent_child(dmmoab->fileset,vtxset);MBERRNM(merr);
3433a4aeca1SVijay Mahadevan   merr = mbiface->unite_meshset(dmmoab->fileset, vtxset);MBERRNM(merr);
3443a4aeca1SVijay Mahadevan 
3453a4aeca1SVijay Mahadevan   if (build_adjacencies) {
3463a4aeca1SVijay Mahadevan     // generate all lower dimensional adjacencies
3473a4aeca1SVijay Mahadevan     merr = mbiface->get_adjacencies( ownedelms, dim-1, true, adj, moab::Interface::UNION );MBERRNM(merr);
3483a4aeca1SVijay Mahadevan     merr = dmmoab->pcomm->get_part_entities(dim2, dim-1);MBERRNM(merr);
3493a4aeca1SVijay Mahadevan     adj.merge(dim2);
3503a4aeca1SVijay Mahadevan 
3513a4aeca1SVijay Mahadevan     /* create face sets */
3523a4aeca1SVijay Mahadevan     merr = mbiface->create_meshset(moab::MESHSET_SET, faceset);MBERRNM(merr);
3533a4aeca1SVijay Mahadevan     merr = mbiface->add_entities(faceset, adj);MBERRNM(merr);
3543a4aeca1SVijay Mahadevan     merr = mbiface->add_parent_child(dmmoab->fileset,faceset);MBERRNM(merr);
3553a4aeca1SVijay Mahadevan     i=dim-1;
3563a4aeca1SVijay Mahadevan     merr = mbiface->tag_set_data(geom_tag, &faceset, 1, &i);MBERRNM(merr);
3573a4aeca1SVijay Mahadevan     merr = mbiface->unite_meshset(dmmoab->fileset, faceset);MBERRNM(merr);
358b8ecf6d3SVijay Mahadevan     PetscInfo2(NULL, "Found %d %d-Dim quantities.\n", adj.size(), dim-1);
3593a4aeca1SVijay Mahadevan 
3603a4aeca1SVijay Mahadevan     if (dim > 2) {
3613a4aeca1SVijay Mahadevan       dim2.clear();
3623a4aeca1SVijay Mahadevan       /* create edge sets, if appropriate i.e., if dim=3 */
3633a4aeca1SVijay Mahadevan       merr = mbiface->create_meshset(moab::MESHSET_SET, edgeset);MBERRNM(merr);
3643a4aeca1SVijay Mahadevan       merr = mbiface->get_adjacencies(adj, dim-1, true, dim2, moab::Interface::UNION );MBERRNM(merr);
3653a4aeca1SVijay Mahadevan       merr = mbiface->add_entities(edgeset, dim2);MBERRNM(merr);
3663a4aeca1SVijay Mahadevan       merr = mbiface->add_parent_child(dmmoab->fileset,edgeset);MBERRNM(merr);
3673a4aeca1SVijay Mahadevan       i=dim-2;
3683a4aeca1SVijay Mahadevan       merr = mbiface->tag_set_data(geom_tag, &edgeset, 1, &i);MBERRNM(merr);
3693a4aeca1SVijay Mahadevan       merr = mbiface->unite_meshset(dmmoab->fileset, edgeset);MBERRNM(merr);
370b8ecf6d3SVijay Mahadevan       PetscInfo2(NULL, "Found %d %d-Dim quantities.\n", adj.size(), dim-2);
3713a4aeca1SVijay Mahadevan     }
3723a4aeca1SVijay Mahadevan   }
3733a4aeca1SVijay Mahadevan 
374a4717116SVijay Mahadevan   /* check the handles */
375a4717116SVijay Mahadevan   merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr);
37651d15aeeSVijay Mahadevan 
377a4717116SVijay Mahadevan   /* resolve the shared entities by exchanging information to adjacent processors */
378a4717116SVijay Mahadevan   merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr);
379ce1fd009SVijay Mahadevan   merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,NULL,&id_tag);MBERRV(mbiface,merr);
38051d15aeeSVijay Mahadevan 
3813a4aeca1SVijay Mahadevan   merr = pcomm->exchange_ghost_cells(dim,0,user_nghost,dim,true,false,&dmmoab->fileset);MBERRV(mbiface,merr);
382c8d5365dSVijay Mahadevan 
383e427d9c9SVijay Mahadevan   /* Reassign global IDs on all entities. */
384e427d9c9SVijay Mahadevan   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,false,true,false);MBERRNM(merr);
385e427d9c9SVijay Mahadevan 
386a4717116SVijay Mahadevan   /* Everything is set up, now just do a tag exchange to update tags
387a4717116SVijay Mahadevan      on all of the ghost vertexes */
388c8d5365dSVijay Mahadevan   merr = mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr);
389c8d5365dSVijay Mahadevan   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, ownedelms);MBERRNM(merr);
390c8d5365dSVijay Mahadevan 
391a4717116SVijay Mahadevan   merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr);
392a4717116SVijay Mahadevan   merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr);
39351d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
39451d15aeeSVijay Mahadevan }
39551d15aeeSVijay Mahadevan 
39651d15aeeSVijay Mahadevan 
397a4717116SVijay Mahadevan #undef __FUNCT__
398a4717116SVijay Mahadevan #define __FUNCT__ "DMMoab_GetReadOptions_Private"
399*2e4e7c01SVijay 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)
40051d15aeeSVijay Mahadevan {
401a4717116SVijay Mahadevan   std::ostringstream str;
40251d15aeeSVijay Mahadevan 
403a4717116SVijay Mahadevan   PetscFunctionBegin;
404e23c60ebSVijay Mahadevan   /* do parallel read unless using only one processor */
405a4717116SVijay Mahadevan   if (numproc > 1) {
406*2e4e7c01SVijay Mahadevan     str << "PARALLEL=" << MoabReadModes[mode] << ";PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;";
4077023aa44SVijay Mahadevan     str << "PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=" << dim << ".0.1;";
408a4717116SVijay Mahadevan     if (by_rank)
409a4717116SVijay Mahadevan       str << "PARTITION_BY_RANK;";
410a4717116SVijay Mahadevan   }
41151d15aeeSVijay Mahadevan 
412*2e4e7c01SVijay Mahadevan   if (strlen(dm_opts)) {
413*2e4e7c01SVijay Mahadevan     str << dm_opts << ";";
414*2e4e7c01SVijay Mahadevan   }
415*2e4e7c01SVijay Mahadevan 
416c8d5365dSVijay Mahadevan   if (dbglevel) {
417c8d5365dSVijay Mahadevan     str << "CPUTIME;DEBUG_IO=" << dbglevel << ";";
418c8d5365dSVijay Mahadevan     if (numproc>1)
419c8d5365dSVijay Mahadevan       str << "DEBUG_PIO=" << dbglevel << ";";
420c8d5365dSVijay Mahadevan   }
42151d15aeeSVijay Mahadevan 
422c8d5365dSVijay Mahadevan   if (extra_opts)
423c8d5365dSVijay Mahadevan     str << extra_opts;
42451d15aeeSVijay Mahadevan 
4257023aa44SVijay Mahadevan   *read_opts = str.str().c_str();
42651d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
42751d15aeeSVijay Mahadevan }
42851d15aeeSVijay Mahadevan 
42951d15aeeSVijay Mahadevan 
43051d15aeeSVijay Mahadevan #undef __FUNCT__
431a4717116SVijay Mahadevan #define __FUNCT__ "DMMoabLoadFromFile"
432aa859218SVijay Mahadevan /*@
433aa859218SVijay Mahadevan   DMMoabLoadFromFile - Creates a DM object by loading the mesh from a user specified file.
434aa859218SVijay Mahadevan 
435aa859218SVijay Mahadevan   Collective on MPI_Comm
436aa859218SVijay Mahadevan 
437aa859218SVijay Mahadevan   Input Parameters:
438aa859218SVijay Mahadevan + comm - The communicator for the DM object
439aa859218SVijay Mahadevan . dim - The spatial dimension
440b8ecf6d3SVijay Mahadevan . filename - The name of the mesh file to be loaded
441b8ecf6d3SVijay Mahadevan . usrreadopts - The options string to read a MOAB mesh.
442aa859218SVijay Mahadevan   Reference (Parallel Mesh Initialization: http://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo)
443aa859218SVijay Mahadevan 
444aa859218SVijay Mahadevan   Output Parameter:
445aa859218SVijay Mahadevan . dm  - The DM object
446aa859218SVijay Mahadevan 
447aa859218SVijay Mahadevan   Level: beginner
448aa859218SVijay Mahadevan 
449aa859218SVijay Mahadevan .keywords: DM, create
450b8ecf6d3SVijay Mahadevan 
451aa859218SVijay Mahadevan .seealso: DMSetType(), DMCreate(), DMMoabCreateBoxMesh()
452aa859218SVijay Mahadevan @*/
453a4717116SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm)
45451d15aeeSVijay Mahadevan {
455a4717116SVijay Mahadevan   moab::ErrorCode merr;
456*2e4e7c01SVijay Mahadevan   PetscInt        nprocs;
457a4717116SVijay Mahadevan   DM_Moab        *dmmoab;
458a4717116SVijay Mahadevan   moab::Interface *mbiface;
459a4717116SVijay Mahadevan   moab::ParallelComm *pcomm;
460a4717116SVijay Mahadevan   moab::Range verts,elems;
4617023aa44SVijay Mahadevan   const char *readopts;
462a4717116SVijay Mahadevan   PetscErrorCode ierr;
46351d15aeeSVijay Mahadevan 
46451d15aeeSVijay Mahadevan   PetscFunctionBegin;
465a4717116SVijay Mahadevan   PetscValidPointer(dm,5);
46651d15aeeSVijay Mahadevan 
467a4717116SVijay Mahadevan   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
468a4717116SVijay Mahadevan   ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr);
46951d15aeeSVijay Mahadevan 
470a4717116SVijay Mahadevan   /* get all the necessary handles from the private DM object */
471a4717116SVijay Mahadevan   dmmoab = (DM_Moab*)(*dm)->data;
472a4717116SVijay Mahadevan   mbiface = dmmoab->mbiface;
473a4717116SVijay Mahadevan   pcomm = dmmoab->pcomm;
474a4717116SVijay Mahadevan   nprocs = pcomm->size();
475aa859218SVijay Mahadevan   /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */
476c8d5365dSVijay Mahadevan   dmmoab->dim = dim;
477a4717116SVijay Mahadevan 
478b5410836SVijay Mahadevan   /* create a file set to associate all entities in current mesh */
479b5410836SVijay Mahadevan   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr);
480b5410836SVijay Mahadevan 
481a4717116SVijay Mahadevan   /* add mesh loading options specific to the DM */
482*2e4e7c01SVijay Mahadevan   ierr = DMMoab_GetReadOptions_Private(dmmoab->partition_by_rank, nprocs, dim, dmmoab->read_mode,
483*2e4e7c01SVijay Mahadevan                                         dmmoab->rw_dbglevel, dmmoab->extra_read_options, usrreadopts, &readopts);CHKERRQ(ierr);
4847023aa44SVijay Mahadevan 
485e5136372SVijay Mahadevan   PetscInfo2(*dm, "Reading file %s with options: %s\n",filename,readopts);
486a4717116SVijay Mahadevan 
487a4717116SVijay Mahadevan   /* Load the mesh from a file. */
4887023aa44SVijay Mahadevan   merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr);
489a4717116SVijay Mahadevan 
4906e40195eSVijay Mahadevan   /* Reassign global IDs on all entities. */
491e5136372SVijay Mahadevan   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr);
492e5136372SVijay Mahadevan 
493e5136372SVijay Mahadevan   /* load the local vertices */
494e5136372SVijay Mahadevan   merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr);
495e5136372SVijay Mahadevan   /* load the local elements */
496e5136372SVijay Mahadevan   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr);
497e5136372SVijay Mahadevan 
498e5136372SVijay Mahadevan   /* Everything is set up, now just do a tag exchange to update tags
499e5136372SVijay Mahadevan      on all of the ghost vertexes */
500e5136372SVijay Mahadevan   merr = pcomm->exchange_tags(dmmoab->ltog_tag,verts);MBERRV(mbiface,merr);
501e5136372SVijay Mahadevan   merr = pcomm->exchange_tags(dmmoab->ltog_tag,elems);MBERRV(mbiface,merr);
502e5136372SVijay Mahadevan 
503e5136372SVijay Mahadevan   merr = pcomm->exchange_ghost_cells(dim,0,1,0,true,true,&dmmoab->fileset);MBERRV(mbiface,merr);
5046e40195eSVijay Mahadevan 
505a4717116SVijay Mahadevan   merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr);
506a4717116SVijay Mahadevan 
507e5136372SVijay Mahadevan   PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size());
50851d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
50951d15aeeSVijay Mahadevan }
51051d15aeeSVijay Mahadevan 
511