xref: /petsc/src/dm/impls/moab/dmmbmg.cxx (revision b117cd09e766aba90b52b0bef2091bbb861995da)
1*b117cd09SVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
2*b117cd09SVijay Mahadevan 
3*b117cd09SVijay Mahadevan #include <petscdmmoab.h>
4*b117cd09SVijay Mahadevan #include <MBTagConventions.hpp>
5*b117cd09SVijay Mahadevan #include <moab/NestedRefine.hpp>
6*b117cd09SVijay Mahadevan 
7*b117cd09SVijay Mahadevan #undef __FUNCT__
8*b117cd09SVijay Mahadevan #define __FUNCT__ "DMMoabGenerateHierarchy"
9*b117cd09SVijay Mahadevan /*@
10*b117cd09SVijay Mahadevan   DMMoabGenerateHierarchy - Generate a multi-level uniform refinement hierarchy
11*b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
12*b117cd09SVijay Mahadevan   provided by the user.
13*b117cd09SVijay Mahadevan 
14*b117cd09SVijay Mahadevan   Collective on MPI_Comm
15*b117cd09SVijay Mahadevan 
16*b117cd09SVijay Mahadevan   Input Parameter:
17*b117cd09SVijay Mahadevan + dmb  - The DMMoab object
18*b117cd09SVijay Mahadevan 
19*b117cd09SVijay Mahadevan   Output Parameter:
20*b117cd09SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
21*b117cd09SVijay Mahadevan . ldegrees  - The degree of refinement at each level in the hierarchy
22*b117cd09SVijay Mahadevan 
23*b117cd09SVijay Mahadevan   Level: beginner
24*b117cd09SVijay Mahadevan 
25*b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
26*b117cd09SVijay Mahadevan @*/
27*b117cd09SVijay Mahadevan PetscErrorCode DMMoabGenerateHierarchy(DM dm,PetscInt nlevels,PetscInt *ldegrees)
28*b117cd09SVijay Mahadevan {
29*b117cd09SVijay Mahadevan   DM_Moab        *dmmoab;
30*b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
31*b117cd09SVijay Mahadevan   moab::ErrorCode merr;
32*b117cd09SVijay Mahadevan   PetscInt *pdegrees,i;
33*b117cd09SVijay Mahadevan 
34*b117cd09SVijay Mahadevan   PetscFunctionBegin;
35*b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36*b117cd09SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
37*b117cd09SVijay Mahadevan 
38*b117cd09SVijay Mahadevan   if (!ldegrees) {
39*b117cd09SVijay Mahadevan     ierr = PetscMalloc1(nlevels,&pdegrees);CHKERRQ(ierr);
40*b117cd09SVijay Mahadevan     for (i=0; i<nlevels; i++) pdegrees[i]=2; /* default = Degree 2 refinement */
41*b117cd09SVijay Mahadevan   }
42*b117cd09SVijay Mahadevan   else pdegrees=ldegrees;
43*b117cd09SVijay Mahadevan 
44*b117cd09SVijay Mahadevan   /* initialize set level refinement data for hierarchy */
45*b117cd09SVijay Mahadevan   dmmoab->nhlevels=nlevels;
46*b117cd09SVijay Mahadevan 
47*b117cd09SVijay Mahadevan   /* Instantiate the nested refinement class */
48*b117cd09SVijay Mahadevan   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), dmmoab->fileset);
49*b117cd09SVijay Mahadevan 
50*b117cd09SVijay Mahadevan   /* generate the mesh hierarchy */
51*b117cd09SVijay Mahadevan   merr = dmmoab->hierarchy->generate_mesh_hierarchy(pdegrees, nlevels, &dmmoab->hsets);MBERRNM(merr);
52*b117cd09SVijay Mahadevan 
53*b117cd09SVijay Mahadevan   if (!ldegrees) {
54*b117cd09SVijay Mahadevan     ierr = PetscFree(pdegrees);CHKERRQ(ierr);
55*b117cd09SVijay Mahadevan   }
56*b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
57*b117cd09SVijay Mahadevan }
58*b117cd09SVijay Mahadevan 
59*b117cd09SVijay Mahadevan #undef __FUNCT__
60*b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefineHierarchy_Moab"
61*b117cd09SVijay Mahadevan /*@
62*b117cd09SVijay Mahadevan   DMRefineHierarchy_Moab - Generate a multi-level uniform refinement hierarchy
63*b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
64*b117cd09SVijay Mahadevan   provided by the user.
65*b117cd09SVijay Mahadevan 
66*b117cd09SVijay Mahadevan   Collective on MPI_Comm
67*b117cd09SVijay Mahadevan 
68*b117cd09SVijay Mahadevan   Input Parameter:
69*b117cd09SVijay Mahadevan . dmb  - The DMMoab object
70*b117cd09SVijay Mahadevan 
71*b117cd09SVijay Mahadevan   Output Parameter:
72*b117cd09SVijay Mahadevan . nlevels   - The number of levels of refinement needed to generate the hierarchy
73*b117cd09SVijay Mahadevan + ldegrees  - The degree of refinement at each level in the hierarchy
74*b117cd09SVijay Mahadevan 
75*b117cd09SVijay Mahadevan   Level: beginner
76*b117cd09SVijay Mahadevan 
77*b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
78*b117cd09SVijay Mahadevan @*/
79*b117cd09SVijay Mahadevan PetscErrorCode  DMRefineHierarchy_Moab(DM dm,PetscInt nlevels,DM dmf[])
80*b117cd09SVijay Mahadevan {
81*b117cd09SVijay Mahadevan   DM_Moab        *dmmoab;
82*b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
83*b117cd09SVijay Mahadevan   moab::ErrorCode merr;
84*b117cd09SVijay Mahadevan   PetscInt        i;
85*b117cd09SVijay Mahadevan 
86*b117cd09SVijay Mahadevan   PetscFunctionBegin;
87*b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
88*b117cd09SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
89*b117cd09SVijay Mahadevan 
90*b117cd09SVijay Mahadevan   for (i=0; i<nlevels; i++) {
91*b117cd09SVijay Mahadevan     dmf[i] = dm;
92*b117cd09SVijay Mahadevan     PetscObjectReference((PetscObject)dm);
93*b117cd09SVijay Mahadevan   }
94*b117cd09SVijay Mahadevan 
95*b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMRefineHierarchy_Moab] :: Placeholder\n");
96*b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
97*b117cd09SVijay Mahadevan }
98*b117cd09SVijay Mahadevan 
99*b117cd09SVijay Mahadevan #undef __FUNCT__
100*b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsenHierarchy_Moab"
101*b117cd09SVijay Mahadevan /*@
102*b117cd09SVijay Mahadevan   DMCoarsenHierarchy_Moab - Generate a multi-level uniform refinement hierarchy
103*b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
104*b117cd09SVijay Mahadevan   provided by the user.
105*b117cd09SVijay Mahadevan 
106*b117cd09SVijay Mahadevan   Collective on MPI_Comm
107*b117cd09SVijay Mahadevan 
108*b117cd09SVijay Mahadevan   Input Parameter:
109*b117cd09SVijay Mahadevan . dmb  - The DMMoab object
110*b117cd09SVijay Mahadevan 
111*b117cd09SVijay Mahadevan   Output Parameter:
112*b117cd09SVijay Mahadevan . nlevels   - The number of levels of refinement needed to generate the hierarchy
113*b117cd09SVijay Mahadevan + ldegrees  - The degree of refinement at each level in the hierarchy
114*b117cd09SVijay Mahadevan 
115*b117cd09SVijay Mahadevan   Level: beginner
116*b117cd09SVijay Mahadevan 
117*b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
118*b117cd09SVijay Mahadevan @*/
119*b117cd09SVijay Mahadevan PetscErrorCode DMCoarsenHierarchy_Moab(DM dm,PetscInt nlevels,DM dmc[])
120*b117cd09SVijay Mahadevan {
121*b117cd09SVijay Mahadevan   DM_Moab        *dmmoab;
122*b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
123*b117cd09SVijay Mahadevan   moab::ErrorCode merr;
124*b117cd09SVijay Mahadevan   PetscInt        i;
125*b117cd09SVijay Mahadevan 
126*b117cd09SVijay Mahadevan   PetscFunctionBegin;
127*b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
128*b117cd09SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
129*b117cd09SVijay Mahadevan 
130*b117cd09SVijay Mahadevan   for (i=0; i<nlevels; i++) {
131*b117cd09SVijay Mahadevan     dmc[i] = dm;
132*b117cd09SVijay Mahadevan     PetscObjectReference((PetscObject)dm);
133*b117cd09SVijay Mahadevan   }
134*b117cd09SVijay Mahadevan 
135*b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMCoarsenHierarchy_Moab] :: Placeholder\n");
136*b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
137*b117cd09SVijay Mahadevan }
138*b117cd09SVijay Mahadevan 
139*b117cd09SVijay Mahadevan 
140*b117cd09SVijay Mahadevan #undef __FUNCT__
141*b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInterpolation_Moab"
142*b117cd09SVijay Mahadevan /*@
143*b117cd09SVijay Mahadevan   DMCreateInterpolation_Moab - Generate a multi-level uniform refinement hierarchy
144*b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
145*b117cd09SVijay Mahadevan   provided by the user.
146*b117cd09SVijay Mahadevan 
147*b117cd09SVijay Mahadevan   Collective on MPI_Comm
148*b117cd09SVijay Mahadevan 
149*b117cd09SVijay Mahadevan   Input Parameter:
150*b117cd09SVijay Mahadevan . dmb  - The DMMoab object
151*b117cd09SVijay Mahadevan 
152*b117cd09SVijay Mahadevan   Output Parameter:
153*b117cd09SVijay Mahadevan . nlevels   - The number of levels of refinement needed to generate the hierarchy
154*b117cd09SVijay Mahadevan + ldegrees  - The degree of refinement at each level in the hierarchy
155*b117cd09SVijay Mahadevan 
156*b117cd09SVijay Mahadevan   Level: beginner
157*b117cd09SVijay Mahadevan 
158*b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
159*b117cd09SVijay Mahadevan @*/
160*b117cd09SVijay Mahadevan PetscErrorCode DMCreateInterpolation_Moab(DM dm1,DM dm2,Mat* interpl,Vec* vec)
161*b117cd09SVijay Mahadevan {
162*b117cd09SVijay Mahadevan   DM_Moab        *dmb1, *dmb2;
163*b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
164*b117cd09SVijay Mahadevan   moab::ErrorCode merr;
165*b117cd09SVijay Mahadevan   Vec             unitv;
166*b117cd09SVijay Mahadevan   PetscInt        dim,dofs_per_element=4;
167*b117cd09SVijay Mahadevan   PetscReal       unitval=1.0;
168*b117cd09SVijay Mahadevan   PetscBool       eonbnd,dbdry[27];
169*b117cd09SVijay Mahadevan   std::vector<int> bndrows;
170*b117cd09SVijay Mahadevan 
171*b117cd09SVijay Mahadevan   PetscFunctionBegin;
172*b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
173*b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
174*b117cd09SVijay Mahadevan   dmb1 = (DM_Moab*)(dm1)->data;
175*b117cd09SVijay Mahadevan   dmb2 = (DM_Moab*)(dm2)->data;
176*b117cd09SVijay Mahadevan 
177*b117cd09SVijay Mahadevan   ierr = MatCreate(PetscObjectComm((PetscObject)dm1), interpl);CHKERRQ(ierr);
178*b117cd09SVijay Mahadevan   ierr = MatSetType(*interpl, MATMPIAIJ);CHKERRQ(ierr);
179*b117cd09SVijay Mahadevan   ierr = MatSetSizes(*interpl, dmb2->nloc, dmb1->nloc, dmb2->n, dmb1->n);CHKERRQ(ierr);
180*b117cd09SVijay Mahadevan 
181*b117cd09SVijay Mahadevan   /* TODO: This is a hack for the rectangular system - decipher NNZ pattern better */
182*b117cd09SVijay Mahadevan   ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr);
183*b117cd09SVijay Mahadevan 
184*b117cd09SVijay Mahadevan   ierr = MatSeqAIJSetPreallocation(*interpl,dmb2->nloc,0);CHKERRQ(ierr);
185*b117cd09SVijay Mahadevan   ierr = MatMPIAIJSetPreallocation(*interpl,dmb2->nloc,0,dmb2->nghost,0);CHKERRQ(ierr);
186*b117cd09SVijay Mahadevan 
187*b117cd09SVijay Mahadevan   /* set up internal matrix data-structures */
188*b117cd09SVijay Mahadevan   ierr = MatSetUp(*interpl);CHKERRQ(ierr);
189*b117cd09SVijay Mahadevan 
190*b117cd09SVijay Mahadevan   ierr = MatZeroEntries(*interpl);CHKERRQ(ierr);
191*b117cd09SVijay Mahadevan 
192*b117cd09SVijay Mahadevan   ierr = DMGetDimension(dm1, &dim);CHKERRQ(ierr);
193*b117cd09SVijay Mahadevan 
194*b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "Creating a %D DIM matrix that is %D X %D and setting diagonal to 1.0", dim, dmb1->nloc, dmb2->nloc);
195*b117cd09SVijay Mahadevan 
196*b117cd09SVijay Mahadevan   double factor = std::pow(2.0,(dmb2->hlevel-dmb1->hlevel)*dmb1->dim*1.0);
197*b117cd09SVijay Mahadevan 
198*b117cd09SVijay Mahadevan     //Loop through the remaining vertices. These vertices appear only on the current refined_level.
199*b117cd09SVijay Mahadevan   for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) {
200*b117cd09SVijay Mahadevan 
201*b117cd09SVijay Mahadevan     const moab::EntityHandle ehandle = *iter;
202*b117cd09SVijay Mahadevan     std::vector<moab::EntityHandle> children;
203*b117cd09SVijay Mahadevan     std::vector<moab::EntityHandle> connp, connc;
204*b117cd09SVijay Mahadevan 
205*b117cd09SVijay Mahadevan     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
206*b117cd09SVijay Mahadevan     merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr);
207*b117cd09SVijay Mahadevan 
208*b117cd09SVijay Mahadevan     /* Get connectivity and coordinates of the parent vertices */
209*b117cd09SVijay Mahadevan     merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr);
210*b117cd09SVijay Mahadevan     for (int ic=0; ic < children.size(); ic++) {
211*b117cd09SVijay Mahadevan       std::vector<moab::EntityHandle> tconnc;
212*b117cd09SVijay Mahadevan       /* Get coordinates of the parent vertices in canonical order */
213*b117cd09SVijay Mahadevan       merr = dmb1->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr);
214*b117cd09SVijay Mahadevan       for (int tc=0; tc<tconnc.size(); tc++) {
215*b117cd09SVijay Mahadevan         connc.push_back(tconnc[tc]);
216*b117cd09SVijay Mahadevan       }
217*b117cd09SVijay Mahadevan     }
218*b117cd09SVijay Mahadevan 
219*b117cd09SVijay Mahadevan     std::vector<double> pcoords(connp.size()*3), ccoords(connc.size()*3), values_phi(connc.size());
220*b117cd09SVijay Mahadevan     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
221*b117cd09SVijay Mahadevan     merr = dmb1->hierarchy->get_coordinates(&connp[0], connp.size(), dmb1->hlevel, &pcoords[0]);MBERRNM(merr);
222*b117cd09SVijay Mahadevan     merr = dmb2->hierarchy->get_coordinates(&connc[0], connc.size(), dmb2->hlevel, &ccoords[0]);MBERRNM(merr);
223*b117cd09SVijay Mahadevan 
224*b117cd09SVijay Mahadevan     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
225*b117cd09SVijay Mahadevan     /* TODO: specific to scalar system - use GetDofs */
226*b117cd09SVijay Mahadevan     ierr = DMMoabGetFieldDofs(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr);
227*b117cd09SVijay Mahadevan     ierr = DMMoabGetFieldDofs(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr);
228*b117cd09SVijay Mahadevan 
229*b117cd09SVijay Mahadevan     // Compute the interoplation weights by determining distance of 1-ring neighbor vertices
230*b117cd09SVijay Mahadevan     // from current vertex
231*b117cd09SVijay Mahadevan     for (int i=0;i<connp.size(); i++) {
232*b117cd09SVijay Mahadevan       double normsum=0.0;
233*b117cd09SVijay Mahadevan       for (int j=0;j<connc.size(); j++) {
234*b117cd09SVijay Mahadevan         int offset=j;
235*b117cd09SVijay Mahadevan         values_phi[offset] = 0.0;
236*b117cd09SVijay Mahadevan         for (int k=0;k<3; k++)
237*b117cd09SVijay Mahadevan           values_phi[offset] += (pcoords[i*3+k]-ccoords[k+j*3])*(pcoords[i*3+k]-ccoords[k+j*3]);
238*b117cd09SVijay Mahadevan         if (values_phi[offset] < 1e-12) {
239*b117cd09SVijay Mahadevan           values_phi[offset] = 1e12;
240*b117cd09SVijay Mahadevan           // PetscPrintf(PETSC_COMM_WORLD, "Found coarse and fine restrictive space: %D, %D\n",dofsp[i],dofsc[j]);
241*b117cd09SVijay Mahadevan         }
242*b117cd09SVijay Mahadevan         else {
243*b117cd09SVijay Mahadevan           values_phi[offset] = 1.0/(values_phi[offset]);
244*b117cd09SVijay Mahadevan           normsum += values_phi[offset];
245*b117cd09SVijay Mahadevan         }
246*b117cd09SVijay Mahadevan       }
247*b117cd09SVijay Mahadevan       //PetscPrintf(PETSC_COMM_WORLD, "\nRow [%D]: ", dofsp[i]);
248*b117cd09SVijay Mahadevan       for (int j=0;j<connc.size(); j++) {
249*b117cd09SVijay Mahadevan         //int offset=i*connc.size()+j;
250*b117cd09SVijay Mahadevan         int offset=j;
251*b117cd09SVijay Mahadevan         if (values_phi[offset] > 1e11)
252*b117cd09SVijay Mahadevan           values_phi[offset] = factor*0.5/connc.size();
253*b117cd09SVijay Mahadevan         else
254*b117cd09SVijay Mahadevan           values_phi[offset] = factor*values_phi[offset]*0.5/(connc.size()*normsum);
255*b117cd09SVijay Mahadevan         //PetscPrintf(PETSC_COMM_WORLD, "%D : %g, ", dofsc[j], values_phi[offset]);
256*b117cd09SVijay Mahadevan       }
257*b117cd09SVijay Mahadevan       ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr);
258*b117cd09SVijay Mahadevan     }
259*b117cd09SVijay Mahadevan 
260*b117cd09SVijay Mahadevan     /* check if element is on the boundary */
261*b117cd09SVijay Mahadevan     //ierr = DMMoabIsEntityOnBoundary(dm1,ehandle,&eonbnd);CHKERRQ(ierr);
262*b117cd09SVijay Mahadevan     ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry);CHKERRQ(ierr);
263*b117cd09SVijay Mahadevan     eonbnd=PETSC_FALSE;
264*b117cd09SVijay Mahadevan     for (int i=0; i< connc.size(); ++i)
265*b117cd09SVijay Mahadevan       if (dbdry[i]) eonbnd=PETSC_TRUE;
266*b117cd09SVijay Mahadevan 
267*b117cd09SVijay Mahadevan     values_phi.clear();
268*b117cd09SVijay Mahadevan     values_phi.resize(connp.size());
269*b117cd09SVijay Mahadevan     /* apply dirichlet boundary conditions */
270*b117cd09SVijay Mahadevan     if (eonbnd) {
271*b117cd09SVijay Mahadevan 
272*b117cd09SVijay Mahadevan       ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
273*b117cd09SVijay Mahadevan       ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
274*b117cd09SVijay Mahadevan       /* get the list of nodes on boundary so that we can enforce dirichlet conditions strongly */
275*b117cd09SVijay Mahadevan       //ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry);CHKERRQ(ierr);
276*b117cd09SVijay Mahadevan       for (int i=0; i < connc.size(); i++) {
277*b117cd09SVijay Mahadevan         if (dmb2->hierarchy->is_boundary_vertex(connc[i])) {  /* dirichlet node */
278*b117cd09SVijay Mahadevan           /* think about strongly imposing dirichlet */
279*b117cd09SVijay Mahadevan           //bndrows.push_back(dofsc[i]);
280*b117cd09SVijay Mahadevan 
281*b117cd09SVijay Mahadevan           ierr = MatSetValues(*interpl, 1, &dofsc[i], connp.size(), &dofsp[0], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr);
282*b117cd09SVijay Mahadevan           //values_phi[0]=1.0;
283*b117cd09SVijay Mahadevan           //ierr = MatSetValues(*interpl, 1, &dofsc[i], 1, &dofsc[i], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr);
284*b117cd09SVijay Mahadevan         }
285*b117cd09SVijay Mahadevan       }
286*b117cd09SVijay Mahadevan 
287*b117cd09SVijay Mahadevan       ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
288*b117cd09SVijay Mahadevan       ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
289*b117cd09SVijay Mahadevan     }
290*b117cd09SVijay Mahadevan 
291*b117cd09SVijay Mahadevan     //get interpolation weights
292*b117cd09SVijay Mahadevan     //ierr = Compute_Quad4_Basis(pcoords, 1, coord, values_phi);CHKERRQ(ierr);
293*b117cd09SVijay Mahadevan     /*for (int j=0;j<dofs_per_element; j++)
294*b117cd09SVijay Mahadevan       std::cout<<"values "<<values_phi[j]<<std::endl;*/
295*b117cd09SVijay Mahadevan 
296*b117cd09SVijay Mahadevan     //get row and column indices, zero weights are ignored
297*b117cd09SVijay Mahadevan     /*
298*b117cd09SVijay Mahadevan     int nz_ind = 0;
299*b117cd09SVijay Mahadevan     idx = dmb2->vowned->index(vhandle);
300*b117cd09SVijay Mahadevan     for (int j=0;j<dofs_per_element; j++){
301*b117cd09SVijay Mahadevan       idy[nz_ind] = dmb1->vowned->index(connectivity[j]);
302*b117cd09SVijay Mahadevan       PetscPrintf(PETSC_COMM_WORLD, "Finding coarse connectivity vertex %D associated with [%D, %D] - set to %D\n", connectivity[j], parent.size(), vhandle, idy[nz_ind]);
303*b117cd09SVijay Mahadevan       //values_phi[nz_ind] = values_phi[j];
304*b117cd09SVijay Mahadevan       nz_ind = nz_ind+1;
305*b117cd09SVijay Mahadevan     }
306*b117cd09SVijay Mahadevan     */
307*b117cd09SVijay Mahadevan 
308*b117cd09SVijay Mahadevan     //PetscPrintf(PETSC_COMM_WORLD, "Setting basis for row %D: %g %g %g %g\n", idx, values_phi[0], values_phi[1], values_phi[2], values_phi[3]);
309*b117cd09SVijay Mahadevan     //PetscPrintf(PETSC_COMM_WORLD, "Setting entries for row %D: %D %D %D %D\n", idx, idy[0], idy[1], idy[2], idy[3]);
310*b117cd09SVijay Mahadevan     //ierr = MatSetValues(*interpl, nz_ind, idy, 1, &idx, values_phi, INSERT_VALUES);CHKERRQ(ierr);
311*b117cd09SVijay Mahadevan     //ierr = MatSetValues(*interpl, connp.size(), dofsp, connc.size(), dofsc, &values_phi[0], INSERT_VALUES);CHKERRQ(ierr);
312*b117cd09SVijay Mahadevan   }
313*b117cd09SVijay Mahadevan 
314*b117cd09SVijay Mahadevan   //PetscPrintf(PETSC_COMM_WORLD, "[Boundary vertices = %D] :: A few: %D %D %D %D \n", bndrows.size(), bndrows[0], bndrows[1], bndrows[2], bndrows[3]);
315*b117cd09SVijay Mahadevan   //ierr = MatZeroRows(*interpl, bndrows.size(), &bndrows[0], 1.0, 0, 0);CHKERRQ(ierr);
316*b117cd09SVijay Mahadevan 
317*b117cd09SVijay Mahadevan   ierr = MatAssemblyBegin(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
318*b117cd09SVijay Mahadevan   ierr = MatAssemblyEnd(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
319*b117cd09SVijay Mahadevan   //MatView(*interpl, 0);
320*b117cd09SVijay Mahadevan 
321*b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInterpolation_Moab] :: Placeholder\n");
322*b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
323*b117cd09SVijay Mahadevan }
324*b117cd09SVijay Mahadevan 
325*b117cd09SVijay Mahadevan #undef __FUNCT__
326*b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInjection_Moab"
327*b117cd09SVijay Mahadevan /*@
328*b117cd09SVijay Mahadevan   DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy
329*b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
330*b117cd09SVijay Mahadevan   provided by the user.
331*b117cd09SVijay Mahadevan 
332*b117cd09SVijay Mahadevan   Collective on MPI_Comm
333*b117cd09SVijay Mahadevan 
334*b117cd09SVijay Mahadevan   Input Parameter:
335*b117cd09SVijay Mahadevan . dmb  - The DMMoab object
336*b117cd09SVijay Mahadevan 
337*b117cd09SVijay Mahadevan   Output Parameter:
338*b117cd09SVijay Mahadevan . nlevels   - The number of levels of refinement needed to generate the hierarchy
339*b117cd09SVijay Mahadevan + ldegrees  - The degree of refinement at each level in the hierarchy
340*b117cd09SVijay Mahadevan 
341*b117cd09SVijay Mahadevan   Level: beginner
342*b117cd09SVijay Mahadevan 
343*b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
344*b117cd09SVijay Mahadevan @*/
345*b117cd09SVijay Mahadevan PetscErrorCode DMCreateInjection_Moab(DM dm1,DM dm2,VecScatter* ctx)
346*b117cd09SVijay Mahadevan {
347*b117cd09SVijay Mahadevan   DM_Moab        *dmmoab;
348*b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
349*b117cd09SVijay Mahadevan   moab::ErrorCode merr;
350*b117cd09SVijay Mahadevan 
351*b117cd09SVijay Mahadevan   PetscFunctionBegin;
352*b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
353*b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
354*b117cd09SVijay Mahadevan   dmmoab = (DM_Moab*)(dm1)->data;
355*b117cd09SVijay Mahadevan 
356*b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n");
357*b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
358*b117cd09SVijay Mahadevan }
359*b117cd09SVijay Mahadevan 
360*b117cd09SVijay Mahadevan 
361*b117cd09SVijay Mahadevan #undef __FUNCT__
362*b117cd09SVijay Mahadevan #define __FUNCT__ "DM_UMR_Moab_Private"
363*b117cd09SVijay Mahadevan PetscErrorCode  DM_UMR_Moab_Private(DM dm,MPI_Comm comm,PetscBool refine,DM *dmref)
364*b117cd09SVijay Mahadevan {
365*b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
366*b117cd09SVijay Mahadevan   PetscInt        M,N,P,i,dim;
367*b117cd09SVijay Mahadevan   DM              dm2;
368*b117cd09SVijay Mahadevan   moab::ErrorCode merr;
369*b117cd09SVijay Mahadevan   DM_Moab        *dmb = (DM_Moab*)dm->data,*dd2;
370*b117cd09SVijay Mahadevan 
371*b117cd09SVijay Mahadevan   PetscFunctionBegin;
372*b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
373*b117cd09SVijay Mahadevan   PetscValidPointer(dmref,4);
374*b117cd09SVijay Mahadevan 
375*b117cd09SVijay Mahadevan   //if (dmb->hlevel+1 > dmb->nhlevels && refine) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Invalid multigrid refinement hierarchy level specified (%D). MOAB UMR max levels = %D\n",dmb->hlevel+1,dmb->nhlevels);
376*b117cd09SVijay Mahadevan   //if (dmb->hlevel-1 < 0 && !refine) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Invalid multigrid coarsen hierarchy level specified (%D).\n",dmb->hlevel-1);
377*b117cd09SVijay Mahadevan 
378*b117cd09SVijay Mahadevan   if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) )  {
379*b117cd09SVijay Mahadevan     dmref = PETSC_NULL;
380*b117cd09SVijay Mahadevan     PetscFunctionReturn(0);
381*b117cd09SVijay Mahadevan   }
382*b117cd09SVijay Mahadevan 
383*b117cd09SVijay Mahadevan   ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr);
384*b117cd09SVijay Mahadevan   dd2 = (DM_Moab*)dm2->data;
385*b117cd09SVijay Mahadevan 
386*b117cd09SVijay Mahadevan   dd2->mbiface = dmb->mbiface;
387*b117cd09SVijay Mahadevan   dd2->pcomm = dmb->pcomm;
388*b117cd09SVijay Mahadevan   dd2->icreatedinstance = PETSC_FALSE;
389*b117cd09SVijay Mahadevan 
390*b117cd09SVijay Mahadevan   /* set the new level based on refinement/coarsening */
391*b117cd09SVijay Mahadevan   if (refine) {
392*b117cd09SVijay Mahadevan     dd2->hlevel=dmb->hlevel+1;
393*b117cd09SVijay Mahadevan   }
394*b117cd09SVijay Mahadevan   else {
395*b117cd09SVijay Mahadevan     dd2->hlevel=dmb->hlevel-1;
396*b117cd09SVijay Mahadevan   }
397*b117cd09SVijay Mahadevan 
398*b117cd09SVijay Mahadevan   /* Copy the multilevel hierarchy pointers in MOAB */
399*b117cd09SVijay Mahadevan   dd2->hierarchy = dmb->hierarchy;
400*b117cd09SVijay Mahadevan   dd2->nhlevels = dmb->nhlevels;
401*b117cd09SVijay Mahadevan   ierr = PetscMalloc1(dd2->nhlevels+1,&dd2->hsets);CHKERRQ(ierr);
402*b117cd09SVijay Mahadevan   for (i=0; i<=dd2->nhlevels; i++) {
403*b117cd09SVijay Mahadevan     dd2->hsets[i]=dmb->hsets[i];
404*b117cd09SVijay Mahadevan   }
405*b117cd09SVijay Mahadevan   dd2->fileset = dd2->hsets[dd2->hlevel];
406*b117cd09SVijay Mahadevan 
407*b117cd09SVijay Mahadevan   /* do the remaining initializations for DMMoab */
408*b117cd09SVijay Mahadevan   dd2->bs = dmb->bs;
409*b117cd09SVijay Mahadevan   dd2->numFields = dmb->numFields;
410*b117cd09SVijay Mahadevan   dd2->rw_dbglevel = dmb->rw_dbglevel;
411*b117cd09SVijay Mahadevan   dd2->partition_by_rank = dmb->partition_by_rank;
412*b117cd09SVijay Mahadevan   ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr);
413*b117cd09SVijay Mahadevan   ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr);
414*b117cd09SVijay Mahadevan   dd2->read_mode = dmb->read_mode;
415*b117cd09SVijay Mahadevan   dd2->write_mode = dmb->write_mode;
416*b117cd09SVijay Mahadevan 
417*b117cd09SVijay Mahadevan   /* set global ID tag handle */
418*b117cd09SVijay Mahadevan   ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr);
419*b117cd09SVijay Mahadevan 
420*b117cd09SVijay Mahadevan   merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr);
421*b117cd09SVijay Mahadevan 
422*b117cd09SVijay Mahadevan   //ierr = DMMoabCreateMoab(PetscObjectComm((PetscObject)dm),dmb->mbiface,dmb->pcomm,&dmb->ltog_tag,NULL,&dm2);CHKERRQ(ierr);
423*b117cd09SVijay Mahadevan   ierr = DMSetOptionsPrefix(dm2,((PetscObject)dm)->prefix);CHKERRQ(ierr);
424*b117cd09SVijay Mahadevan   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
425*b117cd09SVijay Mahadevan   ierr = DMSetDimension(dm2,dim);CHKERRQ(ierr);
426*b117cd09SVijay Mahadevan 
427*b117cd09SVijay Mahadevan   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
428*b117cd09SVijay Mahadevan   dm2->ops->creatematrix = dm->ops->creatematrix;
429*b117cd09SVijay Mahadevan 
430*b117cd09SVijay Mahadevan   /* copy fill information if given */
431*b117cd09SVijay Mahadevan   ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr);
432*b117cd09SVijay Mahadevan 
433*b117cd09SVijay Mahadevan   /* copy vector type information */
434*b117cd09SVijay Mahadevan   ierr = DMSetMatType(dm2,dm->mattype);CHKERRQ(ierr);
435*b117cd09SVijay Mahadevan   ierr = DMSetVecType(dm2,dm->vectype);CHKERRQ(ierr);
436*b117cd09SVijay Mahadevan   ierr = DMMoabSetFieldNames(dm2,dmb->numFields,dmb->fieldNames);CHKERRQ(ierr);
437*b117cd09SVijay Mahadevan 
438*b117cd09SVijay Mahadevan   ierr = DMSetFromOptions(dm2);CHKERRQ(ierr);
439*b117cd09SVijay Mahadevan 
440*b117cd09SVijay Mahadevan   /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */
441*b117cd09SVijay Mahadevan   ierr = DMSetUp(dm2);CHKERRQ(ierr);
442*b117cd09SVijay Mahadevan 
443*b117cd09SVijay Mahadevan   *dmref = dm2;
444*b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
445*b117cd09SVijay Mahadevan }
446*b117cd09SVijay Mahadevan 
447*b117cd09SVijay Mahadevan 
448*b117cd09SVijay Mahadevan #undef __FUNCT__
449*b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefine_Moab"
450*b117cd09SVijay Mahadevan /*@
451*b117cd09SVijay Mahadevan   DMRefine_Moab - Generate a multi-level uniform refinement hierarchy
452*b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
453*b117cd09SVijay Mahadevan   provided by the user.
454*b117cd09SVijay Mahadevan 
455*b117cd09SVijay Mahadevan   Collective on MPI_Comm
456*b117cd09SVijay Mahadevan 
457*b117cd09SVijay Mahadevan   Input Parameter:
458*b117cd09SVijay Mahadevan . dmb  - The DMMoab object
459*b117cd09SVijay Mahadevan 
460*b117cd09SVijay Mahadevan   Output Parameter:
461*b117cd09SVijay Mahadevan . nlevels   - The number of levels of refinement needed to generate the hierarchy
462*b117cd09SVijay Mahadevan + ldegrees  - The degree of refinement at each level in the hierarchy
463*b117cd09SVijay Mahadevan 
464*b117cd09SVijay Mahadevan   Level: beginner
465*b117cd09SVijay Mahadevan 
466*b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
467*b117cd09SVijay Mahadevan @*/
468*b117cd09SVijay Mahadevan PetscErrorCode DMRefine_Moab(DM dm,MPI_Comm comm,DM* dmf)
469*b117cd09SVijay Mahadevan {
470*b117cd09SVijay Mahadevan   DM_Moab        *dmmoab;
471*b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
472*b117cd09SVijay Mahadevan   moab::ErrorCode merr;
473*b117cd09SVijay Mahadevan 
474*b117cd09SVijay Mahadevan   PetscFunctionBegin;
475*b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
476*b117cd09SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
477*b117cd09SVijay Mahadevan 
478*b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMRefine_Moab] :: Placeholder\n");
479*b117cd09SVijay Mahadevan   ierr = DM_UMR_Moab_Private(dm,comm,PETSC_TRUE,dmf);CHKERRQ(ierr);
480*b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMRefine_Moab] :: All done\n");
481*b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
482*b117cd09SVijay Mahadevan }
483*b117cd09SVijay Mahadevan 
484*b117cd09SVijay Mahadevan #undef __FUNCT__
485*b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsen_Moab"
486*b117cd09SVijay Mahadevan /*@
487*b117cd09SVijay Mahadevan   DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy
488*b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
489*b117cd09SVijay Mahadevan   provided by the user.
490*b117cd09SVijay Mahadevan 
491*b117cd09SVijay Mahadevan   Collective on MPI_Comm
492*b117cd09SVijay Mahadevan 
493*b117cd09SVijay Mahadevan   Input Parameter:
494*b117cd09SVijay Mahadevan . dmb  - The DMMoab object
495*b117cd09SVijay Mahadevan 
496*b117cd09SVijay Mahadevan   Output Parameter:
497*b117cd09SVijay Mahadevan . nlevels   - The number of levels of refinement needed to generate the hierarchy
498*b117cd09SVijay Mahadevan + ldegrees  - The degree of refinement at each level in the hierarchy
499*b117cd09SVijay Mahadevan 
500*b117cd09SVijay Mahadevan   Level: beginner
501*b117cd09SVijay Mahadevan 
502*b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
503*b117cd09SVijay Mahadevan @*/
504*b117cd09SVijay Mahadevan PetscErrorCode DMCoarsen_Moab(DM dm,MPI_Comm comm,DM* dmc)
505*b117cd09SVijay Mahadevan {
506*b117cd09SVijay Mahadevan   DM_Moab        *dmmoab;
507*b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
508*b117cd09SVijay Mahadevan   moab::ErrorCode merr;
509*b117cd09SVijay Mahadevan 
510*b117cd09SVijay Mahadevan   PetscFunctionBegin;
511*b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
512*b117cd09SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
513*b117cd09SVijay Mahadevan 
514*b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMCoarsen_Moab] :: Placeholder\n");
515*b117cd09SVijay Mahadevan   ierr = DM_UMR_Moab_Private(dm,comm,PETSC_FALSE,dmc);CHKERRQ(ierr);
516*b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMCoarsen_Moab] :: All done\n");
517*b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
518*b117cd09SVijay Mahadevan }
519