xref: /petsc/src/ksp/pc/impls/wb/wb.c (revision c0decd05c6848b80907752eef350b55c8c90e696)
1174b6946SBarry Smith 
2c6db04a5SJed Brown #include <petscdmda.h>   /*I "petscdmda.h" I*/
3af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h>   /*I "petscksp.h" I*/
482c86c8fSBarry Smith #include <petscctable.h>
57233f9f0SBarry Smith 
68e722e36SBarry Smith typedef struct {
78e722e36SBarry Smith   PCExoticType type;
88e722e36SBarry Smith   Mat          P;            /* the constructed interpolation matrix */
9ace3abfcSBarry Smith   PetscBool    directSolve;  /* use direct LU factorization to construct interpolation */
108e722e36SBarry Smith   KSP          ksp;
118e722e36SBarry Smith } PC_Exotic;
12174b6946SBarry Smith 
136a6fc655SJed Brown const char *const PCExoticTypes[] = {"face","wirebasket","PCExoticType","PC_Exotic",0};
14064c8009SBarry Smith 
15064c8009SBarry Smith 
16064c8009SBarry Smith /*
17aa219208SBarry Smith       DMDAGetWireBasketInterpolation - Gets the interpolation for a wirebasket based coarse space
18064c8009SBarry Smith 
19064c8009SBarry Smith */
20*c0decd05SBarry Smith PetscErrorCode DMDAGetWireBasketInterpolation(PC pc,DM da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P)
21064c8009SBarry Smith {
22064c8009SBarry Smith   PetscErrorCode         ierr;
23064c8009SBarry Smith   PetscInt               dim,i,j,k,m,n,p,dof,Nint,Nface,Nwire,Nsurf,*Iint,*Isurf,cint = 0,csurf = 0,istart,jstart,kstart,*II,N,c = 0;
2428d20b34SBarry Smith   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[26],*globals,Ng,*IIint,*IIsurf,Nt;
25064c8009SBarry Smith   Mat                    Xint, Xsurf,Xint_tmp;
26064c8009SBarry Smith   IS                     isint,issurf,is,row,col;
27064c8009SBarry Smith   ISLocalToGlobalMapping ltg;
28064c8009SBarry Smith   MPI_Comm               comm;
29064c8009SBarry Smith   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
30064c8009SBarry Smith   MatFactorInfo          info;
31064c8009SBarry Smith   PetscScalar            *xsurf,*xint;
328e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
33064c8009SBarry Smith   PetscScalar            tmp;
34064c8009SBarry Smith #endif
35064c8009SBarry Smith   PetscTable             ht;
36064c8009SBarry Smith 
37064c8009SBarry Smith   PetscFunctionBegin;
381321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0,0,0);CHKERRQ(ierr);
39ce94432eSBarry Smith   if (dof != 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only for single field problems");
40ce94432eSBarry Smith   if (dim != 3) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only coded for 3d problems");
41aa219208SBarry Smith   ierr   = DMDAGetCorners(da,0,0,0,&m,&n,&p);CHKERRQ(ierr);
42aa219208SBarry Smith   ierr   = DMDAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr);
43064c8009SBarry Smith   istart = istart ? -1 : 0;
44064c8009SBarry Smith   jstart = jstart ? -1 : 0;
45064c8009SBarry Smith   kstart = kstart ? -1 : 0;
46064c8009SBarry Smith 
47064c8009SBarry Smith   /*
48064c8009SBarry Smith     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge)
49064c8009SBarry Smith     to all the local degrees of freedom (this includes the vertices, edges and faces).
50064c8009SBarry Smith 
51064c8009SBarry Smith     Xint are the subset of the interpolation into the interior
52064c8009SBarry Smith 
53064c8009SBarry Smith     Xface are the interpolation onto faces but not into the interior
54064c8009SBarry Smith 
55064c8009SBarry Smith     Xsurf are the interpolation onto the vertices and edges (the surfbasket)
56064c8009SBarry Smith                                       Xint
57064c8009SBarry Smith     Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain
58064c8009SBarry Smith                                       Xsurf
59064c8009SBarry Smith   */
60064c8009SBarry Smith   N     = (m - istart)*(n - jstart)*(p - kstart);
61064c8009SBarry Smith   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
62064c8009SBarry Smith   Nface = 2*((m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart));
63064c8009SBarry Smith   Nwire = 4*((m-2-istart) + (n-2-jstart) + (p-2-kstart)) + 8;
64064c8009SBarry Smith   Nsurf = Nface + Nwire;
650298fd71SBarry Smith   ierr  = MatCreateSeqDense(MPI_COMM_SELF,Nint,26,NULL,&Xint);CHKERRQ(ierr);
660298fd71SBarry Smith   ierr  = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,26,NULL,&Xsurf);CHKERRQ(ierr);
678c778c55SBarry Smith   ierr  = MatDenseGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
68064c8009SBarry Smith 
69064c8009SBarry Smith   /*
70064c8009SBarry Smith      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
71064c8009SBarry Smith      Xsurf will be all zero (thus making the coarse matrix singular).
72064c8009SBarry Smith   */
73e32f2f54SBarry Smith   if (m-istart < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in X direction must be at least 3");
74e32f2f54SBarry Smith   if (n-jstart < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in Y direction must be at least 3");
75e32f2f54SBarry Smith   if (p-kstart < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in Z direction must be at least 3");
76064c8009SBarry Smith 
77064c8009SBarry Smith   cnt = 0;
782fa5cd67SKarl Rupp 
792fa5cd67SKarl Rupp   xsurf[cnt++] = 1;
802fa5cd67SKarl Rupp   for (i=1; i<m-istart-1; i++) xsurf[cnt++ + Nsurf] = 1;
812fa5cd67SKarl Rupp   xsurf[cnt++ + 2*Nsurf] = 1;
822fa5cd67SKarl Rupp 
832fa5cd67SKarl Rupp   for (j=1; j<n-1-jstart; j++) {
842fa5cd67SKarl Rupp     xsurf[cnt++ + 3*Nsurf] = 1;
852fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 4*Nsurf] = 1;
862fa5cd67SKarl Rupp     xsurf[cnt++ + 5*Nsurf] = 1;
87064c8009SBarry Smith   }
882fa5cd67SKarl Rupp 
892fa5cd67SKarl Rupp   xsurf[cnt++ + 6*Nsurf] = 1;
902fa5cd67SKarl Rupp   for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 7*Nsurf] = 1;
912fa5cd67SKarl Rupp   xsurf[cnt++ + 8*Nsurf] = 1;
922fa5cd67SKarl Rupp 
932fa5cd67SKarl Rupp   for (k=1; k<p-1-kstart; k++) {
942fa5cd67SKarl Rupp     xsurf[cnt++ + 9*Nsurf] = 1;
952fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 10*Nsurf] = 1;
962fa5cd67SKarl Rupp     xsurf[cnt++ + 11*Nsurf] = 1;
972fa5cd67SKarl Rupp 
982fa5cd67SKarl Rupp     for (j=1; j<n-1-jstart; j++) {
992fa5cd67SKarl Rupp       xsurf[cnt++ + 12*Nsurf] = 1;
1002fa5cd67SKarl Rupp       /* these are the interior nodes */
1012fa5cd67SKarl Rupp       xsurf[cnt++ + 13*Nsurf] = 1;
1022fa5cd67SKarl Rupp     }
1032fa5cd67SKarl Rupp 
1042fa5cd67SKarl Rupp     xsurf[cnt++ + 14*Nsurf] = 1;
1052fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 15*Nsurf] = 1;
1062fa5cd67SKarl Rupp     xsurf[cnt++ + 16*Nsurf] = 1;
1072fa5cd67SKarl Rupp   }
1082fa5cd67SKarl Rupp 
1092fa5cd67SKarl Rupp   xsurf[cnt++ + 17*Nsurf] = 1;
1102fa5cd67SKarl Rupp   for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 18*Nsurf] = 1;
1112fa5cd67SKarl Rupp   xsurf[cnt++ + 19*Nsurf] = 1;
1122fa5cd67SKarl Rupp 
1132fa5cd67SKarl Rupp   for (j=1;j<n-1-jstart;j++) {
1142fa5cd67SKarl Rupp     xsurf[cnt++ + 20*Nsurf] = 1;
1152fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 21*Nsurf] = 1;
1162fa5cd67SKarl Rupp     xsurf[cnt++ + 22*Nsurf] = 1;
1172fa5cd67SKarl Rupp   }
1182fa5cd67SKarl Rupp 
1192fa5cd67SKarl Rupp   xsurf[cnt++ + 23*Nsurf] = 1;
1202fa5cd67SKarl Rupp   for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 24*Nsurf] = 1;
1212fa5cd67SKarl Rupp   xsurf[cnt++ + 25*Nsurf] = 1;
1222fa5cd67SKarl Rupp 
123064c8009SBarry Smith 
1248e722e36SBarry Smith   /* interpolations only sum to 1 when using direct solver */
1258e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
126064c8009SBarry Smith   for (i=0; i<Nsurf; i++) {
127064c8009SBarry Smith     tmp = 0.0;
1282fa5cd67SKarl Rupp     for (j=0; j<26; j++) tmp += xsurf[i+j*Nsurf];
12957622a8eSBarry Smith     if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %D value %g",i,(double)PetscAbsScalar(tmp));
130064c8009SBarry Smith   }
131064c8009SBarry Smith #endif
1328c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
133064c8009SBarry Smith   /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/
134064c8009SBarry Smith 
135064c8009SBarry Smith 
136064c8009SBarry Smith   /*
137064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
138064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
1397dae84e0SHong Zhang             (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it
140aa219208SBarry Smith              is NOT the local DMDA ordering.)
141064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
142064c8009SBarry Smith   */
143064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
144dcca6d9dSJed Brown   ierr = PetscMalloc3(N,&II,Nint,&Iint,Nsurf,&Isurf);CHKERRQ(ierr);
145dcca6d9dSJed Brown   ierr = PetscMalloc2(Nint,&IIint,Nsurf,&IIsurf);CHKERRQ(ierr);
146064c8009SBarry Smith   for (k=0; k<p-kstart; k++) {
147064c8009SBarry Smith     for (j=0; j<n-jstart; j++) {
148064c8009SBarry Smith       for (i=0; i<m-istart; i++) {
149064c8009SBarry Smith         II[c++] = i + j*mwidth + k*mwidth*nwidth;
150064c8009SBarry Smith 
151064c8009SBarry Smith         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
152064c8009SBarry Smith           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
153064c8009SBarry Smith           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
154064c8009SBarry Smith         } else {
155064c8009SBarry Smith           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
156064c8009SBarry Smith           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
157064c8009SBarry Smith         }
158064c8009SBarry Smith       }
159064c8009SBarry Smith     }
160064c8009SBarry Smith   }
161e32f2f54SBarry Smith   if (c != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N");
162e32f2f54SBarry Smith   if (cint != Nint) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint");
163e32f2f54SBarry Smith   if (csurf != Nsurf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf");
1641411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltg);CHKERRQ(ierr);
165064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr);
166064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr);
167064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr);
168064c8009SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
16970b3c8c7SBarry Smith   ierr = ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
17070b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint);CHKERRQ(ierr);
17170b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf);CHKERRQ(ierr);
172064c8009SBarry Smith   ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr);
173064c8009SBarry Smith 
1747dae84e0SHong Zhang   ierr = MatCreateSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr);
175064c8009SBarry Smith   A    = *Aholder;
176064c8009SBarry Smith   ierr = PetscFree(Aholder);CHKERRQ(ierr);
177064c8009SBarry Smith 
1787dae84e0SHong Zhang   ierr = MatCreateSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr);
1797dae84e0SHong Zhang   ierr = MatCreateSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr);
1807dae84e0SHong Zhang   ierr = MatCreateSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr);
181064c8009SBarry Smith 
182064c8009SBarry Smith   /*
183064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
184064c8009SBarry Smith   */
18592e4c2edSBarry Smith   ierr = MatMatMult(Ais,Xsurf,MAT_INITIAL_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr);
1868e722e36SBarry Smith   ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr);
1878e722e36SBarry Smith   if (exotic->directSolve) {
1882692d6eeSBarry Smith     ierr = MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr);
189064c8009SBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
1902692d6eeSBarry Smith     ierr = MatGetOrdering(Aii,MATORDERINGND,&row,&col);CHKERRQ(ierr);
191064c8009SBarry Smith     ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr);
192fcfd50ebSBarry Smith     ierr = ISDestroy(&row);CHKERRQ(ierr);
193fcfd50ebSBarry Smith     ierr = ISDestroy(&col);CHKERRQ(ierr);
194064c8009SBarry Smith     ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr);
195064c8009SBarry Smith     ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr);
1966bf464f9SBarry Smith     ierr = MatDestroy(&iAii);CHKERRQ(ierr);
1978e722e36SBarry Smith   } else {
1988e722e36SBarry Smith     Vec         b,x;
1998e722e36SBarry Smith     PetscScalar *xint_tmp;
200064c8009SBarry Smith 
2018c778c55SBarry Smith     ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr);
202778a2246SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,0,&x);CHKERRQ(ierr);
2038c778c55SBarry Smith     ierr = MatDenseGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
204778a2246SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,0,&b);CHKERRQ(ierr);
20523ee1639SBarry Smith     ierr = KSPSetOperators(exotic->ksp,Aii,Aii);CHKERRQ(ierr);
2068e722e36SBarry Smith     for (i=0; i<26; i++) {
2078e722e36SBarry Smith       ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr);
2088e722e36SBarry Smith       ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr);
2098e722e36SBarry Smith       ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr);
210*c0decd05SBarry Smith       ierr = KSPCheckSolve(exotic->ksp,pc,x);CHKERRQ(ierr);
2118e722e36SBarry Smith       ierr = VecResetArray(x);CHKERRQ(ierr);
2128e722e36SBarry Smith       ierr = VecResetArray(b);CHKERRQ(ierr);
2138e722e36SBarry Smith     }
2148c778c55SBarry Smith     ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr);
2158c778c55SBarry Smith     ierr = MatDenseRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
2166bf464f9SBarry Smith     ierr = VecDestroy(&x);CHKERRQ(ierr);
2176bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
2188e722e36SBarry Smith   }
2196bf464f9SBarry Smith   ierr = MatDestroy(&Xint_tmp);CHKERRQ(ierr);
2208e722e36SBarry Smith 
2218e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
2228c778c55SBarry Smith   ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr);
223064c8009SBarry Smith   for (i=0; i<Nint; i++) {
224064c8009SBarry Smith     tmp = 0.0;
2252fa5cd67SKarl Rupp     for (j=0; j<26; j++) tmp += xint[i+j*Nint];
2262fa5cd67SKarl Rupp 
22757622a8eSBarry Smith     if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xint interpolation at i %D value %g",i,(double)PetscAbsScalar(tmp));
228064c8009SBarry Smith   }
2298c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr);
230064c8009SBarry Smith   /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
231064c8009SBarry Smith #endif
232064c8009SBarry Smith 
233064c8009SBarry Smith 
234064c8009SBarry Smith   /*         total vertices             total faces                                  total edges */
235064c8009SBarry Smith   Ntotal = (mp + 1)*(np + 1)*(pp + 1) + mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1) + mp*(np+1)*(pp+1) + np*(mp+1)*(pp+1) +  pp*(mp+1)*(np+1);
236064c8009SBarry Smith 
237064c8009SBarry Smith   /*
238064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
239064c8009SBarry Smith   */
240064c8009SBarry Smith   cnt = 0;
2412fa5cd67SKarl Rupp 
242064c8009SBarry Smith   gl[cnt++] = 0;  { gl[cnt++] = 1;} gl[cnt++] = m-istart-1;
243064c8009SBarry Smith   { gl[cnt++] = mwidth;  { gl[cnt++] = mwidth+1;} gl[cnt++] = mwidth + m-istart-1;}
244064c8009SBarry Smith   gl[cnt++] = mwidth*(n-jstart-1);  { gl[cnt++] = mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*(n-jstart-1) + m-istart-1;
245064c8009SBarry Smith   {
246064c8009SBarry Smith     gl[cnt++] = mwidth*nwidth;  { gl[cnt++] = mwidth*nwidth+1;}  gl[cnt++] = mwidth*nwidth+ m-istart-1;
247064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
248064c8009SBarry Smith     gl[cnt++] = mwidth*nwidth+ mwidth*(n-jstart-1);   { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1) + m-istart-1;
249064c8009SBarry Smith   }
250064c8009SBarry Smith   gl[cnt++] = mwidth*nwidth*(p-kstart-1); { gl[cnt++] = mwidth*nwidth*(p-kstart-1)+1;} gl[cnt++] = mwidth*nwidth*(p-kstart-1) +  m-istart-1;
251064c8009SBarry Smith   { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth;   { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;} gl[cnt++] = mwidth*nwidth*(p-kstart-1)+mwidth+m-istart-1;}
252064c8009SBarry Smith   gl[cnt++] = mwidth*nwidth*(p-kstart-1) +  mwidth*(n-jstart-1);  { gl[cnt++] = mwidth*nwidth*(p-kstart-1)+ mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth*(n-jstart-1) + m-istart-1;
253064c8009SBarry Smith 
254064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
255064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
256064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,26,gl,gl);CHKERRQ(ierr);
257064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
258785e854fSJed Brown   ierr = PetscMalloc1(26*mp*np*pp,&globals);CHKERRQ(ierr);
259ce94432eSBarry Smith   ierr = MPI_Allgather(gl,26,MPIU_INT,globals,26,MPIU_INT,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
260064c8009SBarry Smith 
261064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
2620298fd71SBarry Smith   ierr = MatGetSize(Aglobal,&Nt,NULL);CHKERRQ(ierr);
26312dd4568SBarry Smith   ierr = PetscTableCreate(Ntotal/3,Nt+1,&ht);CHKERRQ(ierr);
264064c8009SBarry Smith   for (i=0; i<26*mp*np*pp; i++) {
265064c8009SBarry Smith     ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr);
266064c8009SBarry Smith   }
267064c8009SBarry Smith   ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr);
268e32f2f54SBarry Smith   if (cnt != Ntotal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hash table size %D not equal to total number coarse grid points %D",cnt,Ntotal);
269064c8009SBarry Smith   ierr = PetscFree(globals);CHKERRQ(ierr);
270064c8009SBarry Smith   for (i=0; i<26; i++) {
271064c8009SBarry Smith     ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr);
272064c8009SBarry Smith     gl[i]--;
273064c8009SBarry Smith   }
2746bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&ht);CHKERRQ(ierr);
275064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
276064c8009SBarry Smith 
277064c8009SBarry Smith   /* construct global interpolation matrix */
2780298fd71SBarry Smith   ierr = MatGetLocalSize(Aglobal,&Ng,NULL);CHKERRQ(ierr);
279064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
280ce94432eSBarry Smith     ierr = MatCreateAIJ(PetscObjectComm((PetscObject)da),Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,NULL,Nint+Nsurf,NULL,P);CHKERRQ(ierr);
281064c8009SBarry Smith   } else {
282064c8009SBarry Smith     ierr = MatZeroEntries(*P);CHKERRQ(ierr);
283064c8009SBarry Smith   }
284064c8009SBarry Smith   ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
2858c778c55SBarry Smith   ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr);
286064c8009SBarry Smith   ierr = MatSetValues(*P,Nint,IIint,26,gl,xint,INSERT_VALUES);CHKERRQ(ierr);
2878c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr);
2888c778c55SBarry Smith   ierr = MatDenseGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
289064c8009SBarry Smith   ierr = MatSetValues(*P,Nsurf,IIsurf,26,gl,xsurf,INSERT_VALUES);CHKERRQ(ierr);
2908c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
291064c8009SBarry Smith   ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
292064c8009SBarry Smith   ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
293064c8009SBarry Smith   ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr);
294064c8009SBarry Smith 
2958e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
296064c8009SBarry Smith   {
297064c8009SBarry Smith     Vec         x,y;
298064c8009SBarry Smith     PetscScalar *yy;
299ce94432eSBarry Smith     ierr = VecCreateMPI(PetscObjectComm((PetscObject)da),Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr);
300ce94432eSBarry Smith     ierr = VecCreateMPI(PetscObjectComm((PetscObject)da),PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr);
301064c8009SBarry Smith     ierr = VecSet(x,1.0);CHKERRQ(ierr);
302064c8009SBarry Smith     ierr = MatMult(*P,x,y);CHKERRQ(ierr);
303064c8009SBarry Smith     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
304064c8009SBarry Smith     for (i=0; i<Ng; i++) {
30557622a8eSBarry Smith       if (PetscAbsScalar(yy[i]-1.0) > 1.e-10) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong p interpolation at i %D value %g",i,(double)PetscAbsScalar(yy[i]));
306064c8009SBarry Smith     }
307064c8009SBarry Smith     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
308064c8009SBarry Smith     ierr = VecDestroy(x);CHKERRQ(ierr);
309064c8009SBarry Smith     ierr = VecDestroy(y);CHKERRQ(ierr);
310064c8009SBarry Smith   }
311064c8009SBarry Smith #endif
312064c8009SBarry Smith 
313fcfd50ebSBarry Smith   ierr = MatDestroy(&Aii);CHKERRQ(ierr);
314fcfd50ebSBarry Smith   ierr = MatDestroy(&Ais);CHKERRQ(ierr);
315fcfd50ebSBarry Smith   ierr = MatDestroy(&Asi);CHKERRQ(ierr);
316fcfd50ebSBarry Smith   ierr = MatDestroy(&A);CHKERRQ(ierr);
317fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
318fcfd50ebSBarry Smith   ierr = ISDestroy(&isint);CHKERRQ(ierr);
319fcfd50ebSBarry Smith   ierr = ISDestroy(&issurf);CHKERRQ(ierr);
320fcfd50ebSBarry Smith   ierr = MatDestroy(&Xint);CHKERRQ(ierr);
321fcfd50ebSBarry Smith   ierr = MatDestroy(&Xsurf);CHKERRQ(ierr);
322064c8009SBarry Smith   PetscFunctionReturn(0);
323064c8009SBarry Smith }
324064c8009SBarry Smith 
325064c8009SBarry Smith /*
326aa219208SBarry Smith       DMDAGetFaceInterpolation - Gets the interpolation for a face based coarse space
327064c8009SBarry Smith 
328064c8009SBarry Smith */
329*c0decd05SBarry Smith PetscErrorCode DMDAGetFaceInterpolation(PC pc,DM da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P)
330064c8009SBarry Smith {
331064c8009SBarry Smith   PetscErrorCode         ierr;
332064c8009SBarry Smith   PetscInt               dim,i,j,k,m,n,p,dof,Nint,Nface,Nwire,Nsurf,*Iint,*Isurf,cint = 0,csurf = 0,istart,jstart,kstart,*II,N,c = 0;
33328d20b34SBarry Smith   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[6],*globals,Ng,*IIint,*IIsurf,Nt;
334064c8009SBarry Smith   Mat                    Xint, Xsurf,Xint_tmp;
335064c8009SBarry Smith   IS                     isint,issurf,is,row,col;
336064c8009SBarry Smith   ISLocalToGlobalMapping ltg;
337064c8009SBarry Smith   MPI_Comm               comm;
338064c8009SBarry Smith   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
339064c8009SBarry Smith   MatFactorInfo          info;
340064c8009SBarry Smith   PetscScalar            *xsurf,*xint;
341064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
342064c8009SBarry Smith   PetscScalar            tmp;
343064c8009SBarry Smith #endif
344064c8009SBarry Smith   PetscTable             ht;
345064c8009SBarry Smith 
346064c8009SBarry Smith   PetscFunctionBegin;
3471321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0,0,0);CHKERRQ(ierr);
348ce94432eSBarry Smith   if (dof != 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only for single field problems");
349ce94432eSBarry Smith   if (dim != 3) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only coded for 3d problems");
350aa219208SBarry Smith   ierr   = DMDAGetCorners(da,0,0,0,&m,&n,&p);CHKERRQ(ierr);
351aa219208SBarry Smith   ierr   = DMDAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr);
352064c8009SBarry Smith   istart = istart ? -1 : 0;
353064c8009SBarry Smith   jstart = jstart ? -1 : 0;
354064c8009SBarry Smith   kstart = kstart ? -1 : 0;
355064c8009SBarry Smith 
356064c8009SBarry Smith   /*
357064c8009SBarry Smith     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge)
358064c8009SBarry Smith     to all the local degrees of freedom (this includes the vertices, edges and faces).
359064c8009SBarry Smith 
360064c8009SBarry Smith     Xint are the subset of the interpolation into the interior
361064c8009SBarry Smith 
362064c8009SBarry Smith     Xface are the interpolation onto faces but not into the interior
363064c8009SBarry Smith 
364064c8009SBarry Smith     Xsurf are the interpolation onto the vertices and edges (the surfbasket)
365064c8009SBarry Smith                                       Xint
366064c8009SBarry Smith     Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain
367064c8009SBarry Smith                                       Xsurf
368064c8009SBarry Smith   */
369064c8009SBarry Smith   N     = (m - istart)*(n - jstart)*(p - kstart);
370064c8009SBarry Smith   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
371064c8009SBarry Smith   Nface = 2*((m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart));
372064c8009SBarry Smith   Nwire = 4*((m-2-istart) + (n-2-jstart) + (p-2-kstart)) + 8;
373064c8009SBarry Smith   Nsurf = Nface + Nwire;
3740298fd71SBarry Smith   ierr  = MatCreateSeqDense(MPI_COMM_SELF,Nint,6,NULL,&Xint);CHKERRQ(ierr);
3750298fd71SBarry Smith   ierr  = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,6,NULL,&Xsurf);CHKERRQ(ierr);
3768c778c55SBarry Smith   ierr  = MatDenseGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
377064c8009SBarry Smith 
378064c8009SBarry Smith   /*
379064c8009SBarry Smith      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
380064c8009SBarry Smith      Xsurf will be all zero (thus making the coarse matrix singular).
381064c8009SBarry Smith   */
382e32f2f54SBarry Smith   if (m-istart < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in X direction must be at least 3");
383e32f2f54SBarry Smith   if (n-jstart < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in Y direction must be at least 3");
384e32f2f54SBarry Smith   if (p-kstart < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in Z direction must be at least 3");
385064c8009SBarry Smith 
386064c8009SBarry Smith   cnt = 0;
3872fa5cd67SKarl Rupp   for (j=1; j<n-1-jstart; j++) {
3882fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 0*Nsurf] = 1;
389064c8009SBarry Smith   }
3902fa5cd67SKarl Rupp 
3912fa5cd67SKarl Rupp   for (k=1; k<p-1-kstart; k++) {
3922fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 1*Nsurf] = 1;
3932fa5cd67SKarl Rupp     for (j=1; j<n-1-jstart; j++) {
3942fa5cd67SKarl Rupp       xsurf[cnt++ + 2*Nsurf] = 1;
3952fa5cd67SKarl Rupp       /* these are the interior nodes */
3962fa5cd67SKarl Rupp       xsurf[cnt++ + 3*Nsurf] = 1;
3972fa5cd67SKarl Rupp     }
3982fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 4*Nsurf] = 1;
3992fa5cd67SKarl Rupp   }
4002fa5cd67SKarl Rupp   for (j=1;j<n-1-jstart;j++) {
4012fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 5*Nsurf] = 1;
4022fa5cd67SKarl Rupp   }
403064c8009SBarry Smith 
404064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
405064c8009SBarry Smith   for (i=0; i<Nsurf; i++) {
406064c8009SBarry Smith     tmp = 0.0;
4072fa5cd67SKarl Rupp     for (j=0; j<6; j++) tmp += xsurf[i+j*Nsurf];
4082fa5cd67SKarl Rupp 
40957622a8eSBarry Smith     if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %D value %g",i,(double)PetscAbsScalar(tmp));
410064c8009SBarry Smith   }
411064c8009SBarry Smith #endif
4128c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
413064c8009SBarry Smith   /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/
414064c8009SBarry Smith 
415064c8009SBarry Smith 
416064c8009SBarry Smith   /*
417064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
418064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
4197dae84e0SHong Zhang             (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it
420aa219208SBarry Smith              is NOT the local DMDA ordering.)
421064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
422064c8009SBarry Smith   */
423064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
424dcca6d9dSJed Brown   ierr = PetscMalloc3(N,&II,Nint,&Iint,Nsurf,&Isurf);CHKERRQ(ierr);
425dcca6d9dSJed Brown   ierr = PetscMalloc2(Nint,&IIint,Nsurf,&IIsurf);CHKERRQ(ierr);
426064c8009SBarry Smith   for (k=0; k<p-kstart; k++) {
427064c8009SBarry Smith     for (j=0; j<n-jstart; j++) {
428064c8009SBarry Smith       for (i=0; i<m-istart; i++) {
429064c8009SBarry Smith         II[c++] = i + j*mwidth + k*mwidth*nwidth;
430064c8009SBarry Smith 
431064c8009SBarry Smith         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
432064c8009SBarry Smith           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
433064c8009SBarry Smith           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
434064c8009SBarry Smith         } else {
435064c8009SBarry Smith           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
436064c8009SBarry Smith           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
437064c8009SBarry Smith         }
438064c8009SBarry Smith       }
439064c8009SBarry Smith     }
440064c8009SBarry Smith   }
441e32f2f54SBarry Smith   if (c != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N");
442e32f2f54SBarry Smith   if (cint != Nint) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint");
443e32f2f54SBarry Smith   if (csurf != Nsurf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf");
4441411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltg);CHKERRQ(ierr);
445064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr);
446064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr);
447064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr);
448064c8009SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
44970b3c8c7SBarry Smith   ierr = ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
45070b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint);CHKERRQ(ierr);
45170b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf);CHKERRQ(ierr);
452064c8009SBarry Smith   ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr);
453064c8009SBarry Smith 
454671e8369SHong Zhang   ierr = ISSort(is);CHKERRQ(ierr);
4557dae84e0SHong Zhang   ierr = MatCreateSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr);
456064c8009SBarry Smith   A    = *Aholder;
457064c8009SBarry Smith   ierr = PetscFree(Aholder);CHKERRQ(ierr);
458064c8009SBarry Smith 
4597dae84e0SHong Zhang   ierr = MatCreateSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr);
4607dae84e0SHong Zhang   ierr = MatCreateSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr);
4617dae84e0SHong Zhang   ierr = MatCreateSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr);
462064c8009SBarry Smith 
463064c8009SBarry Smith   /*
464064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
465064c8009SBarry Smith   */
46692e4c2edSBarry Smith   ierr = MatMatMult(Ais,Xsurf,MAT_INITIAL_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr);
467064c8009SBarry Smith   ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr);
468064c8009SBarry Smith 
4698e722e36SBarry Smith   if (exotic->directSolve) {
4702692d6eeSBarry Smith     ierr = MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr);
471064c8009SBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
4722692d6eeSBarry Smith     ierr = MatGetOrdering(Aii,MATORDERINGND,&row,&col);CHKERRQ(ierr);
473064c8009SBarry Smith     ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr);
474fcfd50ebSBarry Smith     ierr = ISDestroy(&row);CHKERRQ(ierr);
475fcfd50ebSBarry Smith     ierr = ISDestroy(&col);CHKERRQ(ierr);
476064c8009SBarry Smith     ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr);
477064c8009SBarry Smith     ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr);
478fcfd50ebSBarry Smith     ierr = MatDestroy(&iAii);CHKERRQ(ierr);
479064c8009SBarry Smith   } else {
480064c8009SBarry Smith     Vec         b,x;
481064c8009SBarry Smith     PetscScalar *xint_tmp;
482064c8009SBarry Smith 
4838c778c55SBarry Smith     ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr);
484778a2246SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,0,&x);CHKERRQ(ierr);
4858c778c55SBarry Smith     ierr = MatDenseGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
486778a2246SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,0,&b);CHKERRQ(ierr);
48723ee1639SBarry Smith     ierr = KSPSetOperators(exotic->ksp,Aii,Aii);CHKERRQ(ierr);
488064c8009SBarry Smith     for (i=0; i<6; i++) {
489064c8009SBarry Smith       ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr);
490064c8009SBarry Smith       ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr);
4918e722e36SBarry Smith       ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr);
492*c0decd05SBarry Smith       ierr = KSPCheckSolve(exotic->ksp,pc,x);CHKERRQ(ierr);
493064c8009SBarry Smith       ierr = VecResetArray(x);CHKERRQ(ierr);
494064c8009SBarry Smith       ierr = VecResetArray(b);CHKERRQ(ierr);
495064c8009SBarry Smith     }
4968c778c55SBarry Smith     ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr);
4978c778c55SBarry Smith     ierr = MatDenseRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
4986bf464f9SBarry Smith     ierr = VecDestroy(&x);CHKERRQ(ierr);
4996bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
500064c8009SBarry Smith   }
5016bf464f9SBarry Smith   ierr = MatDestroy(&Xint_tmp);CHKERRQ(ierr);
502064c8009SBarry Smith 
503064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
5048c778c55SBarry Smith   ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr);
505064c8009SBarry Smith   for (i=0; i<Nint; i++) {
506064c8009SBarry Smith     tmp = 0.0;
5072fa5cd67SKarl Rupp     for (j=0; j<6; j++) tmp += xint[i+j*Nint];
5082fa5cd67SKarl Rupp 
50957622a8eSBarry Smith     if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xint interpolation at i %D value %g",i,(double)PetscAbsScalar(tmp));
510064c8009SBarry Smith   }
5118c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr);
512064c8009SBarry Smith   /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
513064c8009SBarry Smith #endif
514064c8009SBarry Smith 
515064c8009SBarry Smith 
516064c8009SBarry Smith   /*         total faces    */
517064c8009SBarry Smith   Ntotal =  mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1);
518064c8009SBarry Smith 
519064c8009SBarry Smith   /*
520064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
521064c8009SBarry Smith   */
522064c8009SBarry Smith   cnt = 0;
523064c8009SBarry Smith   { gl[cnt++] = mwidth+1;}
524064c8009SBarry Smith   {
525064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth+1;}
526064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
527064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;}
528064c8009SBarry Smith   }
529064c8009SBarry Smith   { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;}
530064c8009SBarry Smith 
531064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
532064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
533064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,6,gl,gl);CHKERRQ(ierr);
534064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
535785e854fSJed Brown   ierr = PetscMalloc1(6*mp*np*pp,&globals);CHKERRQ(ierr);
536ce94432eSBarry Smith   ierr = MPI_Allgather(gl,6,MPIU_INT,globals,6,MPIU_INT,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
537064c8009SBarry Smith 
538064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
5390298fd71SBarry Smith   ierr = MatGetSize(Aglobal,&Nt,NULL);CHKERRQ(ierr);
54028d20b34SBarry Smith   ierr = PetscTableCreate(Ntotal/3,Nt+1,&ht);CHKERRQ(ierr);
541064c8009SBarry Smith   for (i=0; i<6*mp*np*pp; i++) {
542064c8009SBarry Smith     ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr);
543064c8009SBarry Smith   }
544064c8009SBarry Smith   ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr);
545e32f2f54SBarry Smith   if (cnt != Ntotal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hash table size %D not equal to total number coarse grid points %D",cnt,Ntotal);
546064c8009SBarry Smith   ierr = PetscFree(globals);CHKERRQ(ierr);
547064c8009SBarry Smith   for (i=0; i<6; i++) {
548064c8009SBarry Smith     ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr);
549064c8009SBarry Smith     gl[i]--;
550064c8009SBarry Smith   }
5516bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&ht);CHKERRQ(ierr);
552064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
553064c8009SBarry Smith 
554064c8009SBarry Smith   /* construct global interpolation matrix */
5550298fd71SBarry Smith   ierr = MatGetLocalSize(Aglobal,&Ng,NULL);CHKERRQ(ierr);
556064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
557ce94432eSBarry Smith     ierr = MatCreateAIJ(PetscObjectComm((PetscObject)da),Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,NULL,Nint,NULL,P);CHKERRQ(ierr);
558064c8009SBarry Smith   } else {
559064c8009SBarry Smith     ierr = MatZeroEntries(*P);CHKERRQ(ierr);
560064c8009SBarry Smith   }
561064c8009SBarry Smith   ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
5628c778c55SBarry Smith   ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr);
563064c8009SBarry Smith   ierr = MatSetValues(*P,Nint,IIint,6,gl,xint,INSERT_VALUES);CHKERRQ(ierr);
5648c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr);
5658c778c55SBarry Smith   ierr = MatDenseGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
566064c8009SBarry Smith   ierr = MatSetValues(*P,Nsurf,IIsurf,6,gl,xsurf,INSERT_VALUES);CHKERRQ(ierr);
5678c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
568064c8009SBarry Smith   ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
569064c8009SBarry Smith   ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
570064c8009SBarry Smith   ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr);
571064c8009SBarry Smith 
572064c8009SBarry Smith 
573064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
574064c8009SBarry Smith   {
575064c8009SBarry Smith     Vec         x,y;
576064c8009SBarry Smith     PetscScalar *yy;
577ce94432eSBarry Smith     ierr = VecCreateMPI(PetscObjectComm((PetscObject)da),Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr);
578ce94432eSBarry Smith     ierr = VecCreateMPI(PetscObjectComm((PetscObject)da),PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr);
579064c8009SBarry Smith     ierr = VecSet(x,1.0);CHKERRQ(ierr);
580064c8009SBarry Smith     ierr = MatMult(*P,x,y);CHKERRQ(ierr);
581064c8009SBarry Smith     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
582064c8009SBarry Smith     for (i=0; i<Ng; i++) {
58357622a8eSBarry Smith       if (PetscAbsScalar(yy[i]-1.0) > 1.e-10) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong p interpolation at i %D value %g",i,(double)PetscAbsScalar(yy[i]));
584064c8009SBarry Smith     }
585064c8009SBarry Smith     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
586064c8009SBarry Smith     ierr = VecDestroy(x);CHKERRQ(ierr);
587064c8009SBarry Smith     ierr = VecDestroy(y);CHKERRQ(ierr);
588064c8009SBarry Smith   }
589064c8009SBarry Smith #endif
590064c8009SBarry Smith 
591fcfd50ebSBarry Smith   ierr = MatDestroy(&Aii);CHKERRQ(ierr);
592fcfd50ebSBarry Smith   ierr = MatDestroy(&Ais);CHKERRQ(ierr);
593fcfd50ebSBarry Smith   ierr = MatDestroy(&Asi);CHKERRQ(ierr);
594fcfd50ebSBarry Smith   ierr = MatDestroy(&A);CHKERRQ(ierr);
595fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
596fcfd50ebSBarry Smith   ierr = ISDestroy(&isint);CHKERRQ(ierr);
597fcfd50ebSBarry Smith   ierr = ISDestroy(&issurf);CHKERRQ(ierr);
598fcfd50ebSBarry Smith   ierr = MatDestroy(&Xint);CHKERRQ(ierr);
599fcfd50ebSBarry Smith   ierr = MatDestroy(&Xsurf);CHKERRQ(ierr);
600064c8009SBarry Smith   PetscFunctionReturn(0);
601064c8009SBarry Smith }
602174b6946SBarry Smith 
6037233f9f0SBarry Smith 
6047233f9f0SBarry Smith /*@
6057233f9f0SBarry Smith    PCExoticSetType - Sets the type of coarse grid interpolation to use
6067233f9f0SBarry Smith 
6073f9fe445SBarry Smith    Logically Collective on PC
6087233f9f0SBarry Smith 
6097233f9f0SBarry Smith    Input Parameters:
6107233f9f0SBarry Smith +  pc - the preconditioner context
6117233f9f0SBarry Smith -  type - either PC_EXOTIC_FACE or PC_EXOTIC_WIREBASKET (defaults to face)
6127233f9f0SBarry Smith 
61395452b02SPatrick Sanan    Notes:
61495452b02SPatrick Sanan     The face based interpolation has 1 degree of freedom per face and ignores the
615563e08c6SBarry Smith      edge and vertex values completely in the coarse problem. For any seven point
616563e08c6SBarry Smith      stencil the interpolation of a constant on all faces into the interior is that constant.
617563e08c6SBarry Smith 
618563e08c6SBarry Smith      The wirebasket interpolation has 1 degree of freedom per vertex, per edge and
619563e08c6SBarry Smith      per face. A constant on the subdomain boundary is interpolated as that constant
620563e08c6SBarry Smith      in the interior of the domain.
621563e08c6SBarry Smith 
622563e08c6SBarry Smith      The coarse grid matrix is obtained via the Galerkin computation A_c = R A R^T, hence
623563e08c6SBarry Smith      if A is nonsingular A_c is also nonsingular.
624563e08c6SBarry Smith 
625563e08c6SBarry Smith      Both interpolations are suitable for only scalar problems.
626563e08c6SBarry Smith 
6277233f9f0SBarry Smith    Level: intermediate
6287233f9f0SBarry Smith 
6297233f9f0SBarry Smith 
6307233f9f0SBarry Smith .seealso: PCEXOTIC, PCExoticType()
6317233f9f0SBarry Smith @*/
6327087cfbeSBarry Smith PetscErrorCode  PCExoticSetType(PC pc,PCExoticType type)
6337233f9f0SBarry Smith {
6344ac538c5SBarry Smith   PetscErrorCode ierr;
6357233f9f0SBarry Smith 
6367233f9f0SBarry Smith   PetscFunctionBegin;
6370700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
638c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
6394ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCExoticSetType_C",(PC,PCExoticType),(pc,type));CHKERRQ(ierr);
6407233f9f0SBarry Smith   PetscFunctionReturn(0);
6417233f9f0SBarry Smith }
6427233f9f0SBarry Smith 
643f7a08781SBarry Smith static PetscErrorCode  PCExoticSetType_Exotic(PC pc,PCExoticType type)
6447233f9f0SBarry Smith {
645f3fbd535SBarry Smith   PC_MG     *mg  = (PC_MG*)pc->data;
64631567311SBarry Smith   PC_Exotic *ctx = (PC_Exotic*) mg->innerctx;
6477233f9f0SBarry Smith 
6487233f9f0SBarry Smith   PetscFunctionBegin;
6497233f9f0SBarry Smith   ctx->type = type;
6507233f9f0SBarry Smith   PetscFunctionReturn(0);
6517233f9f0SBarry Smith }
6527233f9f0SBarry Smith 
6537233f9f0SBarry Smith PetscErrorCode PCSetUp_Exotic(PC pc)
654174b6946SBarry Smith {
655174b6946SBarry Smith   PetscErrorCode ierr;
65696bdf778SBarry Smith   Mat            A;
657f3fbd535SBarry Smith   PC_MG          *mg   = (PC_MG*)pc->data;
65831567311SBarry Smith   PC_Exotic      *ex   = (PC_Exotic*) mg->innerctx;
65996bdf778SBarry Smith   MatReuse       reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX;
660174b6946SBarry Smith 
661174b6946SBarry Smith   PetscFunctionBegin;
662ce94432eSBarry Smith   if (!pc->dm) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Need to call PCSetDM() before using this PC");
66323ee1639SBarry Smith   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
6647233f9f0SBarry Smith   if (ex->type == PC_EXOTIC_FACE) {
665*c0decd05SBarry Smith     ierr = DMDAGetFaceInterpolation(pc,pc->dm,ex,A,reuse,&ex->P);CHKERRQ(ierr);
6667233f9f0SBarry Smith   } else if (ex->type == PC_EXOTIC_WIREBASKET) {
667*c0decd05SBarry Smith     ierr = DMDAGetWireBasketInterpolation(pc,pc->dm,ex,A,reuse,&ex->P);CHKERRQ(ierr);
668ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unknown exotic coarse space %d",ex->type);
66996bdf778SBarry Smith   ierr = PCMGSetInterpolation(pc,1,ex->P);CHKERRQ(ierr);
670d0660788SBarry Smith   /* if PC has attached DM we must remove it or the PCMG will use it to compute incorrect sized vectors and interpolations */
6710298fd71SBarry Smith   ierr = PCSetDM(pc,NULL);CHKERRQ(ierr);
6727233f9f0SBarry Smith   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
673174b6946SBarry Smith   PetscFunctionReturn(0);
674174b6946SBarry Smith }
675174b6946SBarry Smith 
6767233f9f0SBarry Smith PetscErrorCode PCDestroy_Exotic(PC pc)
677174b6946SBarry Smith {
678174b6946SBarry Smith   PetscErrorCode ierr;
679f3fbd535SBarry Smith   PC_MG          *mg  = (PC_MG*)pc->data;
68031567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
681174b6946SBarry Smith 
682174b6946SBarry Smith   PetscFunctionBegin;
683fcfd50ebSBarry Smith   ierr = MatDestroy(&ctx->P);CHKERRQ(ierr);
684fcfd50ebSBarry Smith   ierr = KSPDestroy(&ctx->ksp);CHKERRQ(ierr);
6857233f9f0SBarry Smith   ierr = PetscFree(ctx);CHKERRQ(ierr);
6867233f9f0SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
687174b6946SBarry Smith   PetscFunctionReturn(0);
688174b6946SBarry Smith }
689174b6946SBarry Smith 
6907233f9f0SBarry Smith PetscErrorCode PCView_Exotic(PC pc,PetscViewer viewer)
6917233f9f0SBarry Smith {
692f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
6937233f9f0SBarry Smith   PetscErrorCode ierr;
694ace3abfcSBarry Smith   PetscBool      iascii;
69531567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
6967233f9f0SBarry Smith 
6977233f9f0SBarry Smith   PetscFunctionBegin;
698251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
6997233f9f0SBarry Smith   if (iascii) {
7007233f9f0SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    Exotic type = %s\n",PCExoticTypes[ctx->type]);CHKERRQ(ierr);
7018e722e36SBarry Smith     if (ctx->directSolve) {
7028e722e36SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"      Using direct solver to construct interpolation\n");CHKERRQ(ierr);
7038e722e36SBarry Smith     } else {
7048e722e36SBarry Smith       PetscViewer sviewer;
7058e722e36SBarry Smith       PetscMPIInt rank;
7068e722e36SBarry Smith 
7078e722e36SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"      Using iterative solver to construct interpolation\n");CHKERRQ(ierr);
7088e722e36SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
7098e722e36SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* should not need to push this twice? */
7103f08860eSBarry Smith       ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
711ce94432eSBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
7128e722e36SBarry Smith       if (!rank) {
7138e722e36SBarry Smith         ierr = KSPView(ctx->ksp,sviewer);CHKERRQ(ierr);
7148e722e36SBarry Smith       }
7153f08860eSBarry Smith       ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
7168e722e36SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
7178e722e36SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
7188e722e36SBarry Smith     }
7197233f9f0SBarry Smith   }
7207233f9f0SBarry Smith   ierr = PCView_MG(pc,viewer);CHKERRQ(ierr);
7217233f9f0SBarry Smith   PetscFunctionReturn(0);
7227233f9f0SBarry Smith }
7237233f9f0SBarry Smith 
7244416b707SBarry Smith PetscErrorCode PCSetFromOptions_Exotic(PetscOptionItems *PetscOptionsObject,PC pc)
7257233f9f0SBarry Smith {
7267233f9f0SBarry Smith   PetscErrorCode ierr;
727ace3abfcSBarry Smith   PetscBool      flg;
728f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
7297233f9f0SBarry Smith   PCExoticType   mgctype;
73031567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
7317233f9f0SBarry Smith 
7327233f9f0SBarry Smith   PetscFunctionBegin;
733e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Exotic coarse space options");CHKERRQ(ierr);
7347233f9f0SBarry Smith   ierr = PetscOptionsEnum("-pc_exotic_type","face or wirebasket","PCExoticSetType",PCExoticTypes,(PetscEnum)ctx->type,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
7357233f9f0SBarry Smith   if (flg) {
7367233f9f0SBarry Smith     ierr = PCExoticSetType(pc,mgctype);CHKERRQ(ierr);
7377233f9f0SBarry Smith   }
7380298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_exotic_direct_solver","use direct solver to construct interpolation","None",ctx->directSolve,&ctx->directSolve,NULL);CHKERRQ(ierr);
7398e722e36SBarry Smith   if (!ctx->directSolve) {
7408e722e36SBarry Smith     if (!ctx->ksp) {
7418e722e36SBarry Smith       const char *prefix;
7428e722e36SBarry Smith       ierr = KSPCreate(PETSC_COMM_SELF,&ctx->ksp);CHKERRQ(ierr);
743422a814eSBarry Smith       ierr = KSPSetErrorIfNotConverged(ctx->ksp,pc->erroriffailure);CHKERRQ(ierr);
7448e722e36SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
7453bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ctx->ksp);CHKERRQ(ierr);
7468e722e36SBarry Smith       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
7478e722e36SBarry Smith       ierr = KSPSetOptionsPrefix(ctx->ksp,prefix);CHKERRQ(ierr);
7488e722e36SBarry Smith       ierr = KSPAppendOptionsPrefix(ctx->ksp,"exotic_");CHKERRQ(ierr);
7498e722e36SBarry Smith     }
7508e722e36SBarry Smith     ierr = KSPSetFromOptions(ctx->ksp);CHKERRQ(ierr);
7518e722e36SBarry Smith   }
7527233f9f0SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7537233f9f0SBarry Smith   PetscFunctionReturn(0);
7547233f9f0SBarry Smith }
7557233f9f0SBarry Smith 
756174b6946SBarry Smith 
757174b6946SBarry Smith /*MC
7587233f9f0SBarry Smith      PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces
759174b6946SBarry Smith 
7607233f9f0SBarry Smith      This uses the PCMG infrastructure restricted to two levels and the face and wirebasket based coarse
76124c3aa18SBarry Smith    grid spaces.
76224c3aa18SBarry Smith 
76395452b02SPatrick Sanan    Notes:
76495452b02SPatrick Sanan     By default this uses GMRES on the fine grid smoother so this should be used with KSPFGMRES or the smoother changed to not use GMRES
76524c3aa18SBarry Smith 
76696a0c994SBarry Smith    References:
76796a0c994SBarry Smith +  1. - These coarse grid spaces originate in the work of Bramble, Pasciak  and Schatz, "The Construction
76896a0c994SBarry Smith    of Preconditioners for Elliptic Problems by Substructing IV", Mathematics of Computation, volume 53, 1989.
76996a0c994SBarry Smith .  2. - They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith,
7704f02bc6aSBarry Smith    New York University, 1990.
77196a0c994SBarry Smith .  3. - They were then explored in great detail in Dryja, Smith, Widlund, "Schwarz Analysis
7723b65e785SBarry Smith    of Iterative Substructuring Methods for Elliptic Problems in Three Dimensions, SIAM Journal on Numerical
77396a0c994SBarry Smith    Analysis, volume 31. 1994. These were developed in the context of iterative substructuring preconditioners.
77496a0c994SBarry Smith .  4. - They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund.
7753b65e785SBarry Smith    They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example,
7763b65e785SBarry Smith    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Extending theory for domain decomposition algorithms to irregular subdomains. In Ulrich Langer, Marco
7773b65e785SBarry Smith    Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
7783b65e785SBarry Smith    of the 17th International Conference on Domain Decomposition Methods in
77996a0c994SBarry Smith    Science and Engineering, held in Strobl, Austria, 2006, number 60 in
78096a0c994SBarry Smith    Springer Verlag, Lecture Notes in Computational Science and Engineering, 2007.
78196a0c994SBarry Smith .  5. -  Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. A family of energy minimizing coarse spaces for overlapping Schwarz preconditioners. In Ulrich Langer,
7823b65e785SBarry Smith    Marco Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
7833b65e785SBarry Smith    of the 17th International Conference on Domain Decomposition Methods
78496a0c994SBarry Smith    in Science and Engineering, held in Strobl, Austria, 2006, number 60 in
78596a0c994SBarry Smith    Springer Verlag, Lecture Notes in Computational Science and Engineering, 2007
78696a0c994SBarry Smith .  6. - Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Domain decomposition
7873b65e785SBarry Smith    for less regular subdomains: Overlapping Schwarz in two dimensions. SIAM J.
78896a0c994SBarry Smith    Numer. Anal., 46(4), 2008.
78996a0c994SBarry Smith -  7. - Clark R. Dohrmann and Olof B. Widlund. An overlapping Schwarz
7903b65e785SBarry Smith    algorithm for almost incompressible elasticity. Technical Report
79196a0c994SBarry Smith    TR2008 912, Department of Computer Science, Courant Institute
7923b65e785SBarry Smith    of Mathematical Sciences, New York University, May 2008. URL:
7937233f9f0SBarry Smith 
7947233f9f0SBarry Smith    Options Database: The usual PCMG options are supported, such as -mg_levels_pc_type <type> -mg_coarse_pc_type <type>
7957233f9f0SBarry Smith       -pc_mg_type <type>
7967233f9f0SBarry Smith 
79725a35f6fSSatish Balay    Level: advanced
798174b6946SBarry Smith 
7996c699258SBarry Smith .seealso:  PCMG, PCSetDM(), PCExoticType, PCExoticSetType()
800174b6946SBarry Smith M*/
801174b6946SBarry Smith 
8028cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Exotic(PC pc)
803174b6946SBarry Smith {
804174b6946SBarry Smith   PetscErrorCode ierr;
8057233f9f0SBarry Smith   PC_Exotic      *ex;
806f3fbd535SBarry Smith   PC_MG          *mg;
807174b6946SBarry Smith 
808174b6946SBarry Smith   PetscFunctionBegin;
809f91d8e95SBarry Smith   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
8102fa5cd67SKarl Rupp   if (pc->ops->destroy) {
8112fa5cd67SKarl Rupp     ierr =  (*pc->ops->destroy)(pc);CHKERRQ(ierr);
8122fa5cd67SKarl Rupp     pc->data = 0;
8132fa5cd67SKarl Rupp   }
814503cfb0cSBarry Smith   ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr);
815f91d8e95SBarry Smith   ((PetscObject)pc)->type_name = 0;
816f91d8e95SBarry Smith 
817174b6946SBarry Smith   ierr         = PCSetType(pc,PCMG);CHKERRQ(ierr);
8180298fd71SBarry Smith   ierr         = PCMGSetLevels(pc,2,NULL);CHKERRQ(ierr);
8192134b1e4SBarry Smith   ierr         = PCMGSetGalerkin(pc,PC_MG_GALERKIN_PMAT);CHKERRQ(ierr);
820b00a9115SJed Brown   ierr         = PetscNew(&ex);CHKERRQ(ierr); \
8217233f9f0SBarry Smith   ex->type     = PC_EXOTIC_FACE;
822f3fbd535SBarry Smith   mg           = (PC_MG*) pc->data;
82331567311SBarry Smith   mg->innerctx = ex;
8247233f9f0SBarry Smith 
8257233f9f0SBarry Smith 
8267233f9f0SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_Exotic;
8277233f9f0SBarry Smith   pc->ops->view           = PCView_Exotic;
8287233f9f0SBarry Smith   pc->ops->destroy        = PCDestroy_Exotic;
8296c699258SBarry Smith   pc->ops->setup          = PCSetUp_Exotic;
8302fa5cd67SKarl Rupp 
831bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCExoticSetType_C",PCExoticSetType_Exotic);CHKERRQ(ierr);
832174b6946SBarry Smith   PetscFunctionReturn(0);
833174b6946SBarry Smith }
834