xref: /petsc/src/dm/impls/moab/dmmoab.cxx (revision fc4180132621dce802f5cfc3099d116b2d2b65bb)
1032b8ab6SVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I  "petscdm.h"   I*/
2aa768e4cSTim Tautges #include <petsc-private/vecimpl.h> /*I  "petscdm.h"   I*/
31d72bce8STim Tautges 
41d72bce8STim Tautges #include <petscdmmoab.h>
588face26SJed Brown #include <MBTagConventions.hpp>
61cec0304SVijay Mahadevan #include <moab/Skinner.hpp>
7032b8ab6SVijay Mahadevan 
8853cdec3SJed Brown #undef __FUNCT__
9853cdec3SJed Brown #define __FUNCT__ "DMDestroy_Moab"
10853cdec3SJed Brown PetscErrorCode DMDestroy_Moab(DM dm)
11853cdec3SJed Brown {
12853cdec3SJed Brown   PetscErrorCode ierr;
13032b8ab6SVijay Mahadevan   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
14853cdec3SJed Brown 
15853cdec3SJed Brown   PetscFunctionBegin;
16853cdec3SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
17032b8ab6SVijay Mahadevan   if (dmmoab->icreatedinstance) {
18032b8ab6SVijay Mahadevan     delete dmmoab->mbiface;
19853cdec3SJed Brown   }
20032b8ab6SVijay Mahadevan   dmmoab->mbiface = NULL;
21032b8ab6SVijay Mahadevan   dmmoab->pcomm = NULL;
22032b8ab6SVijay Mahadevan   delete dmmoab->vlocal;
23032b8ab6SVijay Mahadevan   delete dmmoab->vowned;
24032b8ab6SVijay Mahadevan   delete dmmoab->vghost;
25032b8ab6SVijay Mahadevan   delete dmmoab->elocal;
26032b8ab6SVijay Mahadevan   delete dmmoab->eghost;
271cec0304SVijay Mahadevan   delete dmmoab->bndyvtx;
281cec0304SVijay Mahadevan   delete dmmoab->bndyfaces;
29*fc418013SVijay Mahadevan 
30*fc418013SVijay Mahadevan   ierr = PetscFree(dmmoab->gsindices);CHKERRQ(ierr);
31032b8ab6SVijay Mahadevan   ierr = VecScatterDestroy(&dmmoab->ltog_sendrecv);CHKERRQ(ierr);
32032b8ab6SVijay Mahadevan   ierr = ISLocalToGlobalMappingDestroy(&dmmoab->ltog_map);CHKERRQ(ierr);
33853cdec3SJed Brown   ierr = PetscFree(dm->data);CHKERRQ(ierr);
34853cdec3SJed Brown   PetscFunctionReturn(0);
35853cdec3SJed Brown }
36853cdec3SJed Brown 
37aa768e4cSTim Tautges #undef __FUNCT__
38032b8ab6SVijay Mahadevan #define __FUNCT__ "DMSetUp_Moab"
39032b8ab6SVijay Mahadevan PetscErrorCode DMSetUp_Moab(DM dm)
40032b8ab6SVijay Mahadevan {
41032b8ab6SVijay Mahadevan   PetscErrorCode          ierr;
42032b8ab6SVijay Mahadevan   moab::ErrorCode         merr;
43032b8ab6SVijay Mahadevan   Vec                     local, global;
44032b8ab6SVijay Mahadevan   IS                      from;
45032b8ab6SVijay Mahadevan   moab::Range::iterator   iter;
46*fc418013SVijay Mahadevan   PetscInt                i,j,bs,gsiz,lsiz;
47032b8ab6SVijay Mahadevan   DM_Moab                *dmmoab = (DM_Moab*)dm->data;
48*fc418013SVijay Mahadevan   PetscInt                totsize,local_min,local_max,global_min;
491cec0304SVijay Mahadevan   PetscSection            section;
50032b8ab6SVijay Mahadevan 
51032b8ab6SVijay Mahadevan   PetscFunctionBegin;
52032b8ab6SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
53032b8ab6SVijay Mahadevan   /* Get the local and shared vertices and cache it */
54032b8ab6SVijay Mahadevan   if (dmmoab->mbiface == PETSC_NULL || dmmoab->pcomm == PETSC_NULL) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ORDER, "Set the MOAB Interface and ParallelComm objects before calling SetUp.");
55032b8ab6SVijay Mahadevan 
561cec0304SVijay Mahadevan   /* Get the entities recursively in the current part of the mesh, if user did not set the local vertices explicitly */
57032b8ab6SVijay Mahadevan   if (dmmoab->vlocal->empty()) {
581cec0304SVijay Mahadevan     merr = dmmoab->mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,*dmmoab->vlocal,true);MBERRNM(merr);
59032b8ab6SVijay Mahadevan 
60032b8ab6SVijay Mahadevan     /* filter based on parallel status */
61*fc418013SVijay Mahadevan     merr = dmmoab->pcomm->filter_pstatus(*dmmoab->vlocal,PSTATUS_NOT_OWNED,PSTATUS_NOT,-1,dmmoab->vowned);MBERRNM(merr);
62032b8ab6SVijay Mahadevan     *dmmoab->vghost = moab::subtract(*dmmoab->vlocal, *dmmoab->vowned);
63032b8ab6SVijay Mahadevan 
64032b8ab6SVijay Mahadevan     dmmoab->nloc = dmmoab->vowned->size();
65032b8ab6SVijay Mahadevan     dmmoab->nghost = dmmoab->vghost->size();
66032b8ab6SVijay Mahadevan     ierr = MPI_Allreduce(&dmmoab->nloc, &dmmoab->n, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
67*fc418013SVijay Mahadevan 
68*fc418013SVijay Mahadevan #if 1
69*fc418013SVijay Mahadevan     if(dmmoab->pcomm->rank() || dmmoab->pcomm->size()==1) {
70*fc418013SVijay Mahadevan       PetscPrintf(PETSC_COMM_SELF, "Vertices: global: %D, local: %D", dmmoab->n, dmmoab->nloc+dmmoab->nghost);
71*fc418013SVijay Mahadevan       dmmoab->vlocal->print(0);
72*fc418013SVijay Mahadevan       PetscPrintf(PETSC_COMM_SELF, "Vertices: owned: %D", dmmoab->nloc);
73*fc418013SVijay Mahadevan       dmmoab->vowned->print(0);
74*fc418013SVijay Mahadevan       PetscPrintf(PETSC_COMM_SELF, "Vertices: ghost: %D", dmmoab->nghost);
75*fc418013SVijay Mahadevan       dmmoab->vghost->print(0);
76*fc418013SVijay Mahadevan     }
77*fc418013SVijay Mahadevan #endif
78032b8ab6SVijay Mahadevan   }
79032b8ab6SVijay Mahadevan 
80032b8ab6SVijay Mahadevan   /* get the information about the local elements in the mesh */
81032b8ab6SVijay Mahadevan   {
82032b8ab6SVijay Mahadevan     dmmoab->eghost->clear();
83*fc418013SVijay Mahadevan 
84*fc418013SVijay Mahadevan     /* first decipher the leading dimension */
85*fc418013SVijay Mahadevan     for (i=3;i>0;i--) {
86*fc418013SVijay Mahadevan       dmmoab->elocal->clear();
87*fc418013SVijay Mahadevan       merr = dmmoab->mbiface->get_entities_by_dimension(dmmoab->fileset, i, *dmmoab->elocal, true);CHKERRQ(merr);
88*fc418013SVijay Mahadevan 
89*fc418013SVijay Mahadevan       /* store the current mesh dimension */
90*fc418013SVijay Mahadevan       if (dmmoab->elocal->size()) {
91*fc418013SVijay Mahadevan         dmmoab->dim=i;
92*fc418013SVijay Mahadevan         break;
93*fc418013SVijay Mahadevan       }
94*fc418013SVijay Mahadevan     }
95*fc418013SVijay Mahadevan 
96032b8ab6SVijay Mahadevan     *dmmoab->eghost = *dmmoab->elocal;
97032b8ab6SVijay Mahadevan     merr = dmmoab->pcomm->filter_pstatus(*dmmoab->elocal,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
98032b8ab6SVijay Mahadevan     *dmmoab->eghost = moab::subtract(*dmmoab->eghost, *dmmoab->elocal);
99032b8ab6SVijay Mahadevan 
100032b8ab6SVijay Mahadevan     dmmoab->neleloc = dmmoab->elocal->size();
101032b8ab6SVijay Mahadevan     ierr = MPI_Allreduce(&dmmoab->neleloc, &dmmoab->nele, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
102032b8ab6SVijay Mahadevan   }
103032b8ab6SVijay Mahadevan 
104032b8ab6SVijay Mahadevan   bs = dmmoab->bs;
105032b8ab6SVijay Mahadevan   if (!dmmoab->ltog_tag) {
106db66d124SVijay Mahadevan     /* Get the global ID tag. The global ID tag is applied to each
107db66d124SVijay Mahadevan        vertex. It acts as an global identifier which MOAB uses to
108db66d124SVijay Mahadevan        assemble the individual pieces of the mesh */
109032b8ab6SVijay Mahadevan     merr = dmmoab->mbiface->tag_get_handle(GLOBAL_ID_TAG_NAME, dmmoab->ltog_tag);MBERRNM(merr);
110032b8ab6SVijay Mahadevan   }
111032b8ab6SVijay Mahadevan 
112032b8ab6SVijay Mahadevan   totsize=dmmoab->vlocal->size();
113*fc418013SVijay Mahadevan   ierr = PetscMalloc(totsize*sizeof(PetscInt), &dmmoab->gsindices);CHKERRQ(ierr);
1141cec0304SVijay Mahadevan   {
115032b8ab6SVijay Mahadevan     /* first get the local indices */
116*fc418013SVijay Mahadevan     merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,*dmmoab->vowned,&dmmoab->gsindices[0]);MBERRNM(merr);
1174a40b570SVijay Mahadevan     /* next get the ghosted indices */
118*fc418013SVijay Mahadevan     if (dmmoab->nghost) {
119*fc418013SVijay Mahadevan       merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,*dmmoab->vghost,&dmmoab->gsindices[dmmoab->nloc]);MBERRNM(merr);
1201cec0304SVijay Mahadevan     }
1216e40195eSVijay Mahadevan 
1226e40195eSVijay Mahadevan     /* find out the local and global minima of GLOBAL_ID */
123*fc418013SVijay Mahadevan     local_min=local_max=dmmoab->gsindices[0];
124*fc418013SVijay Mahadevan     for (i=1; i<totsize; ++i) {
125*fc418013SVijay Mahadevan //      if (dmmoab->pcomm->rank())
126*fc418013SVijay Mahadevan //        PetscPrintf(PETSC_COMM_SELF, "[%D] gsindices[%D] = %D\n", dmmoab->pcomm->rank(), i, dmmoab->gsindices[i]);
127*fc418013SVijay Mahadevan       if(local_min>dmmoab->gsindices[i]) local_min=dmmoab->gsindices[i];
128*fc418013SVijay Mahadevan       if(local_max<dmmoab->gsindices[i]) local_max=dmmoab->gsindices[i];
129*fc418013SVijay Mahadevan     }
1306e40195eSVijay Mahadevan 
1316e40195eSVijay Mahadevan     ierr = MPI_Allreduce(&local_min, &global_min, 1, MPI_INT, MPI_MIN, ((PetscObject)dm)->comm);CHKERRQ(ierr);
132*fc418013SVijay Mahadevan     PetscInfo3(dm, "GLOBAL_ID: Local minima - %D, Local maxima - %D, Global minima - %D.\n", local_min, local_max, global_min);
133*fc418013SVijay Mahadevan //    PetscPrintf(PETSC_COMM_SELF, "[%D] GLOBAL_ID: Local minima - %D, Local maxima - %D, Global minima - %D.\n", dmmoab->pcomm->rank(), local_min, local_max, global_min);
1341cec0304SVijay Mahadevan   }
135032b8ab6SVijay Mahadevan 
1361cec0304SVijay Mahadevan   {
1371cec0304SVijay Mahadevan     ierr = PetscSectionCreate(((PetscObject)dm)->comm, &section);CHKERRQ(ierr);
1381cec0304SVijay Mahadevan     ierr = PetscSectionSetNumFields(section, dmmoab->nfields);CHKERRQ(ierr);
139*fc418013SVijay Mahadevan //    ierr = PetscSectionSetChart(section, dmmoab->gsindices[0], dmmoab->gsindices[dmmoab->nloc-1]+1);CHKERRQ(ierr);
140*fc418013SVijay Mahadevan     ierr = PetscSectionSetChart(section, local_min, local_max+1);CHKERRQ(ierr);
141*fc418013SVijay Mahadevan     for (j=0; j<totsize; ++j) {
142*fc418013SVijay Mahadevan       PetscInt locgid = dmmoab->gsindices[j];
1431cec0304SVijay Mahadevan       for (i=0; i < dmmoab->nfields; ++i) {
1441cec0304SVijay Mahadevan         ierr = PetscSectionSetFieldName(section, i, dmmoab->fields[i]);CHKERRQ(ierr);
1451cec0304SVijay Mahadevan         if (bs>1) {
146*fc418013SVijay Mahadevan           ierr = PetscSectionSetFieldDof(section, locgid, i, (locgid-global_min)*dmmoab->nfields+i);CHKERRQ(ierr);
147*fc418013SVijay Mahadevan           ierr = PetscSectionSetFieldOffset(section, locgid, i, (locgid-global_min)*dmmoab->nfields);
1481cec0304SVijay Mahadevan         }
1491cec0304SVijay Mahadevan         else {
150*fc418013SVijay Mahadevan           ierr = PetscSectionSetFieldDof(section, locgid, i, totsize*i+locgid-global_min);CHKERRQ(ierr);
151*fc418013SVijay Mahadevan           ierr = PetscSectionSetFieldOffset(section, locgid, i, i*totsize);
152*fc418013SVijay Mahadevan           PetscPrintf(PETSC_COMM_SELF, "[%D] Local_GID = %D, FDOF = %D, OFF = %D.\n", dmmoab->pcomm->rank(), locgid, totsize*i+locgid-global_min, totsize );
1531cec0304SVijay Mahadevan         }
1541cec0304SVijay Mahadevan       }
155*fc418013SVijay Mahadevan       ierr = PetscSectionSetDof(section, locgid, dmmoab->nfields);CHKERRQ(ierr);
1561cec0304SVijay Mahadevan     }
1571cec0304SVijay Mahadevan     ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
1581cec0304SVijay Mahadevan     ierr = DMSetDefaultSection(dm, section);CHKERRQ(ierr);
1591cec0304SVijay Mahadevan   }
1601cec0304SVijay Mahadevan 
1611cec0304SVijay Mahadevan   {
162*fc418013SVijay Mahadevan     for (i=0; i<totsize; ++i) {
163*fc418013SVijay Mahadevan       dmmoab->gsindices[i]-=global_min;   /* zero based index needed for IS */
164*fc418013SVijay Mahadevan //      if (dmmoab->pcomm->rank())
165*fc418013SVijay Mahadevan //        PetscPrintf(PETSC_COMM_SELF, "[%D] modified gsindices[%D] = %D\n", dmmoab->pcomm->rank(), i, dmmoab->gsindices[i]);
166*fc418013SVijay Mahadevan     }
167*fc418013SVijay Mahadevan 
168032b8ab6SVijay Mahadevan     /* Create Global to Local Vector Scatter Context */
169032b8ab6SVijay Mahadevan     ierr = DMCreateGlobalVector_Moab(dm, &global);CHKERRQ(ierr);
170032b8ab6SVijay Mahadevan     ierr = DMCreateLocalVector_Moab(dm, &local);CHKERRQ(ierr);
171032b8ab6SVijay Mahadevan 
172032b8ab6SVijay Mahadevan     /* global to local must retrieve ghost points */
173*fc418013SVijay Mahadevan     ierr = ISCreateBlock(((PetscObject)dm)->comm,bs,totsize,&dmmoab->gsindices[0],PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
174032b8ab6SVijay Mahadevan 
175db66d124SVijay Mahadevan     ierr = VecGetLocalSize(global,&gsiz);CHKERRQ(ierr);
176db66d124SVijay Mahadevan     ierr = VecGetLocalSize(local,&lsiz);CHKERRQ(ierr);
177032b8ab6SVijay Mahadevan 
178032b8ab6SVijay Mahadevan     ierr = VecScatterCreate(local,from,global,from,&dmmoab->ltog_sendrecv);CHKERRQ(ierr);
179032b8ab6SVijay Mahadevan     ierr = ISDestroy(&from);CHKERRQ(ierr);
180032b8ab6SVijay Mahadevan     ierr = VecDestroy(&local);CHKERRQ(ierr);
181032b8ab6SVijay Mahadevan     ierr = VecDestroy(&global);CHKERRQ(ierr);
182032b8ab6SVijay Mahadevan   }
183032b8ab6SVijay Mahadevan 
1841cec0304SVijay Mahadevan   /* skin the boundary and store nodes */
1851cec0304SVijay Mahadevan   {
1861cec0304SVijay Mahadevan     // get the skin vertices of those faces and mark them as fixed; we don't want to fix the vertices on a
1871cec0304SVijay Mahadevan     // part boundary, but since we exchanged a layer of ghost faces, those vertices aren't on the skin locally
1881cec0304SVijay Mahadevan     // ok to mark non-owned skin vertices too, I won't move those anyway
1891cec0304SVijay Mahadevan     // use MOAB's skinner class to find the skin
1901cec0304SVijay Mahadevan     moab::Skinner skinner(dmmoab->mbiface);
1911cec0304SVijay Mahadevan     dmmoab->bndyvtx = new moab::Range();
1921cec0304SVijay Mahadevan     dmmoab->bndyfaces = new moab::Range();
1931cec0304SVijay Mahadevan     merr = skinner.find_skin(dmmoab->fileset, *dmmoab->elocal, true, *dmmoab->bndyvtx);MBERRNM(merr); // 'true' param indicates we want vertices back, not faces
1941cec0304SVijay Mahadevan     merr = skinner.find_skin(dmmoab->fileset, *dmmoab->elocal, false, *dmmoab->bndyfaces);MBERRNM(merr); // 'false' param indicates we want faces back, not vertices
1956e40195eSVijay Mahadevan     PetscInfo2(dm, "Found %D boundary vertices and %D faces.\n", dmmoab->bndyvtx->size(), dmmoab->bndyvtx->size());
1961cec0304SVijay Mahadevan   }
197032b8ab6SVijay Mahadevan   PetscFunctionReturn(0);
198032b8ab6SVijay Mahadevan }
199032b8ab6SVijay Mahadevan 
2001cec0304SVijay Mahadevan 
201032b8ab6SVijay Mahadevan #undef __FUNCT__
202aa768e4cSTim Tautges #define __FUNCT__ "DMCreate_Moab"
203853cdec3SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Moab(DM dm)
204aa768e4cSTim Tautges {
205aa768e4cSTim Tautges   PetscErrorCode ierr;
206aa768e4cSTim Tautges 
207aa768e4cSTim Tautges   PetscFunctionBegin;
208aa768e4cSTim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
209032b8ab6SVijay Mahadevan   ierr = PetscNewLog(dm,&dm->data);CHKERRQ(ierr);
210032b8ab6SVijay Mahadevan 
211032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->bs = 1;
212032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->n = 0;
213032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->nloc = 0;
214032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->nele = 0;
215032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->neleloc = 0;
216032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->nghost = 0;
217032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->ltog_map = PETSC_NULL;
218032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->ltog_sendrecv = PETSC_NULL;
219032b8ab6SVijay Mahadevan 
220032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->vlocal = new moab::Range();
221032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->vowned = new moab::Range();
222032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->vghost = new moab::Range();
223032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->elocal = new moab::Range();
224032b8ab6SVijay Mahadevan   ((DM_Moab*)dm->data)->eghost = new moab::Range();
225aa768e4cSTim Tautges 
22697ea90e6SJed Brown   dm->ops->createglobalvector              = DMCreateGlobalVector_Moab;
22797ea90e6SJed Brown   dm->ops->createlocalvector               = DMCreateLocalVector_Moab;
228032b8ab6SVijay Mahadevan   dm->ops->creatematrix                    = DMCreateMatrix_Moab;
229032b8ab6SVijay Mahadevan   dm->ops->setup                           = DMSetUp_Moab;
23097ea90e6SJed Brown   dm->ops->destroy                         = DMDestroy_Moab;
231032b8ab6SVijay Mahadevan   dm->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Moab;
232032b8ab6SVijay Mahadevan   dm->ops->globaltolocalend                = DMGlobalToLocalEnd_Moab;
233032b8ab6SVijay Mahadevan   dm->ops->localtoglobalbegin              = DMLocalToGlobalBegin_Moab;
234032b8ab6SVijay Mahadevan   dm->ops->localtoglobalend                = DMLocalToGlobalEnd_Moab;
235aa768e4cSTim Tautges   PetscFunctionReturn(0);
236aa768e4cSTim Tautges }
237fd349b41STim Tautges 
238fd349b41STim Tautges #undef __FUNCT__
2391d72bce8STim Tautges #define __FUNCT__ "DMMoabCreate"
2401d72bce8STim Tautges /*@
2411d72bce8STim Tautges   DMMoabCreate - Creates a DMMoab object, which encapsulates a moab instance
2421d72bce8STim Tautges 
2431d72bce8STim Tautges   Collective on MPI_Comm
2441d72bce8STim Tautges 
2451d72bce8STim Tautges   Input Parameter:
2461d72bce8STim Tautges . comm - The communicator for the DMMoab object
2471d72bce8STim Tautges 
2481d72bce8STim Tautges   Output Parameter:
249032b8ab6SVijay Mahadevan . dmb  - The DMMoab object
2501d72bce8STim Tautges 
2511d72bce8STim Tautges   Level: beginner
2521d72bce8STim Tautges 
2531d72bce8STim Tautges .keywords: DMMoab, create
2541d72bce8STim Tautges @*/
255032b8ab6SVijay Mahadevan PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *dmb)
2561d72bce8STim Tautges {
2571d72bce8STim Tautges   PetscErrorCode ierr;
2581d72bce8STim Tautges 
2591d72bce8STim Tautges   PetscFunctionBegin;
260032b8ab6SVijay Mahadevan   PetscValidPointer(dmb,2);
261032b8ab6SVijay Mahadevan   ierr = DMCreate(comm, dmb);CHKERRQ(ierr);
262032b8ab6SVijay Mahadevan   ierr = DMSetType(*dmb, DMMOAB);CHKERRQ(ierr);
2631d72bce8STim Tautges   PetscFunctionReturn(0);
2641d72bce8STim Tautges }
2651d72bce8STim Tautges 
2661d72bce8STim Tautges #undef __FUNCT__
267aa768e4cSTim Tautges #define __FUNCT__ "DMMoabCreateMoab"
2681d72bce8STim Tautges /*@
269a4d2169cSTim Tautges   DMMoabCreate - Creates a DMMoab object, optionally from an instance and other data
2701d72bce8STim Tautges 
2711d72bce8STim Tautges   Collective on MPI_Comm
2721d72bce8STim Tautges 
2731d72bce8STim Tautges   Input Parameter:
2741d72bce8STim Tautges . comm - The communicator for the DMMoab object
275032b8ab6SVijay Mahadevan . mbiface - (ptr to) the MOAB Instance; if passed in NULL, MOAB instance is created inside PETSc, and destroyed
276a4d2169cSTim Tautges          along with the DMMoab
277a4d2169cSTim Tautges . pcomm - (ptr to) a ParallelComm; if NULL, creates one internally for the whole communicator
2781d72bce8STim Tautges . ltog_tag - A tag to use to retrieve global id for an entity; if 0, will use GLOBAL_ID_TAG_NAME/tag
2791d72bce8STim Tautges . range - If non-NULL, contains range of entities to which DOFs will be assigned
2801d72bce8STim Tautges 
2811d72bce8STim Tautges   Output Parameter:
282032b8ab6SVijay Mahadevan . dmb  - The DMMoab object
2831d72bce8STim Tautges 
284032b8ab6SVijay Mahadevan   Level: intermediate
2851d72bce8STim Tautges 
2861d72bce8STim Tautges .keywords: DMMoab, create
2871d72bce8STim Tautges @*/
288032b8ab6SVijay Mahadevan PetscErrorCode DMMoabCreateMoab(MPI_Comm comm, moab::Interface *mbiface, moab::ParallelComm *pcomm, moab::Tag *ltog_tag, moab::Range *range, DM *dmb)
2891d72bce8STim Tautges {
2901d72bce8STim Tautges   PetscErrorCode ierr;
291032b8ab6SVijay Mahadevan   moab::ErrorCode merr;
2921cec0304SVijay Mahadevan   moab::EntityHandle partnset;
2931cec0304SVijay Mahadevan   PetscInt rank, nprocs;
294853cdec3SJed Brown   DM_Moab        *dmmoab;
2951d72bce8STim Tautges 
2961d72bce8STim Tautges   PetscFunctionBegin;
297032b8ab6SVijay Mahadevan   PetscValidPointer(dmb,6);
298032b8ab6SVijay Mahadevan   ierr = DMMoabCreate(comm, dmb);CHKERRQ(ierr);
299032b8ab6SVijay Mahadevan   dmmoab = (DM_Moab*)(*dmb)->data;
300a4d2169cSTim Tautges 
301a4d2169cSTim Tautges   if (!mbiface) {
30272ff976dSVijay Mahadevan     dmmoab->mbiface = new moab::Core();
3037d89fc02STim Tautges     dmmoab->icreatedinstance = PETSC_TRUE;
3041d72bce8STim Tautges   }
3051cec0304SVijay Mahadevan   else {
3061cec0304SVijay Mahadevan     dmmoab->mbiface = mbiface;
3077d89fc02STim Tautges     dmmoab->icreatedinstance = PETSC_FALSE;
3081cec0304SVijay Mahadevan   }
3091cec0304SVijay Mahadevan 
3101cec0304SVijay Mahadevan   /* create a fileset to store the hierarchy of entities belonging to current DM */
3111cec0304SVijay Mahadevan   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_ORDERED, dmmoab->fileset);MBERR("Creating file set failed", merr);
3127d89fc02STim Tautges 
313a4d2169cSTim Tautges   if (!pcomm) {
314032b8ab6SVijay Mahadevan     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
315032b8ab6SVijay Mahadevan     ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr);
316032b8ab6SVijay Mahadevan 
317db66d124SVijay Mahadevan     /* Create root sets for each mesh.  Then pass these
318db66d124SVijay Mahadevan        to the load_file functions to be populated. */
31972ff976dSVijay Mahadevan     merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, partnset);
320032b8ab6SVijay Mahadevan     MBERR("Creating partition set failed", merr);
321032b8ab6SVijay Mahadevan 
322db66d124SVijay Mahadevan     /* Create the parallel communicator object with the partition handle associated with MOAB */
32372ff976dSVijay Mahadevan     dmmoab->pcomm = moab::ParallelComm::get_pcomm(dmmoab->mbiface, partnset, &comm);
32472ff976dSVijay Mahadevan   }
32572ff976dSVijay Mahadevan   else {
32672ff976dSVijay Mahadevan     ierr = DMMoabSetParallelComm(*dmb, pcomm);CHKERRQ(ierr);
327032b8ab6SVijay Mahadevan   }
328032b8ab6SVijay Mahadevan 
3294973de03SVijay Mahadevan   /* do the remaining initializations for DMMoab */
3304973de03SVijay Mahadevan   dmmoab->bs       = 1;
3314973de03SVijay Mahadevan 
3324973de03SVijay Mahadevan   /* set global ID tag handle */
333032b8ab6SVijay Mahadevan   if (!ltog_tag) {
3344973de03SVijay Mahadevan     merr = dmmoab->mbiface->tag_get_handle(GLOBAL_ID_TAG_NAME, dmmoab->ltog_tag);MBERRNM(merr);
335032b8ab6SVijay Mahadevan   }
336032b8ab6SVijay Mahadevan   else {
337032b8ab6SVijay Mahadevan     ierr = DMMoabSetLocalToGlobalTag(*dmb, *ltog_tag);CHKERRQ(ierr);
338a4d2169cSTim Tautges   }
339a4d2169cSTim Tautges 
3404973de03SVijay Mahadevan   /* set the local range of entities (vertices) of interest */
341a4d2169cSTim Tautges   if (range) {
3425eb88e9dSVijay Mahadevan     ierr = DMMoabSetLocalVertices(*dmb, range);CHKERRQ(ierr);
343a4d2169cSTim Tautges   }
3441d72bce8STim Tautges   PetscFunctionReturn(0);
3451d72bce8STim Tautges }
3461d72bce8STim Tautges 
3471d72bce8STim Tautges #undef __FUNCT__
3481d72bce8STim Tautges #define __FUNCT__ "DMMoabSetParallelComm"
349aa768e4cSTim Tautges /*@
350aa768e4cSTim Tautges   DMMoabSetParallelComm - Set the ParallelComm used with this DMMoab
351aa768e4cSTim Tautges 
352aa768e4cSTim Tautges   Collective on MPI_Comm
353aa768e4cSTim Tautges 
354aa768e4cSTim Tautges   Input Parameter:
355aa768e4cSTim Tautges . dm    - The DMMoab object being set
356aa768e4cSTim Tautges . pcomm - The ParallelComm being set on the DMMoab
357aa768e4cSTim Tautges 
358aa768e4cSTim Tautges   Level: beginner
359aa768e4cSTim Tautges 
360aa768e4cSTim Tautges .keywords: DMMoab, create
361aa768e4cSTim Tautges @*/
3621d72bce8STim Tautges PetscErrorCode DMMoabSetParallelComm(DM dm,moab::ParallelComm *pcomm)
3631d72bce8STim Tautges {
364032b8ab6SVijay Mahadevan   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
365032b8ab6SVijay Mahadevan 
3661d72bce8STim Tautges   PetscFunctionBegin;
3671d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3681cec0304SVijay Mahadevan   PetscValidPointer(pcomm,2);
369032b8ab6SVijay Mahadevan   dmmoab->pcomm = pcomm;
370032b8ab6SVijay Mahadevan   dmmoab->mbiface = pcomm->get_moab();
371032b8ab6SVijay Mahadevan   dmmoab->icreatedinstance = PETSC_FALSE;
3721d72bce8STim Tautges   PetscFunctionReturn(0);
3731d72bce8STim Tautges }
3741d72bce8STim Tautges 
3751d72bce8STim Tautges 
3761d72bce8STim Tautges #undef __FUNCT__
3771d72bce8STim Tautges #define __FUNCT__ "DMMoabGetParallelComm"
378aa768e4cSTim Tautges /*@
379aa768e4cSTim Tautges   DMMoabGetParallelComm - Get the ParallelComm used with this DMMoab
380aa768e4cSTim Tautges 
381aa768e4cSTim Tautges   Collective on MPI_Comm
382aa768e4cSTim Tautges 
383aa768e4cSTim Tautges   Input Parameter:
384aa768e4cSTim Tautges . dm    - The DMMoab object being set
385aa768e4cSTim Tautges 
386aa768e4cSTim Tautges   Output Parameter:
387aa768e4cSTim Tautges . pcomm - The ParallelComm for the DMMoab
388aa768e4cSTim Tautges 
389aa768e4cSTim Tautges   Level: beginner
390aa768e4cSTim Tautges 
391aa768e4cSTim Tautges .keywords: DMMoab, create
392aa768e4cSTim Tautges @*/
3931d72bce8STim Tautges PetscErrorCode DMMoabGetParallelComm(DM dm,moab::ParallelComm **pcomm)
3941d72bce8STim Tautges {
3951d72bce8STim Tautges   PetscFunctionBegin;
3961d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
397032b8ab6SVijay Mahadevan   *pcomm = ((DM_Moab*)(dm)->data)->pcomm;
3981d72bce8STim Tautges   PetscFunctionReturn(0);
3991d72bce8STim Tautges }
4001d72bce8STim Tautges 
4011d72bce8STim Tautges 
4021d72bce8STim Tautges #undef __FUNCT__
4031d72bce8STim Tautges #define __FUNCT__ "DMMoabSetInterface"
404aa768e4cSTim Tautges /*@
405aa768e4cSTim Tautges   DMMoabSetInterface - Set the MOAB instance used with this DMMoab
406aa768e4cSTim Tautges 
407aa768e4cSTim Tautges   Collective on MPI_Comm
408aa768e4cSTim Tautges 
409aa768e4cSTim Tautges   Input Parameter:
410aa768e4cSTim Tautges . dm      - The DMMoab object being set
411aa768e4cSTim Tautges . mbiface - The MOAB instance being set on this DMMoab
412aa768e4cSTim Tautges 
413aa768e4cSTim Tautges   Level: beginner
414aa768e4cSTim Tautges 
415aa768e4cSTim Tautges .keywords: DMMoab, create
416aa768e4cSTim Tautges @*/
417a4d2169cSTim Tautges PetscErrorCode DMMoabSetInterface(DM dm,moab::Interface *mbiface)
4181d72bce8STim Tautges {
419032b8ab6SVijay Mahadevan   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
420032b8ab6SVijay Mahadevan 
4211d72bce8STim Tautges   PetscFunctionBegin;
4221d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4231cec0304SVijay Mahadevan   PetscValidPointer(mbiface,2);
424032b8ab6SVijay Mahadevan   dmmoab->pcomm = NULL;
425032b8ab6SVijay Mahadevan   dmmoab->mbiface = mbiface;
426032b8ab6SVijay Mahadevan   dmmoab->icreatedinstance = PETSC_FALSE;
4271d72bce8STim Tautges   PetscFunctionReturn(0);
4281d72bce8STim Tautges }
4291d72bce8STim Tautges 
4301d72bce8STim Tautges 
4311d72bce8STim Tautges #undef __FUNCT__
4321d72bce8STim Tautges #define __FUNCT__ "DMMoabGetInterface"
433aa768e4cSTim Tautges /*@
434aa768e4cSTim Tautges   DMMoabGetInterface - Get the MOAB instance used with this DMMoab
435aa768e4cSTim Tautges 
436aa768e4cSTim Tautges   Collective on MPI_Comm
437aa768e4cSTim Tautges 
438aa768e4cSTim Tautges   Input Parameter:
439aa768e4cSTim Tautges . dm      - The DMMoab object being set
440aa768e4cSTim Tautges 
441aa768e4cSTim Tautges   Output Parameter:
442aa768e4cSTim Tautges . mbiface - The MOAB instance set on this DMMoab
443aa768e4cSTim Tautges 
444aa768e4cSTim Tautges   Level: beginner
445aa768e4cSTim Tautges 
446aa768e4cSTim Tautges .keywords: DMMoab, create
447aa768e4cSTim Tautges @*/
448a4d2169cSTim Tautges PetscErrorCode DMMoabGetInterface(DM dm,moab::Interface **mbiface)
4491d72bce8STim Tautges {
4509426e041SSatish Balay   PetscErrorCode   ierr;
451cabb514dSBarry Smith   static PetscBool cite = PETSC_FALSE;
452cabb514dSBarry Smith 
4531d72bce8STim Tautges   PetscFunctionBegin;
4541d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
455cabb514dSBarry Smith   ierr = PetscCitationsRegister("@techreport{tautges_moab:_2004,\n  type = {{SAND2004-1592}},\n  title = {{MOAB:} A Mesh-Oriented Database},  institution = {Sandia National Laboratories},\n  author = {Tautges, T. J. and Meyers, R. and Merkley, K. and Stimpson, C. and Ernst, C.},\n  year = {2004},  note = {Report}\n}\n",&cite);CHKERRQ(ierr);
456a4d2169cSTim Tautges   *mbiface = ((DM_Moab*)dm->data)->mbiface;
4571d72bce8STim Tautges   PetscFunctionReturn(0);
4581d72bce8STim Tautges }
4591d72bce8STim Tautges 
4601d72bce8STim Tautges 
4611d72bce8STim Tautges #undef __FUNCT__
4625eb88e9dSVijay Mahadevan #define __FUNCT__ "DMMoabSetLocalVertices"
463aa768e4cSTim Tautges /*@
4645eb88e9dSVijay Mahadevan   DMMoabSetLocalVertices - Set the entities having DOFs on this DMMoab
465aa768e4cSTim Tautges 
466aa768e4cSTim Tautges   Collective on MPI_Comm
467aa768e4cSTim Tautges 
468aa768e4cSTim Tautges   Input Parameter:
469aa768e4cSTim Tautges . dm    - The DMMoab object being set
470aa768e4cSTim Tautges . range - The entities treated by this DMMoab
471aa768e4cSTim Tautges 
472aa768e4cSTim Tautges   Level: beginner
473aa768e4cSTim Tautges 
474aa768e4cSTim Tautges .keywords: DMMoab, create
475aa768e4cSTim Tautges @*/
4765eb88e9dSVijay Mahadevan PetscErrorCode DMMoabSetLocalVertices(DM dm,moab::Range *range)
4771d72bce8STim Tautges {
478032b8ab6SVijay Mahadevan   moab::ErrorCode merr;
479032b8ab6SVijay Mahadevan   PetscErrorCode  ierr;
480032b8ab6SVijay Mahadevan   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
481032b8ab6SVijay Mahadevan 
4821d72bce8STim Tautges   PetscFunctionBegin;
4831d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
484032b8ab6SVijay Mahadevan   dmmoab->vlocal->clear();
485032b8ab6SVijay Mahadevan   dmmoab->vowned->clear();
486032b8ab6SVijay Mahadevan   dmmoab->vlocal->insert(range->begin(), range->end());
487032b8ab6SVijay Mahadevan   *dmmoab->vowned = *dmmoab->vlocal;
488032b8ab6SVijay Mahadevan   merr = dmmoab->pcomm->filter_pstatus(*dmmoab->vowned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
489032b8ab6SVijay Mahadevan   *dmmoab->vghost = moab::subtract(*range, *dmmoab->vowned);
490032b8ab6SVijay Mahadevan   dmmoab->nloc=dmmoab->vowned->size();
491032b8ab6SVijay Mahadevan   dmmoab->nghost=dmmoab->vghost->size();
492032b8ab6SVijay Mahadevan   ierr = MPI_Allreduce(&dmmoab->nloc, &dmmoab->n, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
4931d72bce8STim Tautges   PetscFunctionReturn(0);
4941d72bce8STim Tautges }
4951d72bce8STim Tautges 
4961d72bce8STim Tautges 
4971d72bce8STim Tautges #undef __FUNCT__
4985eb88e9dSVijay Mahadevan #define __FUNCT__ "DMMoabGetLocalVertices"
499aa768e4cSTim Tautges /*@
5005eb88e9dSVijay Mahadevan   DMMoabGetLocalVertices - Get the entities having DOFs on this DMMoab
501aa768e4cSTim Tautges 
502aa768e4cSTim Tautges   Collective on MPI_Comm
503aa768e4cSTim Tautges 
504aa768e4cSTim Tautges   Input Parameter:
505aa768e4cSTim Tautges . dm    - The DMMoab object being set
506aa768e4cSTim Tautges 
507aa768e4cSTim Tautges   Output Parameter:
5085eb88e9dSVijay Mahadevan . owned - The owned vertex entities in this DMMoab
5095eb88e9dSVijay Mahadevan . ghost - The ghosted entities (non-owned) stored locally in this partition
510aa768e4cSTim Tautges 
511aa768e4cSTim Tautges   Level: beginner
512aa768e4cSTim Tautges 
513aa768e4cSTim Tautges .keywords: DMMoab, create
514aa768e4cSTim Tautges @*/
5151cec0304SVijay Mahadevan PetscErrorCode DMMoabGetLocalVertices(DM dm,moab::Range *owned,moab::Range *ghost)
5161d72bce8STim Tautges {
5171d72bce8STim Tautges   PetscFunctionBegin;
5181d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5191cec0304SVijay Mahadevan   if (owned) *owned = *((DM_Moab*)dm->data)->vowned;
5201cec0304SVijay Mahadevan   if (ghost) *ghost = *((DM_Moab*)dm->data)->vghost;
5211d72bce8STim Tautges   PetscFunctionReturn(0);
5221d72bce8STim Tautges }
5231d72bce8STim Tautges 
5241d72bce8STim Tautges #undef __FUNCT__
5255eb88e9dSVijay Mahadevan #define __FUNCT__ "DMMoabGetLocalElements"
5265eb88e9dSVijay Mahadevan /*@
5275eb88e9dSVijay Mahadevan   DMMoabGetLocalElements - Get the higher-dimensional entities that are locally owned
5285eb88e9dSVijay Mahadevan 
5295eb88e9dSVijay Mahadevan   Collective on MPI_Comm
5305eb88e9dSVijay Mahadevan 
5315eb88e9dSVijay Mahadevan   Input Parameter:
5325eb88e9dSVijay Mahadevan . dm    - The DMMoab object being set
5335eb88e9dSVijay Mahadevan 
5345eb88e9dSVijay Mahadevan   Output Parameter:
5355eb88e9dSVijay Mahadevan . range - The entities owned locally
5365eb88e9dSVijay Mahadevan 
5375eb88e9dSVijay Mahadevan   Level: beginner
5385eb88e9dSVijay Mahadevan 
5395eb88e9dSVijay Mahadevan .keywords: DMMoab, create
5405eb88e9dSVijay Mahadevan @*/
5411cec0304SVijay Mahadevan PetscErrorCode DMMoabGetLocalElements(DM dm,moab::Range *range)
5425eb88e9dSVijay Mahadevan {
5435eb88e9dSVijay Mahadevan   PetscFunctionBegin;
5445eb88e9dSVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5451cec0304SVijay Mahadevan   if (range) *range = *((DM_Moab*)dm->data)->elocal;
5461cec0304SVijay Mahadevan   PetscFunctionReturn(0);
5471cec0304SVijay Mahadevan }
5481cec0304SVijay Mahadevan 
5491cec0304SVijay Mahadevan 
5501cec0304SVijay Mahadevan #undef __FUNCT__
5511cec0304SVijay Mahadevan #define __FUNCT__ "DMMoabSetLocalElements"
5521cec0304SVijay Mahadevan /*@
5531cec0304SVijay Mahadevan   DMMoabSetLocalElements - Set the entities having DOFs on this DMMoab
5541cec0304SVijay Mahadevan 
5551cec0304SVijay Mahadevan   Collective on MPI_Comm
5561cec0304SVijay Mahadevan 
5571cec0304SVijay Mahadevan   Input Parameter:
5581cec0304SVijay Mahadevan . dm    - The DMMoab object being set
5591cec0304SVijay Mahadevan . range - The entities treated by this DMMoab
5601cec0304SVijay Mahadevan 
5611cec0304SVijay Mahadevan   Level: beginner
5621cec0304SVijay Mahadevan 
5631cec0304SVijay Mahadevan .keywords: DMMoab, create
5641cec0304SVijay Mahadevan @*/
5651cec0304SVijay Mahadevan PetscErrorCode DMMoabSetLocalElements(DM dm,moab::Range *range)
5661cec0304SVijay Mahadevan {
5671cec0304SVijay Mahadevan   moab::ErrorCode merr;
5681cec0304SVijay Mahadevan   PetscErrorCode  ierr;
5691cec0304SVijay Mahadevan   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
5701cec0304SVijay Mahadevan 
5711cec0304SVijay Mahadevan   PetscFunctionBegin;
5721cec0304SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5731cec0304SVijay Mahadevan   dmmoab->elocal->clear();
5741cec0304SVijay Mahadevan   dmmoab->eghost->clear();
5751cec0304SVijay Mahadevan   dmmoab->elocal->insert(range->begin(), range->end());
5761cec0304SVijay Mahadevan   PetscInfo2(dm, "Range size = %D; elocal size = %D.\n", range->size(), dmmoab->elocal->size());
5771cec0304SVijay Mahadevan   merr = dmmoab->pcomm->filter_pstatus(*dmmoab->elocal,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
5781cec0304SVijay Mahadevan   *dmmoab->eghost = moab::subtract(*range, *dmmoab->elocal);
5791cec0304SVijay Mahadevan   dmmoab->neleloc=dmmoab->elocal->size();
5801cec0304SVijay Mahadevan   ierr = MPI_Allreduce(&dmmoab->nele, &dmmoab->neleloc, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
5811cec0304SVijay Mahadevan   PetscInfo2(dm, "Created %D local and %D glocal elements.\n", dmmoab->neleloc, dmmoab->nele);
5825eb88e9dSVijay Mahadevan   PetscFunctionReturn(0);
5835eb88e9dSVijay Mahadevan }
5845eb88e9dSVijay Mahadevan 
5855eb88e9dSVijay Mahadevan 
5865eb88e9dSVijay Mahadevan #undef __FUNCT__
5871d72bce8STim Tautges #define __FUNCT__ "DMMoabSetLocalToGlobalTag"
588aa768e4cSTim Tautges /*@
589aa768e4cSTim Tautges   DMMoabSetLocalToGlobalTag - Set the tag used for local to global numbering
590aa768e4cSTim Tautges 
591aa768e4cSTim Tautges   Collective on MPI_Comm
592aa768e4cSTim Tautges 
593aa768e4cSTim Tautges   Input Parameter:
594aa768e4cSTim Tautges . dm      - The DMMoab object being set
595aa768e4cSTim Tautges . ltogtag - The MOAB tag used for local to global ids
596aa768e4cSTim Tautges 
597aa768e4cSTim Tautges   Level: beginner
598aa768e4cSTim Tautges 
599aa768e4cSTim Tautges .keywords: DMMoab, create
600aa768e4cSTim Tautges @*/
6011d72bce8STim Tautges PetscErrorCode DMMoabSetLocalToGlobalTag(DM dm,moab::Tag ltogtag)
6021d72bce8STim Tautges {
6031d72bce8STim Tautges   PetscFunctionBegin;
6041d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6051d72bce8STim Tautges   ((DM_Moab*)dm->data)->ltog_tag = ltogtag;
6061d72bce8STim Tautges   PetscFunctionReturn(0);
6071d72bce8STim Tautges }
6081d72bce8STim Tautges 
6091d72bce8STim Tautges 
6101d72bce8STim Tautges #undef __FUNCT__
6111d72bce8STim Tautges #define __FUNCT__ "DMMoabGetLocalToGlobalTag"
612aa768e4cSTim Tautges /*@
613aa768e4cSTim Tautges   DMMoabGetLocalToGlobalTag - Get the tag used for local to global numbering
614aa768e4cSTim Tautges 
615aa768e4cSTim Tautges   Collective on MPI_Comm
616aa768e4cSTim Tautges 
617aa768e4cSTim Tautges   Input Parameter:
618aa768e4cSTim Tautges . dm      - The DMMoab object being set
619aa768e4cSTim Tautges 
620aa768e4cSTim Tautges   Output Parameter:
621aa768e4cSTim Tautges . ltogtag - The MOAB tag used for local to global ids
622aa768e4cSTim Tautges 
623aa768e4cSTim Tautges   Level: beginner
624aa768e4cSTim Tautges 
625aa768e4cSTim Tautges .keywords: DMMoab, create
626aa768e4cSTim Tautges @*/
6271d72bce8STim Tautges PetscErrorCode DMMoabGetLocalToGlobalTag(DM dm,moab::Tag *ltog_tag)
6281d72bce8STim Tautges {
6291d72bce8STim Tautges   PetscFunctionBegin;
6301d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6311d72bce8STim Tautges   *ltog_tag = ((DM_Moab*)dm->data)->ltog_tag;
6321d72bce8STim Tautges   PetscFunctionReturn(0);
6331d72bce8STim Tautges }
6341d72bce8STim Tautges 
6351d72bce8STim Tautges 
6361d72bce8STim Tautges #undef __FUNCT__
6371d72bce8STim Tautges #define __FUNCT__ "DMMoabSetBlockSize"
638aa768e4cSTim Tautges /*@
639aa768e4cSTim Tautges   DMMoabSetBlockSize - Set the block size used with this DMMoab
640aa768e4cSTim Tautges 
641aa768e4cSTim Tautges   Collective on MPI_Comm
642aa768e4cSTim Tautges 
643aa768e4cSTim Tautges   Input Parameter:
644aa768e4cSTim Tautges . dm - The DMMoab object being set
645aa768e4cSTim Tautges . bs - The block size used with this DMMoab
646aa768e4cSTim Tautges 
647aa768e4cSTim Tautges   Level: beginner
648aa768e4cSTim Tautges 
649aa768e4cSTim Tautges .keywords: DMMoab, create
650aa768e4cSTim Tautges @*/
6511d72bce8STim Tautges PetscErrorCode DMMoabSetBlockSize(DM dm,PetscInt bs)
6521d72bce8STim Tautges {
6531d72bce8STim Tautges   PetscFunctionBegin;
6541d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6551d72bce8STim Tautges   ((DM_Moab*)dm->data)->bs = bs;
6561d72bce8STim Tautges   PetscFunctionReturn(0);
6571d72bce8STim Tautges }
6581d72bce8STim Tautges 
6591d72bce8STim Tautges 
6601d72bce8STim Tautges #undef __FUNCT__
6611d72bce8STim Tautges #define __FUNCT__ "DMMoabGetBlockSize"
662aa768e4cSTim Tautges /*@
663aa768e4cSTim Tautges   DMMoabGetBlockSize - Get the block size used with this DMMoab
664aa768e4cSTim Tautges 
665aa768e4cSTim Tautges   Collective on MPI_Comm
666aa768e4cSTim Tautges 
667aa768e4cSTim Tautges   Input Parameter:
668aa768e4cSTim Tautges . dm - The DMMoab object being set
669aa768e4cSTim Tautges 
670aa768e4cSTim Tautges   Output Parameter:
671aa768e4cSTim Tautges . bs - The block size used with this DMMoab
672aa768e4cSTim Tautges 
673aa768e4cSTim Tautges   Level: beginner
674aa768e4cSTim Tautges 
675aa768e4cSTim Tautges .keywords: DMMoab, create
676aa768e4cSTim Tautges @*/
6771d72bce8STim Tautges PetscErrorCode DMMoabGetBlockSize(DM dm,PetscInt *bs)
6781d72bce8STim Tautges {
6791d72bce8STim Tautges   PetscFunctionBegin;
6801d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6811d72bce8STim Tautges   *bs = ((DM_Moab*)dm->data)->bs;
6821d72bce8STim Tautges   PetscFunctionReturn(0);
6831d72bce8STim Tautges }
6841d72bce8STim Tautges 
6851cec0304SVijay Mahadevan 
6861cec0304SVijay Mahadevan #undef __FUNCT__
6871cec0304SVijay Mahadevan #define __FUNCT__ "DMMoabSetFieldVector"
6881cec0304SVijay Mahadevan PetscErrorCode DMMoabSetFieldVector(DM dm, PetscInt ifield, Vec fvec)
6891cec0304SVijay Mahadevan {
6901cec0304SVijay Mahadevan   DM_Moab        *dmmoab;
6911cec0304SVijay Mahadevan   moab::Tag     vtag,ntag;
6921cec0304SVijay Mahadevan   PetscScalar   *varray;
6931cec0304SVijay Mahadevan   moab::ErrorCode merr;
6941cec0304SVijay Mahadevan   PetscErrorCode  ierr;
6951cec0304SVijay Mahadevan   PetscSection section;
6961cec0304SVijay Mahadevan   PetscInt doff;
6971cec0304SVijay Mahadevan   std::string tag_name;
6981cec0304SVijay Mahadevan   moab::Range::iterator iter;
6991cec0304SVijay Mahadevan 
7001cec0304SVijay Mahadevan   PetscFunctionBegin;
7011cec0304SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7021cec0304SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
7031cec0304SVijay Mahadevan 
7041cec0304SVijay Mahadevan   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
7051cec0304SVijay Mahadevan 
7061cec0304SVijay Mahadevan   /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
7071cec0304SVijay Mahadevan   merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag,
7081cec0304SVijay Mahadevan                                           moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr);
7091cec0304SVijay Mahadevan 
7101cec0304SVijay Mahadevan   ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr);
7111cec0304SVijay Mahadevan 
7121cec0304SVijay Mahadevan   merr = dmmoab->mbiface->tag_get_name(vtag, tag_name);
7131cec0304SVijay Mahadevan   if (!tag_name.length() && merr !=moab::MB_SUCCESS) {
7141cec0304SVijay Mahadevan     ierr = DMMoabVecGetArray(dm,fvec,&varray);CHKERRQ(ierr);
7151cec0304SVijay Mahadevan     for(iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++) {
7161cec0304SVijay Mahadevan       moab::EntityHandle vtx = (*iter);
7171cec0304SVijay Mahadevan 
7181cec0304SVijay Mahadevan       /* get field dof index */
7191cec0304SVijay Mahadevan       ierr = PetscSectionGetFieldOffset(section, vtx, ifield, &doff);
7201cec0304SVijay Mahadevan 
7211cec0304SVijay Mahadevan       /* use the entity handle and the Dof index to set the right value */
7221cec0304SVijay Mahadevan       merr = dmmoab->mbiface->tag_set_data(ntag, &vtx, 1, (const void*)&varray[doff]);MBERRNM(merr);
7231cec0304SVijay Mahadevan     }
7241cec0304SVijay Mahadevan     ierr = DMMoabVecRestoreArray(dm,fvec,&varray);CHKERRQ(ierr);
7251cec0304SVijay Mahadevan   }
7261cec0304SVijay Mahadevan   else {
7271cec0304SVijay Mahadevan     ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&varray);CHKERRQ(ierr);
7281cec0304SVijay Mahadevan     /* we are using a MOAB Vec - directly copy the tag data to new one */
7291cec0304SVijay Mahadevan     merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)varray);MBERRNM(merr);
7301cec0304SVijay Mahadevan     merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)varray);MBERRNM(merr);
7311cec0304SVijay Mahadevan     /* make sure the parallel exchange for ghosts are done appropriately */
7321cec0304SVijay Mahadevan     merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr);
7331cec0304SVijay Mahadevan     ierr = PetscFree(varray);CHKERRQ(ierr);
7341cec0304SVijay Mahadevan   }
7351cec0304SVijay Mahadevan   PetscFunctionReturn(0);
7361cec0304SVijay Mahadevan }
7371cec0304SVijay Mahadevan 
7381cec0304SVijay Mahadevan 
7391cec0304SVijay Mahadevan #undef __FUNCT__
7407023aa44SVijay Mahadevan #define __FUNCT__ "DMMoabGetVertexCoordinates"
7417023aa44SVijay Mahadevan PetscErrorCode DMMoabGetVertexCoordinates(DM dm,PetscInt nconn,const moab::EntityHandle *conn,PetscScalar *vpos)
7427023aa44SVijay Mahadevan {
7437023aa44SVijay Mahadevan   DM_Moab         *dmmoab;
7447023aa44SVijay Mahadevan   PetscErrorCode  ierr;
7457023aa44SVijay Mahadevan   moab::ErrorCode merr;
7467023aa44SVijay Mahadevan 
7477023aa44SVijay Mahadevan   PetscFunctionBegin;
7487023aa44SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7497023aa44SVijay Mahadevan   PetscValidPointer(conn,3);
7507023aa44SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
7517023aa44SVijay Mahadevan 
7527023aa44SVijay Mahadevan   if (!vpos) {
7537023aa44SVijay Mahadevan     ierr = PetscMalloc(sizeof(PetscScalar)*nconn*3, &vpos);CHKERRQ(ierr);
7547023aa44SVijay Mahadevan   }
7557023aa44SVijay Mahadevan 
7567023aa44SVijay Mahadevan   /* Get connectivity information in MOAB canonical ordering */
7577023aa44SVijay Mahadevan   merr = dmmoab->mbiface->get_coords(conn, nconn, vpos);MBERRNM(merr);
7587023aa44SVijay Mahadevan   PetscFunctionReturn(0);
7597023aa44SVijay Mahadevan }
7607023aa44SVijay Mahadevan 
7617023aa44SVijay Mahadevan 
7627023aa44SVijay Mahadevan #undef __FUNCT__
7637023aa44SVijay Mahadevan #define __FUNCT__ "DMMoabGetElementConnectivity"
7647023aa44SVijay Mahadevan PetscErrorCode DMMoabGetElementConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn,const moab::EntityHandle **conn)
7657023aa44SVijay Mahadevan {
7667023aa44SVijay Mahadevan   DM_Moab        *dmmoab;
7677023aa44SVijay Mahadevan   const moab::EntityHandle *connect;
7687023aa44SVijay Mahadevan   moab::ErrorCode merr;
7697023aa44SVijay Mahadevan   PetscInt nnodes;
7707023aa44SVijay Mahadevan 
7717023aa44SVijay Mahadevan   PetscFunctionBegin;
7727023aa44SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7737023aa44SVijay Mahadevan   PetscValidPointer(conn,4);
7747023aa44SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
7757023aa44SVijay Mahadevan 
7767023aa44SVijay Mahadevan   /* Get connectivity information in MOAB canonical ordering */
7777023aa44SVijay Mahadevan   merr = dmmoab->mbiface->get_connectivity(ehandle, connect, nnodes);MBERRNM(merr);
7787023aa44SVijay Mahadevan   if (conn) *conn=connect;
7797023aa44SVijay Mahadevan   if (nconn) *nconn=nnodes;
7807023aa44SVijay Mahadevan   PetscFunctionReturn(0);
7817023aa44SVijay Mahadevan }
7827023aa44SVijay Mahadevan 
7837023aa44SVijay Mahadevan 
7847023aa44SVijay Mahadevan #undef __FUNCT__
7857023aa44SVijay Mahadevan #define __FUNCT__ "DMMoabCheckBoundaryVertices"
7867023aa44SVijay Mahadevan PetscErrorCode DMMoabCheckBoundaryVertices(DM dm,PetscInt nconn,const moab::EntityHandle *cnt,PetscBool* isbdvtx,PetscBool* elem_on_boundary)
7877023aa44SVijay Mahadevan {
7887023aa44SVijay Mahadevan   DM_Moab        *dmmoab;
7897023aa44SVijay Mahadevan   PetscInt       i;
7907023aa44SVijay Mahadevan 
7917023aa44SVijay Mahadevan   PetscFunctionBegin;
7927023aa44SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7937023aa44SVijay Mahadevan   PetscValidPointer(cnt,3);
7947023aa44SVijay Mahadevan   PetscValidPointer(isbdvtx,4);
7957023aa44SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
7967023aa44SVijay Mahadevan 
7977023aa44SVijay Mahadevan   if (elem_on_boundary) *elem_on_boundary = PETSC_FALSE;
7987023aa44SVijay Mahadevan   for (i=0; i < nconn; ++i) {
7997023aa44SVijay Mahadevan     moab::Range::const_iterator giter = dmmoab->bndyvtx->find(cnt[i]);
8007023aa44SVijay Mahadevan     if (giter != dmmoab->bndyvtx->end()) {
8017023aa44SVijay Mahadevan       isbdvtx[i] = PETSC_TRUE;
8027023aa44SVijay Mahadevan       if (elem_on_boundary) *elem_on_boundary = PETSC_TRUE;
8037023aa44SVijay Mahadevan     }
8047023aa44SVijay Mahadevan     else isbdvtx[i] = PETSC_FALSE;
8057023aa44SVijay Mahadevan   }
8067023aa44SVijay Mahadevan   PetscFunctionReturn(0);
8077023aa44SVijay Mahadevan }
8087023aa44SVijay Mahadevan 
8097023aa44SVijay Mahadevan 
8107023aa44SVijay Mahadevan #undef __FUNCT__
8111cec0304SVijay Mahadevan #define __FUNCT__ "DMMoabGetBoundaryEntities"
8121cec0304SVijay Mahadevan PetscErrorCode DMMoabGetBoundaryEntities(DM dm,moab::Range *bdvtx,moab::Range* bdfaces)
8131cec0304SVijay Mahadevan {
8141cec0304SVijay Mahadevan   DM_Moab        *dmmoab;
8151cec0304SVijay Mahadevan 
8161cec0304SVijay Mahadevan   PetscFunctionBegin;
8171cec0304SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
8181cec0304SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
8191cec0304SVijay Mahadevan 
8201cec0304SVijay Mahadevan   if (bdvtx)  *bdvtx = *dmmoab->bndyvtx;
8211cec0304SVijay Mahadevan   if (bdfaces)  *bdfaces = *dmmoab->bndyfaces;
8221cec0304SVijay Mahadevan   PetscFunctionReturn(0);
8231cec0304SVijay Mahadevan }
8241cec0304SVijay Mahadevan 
8251cec0304SVijay Mahadevan 
8261cec0304SVijay Mahadevan #undef __FUNCT__
8271cec0304SVijay Mahadevan #define __FUNCT__ "DMMoabSetFields"
8281cec0304SVijay Mahadevan PetscErrorCode DMMoabSetFields(DM dm,PetscInt nfields,const char** fields)
8291cec0304SVijay Mahadevan {
8301cec0304SVijay Mahadevan   DM_Moab        *dmmoab;
8311cec0304SVijay Mahadevan 
8321cec0304SVijay Mahadevan   PetscFunctionBegin;
8331cec0304SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
8341cec0304SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
8351cec0304SVijay Mahadevan 
8361cec0304SVijay Mahadevan   dmmoab->fields = fields;
8371cec0304SVijay Mahadevan   dmmoab->nfields = nfields;
8381cec0304SVijay Mahadevan   PetscFunctionReturn(0);
8391cec0304SVijay Mahadevan }
8401cec0304SVijay Mahadevan 
8411cec0304SVijay Mahadevan 
8421cec0304SVijay Mahadevan #undef __FUNCT__
8431cec0304SVijay Mahadevan #define __FUNCT__ "DMMoabGetFieldDof"
8441cec0304SVijay Mahadevan PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point,PetscInt field,PetscInt* dof)
8451cec0304SVijay Mahadevan {
8461cec0304SVijay Mahadevan   PetscSection section;
847*fc418013SVijay Mahadevan   PetscInt gid;
8481cec0304SVijay Mahadevan   PetscErrorCode ierr;
849*fc418013SVijay Mahadevan   moab::ErrorCode merr;
850*fc418013SVijay Mahadevan   DM_Moab        *dmmoab;
8511cec0304SVijay Mahadevan 
8521cec0304SVijay Mahadevan   PetscFunctionBegin;
8531cec0304SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
854*fc418013SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
855*fc418013SVijay Mahadevan 
8561cec0304SVijay Mahadevan   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
857*fc418013SVijay Mahadevan 
858*fc418013SVijay Mahadevan   /* first get the global ID for the point */
859*fc418013SVijay Mahadevan   merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&point,1,&gid);MBERRNM(merr);
860*fc418013SVijay Mahadevan 
861*fc418013SVijay Mahadevan   /* get the dof value for the field */
862*fc418013SVijay Mahadevan   ierr = PetscSectionGetFieldDof(section, gid, field, dof);CHKERRQ(ierr);
8631cec0304SVijay Mahadevan   PetscFunctionReturn(0);
8641cec0304SVijay Mahadevan }
8651cec0304SVijay Mahadevan 
8661cec0304SVijay Mahadevan 
8671cec0304SVijay Mahadevan #undef __FUNCT__
8681cec0304SVijay Mahadevan #define __FUNCT__ "DMMoabGetFieldDofs"
8691cec0304SVijay Mahadevan PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof)
8701cec0304SVijay Mahadevan {
871*fc418013SVijay Mahadevan   PetscInt i,gid;
8721cec0304SVijay Mahadevan   PetscSection section;
8731cec0304SVijay Mahadevan   PetscErrorCode  ierr;
874*fc418013SVijay Mahadevan   moab::ErrorCode merr;
875*fc418013SVijay Mahadevan   DM_Moab        *dmmoab;
8761cec0304SVijay Mahadevan 
8771cec0304SVijay Mahadevan   PetscFunctionBegin;
8781cec0304SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
879*fc418013SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
880*fc418013SVijay Mahadevan 
8811cec0304SVijay Mahadevan   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
8821cec0304SVijay Mahadevan   if (!dof) {
8831cec0304SVijay Mahadevan     ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr);
8841cec0304SVijay Mahadevan   }
885*fc418013SVijay Mahadevan 
886*fc418013SVijay Mahadevan   /* first get the local indices */
887*fc418013SVijay Mahadevan   merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr);
888*fc418013SVijay Mahadevan 
8891cec0304SVijay Mahadevan   for (i=0; i<npoints; ++i) {
890*fc418013SVijay Mahadevan     gid=dof[i];
891*fc418013SVijay Mahadevan     ierr = PetscSectionGetFieldDof(section, gid, field, &dof[i]);CHKERRQ(ierr);
8921cec0304SVijay Mahadevan   }
8931cec0304SVijay Mahadevan   PetscFunctionReturn(0);
8941cec0304SVijay Mahadevan }
8951cec0304SVijay Mahadevan 
8961cec0304SVijay Mahadevan 
897*fc418013SVijay Mahadevan #undef __FUNCT__
898*fc418013SVijay Mahadevan #define __FUNCT__ "DMMoab_GetWriteOptions_Private"
899*fc418013SVijay Mahadevan PetscErrorCode DMMoab_GetWriteOptions_Private(PetscInt fsetid, PetscInt numproc, PetscInt dim, MoabWriteMode mode, PetscInt dbglevel, const char* extra_opts, const char** write_opts)
900*fc418013SVijay Mahadevan {
901*fc418013SVijay Mahadevan   std::ostringstream str;
902*fc418013SVijay Mahadevan 
903*fc418013SVijay Mahadevan   PetscFunctionBegin;
904*fc418013SVijay Mahadevan 
905*fc418013SVijay Mahadevan   // do parallel read unless only one processor
906*fc418013SVijay Mahadevan   if (numproc > 1) {
907*fc418013SVijay Mahadevan     str << "PARALLEL=" << mode << ";";
908*fc418013SVijay Mahadevan     if (fsetid>=0) str << "PARALLEL_COMM=" << fsetid << ";";
909*fc418013SVijay Mahadevan   }
910*fc418013SVijay Mahadevan 
911*fc418013SVijay Mahadevan   if (dbglevel)
912*fc418013SVijay Mahadevan     str << "CPUTIME;DEBUG_IO=" << dbglevel << ";";
913*fc418013SVijay Mahadevan 
914*fc418013SVijay Mahadevan   if (extra_opts)
915*fc418013SVijay Mahadevan     str << extra_opts ;
916*fc418013SVijay Mahadevan 
917*fc418013SVijay Mahadevan   *write_opts = str.str().c_str();
918*fc418013SVijay Mahadevan   PetscFunctionReturn(0);
919*fc418013SVijay Mahadevan }
920*fc418013SVijay Mahadevan 
921*fc418013SVijay Mahadevan 
922*fc418013SVijay Mahadevan #undef __FUNCT__
923*fc418013SVijay Mahadevan #define __FUNCT__ "DMMoabOutput"
924*fc418013SVijay Mahadevan PetscErrorCode DMMoabOutput(DM dm,const char* filename,const char* usrwriteopts)
925*fc418013SVijay Mahadevan {
926*fc418013SVijay Mahadevan   DM_Moab        *dmmoab;
927*fc418013SVijay Mahadevan   PetscInt       dbglevel=0;
928*fc418013SVijay Mahadevan   const char *writeopts;
929*fc418013SVijay Mahadevan 
930*fc418013SVijay Mahadevan   PetscErrorCode ierr;
931*fc418013SVijay Mahadevan   moab::ErrorCode merr;
932*fc418013SVijay Mahadevan 
933*fc418013SVijay Mahadevan   PetscFunctionBegin;
934*fc418013SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
935*fc418013SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
936*fc418013SVijay Mahadevan 
937*fc418013SVijay Mahadevan   PetscBarrier((PetscObject)dm);
938*fc418013SVijay Mahadevan 
939*fc418013SVijay Mahadevan   /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */
940*fc418013SVijay Mahadevan   ierr  = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab");
941*fc418013SVijay Mahadevan   ierr  = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr);
942*fc418013SVijay Mahadevan   ierr  = PetscOptionsEnd();
943*fc418013SVijay Mahadevan 
944*fc418013SVijay Mahadevan   /* add mesh loading options specific to the DM */
945*fc418013SVijay Mahadevan   ierr = DMMoab_GetWriteOptions_Private(dmmoab->pcomm->get_id(), dmmoab->pcomm->size(), dmmoab->dim, MOAB_PARWOPTS_WRITE_PART, dbglevel, usrwriteopts, &writeopts);CHKERRQ(ierr);
946*fc418013SVijay Mahadevan   PetscInfo2(dm, "Writing file %s with options: %s\n",filename,writeopts);
947*fc418013SVijay Mahadevan 
948*fc418013SVijay Mahadevan   /* output file, using parallel write */
949*fc418013SVijay Mahadevan   merr = dmmoab->mbiface->write_file(filename, NULL, writeopts, &dmmoab->fileset, 1);MBERRVM(dmmoab->mbiface,"Writing output of DMMoab failed.",merr);
950*fc418013SVijay Mahadevan   PetscFunctionReturn(0);
951*fc418013SVijay Mahadevan }
952*fc418013SVijay Mahadevan 
953