xref: /petsc/src/ksp/pc/impls/wb/wb.c (revision 12dd4568308e9b59cb5ccae9e7b91d823363eff4)
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 
138e722e36SBarry Smith const char *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);
69064c8009SBarry Smith   ierr = MatGetArray(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;
80064c8009SBarry Smith   xsurf[cnt++] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + Nsurf] = 1;} xsurf[cnt++ + 2*Nsurf] = 1;
81064c8009SBarry Smith   for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 3*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 4*Nsurf] = 1;} xsurf[cnt++ + 5*Nsurf] = 1;}
82064c8009SBarry Smith   xsurf[cnt++ + 6*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 7*Nsurf] = 1;} xsurf[cnt++ + 8*Nsurf] = 1;
83064c8009SBarry Smith   for (k=1;k<p-1-kstart;k++) {
84064c8009SBarry Smith     xsurf[cnt++ + 9*Nsurf] = 1;  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 10*Nsurf] = 1;}  xsurf[cnt++ + 11*Nsurf] = 1;
85064c8009SBarry Smith     for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 12*Nsurf] = 1; /* these are the interior nodes */ xsurf[cnt++ + 13*Nsurf] = 1;}
86064c8009SBarry Smith     xsurf[cnt++ + 14*Nsurf] = 1;  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 15*Nsurf] = 1;} xsurf[cnt++ + 16*Nsurf] = 1;
87064c8009SBarry Smith   }
88064c8009SBarry Smith   xsurf[cnt++ + 17*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 18*Nsurf] = 1;} xsurf[cnt++ + 19*Nsurf] = 1;
89064c8009SBarry Smith   for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 20*Nsurf] = 1;  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 21*Nsurf] = 1;} xsurf[cnt++ + 22*Nsurf] = 1;}
90064c8009SBarry Smith   xsurf[cnt++ + 23*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 24*Nsurf] = 1;} xsurf[cnt++ + 25*Nsurf] = 1;
91064c8009SBarry Smith 
928e722e36SBarry Smith   /* interpolations only sum to 1 when using direct solver */
938e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
94064c8009SBarry Smith   for (i=0; i<Nsurf; i++) {
95064c8009SBarry Smith     tmp = 0.0;
96064c8009SBarry Smith     for (j=0; j<26; j++) {
97064c8009SBarry Smith       tmp += xsurf[i+j*Nsurf];
98064c8009SBarry Smith     }
99e32f2f54SBarry 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));
100064c8009SBarry Smith   }
101064c8009SBarry Smith #endif
102064c8009SBarry Smith   ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
103064c8009SBarry Smith   /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/
104064c8009SBarry Smith 
105064c8009SBarry Smith 
106064c8009SBarry Smith   /*
107064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
108064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
109064c8009SBarry Smith             (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it
110aa219208SBarry Smith              is NOT the local DMDA ordering.)
111064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
112064c8009SBarry Smith   */
113064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
114064c8009SBarry Smith   ierr = PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);CHKERRQ(ierr);
115064c8009SBarry Smith   ierr = PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);CHKERRQ(ierr);
116064c8009SBarry Smith   for (k=0; k<p-kstart; k++) {
117064c8009SBarry Smith     for (j=0; j<n-jstart; j++) {
118064c8009SBarry Smith       for (i=0; i<m-istart; i++) {
119064c8009SBarry Smith         II[c++] = i + j*mwidth + k*mwidth*nwidth;
120064c8009SBarry Smith 
121064c8009SBarry Smith         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
122064c8009SBarry Smith           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
123064c8009SBarry Smith           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
124064c8009SBarry Smith         } else {
125064c8009SBarry Smith           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
126064c8009SBarry Smith           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
127064c8009SBarry Smith         }
128064c8009SBarry Smith       }
129064c8009SBarry Smith     }
130064c8009SBarry Smith   }
131e32f2f54SBarry Smith   if (c != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N");
132e32f2f54SBarry Smith   if (cint != Nint) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint");
133e32f2f54SBarry Smith   if (csurf != Nsurf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf");
1341411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltg);CHKERRQ(ierr);
135064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr);
136064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr);
137064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr);
138064c8009SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
13970b3c8c7SBarry Smith   ierr = ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
14070b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint);CHKERRQ(ierr);
14170b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf);CHKERRQ(ierr);
142064c8009SBarry Smith   ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr);
143064c8009SBarry Smith 
144064c8009SBarry Smith   ierr = MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr);
145064c8009SBarry Smith   A    = *Aholder;
146064c8009SBarry Smith   ierr = PetscFree(Aholder);CHKERRQ(ierr);
147064c8009SBarry Smith 
148064c8009SBarry Smith   ierr = MatGetSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr);
149064c8009SBarry Smith   ierr = MatGetSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr);
150064c8009SBarry Smith   ierr = MatGetSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr);
151064c8009SBarry Smith 
152064c8009SBarry Smith   /*
153064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
154064c8009SBarry Smith   */
1558e722e36SBarry Smith   ierr = MatDuplicate(Xint,MAT_DO_NOT_COPY_VALUES,&Xint_tmp);CHKERRQ(ierr);
1568e722e36SBarry Smith   ierr = MatMatMult(Ais,Xsurf,MAT_REUSE_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr);
1578e722e36SBarry Smith   ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr);
1588e722e36SBarry Smith   if (exotic->directSolve) {
1592692d6eeSBarry Smith     ierr = MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr);
160064c8009SBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
1612692d6eeSBarry Smith     ierr = MatGetOrdering(Aii,MATORDERINGND,&row,&col);CHKERRQ(ierr);
162064c8009SBarry Smith     ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr);
163fcfd50ebSBarry Smith     ierr = ISDestroy(&row);CHKERRQ(ierr);
164fcfd50ebSBarry Smith     ierr = ISDestroy(&col);CHKERRQ(ierr);
165064c8009SBarry Smith     ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr);
166064c8009SBarry Smith     ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr);
1676bf464f9SBarry Smith     ierr = MatDestroy(&iAii);CHKERRQ(ierr);
1688e722e36SBarry Smith   } else {
1698e722e36SBarry Smith     Vec         b,x;
1708e722e36SBarry Smith     PetscScalar *xint_tmp;
171064c8009SBarry Smith 
1728e722e36SBarry Smith     ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr);
1738e722e36SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&x);CHKERRQ(ierr);
1748e722e36SBarry Smith     ierr = MatGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
1758e722e36SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&b);CHKERRQ(ierr);
1768e722e36SBarry Smith     ierr = KSPSetOperators(exotic->ksp,Aii,Aii,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1778e722e36SBarry Smith     for (i=0; i<26; i++) {
1788e722e36SBarry Smith       ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr);
1798e722e36SBarry Smith       ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr);
1808e722e36SBarry Smith       ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr);
1818e722e36SBarry Smith       ierr = VecResetArray(x);CHKERRQ(ierr);
1828e722e36SBarry Smith       ierr = VecResetArray(b);CHKERRQ(ierr);
1838e722e36SBarry Smith     }
1848e722e36SBarry Smith     ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr);
1858e722e36SBarry Smith     ierr = MatRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
1866bf464f9SBarry Smith     ierr = VecDestroy(&x);CHKERRQ(ierr);
1876bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
1888e722e36SBarry Smith   }
1896bf464f9SBarry Smith   ierr = MatDestroy(&Xint_tmp);CHKERRQ(ierr);
1908e722e36SBarry Smith 
1918e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
192064c8009SBarry Smith   ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr);
193064c8009SBarry Smith   for (i=0; i<Nint; i++) {
194064c8009SBarry Smith     tmp = 0.0;
195064c8009SBarry Smith     for (j=0; j<26; j++) {
196064c8009SBarry Smith       tmp += xint[i+j*Nint];
197064c8009SBarry Smith     }
198e32f2f54SBarry 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));
199064c8009SBarry Smith   }
200064c8009SBarry Smith   ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr);
201064c8009SBarry Smith   /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
202064c8009SBarry Smith #endif
203064c8009SBarry Smith 
204064c8009SBarry Smith 
205064c8009SBarry Smith   /*         total vertices             total faces                                  total edges */
206064c8009SBarry 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);
207064c8009SBarry Smith 
208064c8009SBarry Smith   /*
209064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
210064c8009SBarry Smith   */
211064c8009SBarry Smith   cnt = 0;
212064c8009SBarry Smith   gl[cnt++] = 0;  { gl[cnt++] = 1;} gl[cnt++] = m-istart-1;
213064c8009SBarry Smith   { gl[cnt++] = mwidth;  { gl[cnt++] = mwidth+1;} gl[cnt++] = mwidth + m-istart-1;}
214064c8009SBarry Smith   gl[cnt++] = mwidth*(n-jstart-1);  { gl[cnt++] = mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*(n-jstart-1) + m-istart-1;
215064c8009SBarry Smith   {
216064c8009SBarry Smith     gl[cnt++] = mwidth*nwidth;  { gl[cnt++] = mwidth*nwidth+1;}  gl[cnt++] = mwidth*nwidth+ m-istart-1;
217064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
218064c8009SBarry 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;
219064c8009SBarry Smith   }
220064c8009SBarry 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;
221064c8009SBarry 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;}
222064c8009SBarry 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;
223064c8009SBarry Smith 
224064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
225064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
226064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,26,gl,gl);CHKERRQ(ierr);
227064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
228064c8009SBarry Smith   ierr = PetscMalloc(26*mp*np*pp*sizeof(PetscInt),&globals);CHKERRQ(ierr);
229064c8009SBarry Smith   ierr = MPI_Allgather(gl,26,MPIU_INT,globals,26,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
230064c8009SBarry Smith 
231064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
23228d20b34SBarry Smith   ierr = MatGetSize(Aglobal,&Nt,PETSC_NULL);CHKERRQ(ierr);
233*12dd4568SBarry Smith   ierr = PetscTableCreate(Ntotal/3,Nt+1,&ht);CHKERRQ(ierr);
234064c8009SBarry Smith   for (i=0; i<26*mp*np*pp; i++){
235064c8009SBarry Smith     ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr);
236064c8009SBarry Smith   }
237064c8009SBarry Smith   ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr);
238e32f2f54SBarry 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);
239064c8009SBarry Smith   ierr = PetscFree(globals);CHKERRQ(ierr);
240064c8009SBarry Smith   for (i=0; i<26; i++) {
241064c8009SBarry Smith     ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr);
242064c8009SBarry Smith     gl[i]--;
243064c8009SBarry Smith   }
2446bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&ht);CHKERRQ(ierr);
245064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
246064c8009SBarry Smith 
247064c8009SBarry Smith   /* construct global interpolation matrix */
24828d20b34SBarry Smith   ierr = MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);CHKERRQ(ierr);
249064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
250064c8009SBarry Smith     ierr = MatCreateMPIAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint+Nsurf,PETSC_NULL,P);CHKERRQ(ierr);
251064c8009SBarry Smith   } else {
252064c8009SBarry Smith     ierr = MatZeroEntries(*P);CHKERRQ(ierr);
253064c8009SBarry Smith   }
254064c8009SBarry Smith   ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
255064c8009SBarry Smith   ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr);
256064c8009SBarry Smith   ierr = MatSetValues(*P,Nint,IIint,26,gl,xint,INSERT_VALUES);CHKERRQ(ierr);
257064c8009SBarry Smith   ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr);
258064c8009SBarry Smith   ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
259064c8009SBarry Smith   ierr = MatSetValues(*P,Nsurf,IIsurf,26,gl,xsurf,INSERT_VALUES);CHKERRQ(ierr);
260064c8009SBarry Smith   ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
261064c8009SBarry Smith   ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
262064c8009SBarry Smith   ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
263064c8009SBarry Smith   ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr);
264064c8009SBarry Smith 
2658e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
266064c8009SBarry Smith   {
267064c8009SBarry Smith     Vec         x,y;
268064c8009SBarry Smith     PetscScalar *yy;
269064c8009SBarry Smith     ierr = VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr);
270064c8009SBarry Smith     ierr = VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr);
271064c8009SBarry Smith     ierr = VecSet(x,1.0);CHKERRQ(ierr);
272064c8009SBarry Smith     ierr = MatMult(*P,x,y);CHKERRQ(ierr);
273064c8009SBarry Smith     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
274064c8009SBarry Smith     for (i=0; i<Ng; i++) {
275e32f2f54SBarry 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]));
276064c8009SBarry Smith     }
277064c8009SBarry Smith     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
278064c8009SBarry Smith     ierr = VecDestroy(x);CHKERRQ(ierr);
279064c8009SBarry Smith     ierr = VecDestroy(y);CHKERRQ(ierr);
280064c8009SBarry Smith   }
281064c8009SBarry Smith #endif
282064c8009SBarry Smith 
283fcfd50ebSBarry Smith   ierr = MatDestroy(&Aii);CHKERRQ(ierr);
284fcfd50ebSBarry Smith   ierr = MatDestroy(&Ais);CHKERRQ(ierr);
285fcfd50ebSBarry Smith   ierr = MatDestroy(&Asi);CHKERRQ(ierr);
286fcfd50ebSBarry Smith   ierr = MatDestroy(&A);CHKERRQ(ierr);
287fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
288fcfd50ebSBarry Smith   ierr = ISDestroy(&isint);CHKERRQ(ierr);
289fcfd50ebSBarry Smith   ierr = ISDestroy(&issurf);CHKERRQ(ierr);
290fcfd50ebSBarry Smith   ierr = MatDestroy(&Xint);CHKERRQ(ierr);
291fcfd50ebSBarry Smith   ierr = MatDestroy(&Xsurf);CHKERRQ(ierr);
292064c8009SBarry Smith   PetscFunctionReturn(0);
293064c8009SBarry Smith }
294064c8009SBarry Smith 
295064c8009SBarry Smith #undef __FUNCT__
296aa219208SBarry Smith #define __FUNCT__ "DMDAGetFaceInterpolation"
297064c8009SBarry Smith /*
298aa219208SBarry Smith       DMDAGetFaceInterpolation - Gets the interpolation for a face based coarse space
299064c8009SBarry Smith 
300064c8009SBarry Smith */
301aa219208SBarry Smith PetscErrorCode DMDAGetFaceInterpolation(DM da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P)
302064c8009SBarry Smith {
303064c8009SBarry Smith   PetscErrorCode         ierr;
304064c8009SBarry 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;
30528d20b34SBarry Smith   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[6],*globals,Ng,*IIint,*IIsurf,Nt;
306064c8009SBarry Smith   Mat                    Xint, Xsurf,Xint_tmp;
307064c8009SBarry Smith   IS                     isint,issurf,is,row,col;
308064c8009SBarry Smith   ISLocalToGlobalMapping ltg;
309064c8009SBarry Smith   MPI_Comm               comm;
310064c8009SBarry Smith   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
311064c8009SBarry Smith   MatFactorInfo          info;
312064c8009SBarry Smith   PetscScalar            *xsurf,*xint;
313064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
314064c8009SBarry Smith   PetscScalar            tmp;
315064c8009SBarry Smith #endif
316064c8009SBarry Smith   PetscTable             ht;
317064c8009SBarry Smith 
318064c8009SBarry Smith   PetscFunctionBegin;
3191321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0,0,0);CHKERRQ(ierr);
320e7e72b3dSBarry Smith   if (dof != 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only for single field problems");
321e7e72b3dSBarry Smith   if (dim != 3) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only coded for 3d problems");
322aa219208SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&m,&n,&p);CHKERRQ(ierr);
323aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr);
324064c8009SBarry Smith   istart = istart ? -1 : 0;
325064c8009SBarry Smith   jstart = jstart ? -1 : 0;
326064c8009SBarry Smith   kstart = kstart ? -1 : 0;
327064c8009SBarry Smith 
328064c8009SBarry Smith   /*
329064c8009SBarry Smith     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge)
330064c8009SBarry Smith     to all the local degrees of freedom (this includes the vertices, edges and faces).
331064c8009SBarry Smith 
332064c8009SBarry Smith     Xint are the subset of the interpolation into the interior
333064c8009SBarry Smith 
334064c8009SBarry Smith     Xface are the interpolation onto faces but not into the interior
335064c8009SBarry Smith 
336064c8009SBarry Smith     Xsurf are the interpolation onto the vertices and edges (the surfbasket)
337064c8009SBarry Smith                                         Xint
338064c8009SBarry Smith     Symbolically one could write P = (  Xface  ) after interchanging the rows to match the natural ordering on the domain
339064c8009SBarry Smith                                         Xsurf
340064c8009SBarry Smith   */
341064c8009SBarry Smith   N     = (m - istart)*(n - jstart)*(p - kstart);
342064c8009SBarry Smith   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
343064c8009SBarry Smith   Nface = 2*( (m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart) );
344064c8009SBarry Smith   Nwire = 4*( (m-2-istart) + (n-2-jstart) + (p-2-kstart) ) + 8;
345064c8009SBarry Smith   Nsurf = Nface + Nwire;
346064c8009SBarry Smith   ierr = MatCreateSeqDense(MPI_COMM_SELF,Nint,6,PETSC_NULL,&Xint);CHKERRQ(ierr);
347064c8009SBarry Smith   ierr = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,6,PETSC_NULL,&Xsurf);CHKERRQ(ierr);
348064c8009SBarry Smith   ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
349064c8009SBarry Smith 
350064c8009SBarry Smith   /*
351064c8009SBarry Smith      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
352064c8009SBarry Smith      Xsurf will be all zero (thus making the coarse matrix singular).
353064c8009SBarry Smith   */
354e32f2f54SBarry 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");
355e32f2f54SBarry 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");
356e32f2f54SBarry 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");
357064c8009SBarry Smith 
358064c8009SBarry Smith   cnt = 0;
359064c8009SBarry Smith   for (j=1;j<n-1-jstart;j++) {  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 0*Nsurf] = 1;} }
360064c8009SBarry Smith    for (k=1;k<p-1-kstart;k++) {
361064c8009SBarry Smith     for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 1*Nsurf] = 1;}
362064c8009SBarry Smith     for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 2*Nsurf] = 1; /* these are the interior nodes */ xsurf[cnt++ + 3*Nsurf] = 1;}
363064c8009SBarry Smith     for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 4*Nsurf] = 1;}
364064c8009SBarry Smith   }
365064c8009SBarry Smith   for (j=1;j<n-1-jstart;j++) {for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 5*Nsurf] = 1;} }
366064c8009SBarry Smith 
367064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
368064c8009SBarry Smith   for (i=0; i<Nsurf; i++) {
369064c8009SBarry Smith     tmp = 0.0;
370064c8009SBarry Smith     for (j=0; j<6; j++) {
371064c8009SBarry Smith       tmp += xsurf[i+j*Nsurf];
372064c8009SBarry Smith     }
373e32f2f54SBarry 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));
374064c8009SBarry Smith   }
375064c8009SBarry Smith #endif
376064c8009SBarry Smith   ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
377064c8009SBarry Smith   /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/
378064c8009SBarry Smith 
379064c8009SBarry Smith 
380064c8009SBarry Smith   /*
381064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
382064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
383064c8009SBarry Smith             (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it
384aa219208SBarry Smith              is NOT the local DMDA ordering.)
385064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
386064c8009SBarry Smith   */
387064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
388064c8009SBarry Smith   ierr = PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);CHKERRQ(ierr);
389064c8009SBarry Smith   ierr = PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);CHKERRQ(ierr);
390064c8009SBarry Smith   for (k=0; k<p-kstart; k++) {
391064c8009SBarry Smith     for (j=0; j<n-jstart; j++) {
392064c8009SBarry Smith       for (i=0; i<m-istart; i++) {
393064c8009SBarry Smith         II[c++] = i + j*mwidth + k*mwidth*nwidth;
394064c8009SBarry Smith 
395064c8009SBarry Smith         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
396064c8009SBarry Smith           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
397064c8009SBarry Smith           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
398064c8009SBarry Smith         } else {
399064c8009SBarry Smith           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
400064c8009SBarry Smith           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
401064c8009SBarry Smith         }
402064c8009SBarry Smith       }
403064c8009SBarry Smith     }
404064c8009SBarry Smith   }
405e32f2f54SBarry Smith   if (c != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N");
406e32f2f54SBarry Smith   if (cint != Nint) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint");
407e32f2f54SBarry Smith   if (csurf != Nsurf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf");
4081411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltg);CHKERRQ(ierr);
409064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr);
410064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr);
411064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr);
412064c8009SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
41370b3c8c7SBarry Smith   ierr = ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
41470b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint);CHKERRQ(ierr);
41570b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf);CHKERRQ(ierr);
416064c8009SBarry Smith   ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr);
417064c8009SBarry Smith 
418064c8009SBarry Smith   ierr = MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr);
419064c8009SBarry Smith   A    = *Aholder;
420064c8009SBarry Smith   ierr = PetscFree(Aholder);CHKERRQ(ierr);
421064c8009SBarry Smith 
422064c8009SBarry Smith   ierr = MatGetSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr);
423064c8009SBarry Smith   ierr = MatGetSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr);
424064c8009SBarry Smith   ierr = MatGetSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr);
425064c8009SBarry Smith 
426064c8009SBarry Smith   /*
427064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
428064c8009SBarry Smith   */
429064c8009SBarry Smith   ierr = MatDuplicate(Xint,MAT_DO_NOT_COPY_VALUES,&Xint_tmp);CHKERRQ(ierr);
430064c8009SBarry Smith   ierr = MatMatMult(Ais,Xsurf,MAT_REUSE_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr);
431064c8009SBarry Smith   ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr);
432064c8009SBarry Smith 
4338e722e36SBarry Smith   if (exotic->directSolve) {
4342692d6eeSBarry Smith     ierr = MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr);
435064c8009SBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
4362692d6eeSBarry Smith     ierr = MatGetOrdering(Aii,MATORDERINGND,&row,&col);CHKERRQ(ierr);
437064c8009SBarry Smith     ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr);
438fcfd50ebSBarry Smith     ierr = ISDestroy(&row);CHKERRQ(ierr);
439fcfd50ebSBarry Smith     ierr = ISDestroy(&col);CHKERRQ(ierr);
440064c8009SBarry Smith     ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr);
441064c8009SBarry Smith     ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr);
442fcfd50ebSBarry Smith     ierr = MatDestroy(&iAii);CHKERRQ(ierr);
443064c8009SBarry Smith   } else {
444064c8009SBarry Smith     Vec         b,x;
445064c8009SBarry Smith     PetscScalar *xint_tmp;
446064c8009SBarry Smith 
447064c8009SBarry Smith     ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr);
448064c8009SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&x);CHKERRQ(ierr);
449064c8009SBarry Smith     ierr = MatGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
450064c8009SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&b);CHKERRQ(ierr);
4518e722e36SBarry Smith     ierr = KSPSetOperators(exotic->ksp,Aii,Aii,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
452064c8009SBarry Smith     for (i=0; i<6; i++) {
453064c8009SBarry Smith       ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr);
454064c8009SBarry Smith       ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr);
4558e722e36SBarry Smith       ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr);
456064c8009SBarry Smith       ierr = VecResetArray(x);CHKERRQ(ierr);
457064c8009SBarry Smith       ierr = VecResetArray(b);CHKERRQ(ierr);
458064c8009SBarry Smith     }
459064c8009SBarry Smith     ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr);
460064c8009SBarry Smith     ierr = MatRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
4616bf464f9SBarry Smith     ierr = VecDestroy(&x);CHKERRQ(ierr);
4626bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
463064c8009SBarry Smith   }
4646bf464f9SBarry Smith   ierr = MatDestroy(&Xint_tmp);CHKERRQ(ierr);
465064c8009SBarry Smith 
466064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
467064c8009SBarry Smith   ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr);
468064c8009SBarry Smith   for (i=0; i<Nint; i++) {
469064c8009SBarry Smith     tmp = 0.0;
470064c8009SBarry Smith     for (j=0; j<6; j++) {
471064c8009SBarry Smith       tmp += xint[i+j*Nint];
472064c8009SBarry Smith     }
473e32f2f54SBarry 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));
474064c8009SBarry Smith   }
475064c8009SBarry Smith   ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr);
476064c8009SBarry Smith   /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
477064c8009SBarry Smith #endif
478064c8009SBarry Smith 
479064c8009SBarry Smith 
480064c8009SBarry Smith   /*         total faces    */
481064c8009SBarry Smith   Ntotal =  mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1);
482064c8009SBarry Smith 
483064c8009SBarry Smith   /*
484064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
485064c8009SBarry Smith   */
486064c8009SBarry Smith   cnt = 0;
487064c8009SBarry Smith   { gl[cnt++] = mwidth+1;}
488064c8009SBarry Smith   {
489064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth+1;}
490064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
491064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;}
492064c8009SBarry Smith   }
493064c8009SBarry Smith   { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;}
494064c8009SBarry Smith 
495064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
496064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
497064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,6,gl,gl);CHKERRQ(ierr);
498064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
499064c8009SBarry Smith   ierr = PetscMalloc(6*mp*np*pp*sizeof(PetscInt),&globals);CHKERRQ(ierr);
500064c8009SBarry Smith   ierr = MPI_Allgather(gl,6,MPIU_INT,globals,6,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
501064c8009SBarry Smith 
502064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
50328d20b34SBarry Smith   ierr = MatGetSize(Aglobal,&Nt,PETSC_NULL);CHKERRQ(ierr);
50428d20b34SBarry Smith   ierr = PetscTableCreate(Ntotal/3,Nt+1,&ht);CHKERRQ(ierr);
505064c8009SBarry Smith   for (i=0; i<6*mp*np*pp; i++){
506064c8009SBarry Smith     ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr);
507064c8009SBarry Smith   }
508064c8009SBarry Smith   ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr);
509e32f2f54SBarry 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);
510064c8009SBarry Smith   ierr = PetscFree(globals);CHKERRQ(ierr);
511064c8009SBarry Smith   for (i=0; i<6; i++) {
512064c8009SBarry Smith     ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr);
513064c8009SBarry Smith     gl[i]--;
514064c8009SBarry Smith   }
5156bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&ht);CHKERRQ(ierr);
516064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
517064c8009SBarry Smith 
518064c8009SBarry Smith   /* construct global interpolation matrix */
51928d20b34SBarry Smith   ierr = MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);CHKERRQ(ierr);
520064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
521064c8009SBarry Smith     ierr = MatCreateMPIAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint,PETSC_NULL,P);CHKERRQ(ierr);
522064c8009SBarry Smith   } else {
523064c8009SBarry Smith     ierr = MatZeroEntries(*P);CHKERRQ(ierr);
524064c8009SBarry Smith   }
525064c8009SBarry Smith   ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
526064c8009SBarry Smith   ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr);
527064c8009SBarry Smith   ierr = MatSetValues(*P,Nint,IIint,6,gl,xint,INSERT_VALUES);CHKERRQ(ierr);
528064c8009SBarry Smith   ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr);
529064c8009SBarry Smith   ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
530064c8009SBarry Smith   ierr = MatSetValues(*P,Nsurf,IIsurf,6,gl,xsurf,INSERT_VALUES);CHKERRQ(ierr);
531064c8009SBarry Smith   ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
532064c8009SBarry Smith   ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
533064c8009SBarry Smith   ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
534064c8009SBarry Smith   ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr);
535064c8009SBarry Smith 
536064c8009SBarry Smith 
537064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
538064c8009SBarry Smith   {
539064c8009SBarry Smith     Vec         x,y;
540064c8009SBarry Smith     PetscScalar *yy;
541064c8009SBarry Smith     ierr = VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr);
542064c8009SBarry Smith     ierr = VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr);
543064c8009SBarry Smith     ierr = VecSet(x,1.0);CHKERRQ(ierr);
544064c8009SBarry Smith     ierr = MatMult(*P,x,y);CHKERRQ(ierr);
545064c8009SBarry Smith     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
546064c8009SBarry Smith     for (i=0; i<Ng; i++) {
547e32f2f54SBarry 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]));
548064c8009SBarry Smith     }
549064c8009SBarry Smith     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
550064c8009SBarry Smith     ierr = VecDestroy(x);CHKERRQ(ierr);
551064c8009SBarry Smith     ierr = VecDestroy(y);CHKERRQ(ierr);
552064c8009SBarry Smith   }
553064c8009SBarry Smith #endif
554064c8009SBarry Smith 
555fcfd50ebSBarry Smith   ierr = MatDestroy(&Aii);CHKERRQ(ierr);
556fcfd50ebSBarry Smith   ierr = MatDestroy(&Ais);CHKERRQ(ierr);
557fcfd50ebSBarry Smith   ierr = MatDestroy(&Asi);CHKERRQ(ierr);
558fcfd50ebSBarry Smith   ierr = MatDestroy(&A);CHKERRQ(ierr);
559fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
560fcfd50ebSBarry Smith   ierr = ISDestroy(&isint);CHKERRQ(ierr);
561fcfd50ebSBarry Smith   ierr = ISDestroy(&issurf);CHKERRQ(ierr);
562fcfd50ebSBarry Smith   ierr = MatDestroy(&Xint);CHKERRQ(ierr);
563fcfd50ebSBarry Smith   ierr = MatDestroy(&Xsurf);CHKERRQ(ierr);
564064c8009SBarry Smith   PetscFunctionReturn(0);
565064c8009SBarry Smith }
566174b6946SBarry Smith 
5677233f9f0SBarry Smith 
568174b6946SBarry Smith #undef __FUNCT__
5697233f9f0SBarry Smith #define __FUNCT__ "PCExoticSetType"
5707233f9f0SBarry Smith /*@
5717233f9f0SBarry Smith    PCExoticSetType - Sets the type of coarse grid interpolation to use
5727233f9f0SBarry Smith 
5733f9fe445SBarry Smith    Logically Collective on PC
5747233f9f0SBarry Smith 
5757233f9f0SBarry Smith    Input Parameters:
5767233f9f0SBarry Smith +  pc - the preconditioner context
5777233f9f0SBarry Smith -  type - either PC_EXOTIC_FACE or PC_EXOTIC_WIREBASKET (defaults to face)
5787233f9f0SBarry Smith 
579563e08c6SBarry Smith    Notes: The face based interpolation has 1 degree of freedom per face and ignores the
580563e08c6SBarry Smith      edge and vertex values completely in the coarse problem. For any seven point
581563e08c6SBarry Smith      stencil the interpolation of a constant on all faces into the interior is that constant.
582563e08c6SBarry Smith 
583563e08c6SBarry Smith      The wirebasket interpolation has 1 degree of freedom per vertex, per edge and
584563e08c6SBarry Smith      per face. A constant on the subdomain boundary is interpolated as that constant
585563e08c6SBarry Smith      in the interior of the domain.
586563e08c6SBarry Smith 
587563e08c6SBarry Smith      The coarse grid matrix is obtained via the Galerkin computation A_c = R A R^T, hence
588563e08c6SBarry Smith      if A is nonsingular A_c is also nonsingular.
589563e08c6SBarry Smith 
590563e08c6SBarry Smith      Both interpolations are suitable for only scalar problems.
591563e08c6SBarry Smith 
5927233f9f0SBarry Smith    Level: intermediate
5937233f9f0SBarry Smith 
5947233f9f0SBarry Smith 
5957233f9f0SBarry Smith .seealso: PCEXOTIC, PCExoticType()
5967233f9f0SBarry Smith @*/
5977087cfbeSBarry Smith PetscErrorCode  PCExoticSetType(PC pc,PCExoticType type)
5987233f9f0SBarry Smith {
5994ac538c5SBarry Smith   PetscErrorCode ierr;
6007233f9f0SBarry Smith 
6017233f9f0SBarry Smith   PetscFunctionBegin;
6020700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
603c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
6044ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCExoticSetType_C",(PC,PCExoticType),(pc,type));CHKERRQ(ierr);
6057233f9f0SBarry Smith   PetscFunctionReturn(0);
6067233f9f0SBarry Smith }
6077233f9f0SBarry Smith 
6087233f9f0SBarry Smith #undef __FUNCT__
6097233f9f0SBarry Smith #define __FUNCT__ "PCExoticSetType_Exotic"
6107087cfbeSBarry Smith PetscErrorCode  PCExoticSetType_Exotic(PC pc,PCExoticType type)
6117233f9f0SBarry Smith {
612f3fbd535SBarry Smith   PC_MG     *mg = (PC_MG*)pc->data;
61331567311SBarry Smith   PC_Exotic *ctx = (PC_Exotic*) mg->innerctx;
6147233f9f0SBarry Smith 
6157233f9f0SBarry Smith   PetscFunctionBegin;
6167233f9f0SBarry Smith   ctx->type = type;
6177233f9f0SBarry Smith   PetscFunctionReturn(0);
6187233f9f0SBarry Smith }
6197233f9f0SBarry Smith 
6207233f9f0SBarry Smith #undef __FUNCT__
6217233f9f0SBarry Smith #define __FUNCT__ "PCSetUp_Exotic"
6227233f9f0SBarry Smith PetscErrorCode PCSetUp_Exotic(PC pc)
623174b6946SBarry Smith {
624174b6946SBarry Smith   PetscErrorCode ierr;
62596bdf778SBarry Smith   Mat            A;
626f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
62731567311SBarry Smith   PC_Exotic      *ex = (PC_Exotic*) mg->innerctx;
62896bdf778SBarry Smith   MatReuse       reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX;
629174b6946SBarry Smith 
630174b6946SBarry Smith   PetscFunctionBegin;
631e7e72b3dSBarry Smith   if (!pc->dm) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to call PCSetDM() before using this PC");
632174b6946SBarry Smith   ierr = PCGetOperators(pc,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr);
6337233f9f0SBarry Smith   if (ex->type == PC_EXOTIC_FACE) {
634aa219208SBarry Smith     ierr = DMDAGetFaceInterpolation(pc->dm,ex,A,reuse,&ex->P);CHKERRQ(ierr);
6357233f9f0SBarry Smith   } else if (ex->type == PC_EXOTIC_WIREBASKET) {
636aa219208SBarry Smith     ierr = DMDAGetWireBasketInterpolation(pc->dm,ex,A,reuse,&ex->P);CHKERRQ(ierr);
63765e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_PLIB,"Unknown exotic coarse space %d",ex->type);
63896bdf778SBarry Smith   ierr = PCMGSetInterpolation(pc,1,ex->P);CHKERRQ(ierr);
639d0660788SBarry Smith   /* if PC has attached DM we must remove it or the PCMG will use it to compute incorrect sized vectors and interpolations */
640d0660788SBarry Smith   ierr = PCSetDM(pc,PETSC_NULL);CHKERRQ(ierr);
6417233f9f0SBarry Smith   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
642174b6946SBarry Smith   PetscFunctionReturn(0);
643174b6946SBarry Smith }
644174b6946SBarry Smith 
645174b6946SBarry Smith #undef __FUNCT__
6467233f9f0SBarry Smith #define __FUNCT__ "PCDestroy_Exotic"
6477233f9f0SBarry Smith PetscErrorCode PCDestroy_Exotic(PC pc)
648174b6946SBarry Smith {
649174b6946SBarry Smith   PetscErrorCode ierr;
650f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
65131567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
652174b6946SBarry Smith 
653174b6946SBarry Smith   PetscFunctionBegin;
654fcfd50ebSBarry Smith   ierr = MatDestroy(&ctx->P);CHKERRQ(ierr);
655fcfd50ebSBarry Smith   ierr = KSPDestroy(&ctx->ksp);CHKERRQ(ierr);
6567233f9f0SBarry Smith   ierr = PetscFree(ctx);CHKERRQ(ierr);
6577233f9f0SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
658174b6946SBarry Smith   PetscFunctionReturn(0);
659174b6946SBarry Smith }
660174b6946SBarry Smith 
661174b6946SBarry Smith #undef __FUNCT__
6627233f9f0SBarry Smith #define __FUNCT__ "PCView_Exotic"
6637233f9f0SBarry Smith PetscErrorCode PCView_Exotic(PC pc,PetscViewer viewer)
6647233f9f0SBarry Smith {
665f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
6667233f9f0SBarry Smith   PetscErrorCode ierr;
667ace3abfcSBarry Smith   PetscBool      iascii;
66831567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
6697233f9f0SBarry Smith 
6707233f9f0SBarry Smith   PetscFunctionBegin;
6712692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
6727233f9f0SBarry Smith   if (iascii) {
6737233f9f0SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    Exotic type = %s\n",PCExoticTypes[ctx->type]);CHKERRQ(ierr);
6748e722e36SBarry Smith     if (ctx->directSolve) {
6758e722e36SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"      Using direct solver to construct interpolation\n");CHKERRQ(ierr);
6768e722e36SBarry Smith     } else {
6778e722e36SBarry Smith       PetscViewer sviewer;
6788e722e36SBarry Smith       PetscMPIInt rank;
6798e722e36SBarry Smith 
6808e722e36SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"      Using iterative solver to construct interpolation\n");CHKERRQ(ierr);
6818e722e36SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
6828e722e36SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* should not need to push this twice? */
6838e722e36SBarry Smith       ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
6848e722e36SBarry Smith       ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
6858e722e36SBarry Smith       if (!rank) {
6868e722e36SBarry Smith 	ierr = KSPView(ctx->ksp,sviewer);CHKERRQ(ierr);
6878e722e36SBarry Smith       }
6888e722e36SBarry Smith       ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
6898e722e36SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
6908e722e36SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
6918e722e36SBarry Smith     }
6927233f9f0SBarry Smith   }
6937233f9f0SBarry Smith   ierr = PCView_MG(pc,viewer);CHKERRQ(ierr);
6947233f9f0SBarry Smith   PetscFunctionReturn(0);
6957233f9f0SBarry Smith }
6967233f9f0SBarry Smith 
6977233f9f0SBarry Smith #undef __FUNCT__
6987233f9f0SBarry Smith #define __FUNCT__ "PCSetFromOptions_Exotic"
6997233f9f0SBarry Smith PetscErrorCode PCSetFromOptions_Exotic(PC pc)
7007233f9f0SBarry Smith {
7017233f9f0SBarry Smith   PetscErrorCode ierr;
702ace3abfcSBarry Smith   PetscBool      flg;
703f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
7047233f9f0SBarry Smith   PCExoticType   mgctype;
70531567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
7067233f9f0SBarry Smith 
7077233f9f0SBarry Smith   PetscFunctionBegin;
7087233f9f0SBarry Smith   ierr = PetscOptionsHead("Exotic coarse space options");CHKERRQ(ierr);
7097233f9f0SBarry Smith     ierr = PetscOptionsEnum("-pc_exotic_type","face or wirebasket","PCExoticSetType",PCExoticTypes,(PetscEnum)ctx->type,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
7107233f9f0SBarry Smith     if (flg) {
7117233f9f0SBarry Smith       ierr = PCExoticSetType(pc,mgctype);CHKERRQ(ierr);
7127233f9f0SBarry Smith     }
713acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_exotic_direct_solver","use direct solver to construct interpolation","None",ctx->directSolve,&ctx->directSolve,PETSC_NULL);CHKERRQ(ierr);
7148e722e36SBarry Smith     if (!ctx->directSolve) {
7158e722e36SBarry Smith       if (!ctx->ksp) {
7168e722e36SBarry Smith         const char *prefix;
7178e722e36SBarry Smith         ierr = KSPCreate(PETSC_COMM_SELF,&ctx->ksp);CHKERRQ(ierr);
7188e722e36SBarry Smith         ierr = PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
7198e722e36SBarry Smith         ierr = PetscLogObjectParent(pc,ctx->ksp);CHKERRQ(ierr);
7208e722e36SBarry Smith         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
7218e722e36SBarry Smith         ierr = KSPSetOptionsPrefix(ctx->ksp,prefix);CHKERRQ(ierr);
7228e722e36SBarry Smith         ierr = KSPAppendOptionsPrefix(ctx->ksp,"exotic_");CHKERRQ(ierr);
7238e722e36SBarry Smith       }
7248e722e36SBarry Smith       ierr = KSPSetFromOptions(ctx->ksp);CHKERRQ(ierr);
7258e722e36SBarry Smith     }
7267233f9f0SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7277233f9f0SBarry Smith   PetscFunctionReturn(0);
7287233f9f0SBarry Smith }
7297233f9f0SBarry Smith 
730174b6946SBarry Smith 
731174b6946SBarry Smith /*MC
7327233f9f0SBarry Smith      PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces
733174b6946SBarry Smith 
7347233f9f0SBarry Smith      This uses the PCMG infrastructure restricted to two levels and the face and wirebasket based coarse
73524c3aa18SBarry Smith    grid spaces.
73624c3aa18SBarry Smith 
73724c3aa18SBarry 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
73824c3aa18SBarry Smith 
73924c3aa18SBarry Smith    References: These coarse grid spaces originate in the work of Bramble, Pasciak  and Schatz, "The Construction
7403b65e785SBarry Smith    of Preconditioners for Elliptic Problems by Substructing IV", Mathematics of Computation, volume 53 pages 1--24, 1989.
7413b65e785SBarry Smith    They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith,
7423b65e785SBarry Smith    New York University, 1990. They were then explored in great detail in Dryja, Smith, Widlund, "Schwarz Analysis
7433b65e785SBarry Smith    of Iterative Substructuring Methods for Elliptic Problems in Three Dimensions, SIAM Journal on Numerical
7443b65e785SBarry Smith    Analysis, volume 31. pages 1662-1694, 1994. These were developed in the context of iterative substructuring preconditioners.
7453b65e785SBarry Smith    They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund.
7463b65e785SBarry Smith    They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example,
7473b65e785SBarry Smith    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Extending theory for domain decomposition algorithms to irregular subdomains. In Ulrich Langer, Marco
7483b65e785SBarry Smith    Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
7493b65e785SBarry Smith    of the 17th International Conference on Domain Decomposition Methods in
7503b65e785SBarry Smith    Science and Engineering, held in Strobl, Austria, July 3-7, 2006, number 60 in
7513b65e785SBarry Smith    Springer-Verlag, Lecture Notes in Computational Science and Engineering, pages 255-261, 2007.
7523b65e785SBarry Smith    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. A family of energy min-
7533b65e785SBarry Smith    imizing coarse spaces for overlapping Schwarz preconditioners. In Ulrich Langer,
7543b65e785SBarry Smith    Marco Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
7553b65e785SBarry Smith    of the 17th International Conference on Domain Decomposition Methods
7563b65e785SBarry Smith    in Science and Engineering, held in Strobl, Austria, July 3-7, 2006, number 60 in
7573b65e785SBarry Smith    Springer-Verlag, Lecture Notes in Computational Science and Engineering, pages 247-254, 2007
7583b65e785SBarry Smith    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Domain decomposition
7593b65e785SBarry Smith    for less regular subdomains: Overlapping Schwarz in two dimensions. SIAM J.
7603b65e785SBarry Smith    Numer. Anal., 46(4):2153-2168, 2008.
7613b65e785SBarry Smith    Clark R. Dohrmann and Olof B. Widlund. An overlapping Schwarz
7623b65e785SBarry Smith    algorithm for almost incompressible elasticity. Technical Report
7633b65e785SBarry Smith    TR2008-912, Department of Computer Science, Courant Institute
7643b65e785SBarry Smith    of Mathematical Sciences, New York University, May 2008. URL:
7653b65e785SBarry Smith    http://cs.nyu.edu/csweb/Research/TechReports/TR2008-912/TR2008-912.pdf
7667233f9f0SBarry Smith 
7677233f9f0SBarry Smith    Options Database: The usual PCMG options are supported, such as -mg_levels_pc_type <type> -mg_coarse_pc_type <type>
7687233f9f0SBarry Smith       -pc_mg_type <type>
7697233f9f0SBarry Smith 
77025a35f6fSSatish Balay    Level: advanced
771174b6946SBarry Smith 
7726c699258SBarry Smith .seealso:  PCMG, PCSetDM(), PCExoticType, PCExoticSetType()
773174b6946SBarry Smith M*/
774174b6946SBarry Smith 
775174b6946SBarry Smith EXTERN_C_BEGIN
776174b6946SBarry Smith #undef __FUNCT__
7777233f9f0SBarry Smith #define __FUNCT__ "PCCreate_Exotic"
7787087cfbeSBarry Smith PetscErrorCode  PCCreate_Exotic(PC pc)
779174b6946SBarry Smith {
780174b6946SBarry Smith   PetscErrorCode ierr;
7817233f9f0SBarry Smith   PC_Exotic      *ex;
782f3fbd535SBarry Smith   PC_MG          *mg;
783174b6946SBarry Smith 
784174b6946SBarry Smith   PetscFunctionBegin;
785f91d8e95SBarry Smith   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
786f91d8e95SBarry Smith   if (pc->ops->destroy) { ierr =  (*pc->ops->destroy)(pc);CHKERRQ(ierr); pc->data = 0;}
787503cfb0cSBarry Smith   ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr);
788f91d8e95SBarry Smith   ((PetscObject)pc)->type_name = 0;
789f91d8e95SBarry Smith 
790174b6946SBarry Smith   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr);
791174b6946SBarry Smith   ierr = PCMGSetLevels(pc,2,PETSC_NULL);CHKERRQ(ierr);
792c91913d3SJed Brown   ierr = PCMGSetGalerkin(pc,PETSC_TRUE);CHKERRQ(ierr);
79396bdf778SBarry Smith   ierr = PetscNew(PC_Exotic,&ex);CHKERRQ(ierr);\
7947233f9f0SBarry Smith   ex->type = PC_EXOTIC_FACE;
795f3fbd535SBarry Smith   mg = (PC_MG*) pc->data;
79631567311SBarry Smith   mg->innerctx = ex;
7977233f9f0SBarry Smith 
7987233f9f0SBarry Smith 
7997233f9f0SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_Exotic;
8007233f9f0SBarry Smith   pc->ops->view           = PCView_Exotic;
8017233f9f0SBarry Smith   pc->ops->destroy        = PCDestroy_Exotic;
8026c699258SBarry Smith   pc->ops->setup          = PCSetUp_Exotic;
8037233f9f0SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCExoticSetType_C","PCExoticSetType_Exotic",PCExoticSetType_Exotic);CHKERRQ(ierr);
804174b6946SBarry Smith   PetscFunctionReturn(0);
805174b6946SBarry Smith }
806174b6946SBarry Smith EXTERN_C_END
807