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__ 56*c3dd00c7SVijay Mahadevan #define __FUNCT__ "DMMoab_SetTensorElementConnectivity_Private" 57*c3dd00c7SVijay 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) 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__ 89*c3dd00c7SVijay Mahadevan #define __FUNCT__ "DMMoab_SetSimplexElementConnectivity_Private" 90*c3dd00c7SVijay 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) 91*c3dd00c7SVijay Mahadevan { 92*c3dd00c7SVijay Mahadevan std::vector<int> subent_conn(pow(2,dim)); /* only linear edge, quad, hex supported now */ 93*c3dd00c7SVijay Mahadevan 94*c3dd00c7SVijay Mahadevan moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data()); 95*c3dd00c7SVijay Mahadevan 96*c3dd00c7SVijay Mahadevan switch(dim) { 97*c3dd00c7SVijay Mahadevan case 1: 98*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i; 99*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1); 100*c3dd00c7SVijay Mahadevan break; 101*c3dd00c7SVijay Mahadevan case 2: 102*c3dd00c7SVijay Mahadevan if (subelem) { /* 1 2 3 of a QUAD */ 103*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1); 104*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1); 105*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1); 106*c3dd00c7SVijay Mahadevan } 107*c3dd00c7SVijay Mahadevan else { /* 3 4 1 of a QUAD */ 108*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+(i+1)+(j+1)*(nele+1); 109*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+i+(j+1)*(nele+1); 110*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+i+j*(nele+1); 111*c3dd00c7SVijay Mahadevan } 112*c3dd00c7SVijay Mahadevan break; 113*c3dd00c7SVijay Mahadevan case 3: 114*c3dd00c7SVijay Mahadevan default: 115*c3dd00c7SVijay Mahadevan switch(subelem) { 116*c3dd00c7SVijay Mahadevan case 0: /* 0 1 2 5 of a HEX */ 117*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k); 118*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k); 119*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k); 120*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1)); 121*c3dd00c7SVijay Mahadevan break; 122*c3dd00c7SVijay Mahadevan case 1: /* 0 2 7 5 of a HEX */ 123*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k); 124*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k); 125*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1)); 126*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1)); 127*c3dd00c7SVijay Mahadevan break; 128*c3dd00c7SVijay Mahadevan case 2: /* 0 2 3 7 of a HEX */ 129*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k); 130*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k); 131*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k); 132*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1)); 133*c3dd00c7SVijay Mahadevan break; 134*c3dd00c7SVijay Mahadevan case 3: /* 0 5 7 4 of a HEX */ 135*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k); 136*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1)); 137*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1)); 138*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1)); 139*c3dd00c7SVijay Mahadevan break; 140*c3dd00c7SVijay Mahadevan case 4: /* 2 7 5 6 of a HEX */ 141*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k); 142*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1)); 143*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1)); 144*c3dd00c7SVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1)); 145*c3dd00c7SVijay Mahadevan break; 146*c3dd00c7SVijay Mahadevan } 147*c3dd00c7SVijay Mahadevan break; 148*c3dd00c7SVijay Mahadevan } 149*c3dd00c7SVijay Mahadevan } 150*c3dd00c7SVijay Mahadevan 151*c3dd00c7SVijay Mahadevan #undef __FUNCT__ 152*c3dd00c7SVijay Mahadevan #define __FUNCT__ "DMMoab_SetElementConnectivity_Private" 153*c3dd00c7SVijay 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) 154*c3dd00c7SVijay Mahadevan { 155*c3dd00c7SVijay Mahadevan PetscInt m,subelem; 156*c3dd00c7SVijay Mahadevan if (useSimplex) { 157*c3dd00c7SVijay Mahadevan switch (dim) { 158*c3dd00c7SVijay Mahadevan case 1: 159*c3dd00c7SVijay Mahadevan subelem=1; 160*c3dd00c7SVijay Mahadevan DMMoab_SetSimplexElementConnectivity_Private(dim, 0, etype, (*ecount)*vpere, nele, i, j, k, vfirst, connectivity); 161*c3dd00c7SVijay Mahadevan break; 162*c3dd00c7SVijay Mahadevan case 2: 163*c3dd00c7SVijay Mahadevan subelem=2; 164*c3dd00c7SVijay Mahadevan for (m=0; m<subelem; m++) 165*c3dd00c7SVijay Mahadevan DMMoab_SetSimplexElementConnectivity_Private(dim, m, etype, (*ecount+m)*vpere, nele, i, j, k, vfirst, connectivity); 166*c3dd00c7SVijay Mahadevan break; 167*c3dd00c7SVijay Mahadevan case 3: 168*c3dd00c7SVijay Mahadevan subelem=5; 169*c3dd00c7SVijay Mahadevan for (m=0; m<subelem; m++) 170*c3dd00c7SVijay Mahadevan DMMoab_SetSimplexElementConnectivity_Private(dim, m, etype, (*ecount+m)*vpere, nele, i, j, k, vfirst, connectivity); 171*c3dd00c7SVijay Mahadevan break; 172*c3dd00c7SVijay Mahadevan } 173*c3dd00c7SVijay Mahadevan } 174*c3dd00c7SVijay Mahadevan else { 175*c3dd00c7SVijay Mahadevan subelem=1; 176*c3dd00c7SVijay Mahadevan DMMoab_SetTensorElementConnectivity_Private(dim, etype, (*ecount)*vpere, nele, i, j, k, vfirst, connectivity); 177*c3dd00c7SVijay Mahadevan } 178*c3dd00c7SVijay Mahadevan *ecount+=subelem; 179*c3dd00c7SVijay Mahadevan } 180*c3dd00c7SVijay Mahadevan 181*c3dd00c7SVijay Mahadevan 182*c3dd00c7SVijay Mahadevan #undef __FUNCT__ 18351d15aeeSVijay Mahadevan #define __FUNCT__ "DMMoabCreateBoxMesh" 184aa859218SVijay Mahadevan /*@ 185aa859218SVijay Mahadevan DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with user specified bounds. 186aa859218SVijay Mahadevan 187aa859218SVijay Mahadevan Collective on MPI_Comm 188aa859218SVijay Mahadevan 189aa859218SVijay Mahadevan Input Parameters: 190aa859218SVijay Mahadevan + comm - The communicator for the DM object 191aa859218SVijay Mahadevan . dim - The spatial dimension 192b8ecf6d3SVijay 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 193b8ecf6d3SVijay Mahadevan . nele - The number of discrete elements in each direction 194b8ecf6d3SVijay Mahadevan . user_nghost - The number of ghosted layers needed in the partitioned mesh 195aa859218SVijay Mahadevan 196aa859218SVijay Mahadevan Output Parameter: 197aa859218SVijay Mahadevan . dm - The DM object 198aa859218SVijay Mahadevan 199aa859218SVijay Mahadevan Level: beginner 200aa859218SVijay Mahadevan 201aa859218SVijay Mahadevan .keywords: DM, create 202aa859218SVijay Mahadevan .seealso: DMSetType(), DMCreate(), DMMoabLoadFromFile() 203aa859218SVijay Mahadevan @*/ 204*c3dd00c7SVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool useSimplex, const PetscReal* bounds, PetscInt nele, PetscInt user_nghost, DM *dm) 20551d15aeeSVijay Mahadevan { 20651d15aeeSVijay Mahadevan PetscErrorCode ierr; 20751d15aeeSVijay Mahadevan moab::ErrorCode merr; 2083a4aeca1SVijay Mahadevan PetscInt i,j,k,n,nprocs; 20951d15aeeSVijay Mahadevan DM_Moab *dmmoab; 21051d15aeeSVijay Mahadevan moab::Interface *mbiface; 21151d15aeeSVijay Mahadevan moab::ParallelComm *pcomm; 212a4717116SVijay Mahadevan moab::ReadUtilIface* readMeshIface; 213a4717116SVijay Mahadevan 21451d15aeeSVijay Mahadevan moab::Tag id_tag=PETSC_NULL; 215a4717116SVijay Mahadevan moab::Range ownedvtx,ownedelms; 2163a4aeca1SVijay Mahadevan moab::EntityHandle vfirst,efirst,regionset,faceset,edgeset,vtxset; 217a4717116SVijay Mahadevan std::vector<double*> vcoords; 218a4717116SVijay Mahadevan moab::EntityHandle *connectivity = 0; 21951d15aeeSVijay Mahadevan moab::EntityType etype; 220c8d5365dSVijay Mahadevan PetscInt ise[6]; 22166f88a78SVijay Mahadevan PetscReal xse[6],defbounds[6]; 2223a4aeca1SVijay Mahadevan /* TODO: Fix nghost > 0 - now relying on exchange_ghost_cells */ 2233a4aeca1SVijay Mahadevan const PetscInt nghost=0; 2243a4aeca1SVijay Mahadevan 2253a4aeca1SVijay Mahadevan moab::Tag geom_tag; 2263a4aeca1SVijay Mahadevan 2273a4aeca1SVijay Mahadevan moab::Range adj,dim3,dim2; 2283a4aeca1SVijay Mahadevan bool build_adjacencies=false; 22951d15aeeSVijay Mahadevan 230a4717116SVijay Mahadevan const PetscInt npts=nele+1; /* Number of points in every dimension */ 231e5136372SVijay Mahadevan PetscInt vpere,locnele,locnpts,ghnele,ghnpts; /* Number of verts/element, vertices, elements owned by this process */ 23251d15aeeSVijay Mahadevan 23351d15aeeSVijay Mahadevan PetscFunctionBegin; 234e5136372SVijay Mahadevan if(dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid dimension argument for mesh: dim=[1,3].\n"); 235e5136372SVijay Mahadevan 236c8d5365dSVijay Mahadevan ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr); 237e5136372SVijay Mahadevan /* total number of vertices in all dimensions */ 238e5136372SVijay Mahadevan n=pow(npts,dim); 239e5136372SVijay Mahadevan 240e5136372SVijay Mahadevan /* do some error checking */ 241e5136372SVijay Mahadevan if(n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2.\n"); 242e5136372SVijay 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"); 243e5136372SVijay Mahadevan if(nghost < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of ghost layers cannot be negative.\n"); 244e5136372SVijay Mahadevan 245a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 24651d15aeeSVijay Mahadevan ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 24751d15aeeSVijay Mahadevan 248a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 24951d15aeeSVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 25051d15aeeSVijay Mahadevan mbiface = dmmoab->mbiface; 25151d15aeeSVijay Mahadevan pcomm = dmmoab->pcomm; 25251d15aeeSVijay Mahadevan id_tag = dmmoab->ltog_tag; 25351d15aeeSVijay Mahadevan nprocs = pcomm->size(); 254c8d5365dSVijay Mahadevan dmmoab->dim = dim; 25551d15aeeSVijay Mahadevan 256b5410836SVijay Mahadevan /* create a file set to associate all entities in current mesh */ 257b5410836SVijay Mahadevan merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr); 258b5410836SVijay Mahadevan 259a4717116SVijay Mahadevan /* No errors yet; proceed with building the mesh */ 260a4717116SVijay Mahadevan merr = mbiface->query_interface(readMeshIface);MBERRNM(merr); 26151d15aeeSVijay Mahadevan 26251d15aeeSVijay Mahadevan ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr); 263a4717116SVijay Mahadevan 264a4717116SVijay Mahadevan /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */ 26551d15aeeSVijay Mahadevan ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr); 26651d15aeeSVijay Mahadevan 267a4717116SVijay Mahadevan /* set some variables based on dimension */ 26851d15aeeSVijay Mahadevan switch(dim) { 26951d15aeeSVijay Mahadevan case 1: 27051d15aeeSVijay Mahadevan vpere = 2; 271a4717116SVijay Mahadevan locnele = (ise[1]-ise[0]); 272a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1); 273c8d5365dSVijay Mahadevan ghnele = (nghost > 0 ? (ise[0] > nghost ? 1 : 0) + (ise[1] < nele - nghost ? 1 : 0) : 0 ); 274c8d5365dSVijay Mahadevan ghnpts = (nghost > 0 ? (ise[0] > 0 ? 1 : 0) + (ise[1] < nele ? 1 : 0) : 0); 27551d15aeeSVijay Mahadevan etype = moab::MBEDGE; 27651d15aeeSVijay Mahadevan break; 27751d15aeeSVijay Mahadevan case 2: 278*c3dd00c7SVijay Mahadevan locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1); 279*c3dd00c7SVijay Mahadevan ghnpts = (nghost > 0 ? (ise[2] > 0 ? npts : 0) + (ise[3] < nele ? npts : 0) : 0); 280*c3dd00c7SVijay Mahadevan if (useSimplex) { 281*c3dd00c7SVijay Mahadevan vpere = 3; 282*c3dd00c7SVijay Mahadevan locnele = 2*(ise[1]-ise[0])*(ise[3]-ise[2]); 283*c3dd00c7SVijay Mahadevan ghnele = 2*(nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0); 284*c3dd00c7SVijay Mahadevan etype = moab::MBTRI; 285*c3dd00c7SVijay Mahadevan } 286*c3dd00c7SVijay Mahadevan else { 28751d15aeeSVijay Mahadevan vpere = 4; 288a4717116SVijay Mahadevan locnele = (ise[1]-ise[0])*(ise[3]-ise[2]); 289c8d5365dSVijay Mahadevan ghnele = (nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0); 29051d15aeeSVijay Mahadevan etype = moab::MBQUAD; 291*c3dd00c7SVijay Mahadevan } 29251d15aeeSVijay Mahadevan break; 29351d15aeeSVijay Mahadevan case 3: 294*c3dd00c7SVijay Mahadevan locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1); 295*c3dd00c7SVijay Mahadevan ghnpts = (nghost > 0 ? (ise[4] > 0 ? npts*npts : 0) + (ise[5] < nele ? npts*npts : 0) : 0); 296*c3dd00c7SVijay Mahadevan if (useSimplex) { 297*c3dd00c7SVijay Mahadevan vpere = 4; 298*c3dd00c7SVijay Mahadevan locnele = 5*(ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]); 299*c3dd00c7SVijay Mahadevan ghnele = 5*(nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0); 300*c3dd00c7SVijay Mahadevan etype = moab::MBTET; 301*c3dd00c7SVijay Mahadevan } 302*c3dd00c7SVijay Mahadevan else { 30351d15aeeSVijay Mahadevan vpere = 8; 304a4717116SVijay Mahadevan locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]); 305c8d5365dSVijay Mahadevan ghnele = (nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0); 30651d15aeeSVijay Mahadevan etype = moab::MBHEX; 307*c3dd00c7SVijay Mahadevan } 30851d15aeeSVijay Mahadevan break; 30951d15aeeSVijay Mahadevan } 31051d15aeeSVijay Mahadevan 31151d15aeeSVijay Mahadevan /* we have a domain of size [1,1,1] - now compute local co-ordinate box */ 31251d15aeeSVijay Mahadevan ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr); 313c8d5365dSVijay Mahadevan for(i=0; i<6; ++i) { 31451d15aeeSVijay Mahadevan xse[i]=(PetscReal)ise[i]/nele; 31551d15aeeSVijay Mahadevan } 31651d15aeeSVijay Mahadevan 317a4717116SVijay Mahadevan /* Create vertexes and set the coodinate of each vertex */ 318c8d5365dSVijay Mahadevan merr = readMeshIface->get_node_coords(3,locnpts+ghnpts,0,vfirst,vcoords,n);MBERRNM(merr); 31951d15aeeSVijay Mahadevan 320a4717116SVijay Mahadevan /* Compute the co-ordinates of vertices and global IDs */ 321c8d5365dSVijay Mahadevan std::vector<int> vgid(locnpts+ghnpts); 32251d15aeeSVijay Mahadevan int vcount=0; 32366f88a78SVijay Mahadevan 32466f88a78SVijay Mahadevan if (!bounds) { /* default box mesh is defined on a unit-cube */ 32566f88a78SVijay Mahadevan defbounds[0]=0.0; defbounds[1]=1.0; 32666f88a78SVijay Mahadevan defbounds[2]=0.0; defbounds[3]=1.0; 32766f88a78SVijay Mahadevan defbounds[4]=0.0; defbounds[5]=1.0; 32866f88a78SVijay Mahadevan bounds=defbounds; 32966f88a78SVijay Mahadevan } 33066f88a78SVijay Mahadevan else { 33166f88a78SVijay Mahadevan /* validate the bounds data */ 33266f88a78SVijay 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]); 33366f88a78SVijay 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]); 33466f88a78SVijay 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]); 33566f88a78SVijay Mahadevan } 33666f88a78SVijay Mahadevan 33766f88a78SVijay Mahadevan const double hx=(bounds[1]-bounds[0])/nele; 33866f88a78SVijay Mahadevan const double hy=(dim > 1 ? (bounds[3]-bounds[2])/nele : 0.0); 33966f88a78SVijay Mahadevan const double hz=(dim > 2 ? (bounds[5]-bounds[4])/nele : 0.0); 340e5136372SVijay Mahadevan 341c8d5365dSVijay Mahadevan /* create all the owned vertices */ 342c8d5365dSVijay Mahadevan for (k = ise[4]; k <= ise[5]; k++) { 343c8d5365dSVijay Mahadevan for (j = ise[2]; j <= ise[3]; j++) { 344c8d5365dSVijay Mahadevan for (i = ise[0]; i <= ise[1]; i++, vcount++) { 345aa859218SVijay Mahadevan DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords); 346c8d5365dSVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 347c8d5365dSVijay Mahadevan } 348c8d5365dSVijay Mahadevan } 349c8d5365dSVijay Mahadevan } 350c8d5365dSVijay Mahadevan 351e5136372SVijay Mahadevan /* create ghosted vertices requested by user - below the current plane */ 352e5136372SVijay Mahadevan if (ise[2*dim-2] > 0) { 353e5136372SVijay Mahadevan for (k = (dim==3?ise[4]-nghost:ise[4]); k <= (dim==3?ise[4]-1:ise[5]); k++) { 354e5136372SVijay Mahadevan for (j = (dim==2?ise[2]-nghost:ise[2]); j <= (dim==2?ise[2]-1:ise[3]); j++) { 355e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[0]-nghost); i <= (dim>1?ise[1]:ise[0]-1); i++, vcount++) { 356aa859218SVijay Mahadevan DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords); 357e5136372SVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 358e5136372SVijay Mahadevan } 359e5136372SVijay Mahadevan } 360e5136372SVijay Mahadevan } 361e5136372SVijay Mahadevan } 362e5136372SVijay Mahadevan 363e5136372SVijay Mahadevan /* create ghosted vertices requested by user - above the current plane */ 364e5136372SVijay Mahadevan if (ise[2*dim-1] < nele) { 365e5136372SVijay Mahadevan for (k = (dim==3?ise[5]+1:ise[4]); k <= (dim==3?ise[5]+nghost:ise[5]); k++) { 366e5136372SVijay Mahadevan for (j = (dim==2?ise[3]+1:ise[2]); j <= (dim==2?ise[3]+nghost:ise[3]); j++) { 367e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[1]+1); i <= (dim>1?ise[1]:ise[1]+nghost); i++, vcount++) { 368aa859218SVijay Mahadevan DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords); 369e5136372SVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 370e5136372SVijay Mahadevan } 37151d15aeeSVijay Mahadevan } 37251d15aeeSVijay Mahadevan } 37351d15aeeSVijay Mahadevan } 37451d15aeeSVijay Mahadevan 375c8d5365dSVijay Mahadevan if (locnpts+ghnpts != vcount) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,vcount); 376c8d5365dSVijay Mahadevan 37751d15aeeSVijay Mahadevan merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr); 37851d15aeeSVijay Mahadevan 379a4717116SVijay Mahadevan /* The global ID tag is applied to each owned 380a4717116SVijay Mahadevan vertex. It acts as an global identifier which MOAB uses to 381a4717116SVijay Mahadevan assemble the individual pieces of the mesh: 382a4717116SVijay Mahadevan Set the global ID indices */ 38351d15aeeSVijay Mahadevan merr = mbiface->tag_set_data(id_tag,ownedvtx,vgid.data());MBERRNM(merr); 38451d15aeeSVijay Mahadevan 385a4717116SVijay Mahadevan /* Create elements between mesh points using the ReadUtilInterface 386a4717116SVijay Mahadevan get the reference to element connectivities for all local elements from the ReadUtilInterface */ 387e5136372SVijay Mahadevan merr = readMeshIface->get_element_connect (locnele+ghnele,vpere,etype,1,efirst,connectivity);MBERRNM(merr); 38851d15aeeSVijay Mahadevan 389a4717116SVijay Mahadevan /* offset appropriately so that only local ID and not global ID numbers are set for connectivity array */ 390e5136372SVijay Mahadevan vfirst-=vgid[0]-1; 39151d15aeeSVijay Mahadevan 392a4717116SVijay Mahadevan /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k): 393a4717116SVijay Mahadevan and then set the connectivity for each element appropriately */ 39451d15aeeSVijay Mahadevan int ecount=0; 39551d15aeeSVijay Mahadevan 396e5136372SVijay Mahadevan /* create ghosted elements requested by user - below the current plane */ 397e5136372SVijay Mahadevan if (ise[2*dim-2] >= nghost) { 398e5136372SVijay Mahadevan for (k = (dim==3?ise[4]-nghost:ise[4]); k < (dim==3?ise[4]:std::max(ise[5],1)); k++) { 399e5136372SVijay Mahadevan for (j = (dim==2?ise[2]-nghost:ise[2]); j < (dim==2?ise[2]:std::max(ise[3],1)); j++) { 400*c3dd00c7SVijay Mahadevan for (i = (dim>1?ise[0]:ise[0]-nghost); i < (dim>1?std::max(ise[1],1):ise[0]); i++) { 401*c3dd00c7SVijay Mahadevan DMMoab_SetElementConnectivity_Private(useSimplex, dim, etype, &ecount, vpere, nele, i, j, k, vfirst, connectivity); 40251d15aeeSVijay Mahadevan } 40351d15aeeSVijay Mahadevan } 40451d15aeeSVijay Mahadevan } 40551d15aeeSVijay Mahadevan } 406e5136372SVijay Mahadevan 407e5136372SVijay Mahadevan /* create owned elements requested by user */ 408e5136372SVijay Mahadevan for (k = ise[4]; k < std::max(ise[5],1); k++) { 409e5136372SVijay Mahadevan for (j = ise[2]; j < std::max(ise[3],1); j++) { 410*c3dd00c7SVijay Mahadevan for (i = ise[0]; i < std::max(ise[1],1); i++) { 411*c3dd00c7SVijay Mahadevan DMMoab_SetElementConnectivity_Private(useSimplex, dim, etype, &ecount, vpere, nele, i, j, k, vfirst, connectivity); 412e5136372SVijay Mahadevan } 413e5136372SVijay Mahadevan } 414e5136372SVijay Mahadevan } 415e5136372SVijay Mahadevan 416e5136372SVijay Mahadevan /* create ghosted elements requested by user - above the current plane */ 417e5136372SVijay Mahadevan if (ise[2*dim-1] <= nele-nghost) { 418e5136372SVijay Mahadevan for (k = (dim==3?ise[5]:ise[4]); k < (dim==3?ise[5]+nghost:std::max(ise[5],1)); k++) { 419e5136372SVijay Mahadevan for (j = (dim==2?ise[3]:ise[2]); j < (dim==2?ise[3]+nghost:std::max(ise[3],1)); j++) { 420*c3dd00c7SVijay Mahadevan for (i = (dim>1?ise[0]:ise[1]); i < (dim>1?std::max(ise[1],1):ise[1]+nghost); i++) { 421*c3dd00c7SVijay Mahadevan DMMoab_SetElementConnectivity_Private(useSimplex, dim, etype, &ecount, vpere, nele, i, j, k, vfirst, connectivity); 422e5136372SVijay Mahadevan } 423e5136372SVijay Mahadevan } 424e5136372SVijay Mahadevan } 425e5136372SVijay Mahadevan } 426e5136372SVijay Mahadevan 427e5136372SVijay Mahadevan merr = readMeshIface->update_adjacencies(efirst,locnele+ghnele,vpere,connectivity);MBERRNM(merr); 42851d15aeeSVijay Mahadevan 429a4717116SVijay 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. 430a4717116SVijay Mahadevan first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */ 43151d15aeeSVijay Mahadevan merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr); 432a4717116SVijay Mahadevan 433e5136372SVijay Mahadevan if (locnele+ghnele != (int) ownedelms.size()) 434e5136372SVijay Mahadevan SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele+ghnele,ownedelms.size()); 435e5136372SVijay Mahadevan else if(locnpts+ghnpts != (int) ownedvtx.size()) 436e5136372SVijay Mahadevan SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,ownedvtx.size()); 43751d15aeeSVijay Mahadevan else 438e427d9c9SVijay Mahadevan PetscInfo2(NULL, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size()); 43951d15aeeSVijay Mahadevan 4403a4aeca1SVijay Mahadevan /* lets create some sets */ 4413a4aeca1SVijay 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); 4423a4aeca1SVijay Mahadevan 4433a4aeca1SVijay Mahadevan merr = mbiface->create_meshset(moab::MESHSET_SET, regionset);MBERRNM(merr); 4443a4aeca1SVijay Mahadevan merr = mbiface->add_entities(regionset, ownedelms);MBERRNM(merr); 4453a4aeca1SVijay Mahadevan merr = mbiface->tag_set_data(geom_tag, ®ionset, 1, &dmmoab->dim);MBERRNM(merr); 446b5410836SVijay Mahadevan merr = mbiface->add_parent_child(dmmoab->fileset,regionset);MBERRNM(merr); 4473a4aeca1SVijay Mahadevan merr = mbiface->unite_meshset(dmmoab->fileset, regionset);MBERRNM(merr); 4483a4aeca1SVijay Mahadevan 4493a4aeca1SVijay Mahadevan merr = mbiface->create_meshset(moab::MESHSET_SET, vtxset);MBERRNM(merr); 4503a4aeca1SVijay Mahadevan merr = mbiface->add_entities(vtxset, ownedvtx);MBERRNM(merr); 4513a4aeca1SVijay Mahadevan merr = mbiface->add_parent_child(dmmoab->fileset,vtxset);MBERRNM(merr); 4523a4aeca1SVijay Mahadevan merr = mbiface->unite_meshset(dmmoab->fileset, vtxset);MBERRNM(merr); 4533a4aeca1SVijay Mahadevan 4543a4aeca1SVijay Mahadevan if (build_adjacencies) { 4553a4aeca1SVijay Mahadevan // generate all lower dimensional adjacencies 4563a4aeca1SVijay Mahadevan merr = mbiface->get_adjacencies( ownedelms, dim-1, true, adj, moab::Interface::UNION );MBERRNM(merr); 4573a4aeca1SVijay Mahadevan merr = dmmoab->pcomm->get_part_entities(dim2, dim-1);MBERRNM(merr); 4583a4aeca1SVijay Mahadevan adj.merge(dim2); 4593a4aeca1SVijay Mahadevan 4603a4aeca1SVijay Mahadevan /* create face sets */ 4613a4aeca1SVijay Mahadevan merr = mbiface->create_meshset(moab::MESHSET_SET, faceset);MBERRNM(merr); 4623a4aeca1SVijay Mahadevan merr = mbiface->add_entities(faceset, adj);MBERRNM(merr); 4633a4aeca1SVijay Mahadevan merr = mbiface->add_parent_child(dmmoab->fileset,faceset);MBERRNM(merr); 4643a4aeca1SVijay Mahadevan i=dim-1; 4653a4aeca1SVijay Mahadevan merr = mbiface->tag_set_data(geom_tag, &faceset, 1, &i);MBERRNM(merr); 4663a4aeca1SVijay Mahadevan merr = mbiface->unite_meshset(dmmoab->fileset, faceset);MBERRNM(merr); 467b8ecf6d3SVijay Mahadevan PetscInfo2(NULL, "Found %d %d-Dim quantities.\n", adj.size(), dim-1); 4683a4aeca1SVijay Mahadevan 4693a4aeca1SVijay Mahadevan if (dim > 2) { 4703a4aeca1SVijay Mahadevan dim2.clear(); 4713a4aeca1SVijay Mahadevan /* create edge sets, if appropriate i.e., if dim=3 */ 4723a4aeca1SVijay Mahadevan merr = mbiface->create_meshset(moab::MESHSET_SET, edgeset);MBERRNM(merr); 4733a4aeca1SVijay Mahadevan merr = mbiface->get_adjacencies(adj, dim-1, true, dim2, moab::Interface::UNION );MBERRNM(merr); 4743a4aeca1SVijay Mahadevan merr = mbiface->add_entities(edgeset, dim2);MBERRNM(merr); 4753a4aeca1SVijay Mahadevan merr = mbiface->add_parent_child(dmmoab->fileset,edgeset);MBERRNM(merr); 4763a4aeca1SVijay Mahadevan i=dim-2; 4773a4aeca1SVijay Mahadevan merr = mbiface->tag_set_data(geom_tag, &edgeset, 1, &i);MBERRNM(merr); 4783a4aeca1SVijay Mahadevan merr = mbiface->unite_meshset(dmmoab->fileset, edgeset);MBERRNM(merr); 479b8ecf6d3SVijay Mahadevan PetscInfo2(NULL, "Found %d %d-Dim quantities.\n", adj.size(), dim-2); 4803a4aeca1SVijay Mahadevan } 4813a4aeca1SVijay Mahadevan } 4823a4aeca1SVijay Mahadevan 483a4717116SVijay Mahadevan /* check the handles */ 484a4717116SVijay Mahadevan merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr); 48551d15aeeSVijay Mahadevan 486a4717116SVijay Mahadevan /* resolve the shared entities by exchanging information to adjacent processors */ 487a4717116SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr); 488ce1fd009SVijay Mahadevan merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,NULL,&id_tag);MBERRV(mbiface,merr); 48951d15aeeSVijay Mahadevan 4903a4aeca1SVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim,0,user_nghost,dim,true,false,&dmmoab->fileset);MBERRV(mbiface,merr); 491c8d5365dSVijay Mahadevan 492e427d9c9SVijay Mahadevan /* Reassign global IDs on all entities. */ 493e427d9c9SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,false,true,false);MBERRNM(merr); 494e427d9c9SVijay Mahadevan 495a4717116SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 496a4717116SVijay Mahadevan on all of the ghost vertexes */ 497c8d5365dSVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr); 498c8d5365dSVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, ownedelms);MBERRNM(merr); 499c8d5365dSVijay Mahadevan 500a4717116SVijay Mahadevan merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr); 501a4717116SVijay Mahadevan merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr); 50251d15aeeSVijay Mahadevan PetscFunctionReturn(0); 50351d15aeeSVijay Mahadevan } 50451d15aeeSVijay Mahadevan 50551d15aeeSVijay Mahadevan 506a4717116SVijay Mahadevan #undef __FUNCT__ 507a4717116SVijay Mahadevan #define __FUNCT__ "DMMoab_GetReadOptions_Private" 5082e4e7c01SVijay 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) 50951d15aeeSVijay Mahadevan { 510f90c3b0eSVijay Mahadevan char ropts[PETSC_MAX_PATH_LEN]; 51161eb6e27SVijay Mahadevan char ropts_par[PETSC_MAX_PATH_LEN]; 51261eb6e27SVijay Mahadevan char ropts_dbg[PETSC_MAX_PATH_LEN]; 513f90c3b0eSVijay Mahadevan PetscErrorCode ierr; 51451d15aeeSVijay Mahadevan 515a4717116SVijay Mahadevan PetscFunctionBegin; 51661eb6e27SVijay Mahadevan ierr = PetscMemzero(&ropts,PETSC_MAX_PATH_LEN);CHKERRQ(ierr); 51761eb6e27SVijay Mahadevan ierr = PetscMemzero(&ropts_par,PETSC_MAX_PATH_LEN);CHKERRQ(ierr); 51861eb6e27SVijay Mahadevan ierr = PetscMemzero(&ropts_dbg,PETSC_MAX_PATH_LEN);CHKERRQ(ierr); 51961eb6e27SVijay Mahadevan 520e23c60ebSVijay Mahadevan /* do parallel read unless using only one processor */ 521a4717116SVijay Mahadevan if (numproc > 1) { 52261eb6e27SVijay Mahadevan ierr = PetscSNPrintf(ropts_par, sizeof(ropts_par), "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); 5232e4e7c01SVijay Mahadevan } 5242e4e7c01SVijay Mahadevan 525c8d5365dSVijay Mahadevan if (dbglevel) { 526f90c3b0eSVijay Mahadevan if (numproc>1) { 52761eb6e27SVijay Mahadevan ierr = PetscSNPrintf(ropts_dbg, sizeof(ropts_dbg), "%sCPUTIME;DEBUG_IO=%d;DEBUG_PIO=%d;",dbglevel,dbglevel);CHKERRQ(ierr); 52861eb6e27SVijay Mahadevan } 52961eb6e27SVijay Mahadevan else { 53061eb6e27SVijay Mahadevan ierr = PetscSNPrintf(ropts_dbg, sizeof(ropts_dbg), "%sCPUTIME;DEBUG_IO=%d;",dbglevel);CHKERRQ(ierr); 531f90c3b0eSVijay Mahadevan } 532c8d5365dSVijay Mahadevan } 53351d15aeeSVijay Mahadevan 53461eb6e27SVijay Mahadevan ierr = PetscSNPrintf(ropts, sizeof(ropts), "%s%s%s%s",ropts_par,ropts_dbg,(extra_opts?extra_opts:""),(dm_opts?dm_opts:""));CHKERRQ(ierr); 535f90c3b0eSVijay Mahadevan *read_opts = ropts; 53651d15aeeSVijay Mahadevan PetscFunctionReturn(0); 53751d15aeeSVijay Mahadevan } 53851d15aeeSVijay Mahadevan 53951d15aeeSVijay Mahadevan 54051d15aeeSVijay Mahadevan #undef __FUNCT__ 541a4717116SVijay Mahadevan #define __FUNCT__ "DMMoabLoadFromFile" 542aa859218SVijay Mahadevan /*@ 543aa859218SVijay Mahadevan DMMoabLoadFromFile - Creates a DM object by loading the mesh from a user specified file. 544aa859218SVijay Mahadevan 545aa859218SVijay Mahadevan Collective on MPI_Comm 546aa859218SVijay Mahadevan 547aa859218SVijay Mahadevan Input Parameters: 548aa859218SVijay Mahadevan + comm - The communicator for the DM object 549aa859218SVijay Mahadevan . dim - The spatial dimension 550b8ecf6d3SVijay Mahadevan . filename - The name of the mesh file to be loaded 551b8ecf6d3SVijay Mahadevan . usrreadopts - The options string to read a MOAB mesh. 552aa859218SVijay Mahadevan Reference (Parallel Mesh Initialization: http://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo) 553aa859218SVijay Mahadevan 554aa859218SVijay Mahadevan Output Parameter: 555aa859218SVijay Mahadevan . dm - The DM object 556aa859218SVijay Mahadevan 557aa859218SVijay Mahadevan Level: beginner 558aa859218SVijay Mahadevan 559aa859218SVijay Mahadevan .keywords: DM, create 560b8ecf6d3SVijay Mahadevan 561aa859218SVijay Mahadevan .seealso: DMSetType(), DMCreate(), DMMoabCreateBoxMesh() 562aa859218SVijay Mahadevan @*/ 563a4717116SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm) 56451d15aeeSVijay Mahadevan { 565a4717116SVijay Mahadevan moab::ErrorCode merr; 5662e4e7c01SVijay Mahadevan PetscInt nprocs; 567a4717116SVijay Mahadevan DM_Moab *dmmoab; 568a4717116SVijay Mahadevan moab::Interface *mbiface; 569a4717116SVijay Mahadevan moab::ParallelComm *pcomm; 570a4717116SVijay Mahadevan moab::Range verts,elems; 5717023aa44SVijay Mahadevan const char *readopts; 572a4717116SVijay Mahadevan PetscErrorCode ierr; 57351d15aeeSVijay Mahadevan 57451d15aeeSVijay Mahadevan PetscFunctionBegin; 575a4717116SVijay Mahadevan PetscValidPointer(dm,5); 57651d15aeeSVijay Mahadevan 577a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 578a4717116SVijay Mahadevan ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 57951d15aeeSVijay Mahadevan 580a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 581a4717116SVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 582a4717116SVijay Mahadevan mbiface = dmmoab->mbiface; 583a4717116SVijay Mahadevan pcomm = dmmoab->pcomm; 584a4717116SVijay Mahadevan nprocs = pcomm->size(); 585aa859218SVijay Mahadevan /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */ 586c8d5365dSVijay Mahadevan dmmoab->dim = dim; 587a4717116SVijay Mahadevan 588b5410836SVijay Mahadevan /* create a file set to associate all entities in current mesh */ 589b5410836SVijay Mahadevan merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr); 590b5410836SVijay Mahadevan 591a4717116SVijay Mahadevan /* add mesh loading options specific to the DM */ 5922e4e7c01SVijay Mahadevan ierr = DMMoab_GetReadOptions_Private(dmmoab->partition_by_rank, nprocs, dim, dmmoab->read_mode, 5932e4e7c01SVijay Mahadevan dmmoab->rw_dbglevel, dmmoab->extra_read_options, usrreadopts, &readopts);CHKERRQ(ierr); 5947023aa44SVijay Mahadevan 595e5136372SVijay Mahadevan PetscInfo2(*dm, "Reading file %s with options: %s\n",filename,readopts); 596a4717116SVijay Mahadevan 597a4717116SVijay Mahadevan /* Load the mesh from a file. */ 5987023aa44SVijay Mahadevan merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr); 599a4717116SVijay Mahadevan 6006e40195eSVijay Mahadevan /* Reassign global IDs on all entities. */ 601e5136372SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr); 602e5136372SVijay Mahadevan 603e5136372SVijay Mahadevan /* load the local vertices */ 604e5136372SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr); 605e5136372SVijay Mahadevan /* load the local elements */ 606e5136372SVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr); 607e5136372SVijay Mahadevan 608e5136372SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 609e5136372SVijay Mahadevan on all of the ghost vertexes */ 610e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag,verts);MBERRV(mbiface,merr); 611e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag,elems);MBERRV(mbiface,merr); 612e5136372SVijay Mahadevan 613e5136372SVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim,0,1,0,true,true,&dmmoab->fileset);MBERRV(mbiface,merr); 6146e40195eSVijay Mahadevan 615a4717116SVijay Mahadevan merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr); 616a4717116SVijay Mahadevan 617e5136372SVijay Mahadevan PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size()); 61851d15aeeSVijay Mahadevan PetscFunctionReturn(0); 61951d15aeeSVijay Mahadevan } 62051d15aeeSVijay Mahadevan 621