xref: /petsc/src/ksp/pc/impls/wb/wb.c (revision 2fa5cd679192b9b390e47ae2d0650965e6b1d9fa)
1174b6946SBarry Smith 
2c6db04a5SJed Brown #include <petscpcmg.h>   /*I "petscpcmg.h" I*/
3c6db04a5SJed Brown #include <petscdmda.h>   /*I "petscdmda.h" I*/
4c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.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 #undef __FUNCT__
17aa219208SBarry Smith #define __FUNCT__ "DMDAGetWireBasketInterpolation"
18064c8009SBarry Smith /*
19aa219208SBarry Smith       DMDAGetWireBasketInterpolation - Gets the interpolation for a wirebasket based coarse space
20064c8009SBarry Smith 
21064c8009SBarry Smith */
22aa219208SBarry Smith PetscErrorCode DMDAGetWireBasketInterpolation(DM da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P)
23064c8009SBarry Smith {
24064c8009SBarry Smith   PetscErrorCode         ierr;
25064c8009SBarry 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;
2628d20b34SBarry Smith   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[26],*globals,Ng,*IIint,*IIsurf,Nt;
27064c8009SBarry Smith   Mat                    Xint, Xsurf,Xint_tmp;
28064c8009SBarry Smith   IS                     isint,issurf,is,row,col;
29064c8009SBarry Smith   ISLocalToGlobalMapping ltg;
30064c8009SBarry Smith   MPI_Comm               comm;
31064c8009SBarry Smith   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
32064c8009SBarry Smith   MatFactorInfo          info;
33064c8009SBarry Smith   PetscScalar            *xsurf,*xint;
348e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
35064c8009SBarry Smith   PetscScalar            tmp;
36064c8009SBarry Smith #endif
37064c8009SBarry Smith   PetscTable             ht;
38064c8009SBarry Smith 
39064c8009SBarry Smith   PetscFunctionBegin;
401321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0,0,0);CHKERRQ(ierr);
41e7e72b3dSBarry Smith   if (dof != 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only for single field problems");
42e7e72b3dSBarry Smith   if (dim != 3) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only coded for 3d problems");
43aa219208SBarry Smith   ierr   = DMDAGetCorners(da,0,0,0,&m,&n,&p);CHKERRQ(ierr);
44aa219208SBarry Smith   ierr   = DMDAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr);
45064c8009SBarry Smith   istart = istart ? -1 : 0;
46064c8009SBarry Smith   jstart = jstart ? -1 : 0;
47064c8009SBarry Smith   kstart = kstart ? -1 : 0;
48064c8009SBarry Smith 
49064c8009SBarry Smith   /*
50064c8009SBarry Smith     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge)
51064c8009SBarry Smith     to all the local degrees of freedom (this includes the vertices, edges and faces).
52064c8009SBarry Smith 
53064c8009SBarry Smith     Xint are the subset of the interpolation into the interior
54064c8009SBarry Smith 
55064c8009SBarry Smith     Xface are the interpolation onto faces but not into the interior
56064c8009SBarry Smith 
57064c8009SBarry Smith     Xsurf are the interpolation onto the vertices and edges (the surfbasket)
58064c8009SBarry Smith                                       Xint
59064c8009SBarry Smith     Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain
60064c8009SBarry Smith                                       Xsurf
61064c8009SBarry Smith   */
62064c8009SBarry Smith   N     = (m - istart)*(n - jstart)*(p - kstart);
63064c8009SBarry Smith   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
64064c8009SBarry Smith   Nface = 2*((m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart));
65064c8009SBarry Smith   Nwire = 4*((m-2-istart) + (n-2-jstart) + (p-2-kstart)) + 8;
66064c8009SBarry Smith   Nsurf = Nface + Nwire;
67064c8009SBarry Smith   ierr  = MatCreateSeqDense(MPI_COMM_SELF,Nint,26,PETSC_NULL,&Xint);CHKERRQ(ierr);
68064c8009SBarry Smith   ierr  = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,26,PETSC_NULL,&Xsurf);CHKERRQ(ierr);
698c778c55SBarry Smith   ierr  = MatDenseGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
70064c8009SBarry Smith 
71064c8009SBarry Smith   /*
72064c8009SBarry Smith      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
73064c8009SBarry Smith      Xsurf will be all zero (thus making the coarse matrix singular).
74064c8009SBarry Smith   */
75e32f2f54SBarry 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");
76e32f2f54SBarry 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");
77e32f2f54SBarry 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");
78064c8009SBarry Smith 
79064c8009SBarry Smith   cnt = 0;
80*2fa5cd67SKarl Rupp 
81*2fa5cd67SKarl Rupp   xsurf[cnt++] = 1;
82*2fa5cd67SKarl Rupp   for (i=1; i<m-istart-1; i++) xsurf[cnt++ + Nsurf] = 1;
83*2fa5cd67SKarl Rupp   xsurf[cnt++ + 2*Nsurf] = 1;
84*2fa5cd67SKarl Rupp 
85*2fa5cd67SKarl Rupp   for (j=1; j<n-1-jstart; j++) {
86*2fa5cd67SKarl Rupp     xsurf[cnt++ + 3*Nsurf] = 1;
87*2fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 4*Nsurf] = 1;
88*2fa5cd67SKarl Rupp     xsurf[cnt++ + 5*Nsurf] = 1;
89064c8009SBarry Smith   }
90*2fa5cd67SKarl Rupp 
91*2fa5cd67SKarl Rupp   xsurf[cnt++ + 6*Nsurf] = 1;
92*2fa5cd67SKarl Rupp   for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 7*Nsurf] = 1;
93*2fa5cd67SKarl Rupp   xsurf[cnt++ + 8*Nsurf] = 1;
94*2fa5cd67SKarl Rupp 
95*2fa5cd67SKarl Rupp   for (k=1; k<p-1-kstart; k++) {
96*2fa5cd67SKarl Rupp     xsurf[cnt++ + 9*Nsurf] = 1;
97*2fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 10*Nsurf] = 1;
98*2fa5cd67SKarl Rupp     xsurf[cnt++ + 11*Nsurf] = 1;
99*2fa5cd67SKarl Rupp 
100*2fa5cd67SKarl Rupp     for (j=1; j<n-1-jstart; j++) {
101*2fa5cd67SKarl Rupp       xsurf[cnt++ + 12*Nsurf] = 1;
102*2fa5cd67SKarl Rupp       /* these are the interior nodes */
103*2fa5cd67SKarl Rupp       xsurf[cnt++ + 13*Nsurf] = 1;
104*2fa5cd67SKarl Rupp     }
105*2fa5cd67SKarl Rupp 
106*2fa5cd67SKarl Rupp     xsurf[cnt++ + 14*Nsurf] = 1;
107*2fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 15*Nsurf] = 1;
108*2fa5cd67SKarl Rupp     xsurf[cnt++ + 16*Nsurf] = 1;
109*2fa5cd67SKarl Rupp   }
110*2fa5cd67SKarl Rupp 
111*2fa5cd67SKarl Rupp   xsurf[cnt++ + 17*Nsurf] = 1;
112*2fa5cd67SKarl Rupp   for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 18*Nsurf] = 1;
113*2fa5cd67SKarl Rupp   xsurf[cnt++ + 19*Nsurf] = 1;
114*2fa5cd67SKarl Rupp 
115*2fa5cd67SKarl Rupp   for (j=1;j<n-1-jstart;j++) {
116*2fa5cd67SKarl Rupp     xsurf[cnt++ + 20*Nsurf] = 1;
117*2fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 21*Nsurf] = 1;
118*2fa5cd67SKarl Rupp     xsurf[cnt++ + 22*Nsurf] = 1;
119*2fa5cd67SKarl Rupp   }
120*2fa5cd67SKarl Rupp 
121*2fa5cd67SKarl Rupp   xsurf[cnt++ + 23*Nsurf] = 1;
122*2fa5cd67SKarl Rupp   for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 24*Nsurf] = 1;
123*2fa5cd67SKarl Rupp   xsurf[cnt++ + 25*Nsurf] = 1;
124*2fa5cd67SKarl Rupp 
125064c8009SBarry Smith 
1268e722e36SBarry Smith   /* interpolations only sum to 1 when using direct solver */
1278e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
128064c8009SBarry Smith   for (i=0; i<Nsurf; i++) {
129064c8009SBarry Smith     tmp = 0.0;
130*2fa5cd67SKarl Rupp     for (j=0; j<26; j++) tmp += xsurf[i+j*Nsurf];
131e32f2f54SBarry 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,PetscAbsScalar(tmp));
132064c8009SBarry Smith   }
133064c8009SBarry Smith #endif
1348c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
135064c8009SBarry Smith   /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/
136064c8009SBarry Smith 
137064c8009SBarry Smith 
138064c8009SBarry Smith   /*
139064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
140064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
141064c8009SBarry Smith             (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it
142aa219208SBarry Smith              is NOT the local DMDA ordering.)
143064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
144064c8009SBarry Smith   */
145064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
146064c8009SBarry Smith   ierr = PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);CHKERRQ(ierr);
147064c8009SBarry Smith   ierr = PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);CHKERRQ(ierr);
148064c8009SBarry Smith   for (k=0; k<p-kstart; k++) {
149064c8009SBarry Smith     for (j=0; j<n-jstart; j++) {
150064c8009SBarry Smith       for (i=0; i<m-istart; i++) {
151064c8009SBarry Smith         II[c++] = i + j*mwidth + k*mwidth*nwidth;
152064c8009SBarry Smith 
153064c8009SBarry Smith         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
154064c8009SBarry Smith           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
155064c8009SBarry Smith           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
156064c8009SBarry Smith         } else {
157064c8009SBarry Smith           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
158064c8009SBarry Smith           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
159064c8009SBarry Smith         }
160064c8009SBarry Smith       }
161064c8009SBarry Smith     }
162064c8009SBarry Smith   }
163e32f2f54SBarry Smith   if (c != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N");
164e32f2f54SBarry Smith   if (cint != Nint) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint");
165e32f2f54SBarry Smith   if (csurf != Nsurf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf");
1661411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltg);CHKERRQ(ierr);
167064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr);
168064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr);
169064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr);
170064c8009SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
17170b3c8c7SBarry Smith   ierr = ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
17270b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint);CHKERRQ(ierr);
17370b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf);CHKERRQ(ierr);
174064c8009SBarry Smith   ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr);
175064c8009SBarry Smith 
176064c8009SBarry Smith   ierr = MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr);
177064c8009SBarry Smith   A    = *Aholder;
178064c8009SBarry Smith   ierr = PetscFree(Aholder);CHKERRQ(ierr);
179064c8009SBarry Smith 
180064c8009SBarry Smith   ierr = MatGetSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr);
181064c8009SBarry Smith   ierr = MatGetSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr);
182064c8009SBarry Smith   ierr = MatGetSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr);
183064c8009SBarry Smith 
184064c8009SBarry Smith   /*
185064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
186064c8009SBarry Smith   */
18792e4c2edSBarry Smith   ierr = MatMatMult(Ais,Xsurf,MAT_INITIAL_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr);
1888e722e36SBarry Smith   ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr);
1898e722e36SBarry Smith   if (exotic->directSolve) {
1902692d6eeSBarry Smith     ierr = MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr);
191064c8009SBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
1922692d6eeSBarry Smith     ierr = MatGetOrdering(Aii,MATORDERINGND,&row,&col);CHKERRQ(ierr);
193064c8009SBarry Smith     ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr);
194fcfd50ebSBarry Smith     ierr = ISDestroy(&row);CHKERRQ(ierr);
195fcfd50ebSBarry Smith     ierr = ISDestroy(&col);CHKERRQ(ierr);
196064c8009SBarry Smith     ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr);
197064c8009SBarry Smith     ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr);
1986bf464f9SBarry Smith     ierr = MatDestroy(&iAii);CHKERRQ(ierr);
1998e722e36SBarry Smith   } else {
2008e722e36SBarry Smith     Vec         b,x;
2018e722e36SBarry Smith     PetscScalar *xint_tmp;
202064c8009SBarry Smith 
2038c778c55SBarry Smith     ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr);
204778a2246SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,0,&x);CHKERRQ(ierr);
2058c778c55SBarry Smith     ierr = MatDenseGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
206778a2246SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,0,&b);CHKERRQ(ierr);
2078e722e36SBarry Smith     ierr = KSPSetOperators(exotic->ksp,Aii,Aii,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
2088e722e36SBarry Smith     for (i=0; i<26; i++) {
2098e722e36SBarry Smith       ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr);
2108e722e36SBarry Smith       ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr);
2118e722e36SBarry Smith       ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr);
2128e722e36SBarry Smith       ierr = VecResetArray(x);CHKERRQ(ierr);
2138e722e36SBarry Smith       ierr = VecResetArray(b);CHKERRQ(ierr);
2148e722e36SBarry Smith     }
2158c778c55SBarry Smith     ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr);
2168c778c55SBarry Smith     ierr = MatDenseRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
2176bf464f9SBarry Smith     ierr = VecDestroy(&x);CHKERRQ(ierr);
2186bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
2198e722e36SBarry Smith   }
2206bf464f9SBarry Smith   ierr = MatDestroy(&Xint_tmp);CHKERRQ(ierr);
2218e722e36SBarry Smith 
2228e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
2238c778c55SBarry Smith   ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr);
224064c8009SBarry Smith   for (i=0; i<Nint; i++) {
225064c8009SBarry Smith     tmp = 0.0;
226*2fa5cd67SKarl Rupp     for (j=0; j<26; j++) tmp += xint[i+j*Nint];
227*2fa5cd67SKarl Rupp 
228e32f2f54SBarry 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,PetscAbsScalar(tmp));
229064c8009SBarry Smith   }
2308c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr);
231064c8009SBarry Smith   /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
232064c8009SBarry Smith #endif
233064c8009SBarry Smith 
234064c8009SBarry Smith 
235064c8009SBarry Smith   /*         total vertices             total faces                                  total edges */
236064c8009SBarry 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);
237064c8009SBarry Smith 
238064c8009SBarry Smith   /*
239064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
240064c8009SBarry Smith   */
241064c8009SBarry Smith   cnt = 0;
242*2fa5cd67SKarl Rupp 
243064c8009SBarry Smith   gl[cnt++] = 0;  { gl[cnt++] = 1;} gl[cnt++] = m-istart-1;
244064c8009SBarry Smith   { gl[cnt++] = mwidth;  { gl[cnt++] = mwidth+1;} gl[cnt++] = mwidth + m-istart-1;}
245064c8009SBarry Smith   gl[cnt++] = mwidth*(n-jstart-1);  { gl[cnt++] = mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*(n-jstart-1) + m-istart-1;
246064c8009SBarry Smith   {
247064c8009SBarry Smith     gl[cnt++] = mwidth*nwidth;  { gl[cnt++] = mwidth*nwidth+1;}  gl[cnt++] = mwidth*nwidth+ m-istart-1;
248064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
249064c8009SBarry 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;
250064c8009SBarry Smith   }
251064c8009SBarry 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;
252064c8009SBarry 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;}
253064c8009SBarry 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;
254064c8009SBarry Smith 
255064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
256064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
257064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,26,gl,gl);CHKERRQ(ierr);
258064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
259064c8009SBarry Smith   ierr = PetscMalloc(26*mp*np*pp*sizeof(PetscInt),&globals);CHKERRQ(ierr);
260064c8009SBarry Smith   ierr = MPI_Allgather(gl,26,MPIU_INT,globals,26,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
261064c8009SBarry Smith 
262064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
26328d20b34SBarry Smith   ierr = MatGetSize(Aglobal,&Nt,PETSC_NULL);CHKERRQ(ierr);
26412dd4568SBarry Smith   ierr = PetscTableCreate(Ntotal/3,Nt+1,&ht);CHKERRQ(ierr);
265064c8009SBarry Smith   for (i=0; i<26*mp*np*pp; i++) {
266064c8009SBarry Smith     ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr);
267064c8009SBarry Smith   }
268064c8009SBarry Smith   ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr);
269e32f2f54SBarry 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);
270064c8009SBarry Smith   ierr = PetscFree(globals);CHKERRQ(ierr);
271064c8009SBarry Smith   for (i=0; i<26; i++) {
272064c8009SBarry Smith     ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr);
273064c8009SBarry Smith     gl[i]--;
274064c8009SBarry Smith   }
2756bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&ht);CHKERRQ(ierr);
276064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
277064c8009SBarry Smith 
278064c8009SBarry Smith   /* construct global interpolation matrix */
27928d20b34SBarry Smith   ierr = MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);CHKERRQ(ierr);
280064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
28169b1f4b7SBarry Smith     ierr = MatCreateAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint+Nsurf,PETSC_NULL,P);CHKERRQ(ierr);
282064c8009SBarry Smith   } else {
283064c8009SBarry Smith     ierr = MatZeroEntries(*P);CHKERRQ(ierr);
284064c8009SBarry Smith   }
285064c8009SBarry Smith   ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
2868c778c55SBarry Smith   ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr);
287064c8009SBarry Smith   ierr = MatSetValues(*P,Nint,IIint,26,gl,xint,INSERT_VALUES);CHKERRQ(ierr);
2888c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr);
2898c778c55SBarry Smith   ierr = MatDenseGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
290064c8009SBarry Smith   ierr = MatSetValues(*P,Nsurf,IIsurf,26,gl,xsurf,INSERT_VALUES);CHKERRQ(ierr);
2918c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
292064c8009SBarry Smith   ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
293064c8009SBarry Smith   ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
294064c8009SBarry Smith   ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr);
295064c8009SBarry Smith 
2968e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
297064c8009SBarry Smith   {
298064c8009SBarry Smith     Vec         x,y;
299064c8009SBarry Smith     PetscScalar *yy;
300064c8009SBarry Smith     ierr = VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr);
301064c8009SBarry Smith     ierr = VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr);
302064c8009SBarry Smith     ierr = VecSet(x,1.0);CHKERRQ(ierr);
303064c8009SBarry Smith     ierr = MatMult(*P,x,y);CHKERRQ(ierr);
304064c8009SBarry Smith     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
305064c8009SBarry Smith     for (i=0; i<Ng; i++) {
306e32f2f54SBarry 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,PetscAbsScalar(yy[i]));
307064c8009SBarry Smith     }
308064c8009SBarry Smith     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
309064c8009SBarry Smith     ierr = VecDestroy(x);CHKERRQ(ierr);
310064c8009SBarry Smith     ierr = VecDestroy(y);CHKERRQ(ierr);
311064c8009SBarry Smith   }
312064c8009SBarry Smith #endif
313064c8009SBarry Smith 
314fcfd50ebSBarry Smith   ierr = MatDestroy(&Aii);CHKERRQ(ierr);
315fcfd50ebSBarry Smith   ierr = MatDestroy(&Ais);CHKERRQ(ierr);
316fcfd50ebSBarry Smith   ierr = MatDestroy(&Asi);CHKERRQ(ierr);
317fcfd50ebSBarry Smith   ierr = MatDestroy(&A);CHKERRQ(ierr);
318fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
319fcfd50ebSBarry Smith   ierr = ISDestroy(&isint);CHKERRQ(ierr);
320fcfd50ebSBarry Smith   ierr = ISDestroy(&issurf);CHKERRQ(ierr);
321fcfd50ebSBarry Smith   ierr = MatDestroy(&Xint);CHKERRQ(ierr);
322fcfd50ebSBarry Smith   ierr = MatDestroy(&Xsurf);CHKERRQ(ierr);
323064c8009SBarry Smith   PetscFunctionReturn(0);
324064c8009SBarry Smith }
325064c8009SBarry Smith 
326064c8009SBarry Smith #undef __FUNCT__
327aa219208SBarry Smith #define __FUNCT__ "DMDAGetFaceInterpolation"
328064c8009SBarry Smith /*
329aa219208SBarry Smith       DMDAGetFaceInterpolation - Gets the interpolation for a face based coarse space
330064c8009SBarry Smith 
331064c8009SBarry Smith */
332aa219208SBarry Smith PetscErrorCode DMDAGetFaceInterpolation(DM da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P)
333064c8009SBarry Smith {
334064c8009SBarry Smith   PetscErrorCode         ierr;
335064c8009SBarry 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;
33628d20b34SBarry Smith   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[6],*globals,Ng,*IIint,*IIsurf,Nt;
337064c8009SBarry Smith   Mat                    Xint, Xsurf,Xint_tmp;
338064c8009SBarry Smith   IS                     isint,issurf,is,row,col;
339064c8009SBarry Smith   ISLocalToGlobalMapping ltg;
340064c8009SBarry Smith   MPI_Comm               comm;
341064c8009SBarry Smith   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
342064c8009SBarry Smith   MatFactorInfo          info;
343064c8009SBarry Smith   PetscScalar            *xsurf,*xint;
344064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
345064c8009SBarry Smith   PetscScalar            tmp;
346064c8009SBarry Smith #endif
347064c8009SBarry Smith   PetscTable             ht;
348064c8009SBarry Smith 
349064c8009SBarry Smith   PetscFunctionBegin;
3501321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0,0,0);CHKERRQ(ierr);
351e7e72b3dSBarry Smith   if (dof != 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only for single field problems");
352e7e72b3dSBarry Smith   if (dim != 3) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only coded for 3d problems");
353aa219208SBarry Smith   ierr   = DMDAGetCorners(da,0,0,0,&m,&n,&p);CHKERRQ(ierr);
354aa219208SBarry Smith   ierr   = DMDAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr);
355064c8009SBarry Smith   istart = istart ? -1 : 0;
356064c8009SBarry Smith   jstart = jstart ? -1 : 0;
357064c8009SBarry Smith   kstart = kstart ? -1 : 0;
358064c8009SBarry Smith 
359064c8009SBarry Smith   /*
360064c8009SBarry Smith     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge)
361064c8009SBarry Smith     to all the local degrees of freedom (this includes the vertices, edges and faces).
362064c8009SBarry Smith 
363064c8009SBarry Smith     Xint are the subset of the interpolation into the interior
364064c8009SBarry Smith 
365064c8009SBarry Smith     Xface are the interpolation onto faces but not into the interior
366064c8009SBarry Smith 
367064c8009SBarry Smith     Xsurf are the interpolation onto the vertices and edges (the surfbasket)
368064c8009SBarry Smith                                       Xint
369064c8009SBarry Smith     Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain
370064c8009SBarry Smith                                       Xsurf
371064c8009SBarry Smith   */
372064c8009SBarry Smith   N     = (m - istart)*(n - jstart)*(p - kstart);
373064c8009SBarry Smith   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
374064c8009SBarry Smith   Nface = 2*((m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart));
375064c8009SBarry Smith   Nwire = 4*((m-2-istart) + (n-2-jstart) + (p-2-kstart)) + 8;
376064c8009SBarry Smith   Nsurf = Nface + Nwire;
377064c8009SBarry Smith   ierr  = MatCreateSeqDense(MPI_COMM_SELF,Nint,6,PETSC_NULL,&Xint);CHKERRQ(ierr);
378064c8009SBarry Smith   ierr  = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,6,PETSC_NULL,&Xsurf);CHKERRQ(ierr);
3798c778c55SBarry Smith   ierr  = MatDenseGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
380064c8009SBarry Smith 
381064c8009SBarry Smith   /*
382064c8009SBarry Smith      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
383064c8009SBarry Smith      Xsurf will be all zero (thus making the coarse matrix singular).
384064c8009SBarry Smith   */
385e32f2f54SBarry 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");
386e32f2f54SBarry 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");
387e32f2f54SBarry 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");
388064c8009SBarry Smith 
389064c8009SBarry Smith   cnt = 0;
390*2fa5cd67SKarl Rupp   for (j=1; j<n-1-jstart; j++) {
391*2fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 0*Nsurf] = 1;
392064c8009SBarry Smith   }
393*2fa5cd67SKarl Rupp 
394*2fa5cd67SKarl Rupp   for (k=1; k<p-1-kstart; k++) {
395*2fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 1*Nsurf] = 1;
396*2fa5cd67SKarl Rupp     for (j=1; j<n-1-jstart; j++) {
397*2fa5cd67SKarl Rupp       xsurf[cnt++ + 2*Nsurf] = 1;
398*2fa5cd67SKarl Rupp       /* these are the interior nodes */
399*2fa5cd67SKarl Rupp       xsurf[cnt++ + 3*Nsurf] = 1;
400*2fa5cd67SKarl Rupp     }
401*2fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 4*Nsurf] = 1;
402*2fa5cd67SKarl Rupp   }
403*2fa5cd67SKarl Rupp   for (j=1;j<n-1-jstart;j++) {
404*2fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 5*Nsurf] = 1;
405*2fa5cd67SKarl Rupp   }
406064c8009SBarry Smith 
407064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
408064c8009SBarry Smith   for (i=0; i<Nsurf; i++) {
409064c8009SBarry Smith     tmp = 0.0;
410*2fa5cd67SKarl Rupp     for (j=0; j<6; j++) tmp += xsurf[i+j*Nsurf];
411*2fa5cd67SKarl Rupp 
412e32f2f54SBarry 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,PetscAbsScalar(tmp));
413064c8009SBarry Smith   }
414064c8009SBarry Smith #endif
4158c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
416064c8009SBarry Smith   /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/
417064c8009SBarry Smith 
418064c8009SBarry Smith 
419064c8009SBarry Smith   /*
420064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
421064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
422064c8009SBarry Smith             (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it
423aa219208SBarry Smith              is NOT the local DMDA ordering.)
424064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
425064c8009SBarry Smith   */
426064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
427064c8009SBarry Smith   ierr = PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);CHKERRQ(ierr);
428064c8009SBarry Smith   ierr = PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);CHKERRQ(ierr);
429064c8009SBarry Smith   for (k=0; k<p-kstart; k++) {
430064c8009SBarry Smith     for (j=0; j<n-jstart; j++) {
431064c8009SBarry Smith       for (i=0; i<m-istart; i++) {
432064c8009SBarry Smith         II[c++] = i + j*mwidth + k*mwidth*nwidth;
433064c8009SBarry Smith 
434064c8009SBarry Smith         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
435064c8009SBarry Smith           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
436064c8009SBarry Smith           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
437064c8009SBarry Smith         } else {
438064c8009SBarry Smith           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
439064c8009SBarry Smith           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
440064c8009SBarry Smith         }
441064c8009SBarry Smith       }
442064c8009SBarry Smith     }
443064c8009SBarry Smith   }
444e32f2f54SBarry Smith   if (c != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N");
445e32f2f54SBarry Smith   if (cint != Nint) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint");
446e32f2f54SBarry Smith   if (csurf != Nsurf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf");
4471411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltg);CHKERRQ(ierr);
448064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr);
449064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr);
450064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr);
451064c8009SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
45270b3c8c7SBarry Smith   ierr = ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
45370b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint);CHKERRQ(ierr);
45470b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf);CHKERRQ(ierr);
455064c8009SBarry Smith   ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr);
456064c8009SBarry Smith 
457671e8369SHong Zhang   ierr = ISSort(is);CHKERRQ(ierr);
458064c8009SBarry Smith   ierr = MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr);
459064c8009SBarry Smith   A    = *Aholder;
460064c8009SBarry Smith   ierr = PetscFree(Aholder);CHKERRQ(ierr);
461064c8009SBarry Smith 
462064c8009SBarry Smith   ierr = MatGetSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr);
463064c8009SBarry Smith   ierr = MatGetSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr);
464064c8009SBarry Smith   ierr = MatGetSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr);
465064c8009SBarry Smith 
466064c8009SBarry Smith   /*
467064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
468064c8009SBarry Smith   */
46992e4c2edSBarry Smith   ierr = MatMatMult(Ais,Xsurf,MAT_INITIAL_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr);
470064c8009SBarry Smith   ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr);
471064c8009SBarry Smith 
4728e722e36SBarry Smith   if (exotic->directSolve) {
4732692d6eeSBarry Smith     ierr = MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr);
474064c8009SBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
4752692d6eeSBarry Smith     ierr = MatGetOrdering(Aii,MATORDERINGND,&row,&col);CHKERRQ(ierr);
476064c8009SBarry Smith     ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr);
477fcfd50ebSBarry Smith     ierr = ISDestroy(&row);CHKERRQ(ierr);
478fcfd50ebSBarry Smith     ierr = ISDestroy(&col);CHKERRQ(ierr);
479064c8009SBarry Smith     ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr);
480064c8009SBarry Smith     ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr);
481fcfd50ebSBarry Smith     ierr = MatDestroy(&iAii);CHKERRQ(ierr);
482064c8009SBarry Smith   } else {
483064c8009SBarry Smith     Vec         b,x;
484064c8009SBarry Smith     PetscScalar *xint_tmp;
485064c8009SBarry Smith 
4868c778c55SBarry Smith     ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr);
487778a2246SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,0,&x);CHKERRQ(ierr);
4888c778c55SBarry Smith     ierr = MatDenseGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
489778a2246SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,0,&b);CHKERRQ(ierr);
4908e722e36SBarry Smith     ierr = KSPSetOperators(exotic->ksp,Aii,Aii,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
491064c8009SBarry Smith     for (i=0; i<6; i++) {
492064c8009SBarry Smith       ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr);
493064c8009SBarry Smith       ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr);
4948e722e36SBarry Smith       ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr);
495064c8009SBarry Smith       ierr = VecResetArray(x);CHKERRQ(ierr);
496064c8009SBarry Smith       ierr = VecResetArray(b);CHKERRQ(ierr);
497064c8009SBarry Smith     }
4988c778c55SBarry Smith     ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr);
4998c778c55SBarry Smith     ierr = MatDenseRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
5006bf464f9SBarry Smith     ierr = VecDestroy(&x);CHKERRQ(ierr);
5016bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
502064c8009SBarry Smith   }
5036bf464f9SBarry Smith   ierr = MatDestroy(&Xint_tmp);CHKERRQ(ierr);
504064c8009SBarry Smith 
505064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
5068c778c55SBarry Smith   ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr);
507064c8009SBarry Smith   for (i=0; i<Nint; i++) {
508064c8009SBarry Smith     tmp = 0.0;
509*2fa5cd67SKarl Rupp     for (j=0; j<6; j++) tmp += xint[i+j*Nint];
510*2fa5cd67SKarl Rupp 
511e32f2f54SBarry 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,PetscAbsScalar(tmp));
512064c8009SBarry Smith   }
5138c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr);
514064c8009SBarry Smith   /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
515064c8009SBarry Smith #endif
516064c8009SBarry Smith 
517064c8009SBarry Smith 
518064c8009SBarry Smith   /*         total faces    */
519064c8009SBarry Smith   Ntotal =  mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1);
520064c8009SBarry Smith 
521064c8009SBarry Smith   /*
522064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
523064c8009SBarry Smith   */
524064c8009SBarry Smith   cnt = 0;
525064c8009SBarry Smith   { gl[cnt++] = mwidth+1;}
526064c8009SBarry Smith   {
527064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth+1;}
528064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
529064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;}
530064c8009SBarry Smith   }
531064c8009SBarry Smith   { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;}
532064c8009SBarry Smith 
533064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
534064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
535064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,6,gl,gl);CHKERRQ(ierr);
536064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
537064c8009SBarry Smith   ierr = PetscMalloc(6*mp*np*pp*sizeof(PetscInt),&globals);CHKERRQ(ierr);
538064c8009SBarry Smith   ierr = MPI_Allgather(gl,6,MPIU_INT,globals,6,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
539064c8009SBarry Smith 
540064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
54128d20b34SBarry Smith   ierr = MatGetSize(Aglobal,&Nt,PETSC_NULL);CHKERRQ(ierr);
54228d20b34SBarry Smith   ierr = PetscTableCreate(Ntotal/3,Nt+1,&ht);CHKERRQ(ierr);
543064c8009SBarry Smith   for (i=0; i<6*mp*np*pp; i++) {
544064c8009SBarry Smith     ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr);
545064c8009SBarry Smith   }
546064c8009SBarry Smith   ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr);
547e32f2f54SBarry 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);
548064c8009SBarry Smith   ierr = PetscFree(globals);CHKERRQ(ierr);
549064c8009SBarry Smith   for (i=0; i<6; i++) {
550064c8009SBarry Smith     ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr);
551064c8009SBarry Smith     gl[i]--;
552064c8009SBarry Smith   }
5536bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&ht);CHKERRQ(ierr);
554064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
555064c8009SBarry Smith 
556064c8009SBarry Smith   /* construct global interpolation matrix */
55728d20b34SBarry Smith   ierr = MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);CHKERRQ(ierr);
558064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
55969b1f4b7SBarry Smith     ierr = MatCreateAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint,PETSC_NULL,P);CHKERRQ(ierr);
560064c8009SBarry Smith   } else {
561064c8009SBarry Smith     ierr = MatZeroEntries(*P);CHKERRQ(ierr);
562064c8009SBarry Smith   }
563064c8009SBarry Smith   ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
5648c778c55SBarry Smith   ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr);
565064c8009SBarry Smith   ierr = MatSetValues(*P,Nint,IIint,6,gl,xint,INSERT_VALUES);CHKERRQ(ierr);
5668c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr);
5678c778c55SBarry Smith   ierr = MatDenseGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
568064c8009SBarry Smith   ierr = MatSetValues(*P,Nsurf,IIsurf,6,gl,xsurf,INSERT_VALUES);CHKERRQ(ierr);
5698c778c55SBarry Smith   ierr = MatDenseRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
570064c8009SBarry Smith   ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
571064c8009SBarry Smith   ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
572064c8009SBarry Smith   ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr);
573064c8009SBarry Smith 
574064c8009SBarry Smith 
575064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
576064c8009SBarry Smith   {
577064c8009SBarry Smith     Vec         x,y;
578064c8009SBarry Smith     PetscScalar *yy;
579064c8009SBarry Smith     ierr = VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr);
580064c8009SBarry Smith     ierr = VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr);
581064c8009SBarry Smith     ierr = VecSet(x,1.0);CHKERRQ(ierr);
582064c8009SBarry Smith     ierr = MatMult(*P,x,y);CHKERRQ(ierr);
583064c8009SBarry Smith     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
584064c8009SBarry Smith     for (i=0; i<Ng; i++) {
585e32f2f54SBarry 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,PetscAbsScalar(yy[i]));
586064c8009SBarry Smith     }
587064c8009SBarry Smith     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
588064c8009SBarry Smith     ierr = VecDestroy(x);CHKERRQ(ierr);
589064c8009SBarry Smith     ierr = VecDestroy(y);CHKERRQ(ierr);
590064c8009SBarry Smith   }
591064c8009SBarry Smith #endif
592064c8009SBarry Smith 
593fcfd50ebSBarry Smith   ierr = MatDestroy(&Aii);CHKERRQ(ierr);
594fcfd50ebSBarry Smith   ierr = MatDestroy(&Ais);CHKERRQ(ierr);
595fcfd50ebSBarry Smith   ierr = MatDestroy(&Asi);CHKERRQ(ierr);
596fcfd50ebSBarry Smith   ierr = MatDestroy(&A);CHKERRQ(ierr);
597fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
598fcfd50ebSBarry Smith   ierr = ISDestroy(&isint);CHKERRQ(ierr);
599fcfd50ebSBarry Smith   ierr = ISDestroy(&issurf);CHKERRQ(ierr);
600fcfd50ebSBarry Smith   ierr = MatDestroy(&Xint);CHKERRQ(ierr);
601fcfd50ebSBarry Smith   ierr = MatDestroy(&Xsurf);CHKERRQ(ierr);
602064c8009SBarry Smith   PetscFunctionReturn(0);
603064c8009SBarry Smith }
604174b6946SBarry Smith 
6057233f9f0SBarry Smith 
606174b6946SBarry Smith #undef __FUNCT__
6077233f9f0SBarry Smith #define __FUNCT__ "PCExoticSetType"
6087233f9f0SBarry Smith /*@
6097233f9f0SBarry Smith    PCExoticSetType - Sets the type of coarse grid interpolation to use
6107233f9f0SBarry Smith 
6113f9fe445SBarry Smith    Logically Collective on PC
6127233f9f0SBarry Smith 
6137233f9f0SBarry Smith    Input Parameters:
6147233f9f0SBarry Smith +  pc - the preconditioner context
6157233f9f0SBarry Smith -  type - either PC_EXOTIC_FACE or PC_EXOTIC_WIREBASKET (defaults to face)
6167233f9f0SBarry Smith 
617563e08c6SBarry Smith    Notes: The face based interpolation has 1 degree of freedom per face and ignores the
618563e08c6SBarry Smith      edge and vertex values completely in the coarse problem. For any seven point
619563e08c6SBarry Smith      stencil the interpolation of a constant on all faces into the interior is that constant.
620563e08c6SBarry Smith 
621563e08c6SBarry Smith      The wirebasket interpolation has 1 degree of freedom per vertex, per edge and
622563e08c6SBarry Smith      per face. A constant on the subdomain boundary is interpolated as that constant
623563e08c6SBarry Smith      in the interior of the domain.
624563e08c6SBarry Smith 
625563e08c6SBarry Smith      The coarse grid matrix is obtained via the Galerkin computation A_c = R A R^T, hence
626563e08c6SBarry Smith      if A is nonsingular A_c is also nonsingular.
627563e08c6SBarry Smith 
628563e08c6SBarry Smith      Both interpolations are suitable for only scalar problems.
629563e08c6SBarry Smith 
6307233f9f0SBarry Smith    Level: intermediate
6317233f9f0SBarry Smith 
6327233f9f0SBarry Smith 
6337233f9f0SBarry Smith .seealso: PCEXOTIC, PCExoticType()
6347233f9f0SBarry Smith @*/
6357087cfbeSBarry Smith PetscErrorCode  PCExoticSetType(PC pc,PCExoticType type)
6367233f9f0SBarry Smith {
6374ac538c5SBarry Smith   PetscErrorCode ierr;
6387233f9f0SBarry Smith 
6397233f9f0SBarry Smith   PetscFunctionBegin;
6400700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
641c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
6424ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCExoticSetType_C",(PC,PCExoticType),(pc,type));CHKERRQ(ierr);
6437233f9f0SBarry Smith   PetscFunctionReturn(0);
6447233f9f0SBarry Smith }
6457233f9f0SBarry Smith 
6467233f9f0SBarry Smith #undef __FUNCT__
6477233f9f0SBarry Smith #define __FUNCT__ "PCExoticSetType_Exotic"
6487087cfbeSBarry Smith PetscErrorCode  PCExoticSetType_Exotic(PC pc,PCExoticType type)
6497233f9f0SBarry Smith {
650f3fbd535SBarry Smith   PC_MG     *mg  = (PC_MG*)pc->data;
65131567311SBarry Smith   PC_Exotic *ctx = (PC_Exotic*) mg->innerctx;
6527233f9f0SBarry Smith 
6537233f9f0SBarry Smith   PetscFunctionBegin;
6547233f9f0SBarry Smith   ctx->type = type;
6557233f9f0SBarry Smith   PetscFunctionReturn(0);
6567233f9f0SBarry Smith }
6577233f9f0SBarry Smith 
6587233f9f0SBarry Smith #undef __FUNCT__
6597233f9f0SBarry Smith #define __FUNCT__ "PCSetUp_Exotic"
6607233f9f0SBarry Smith PetscErrorCode PCSetUp_Exotic(PC pc)
661174b6946SBarry Smith {
662174b6946SBarry Smith   PetscErrorCode ierr;
66396bdf778SBarry Smith   Mat            A;
664f3fbd535SBarry Smith   PC_MG          *mg   = (PC_MG*)pc->data;
66531567311SBarry Smith   PC_Exotic      *ex   = (PC_Exotic*) mg->innerctx;
66696bdf778SBarry Smith   MatReuse       reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX;
667174b6946SBarry Smith 
668174b6946SBarry Smith   PetscFunctionBegin;
669e7e72b3dSBarry Smith   if (!pc->dm) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to call PCSetDM() before using this PC");
670174b6946SBarry Smith   ierr = PCGetOperators(pc,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr);
6717233f9f0SBarry Smith   if (ex->type == PC_EXOTIC_FACE) {
672aa219208SBarry Smith     ierr = DMDAGetFaceInterpolation(pc->dm,ex,A,reuse,&ex->P);CHKERRQ(ierr);
6737233f9f0SBarry Smith   } else if (ex->type == PC_EXOTIC_WIREBASKET) {
674aa219208SBarry Smith     ierr = DMDAGetWireBasketInterpolation(pc->dm,ex,A,reuse,&ex->P);CHKERRQ(ierr);
67565e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_PLIB,"Unknown exotic coarse space %d",ex->type);
67696bdf778SBarry Smith   ierr = PCMGSetInterpolation(pc,1,ex->P);CHKERRQ(ierr);
677d0660788SBarry Smith   /* if PC has attached DM we must remove it or the PCMG will use it to compute incorrect sized vectors and interpolations */
678d0660788SBarry Smith   ierr = PCSetDM(pc,PETSC_NULL);CHKERRQ(ierr);
6797233f9f0SBarry Smith   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
680174b6946SBarry Smith   PetscFunctionReturn(0);
681174b6946SBarry Smith }
682174b6946SBarry Smith 
683174b6946SBarry Smith #undef __FUNCT__
6847233f9f0SBarry Smith #define __FUNCT__ "PCDestroy_Exotic"
6857233f9f0SBarry Smith PetscErrorCode PCDestroy_Exotic(PC pc)
686174b6946SBarry Smith {
687174b6946SBarry Smith   PetscErrorCode ierr;
688f3fbd535SBarry Smith   PC_MG          *mg  = (PC_MG*)pc->data;
68931567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
690174b6946SBarry Smith 
691174b6946SBarry Smith   PetscFunctionBegin;
692fcfd50ebSBarry Smith   ierr = MatDestroy(&ctx->P);CHKERRQ(ierr);
693fcfd50ebSBarry Smith   ierr = KSPDestroy(&ctx->ksp);CHKERRQ(ierr);
6947233f9f0SBarry Smith   ierr = PetscFree(ctx);CHKERRQ(ierr);
6957233f9f0SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
696174b6946SBarry Smith   PetscFunctionReturn(0);
697174b6946SBarry Smith }
698174b6946SBarry Smith 
699174b6946SBarry Smith #undef __FUNCT__
7007233f9f0SBarry Smith #define __FUNCT__ "PCView_Exotic"
7017233f9f0SBarry Smith PetscErrorCode PCView_Exotic(PC pc,PetscViewer viewer)
7027233f9f0SBarry Smith {
703f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
7047233f9f0SBarry Smith   PetscErrorCode ierr;
705ace3abfcSBarry Smith   PetscBool      iascii;
70631567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
7077233f9f0SBarry Smith 
7087233f9f0SBarry Smith   PetscFunctionBegin;
709251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
7107233f9f0SBarry Smith   if (iascii) {
7117233f9f0SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    Exotic type = %s\n",PCExoticTypes[ctx->type]);CHKERRQ(ierr);
7128e722e36SBarry Smith     if (ctx->directSolve) {
7138e722e36SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"      Using direct solver to construct interpolation\n");CHKERRQ(ierr);
7148e722e36SBarry Smith     } else {
7158e722e36SBarry Smith       PetscViewer sviewer;
7168e722e36SBarry Smith       PetscMPIInt rank;
7178e722e36SBarry Smith 
7188e722e36SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"      Using iterative solver to construct interpolation\n");CHKERRQ(ierr);
7198e722e36SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
7208e722e36SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* should not need to push this twice? */
7218e722e36SBarry Smith       ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
7228e722e36SBarry Smith       ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
7238e722e36SBarry Smith       if (!rank) {
7248e722e36SBarry Smith         ierr = KSPView(ctx->ksp,sviewer);CHKERRQ(ierr);
7258e722e36SBarry Smith       }
7268e722e36SBarry Smith       ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
7278e722e36SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
7288e722e36SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
7298e722e36SBarry Smith     }
7307233f9f0SBarry Smith   }
7317233f9f0SBarry Smith   ierr = PCView_MG(pc,viewer);CHKERRQ(ierr);
7327233f9f0SBarry Smith   PetscFunctionReturn(0);
7337233f9f0SBarry Smith }
7347233f9f0SBarry Smith 
7357233f9f0SBarry Smith #undef __FUNCT__
7367233f9f0SBarry Smith #define __FUNCT__ "PCSetFromOptions_Exotic"
7377233f9f0SBarry Smith PetscErrorCode PCSetFromOptions_Exotic(PC pc)
7387233f9f0SBarry Smith {
7397233f9f0SBarry Smith   PetscErrorCode ierr;
740ace3abfcSBarry Smith   PetscBool      flg;
741f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
7427233f9f0SBarry Smith   PCExoticType   mgctype;
74331567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
7447233f9f0SBarry Smith 
7457233f9f0SBarry Smith   PetscFunctionBegin;
7467233f9f0SBarry Smith   ierr = PetscOptionsHead("Exotic coarse space options");CHKERRQ(ierr);
7477233f9f0SBarry Smith   ierr = PetscOptionsEnum("-pc_exotic_type","face or wirebasket","PCExoticSetType",PCExoticTypes,(PetscEnum)ctx->type,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
7487233f9f0SBarry Smith   if (flg) {
7497233f9f0SBarry Smith     ierr = PCExoticSetType(pc,mgctype);CHKERRQ(ierr);
7507233f9f0SBarry Smith   }
751acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_exotic_direct_solver","use direct solver to construct interpolation","None",ctx->directSolve,&ctx->directSolve,PETSC_NULL);CHKERRQ(ierr);
7528e722e36SBarry Smith   if (!ctx->directSolve) {
7538e722e36SBarry Smith     if (!ctx->ksp) {
7548e722e36SBarry Smith       const char *prefix;
7558e722e36SBarry Smith       ierr = KSPCreate(PETSC_COMM_SELF,&ctx->ksp);CHKERRQ(ierr);
7568e722e36SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
7578e722e36SBarry Smith       ierr = PetscLogObjectParent(pc,ctx->ksp);CHKERRQ(ierr);
7588e722e36SBarry Smith       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
7598e722e36SBarry Smith       ierr = KSPSetOptionsPrefix(ctx->ksp,prefix);CHKERRQ(ierr);
7608e722e36SBarry Smith       ierr = KSPAppendOptionsPrefix(ctx->ksp,"exotic_");CHKERRQ(ierr);
7618e722e36SBarry Smith     }
7628e722e36SBarry Smith     ierr = KSPSetFromOptions(ctx->ksp);CHKERRQ(ierr);
7638e722e36SBarry Smith   }
7647233f9f0SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7657233f9f0SBarry Smith   PetscFunctionReturn(0);
7667233f9f0SBarry Smith }
7677233f9f0SBarry Smith 
768174b6946SBarry Smith 
769174b6946SBarry Smith /*MC
7707233f9f0SBarry Smith      PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces
771174b6946SBarry Smith 
7727233f9f0SBarry Smith      This uses the PCMG infrastructure restricted to two levels and the face and wirebasket based coarse
77324c3aa18SBarry Smith    grid spaces.
77424c3aa18SBarry Smith 
77524c3aa18SBarry Smith    Notes: 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
77624c3aa18SBarry Smith 
77724c3aa18SBarry Smith    References: These coarse grid spaces originate in the work of Bramble, Pasciak  and Schatz, "The Construction
7783b65e785SBarry Smith    of Preconditioners for Elliptic Problems by Substructing IV", Mathematics of Computation, volume 53 pages 1--24, 1989.
7793b65e785SBarry Smith    They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith,
7803b65e785SBarry Smith    New York University, 1990. They were then explored in great detail in Dryja, Smith, Widlund, "Schwarz Analysis
7813b65e785SBarry Smith    of Iterative Substructuring Methods for Elliptic Problems in Three Dimensions, SIAM Journal on Numerical
7823b65e785SBarry Smith    Analysis, volume 31. pages 1662-1694, 1994. These were developed in the context of iterative substructuring preconditioners.
7833b65e785SBarry Smith    They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund.
7843b65e785SBarry Smith    They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example,
7853b65e785SBarry Smith    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Extending theory for domain decomposition algorithms to irregular subdomains. In Ulrich Langer, Marco
7863b65e785SBarry Smith    Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
7873b65e785SBarry Smith    of the 17th International Conference on Domain Decomposition Methods in
7883b65e785SBarry Smith    Science and Engineering, held in Strobl, Austria, July 3-7, 2006, number 60 in
7893b65e785SBarry Smith    Springer-Verlag, Lecture Notes in Computational Science and Engineering, pages 255-261, 2007.
7903b65e785SBarry Smith    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. A family of energy min-
7913b65e785SBarry Smith    imizing coarse spaces for overlapping Schwarz preconditioners. In Ulrich Langer,
7923b65e785SBarry Smith    Marco Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
7933b65e785SBarry Smith    of the 17th International Conference on Domain Decomposition Methods
7943b65e785SBarry Smith    in Science and Engineering, held in Strobl, Austria, July 3-7, 2006, number 60 in
7953b65e785SBarry Smith    Springer-Verlag, Lecture Notes in Computational Science and Engineering, pages 247-254, 2007
7963b65e785SBarry Smith    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Domain decomposition
7973b65e785SBarry Smith    for less regular subdomains: Overlapping Schwarz in two dimensions. SIAM J.
7983b65e785SBarry Smith    Numer. Anal., 46(4):2153-2168, 2008.
7993b65e785SBarry Smith    Clark R. Dohrmann and Olof B. Widlund. An overlapping Schwarz
8003b65e785SBarry Smith    algorithm for almost incompressible elasticity. Technical Report
8013b65e785SBarry Smith    TR2008-912, Department of Computer Science, Courant Institute
8023b65e785SBarry Smith    of Mathematical Sciences, New York University, May 2008. URL:
8033b65e785SBarry Smith    http://cs.nyu.edu/csweb/Research/TechReports/TR2008-912/TR2008-912.pdf
8047233f9f0SBarry Smith 
8057233f9f0SBarry Smith    Options Database: The usual PCMG options are supported, such as -mg_levels_pc_type <type> -mg_coarse_pc_type <type>
8067233f9f0SBarry Smith       -pc_mg_type <type>
8077233f9f0SBarry Smith 
80825a35f6fSSatish Balay    Level: advanced
809174b6946SBarry Smith 
8106c699258SBarry Smith .seealso:  PCMG, PCSetDM(), PCExoticType, PCExoticSetType()
811174b6946SBarry Smith M*/
812174b6946SBarry Smith 
813174b6946SBarry Smith EXTERN_C_BEGIN
814174b6946SBarry Smith #undef __FUNCT__
8157233f9f0SBarry Smith #define __FUNCT__ "PCCreate_Exotic"
8167087cfbeSBarry Smith PetscErrorCode  PCCreate_Exotic(PC pc)
817174b6946SBarry Smith {
818174b6946SBarry Smith   PetscErrorCode ierr;
8197233f9f0SBarry Smith   PC_Exotic      *ex;
820f3fbd535SBarry Smith   PC_MG          *mg;
821174b6946SBarry Smith 
822174b6946SBarry Smith   PetscFunctionBegin;
823f91d8e95SBarry Smith   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
824*2fa5cd67SKarl Rupp   if (pc->ops->destroy) {
825*2fa5cd67SKarl Rupp     ierr =  (*pc->ops->destroy)(pc);CHKERRQ(ierr);
826*2fa5cd67SKarl Rupp     pc->data = 0;
827*2fa5cd67SKarl Rupp   }
828503cfb0cSBarry Smith   ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr);
829f91d8e95SBarry Smith   ((PetscObject)pc)->type_name = 0;
830f91d8e95SBarry Smith 
831174b6946SBarry Smith   ierr         = PCSetType(pc,PCMG);CHKERRQ(ierr);
832174b6946SBarry Smith   ierr         = PCMGSetLevels(pc,2,PETSC_NULL);CHKERRQ(ierr);
833c91913d3SJed Brown   ierr         = PCMGSetGalerkin(pc,PETSC_TRUE);CHKERRQ(ierr);
83496bdf778SBarry Smith   ierr         = PetscNew(PC_Exotic,&ex);CHKERRQ(ierr); \
8357233f9f0SBarry Smith   ex->type     = PC_EXOTIC_FACE;
836f3fbd535SBarry Smith   mg           = (PC_MG*) pc->data;
83731567311SBarry Smith   mg->innerctx = ex;
8387233f9f0SBarry Smith 
8397233f9f0SBarry Smith 
8407233f9f0SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_Exotic;
8417233f9f0SBarry Smith   pc->ops->view           = PCView_Exotic;
8427233f9f0SBarry Smith   pc->ops->destroy        = PCDestroy_Exotic;
8436c699258SBarry Smith   pc->ops->setup          = PCSetUp_Exotic;
844*2fa5cd67SKarl Rupp 
8457233f9f0SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCExoticSetType_C","PCExoticSetType_Exotic",PCExoticSetType_Exotic);CHKERRQ(ierr);
846174b6946SBarry Smith   PetscFunctionReturn(0);
847174b6946SBarry Smith }
848174b6946SBarry Smith EXTERN_C_END
849