xref: /petsc/src/ksp/pc/impls/wb/wb.c (revision acfcf0e5ba359b58e6a8a7af3f239cd7334278e8)
1174b6946SBarry Smith #define PETSCKSP_DLL
2174b6946SBarry Smith 
3174b6946SBarry Smith 
4174b6946SBarry Smith #include "petscpc.h"   /*I "petscpc.h" I*/
5a9f2baa8SLisandro Dalcin #include "petscmg.h"   /*I "petscmg.h" I*/
6174b6946SBarry Smith #include "petscda.h"   /*I "petscda.h" I*/
77233f9f0SBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"
87233f9f0SBarry Smith 
98e722e36SBarry Smith typedef struct {
108e722e36SBarry Smith   PCExoticType type;
118e722e36SBarry Smith   Mat          P;            /* the constructed interpolation matrix */
12ace3abfcSBarry Smith   PetscBool    directSolve;  /* use direct LU factorization to construct interpolation */
138e722e36SBarry Smith   KSP          ksp;
148e722e36SBarry Smith } PC_Exotic;
15174b6946SBarry Smith 
168e722e36SBarry Smith const char *PCExoticTypes[] = {"face","wirebasket","PCExoticType","PC_Exotic",0};
17064c8009SBarry Smith 
18064c8009SBarry Smith 
19064c8009SBarry Smith #undef __FUNCT__
20064c8009SBarry Smith #define __FUNCT__ "DAGetWireBasketInterpolation"
21064c8009SBarry Smith /*
22064c8009SBarry Smith       DAGetWireBasketInterpolation - Gets the interpolation for a wirebasket based coarse space
23064c8009SBarry Smith 
24064c8009SBarry Smith */
256c699258SBarry Smith PetscErrorCode DAGetWireBasketInterpolation(DA da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P)
26064c8009SBarry Smith {
27064c8009SBarry Smith   PetscErrorCode         ierr;
28064c8009SBarry 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;
29064c8009SBarry Smith   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[26],*globals,Ng,*IIint,*IIsurf;
30064c8009SBarry Smith   Mat                    Xint, Xsurf,Xint_tmp;
31064c8009SBarry Smith   IS                     isint,issurf,is,row,col;
32064c8009SBarry Smith   ISLocalToGlobalMapping ltg;
33064c8009SBarry Smith   MPI_Comm               comm;
34064c8009SBarry Smith   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
35064c8009SBarry Smith   MatFactorInfo          info;
36064c8009SBarry Smith   PetscScalar            *xsurf,*xint;
378e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
38064c8009SBarry Smith   PetscScalar            tmp;
39064c8009SBarry Smith #endif
40064c8009SBarry Smith   PetscTable             ht;
41064c8009SBarry Smith 
42064c8009SBarry Smith   PetscFunctionBegin;
43064c8009SBarry Smith   ierr = DAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0);CHKERRQ(ierr);
44e7e72b3dSBarry Smith   if (dof != 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only for single field problems");
45e7e72b3dSBarry Smith   if (dim != 3) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only coded for 3d problems");
46064c8009SBarry Smith   ierr = DAGetCorners(da,0,0,0,&m,&n,&p);CHKERRQ(ierr);
47064c8009SBarry Smith   ierr = DAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr);
48064c8009SBarry Smith   istart = istart ? -1 : 0;
49064c8009SBarry Smith   jstart = jstart ? -1 : 0;
50064c8009SBarry Smith   kstart = kstart ? -1 : 0;
51064c8009SBarry Smith 
52064c8009SBarry Smith   /*
53064c8009SBarry Smith     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge)
54064c8009SBarry Smith     to all the local degrees of freedom (this includes the vertices, edges and faces).
55064c8009SBarry Smith 
56064c8009SBarry Smith     Xint are the subset of the interpolation into the interior
57064c8009SBarry Smith 
58064c8009SBarry Smith     Xface are the interpolation onto faces but not into the interior
59064c8009SBarry Smith 
60064c8009SBarry Smith     Xsurf are the interpolation onto the vertices and edges (the surfbasket)
61064c8009SBarry Smith                                         Xint
62064c8009SBarry Smith     Symbolically one could write P = (  Xface  ) after interchanging the rows to match the natural ordering on the domain
63064c8009SBarry Smith                                         Xsurf
64064c8009SBarry Smith   */
65064c8009SBarry Smith   N     = (m - istart)*(n - jstart)*(p - kstart);
66064c8009SBarry Smith   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
67064c8009SBarry Smith   Nface = 2*( (m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart) );
68064c8009SBarry Smith   Nwire = 4*( (m-2-istart) + (n-2-jstart) + (p-2-kstart) ) + 8;
69064c8009SBarry Smith   Nsurf = Nface + Nwire;
70064c8009SBarry Smith   ierr = MatCreateSeqDense(MPI_COMM_SELF,Nint,26,PETSC_NULL,&Xint);CHKERRQ(ierr);
71064c8009SBarry Smith   ierr = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,26,PETSC_NULL,&Xsurf);CHKERRQ(ierr);
72064c8009SBarry Smith   ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
73064c8009SBarry Smith 
74064c8009SBarry Smith   /*
75064c8009SBarry Smith      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
76064c8009SBarry Smith      Xsurf will be all zero (thus making the coarse matrix singular).
77064c8009SBarry Smith   */
78e32f2f54SBarry 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");
79e32f2f54SBarry 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");
80e32f2f54SBarry 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");
81064c8009SBarry Smith 
82064c8009SBarry Smith   cnt = 0;
83064c8009SBarry Smith   xsurf[cnt++] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + Nsurf] = 1;} xsurf[cnt++ + 2*Nsurf] = 1;
84064c8009SBarry 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;}
85064c8009SBarry Smith   xsurf[cnt++ + 6*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 7*Nsurf] = 1;} xsurf[cnt++ + 8*Nsurf] = 1;
86064c8009SBarry Smith   for (k=1;k<p-1-kstart;k++) {
87064c8009SBarry Smith     xsurf[cnt++ + 9*Nsurf] = 1;  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 10*Nsurf] = 1;}  xsurf[cnt++ + 11*Nsurf] = 1;
88064c8009SBarry Smith     for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 12*Nsurf] = 1; /* these are the interior nodes */ xsurf[cnt++ + 13*Nsurf] = 1;}
89064c8009SBarry Smith     xsurf[cnt++ + 14*Nsurf] = 1;  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 15*Nsurf] = 1;} xsurf[cnt++ + 16*Nsurf] = 1;
90064c8009SBarry Smith   }
91064c8009SBarry Smith   xsurf[cnt++ + 17*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 18*Nsurf] = 1;} xsurf[cnt++ + 19*Nsurf] = 1;
92064c8009SBarry 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;}
93064c8009SBarry Smith   xsurf[cnt++ + 23*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 24*Nsurf] = 1;} xsurf[cnt++ + 25*Nsurf] = 1;
94064c8009SBarry Smith 
958e722e36SBarry Smith   /* interpolations only sum to 1 when using direct solver */
968e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
97064c8009SBarry Smith   for (i=0; i<Nsurf; i++) {
98064c8009SBarry Smith     tmp = 0.0;
99064c8009SBarry Smith     for (j=0; j<26; j++) {
100064c8009SBarry Smith       tmp += xsurf[i+j*Nsurf];
101064c8009SBarry Smith     }
102e32f2f54SBarry 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));
103064c8009SBarry Smith   }
104064c8009SBarry Smith #endif
105064c8009SBarry Smith   ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
106064c8009SBarry Smith   /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/
107064c8009SBarry Smith 
108064c8009SBarry Smith 
109064c8009SBarry Smith   /*
110064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
111064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
112064c8009SBarry Smith             (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it
113064c8009SBarry Smith              is NOT the local DA ordering.)
114064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
115064c8009SBarry Smith   */
116064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
117064c8009SBarry Smith   ierr = PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);CHKERRQ(ierr);
118064c8009SBarry Smith   ierr = PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);CHKERRQ(ierr);
119064c8009SBarry Smith   for (k=0; k<p-kstart; k++) {
120064c8009SBarry Smith     for (j=0; j<n-jstart; j++) {
121064c8009SBarry Smith       for (i=0; i<m-istart; i++) {
122064c8009SBarry Smith         II[c++] = i + j*mwidth + k*mwidth*nwidth;
123064c8009SBarry Smith 
124064c8009SBarry Smith         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
125064c8009SBarry Smith           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
126064c8009SBarry Smith           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
127064c8009SBarry Smith         } else {
128064c8009SBarry Smith           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
129064c8009SBarry Smith           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
130064c8009SBarry Smith         }
131064c8009SBarry Smith       }
132064c8009SBarry Smith     }
133064c8009SBarry Smith   }
134e32f2f54SBarry Smith   if (c != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N");
135e32f2f54SBarry Smith   if (cint != Nint) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint");
136e32f2f54SBarry Smith   if (csurf != Nsurf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf");
137064c8009SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltg);CHKERRQ(ierr);
138064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr);
139064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr);
140064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr);
141064c8009SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
14270b3c8c7SBarry Smith   ierr = ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
14370b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint);CHKERRQ(ierr);
14470b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf);CHKERRQ(ierr);
145064c8009SBarry Smith   ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr);
146064c8009SBarry Smith 
147064c8009SBarry Smith   ierr = MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr);
148064c8009SBarry Smith   A    = *Aholder;
149064c8009SBarry Smith   ierr = PetscFree(Aholder);CHKERRQ(ierr);
150064c8009SBarry Smith 
151064c8009SBarry Smith   ierr = MatGetSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr);
152064c8009SBarry Smith   ierr = MatGetSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr);
153064c8009SBarry Smith   ierr = MatGetSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr);
154064c8009SBarry Smith 
155064c8009SBarry Smith   /*
156064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
157064c8009SBarry Smith   */
1588e722e36SBarry Smith   ierr = MatDuplicate(Xint,MAT_DO_NOT_COPY_VALUES,&Xint_tmp);CHKERRQ(ierr);
1598e722e36SBarry Smith   ierr = MatMatMult(Ais,Xsurf,MAT_REUSE_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr);
1608e722e36SBarry Smith   ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr);
1618e722e36SBarry Smith   if (exotic->directSolve) {
1622692d6eeSBarry Smith     ierr = MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr);
163064c8009SBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
1642692d6eeSBarry Smith     ierr = MatGetOrdering(Aii,MATORDERINGND,&row,&col);CHKERRQ(ierr);
165064c8009SBarry Smith     ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr);
166064c8009SBarry Smith     ierr = ISDestroy(row);CHKERRQ(ierr);
167064c8009SBarry Smith     ierr = ISDestroy(col);CHKERRQ(ierr);
168064c8009SBarry Smith     ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr);
169064c8009SBarry Smith     ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr);
170064c8009SBarry Smith     ierr = MatDestroy(iAii);CHKERRQ(ierr);
1718e722e36SBarry Smith   } else {
1728e722e36SBarry Smith     Vec         b,x;
1738e722e36SBarry Smith     PetscScalar *xint_tmp;
174064c8009SBarry Smith 
1758e722e36SBarry Smith     ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr);
1768e722e36SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&x);CHKERRQ(ierr);
1778e722e36SBarry Smith     ierr = MatGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
1788e722e36SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&b);CHKERRQ(ierr);
1798e722e36SBarry Smith     ierr = KSPSetOperators(exotic->ksp,Aii,Aii,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1808e722e36SBarry Smith     for (i=0; i<26; i++) {
1818e722e36SBarry Smith       ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr);
1828e722e36SBarry Smith       ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr);
1838e722e36SBarry Smith       ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr);
1848e722e36SBarry Smith       ierr = VecResetArray(x);CHKERRQ(ierr);
1858e722e36SBarry Smith       ierr = VecResetArray(b);CHKERRQ(ierr);
1868e722e36SBarry Smith     }
1878e722e36SBarry Smith     ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr);
1888e722e36SBarry Smith     ierr = MatRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
1898e722e36SBarry Smith     ierr = VecDestroy(x);CHKERRQ(ierr);
1908e722e36SBarry Smith     ierr = VecDestroy(b);CHKERRQ(ierr);
1918e722e36SBarry Smith   }
1928e722e36SBarry Smith   ierr = MatDestroy(Xint_tmp);CHKERRQ(ierr);
1938e722e36SBarry Smith 
1948e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
195064c8009SBarry Smith   ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr);
196064c8009SBarry Smith   for (i=0; i<Nint; i++) {
197064c8009SBarry Smith     tmp = 0.0;
198064c8009SBarry Smith     for (j=0; j<26; j++) {
199064c8009SBarry Smith       tmp += xint[i+j*Nint];
200064c8009SBarry Smith     }
201e32f2f54SBarry 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));
202064c8009SBarry Smith   }
203064c8009SBarry Smith   ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr);
204064c8009SBarry Smith   /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
205064c8009SBarry Smith #endif
206064c8009SBarry Smith 
207064c8009SBarry Smith 
208064c8009SBarry Smith   /*         total vertices             total faces                                  total edges */
209064c8009SBarry 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);
210064c8009SBarry Smith 
211064c8009SBarry Smith   /*
212064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
213064c8009SBarry Smith   */
214064c8009SBarry Smith   cnt = 0;
215064c8009SBarry Smith   gl[cnt++] = 0;  { gl[cnt++] = 1;} gl[cnt++] = m-istart-1;
216064c8009SBarry Smith   { gl[cnt++] = mwidth;  { gl[cnt++] = mwidth+1;} gl[cnt++] = mwidth + m-istart-1;}
217064c8009SBarry Smith   gl[cnt++] = mwidth*(n-jstart-1);  { gl[cnt++] = mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*(n-jstart-1) + m-istart-1;
218064c8009SBarry Smith   {
219064c8009SBarry Smith     gl[cnt++] = mwidth*nwidth;  { gl[cnt++] = mwidth*nwidth+1;}  gl[cnt++] = mwidth*nwidth+ m-istart-1;
220064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
221064c8009SBarry 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;
222064c8009SBarry Smith   }
223064c8009SBarry 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;
224064c8009SBarry 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;}
225064c8009SBarry 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;
226064c8009SBarry Smith 
227064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
228064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
229064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,26,gl,gl);CHKERRQ(ierr);
230064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
231064c8009SBarry Smith   ierr = PetscMalloc(26*mp*np*pp*sizeof(PetscInt),&globals);CHKERRQ(ierr);
232064c8009SBarry Smith   ierr = MPI_Allgather(gl,26,MPIU_INT,globals,26,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
233064c8009SBarry Smith 
234064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
235064c8009SBarry Smith   ierr = PetscTableCreate(Ntotal/3,&ht);CHKERRQ(ierr);
236064c8009SBarry Smith   for (i=0; i<26*mp*np*pp; i++){
237064c8009SBarry Smith     ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr);
238064c8009SBarry Smith   }
239064c8009SBarry Smith   ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr);
240e32f2f54SBarry 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);
241064c8009SBarry Smith   ierr = PetscFree(globals);CHKERRQ(ierr);
242064c8009SBarry Smith   for (i=0; i<26; i++) {
243064c8009SBarry Smith     ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr);
244064c8009SBarry Smith     gl[i]--;
245064c8009SBarry Smith   }
246064c8009SBarry Smith   ierr = PetscTableDestroy(ht);CHKERRQ(ierr);
247064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
248064c8009SBarry Smith 
249064c8009SBarry Smith   /* construct global interpolation matrix */
250064c8009SBarry Smith   ierr = MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);CHKERRQ(ierr);
251064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
252064c8009SBarry Smith     ierr = MatCreateMPIAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint+Nsurf,PETSC_NULL,P);CHKERRQ(ierr);
253064c8009SBarry Smith   } else {
254064c8009SBarry Smith     ierr = MatZeroEntries(*P);CHKERRQ(ierr);
255064c8009SBarry Smith   }
256064c8009SBarry Smith   ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
257064c8009SBarry Smith   ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr);
258064c8009SBarry Smith   ierr = MatSetValues(*P,Nint,IIint,26,gl,xint,INSERT_VALUES);CHKERRQ(ierr);
259064c8009SBarry Smith   ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr);
260064c8009SBarry Smith   ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
261064c8009SBarry Smith   ierr = MatSetValues(*P,Nsurf,IIsurf,26,gl,xsurf,INSERT_VALUES);CHKERRQ(ierr);
262064c8009SBarry Smith   ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
263064c8009SBarry Smith   ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
264064c8009SBarry Smith   ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
265064c8009SBarry Smith   ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr);
266064c8009SBarry Smith 
2678e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
268064c8009SBarry Smith   {
269064c8009SBarry Smith     Vec         x,y;
270064c8009SBarry Smith     PetscScalar *yy;
271064c8009SBarry Smith     ierr = VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr);
272064c8009SBarry Smith     ierr = VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr);
273064c8009SBarry Smith     ierr = VecSet(x,1.0);CHKERRQ(ierr);
274064c8009SBarry Smith     ierr = MatMult(*P,x,y);CHKERRQ(ierr);
275064c8009SBarry Smith     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
276064c8009SBarry Smith     for (i=0; i<Ng; i++) {
277e32f2f54SBarry 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]));
278064c8009SBarry Smith     }
279064c8009SBarry Smith     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
280064c8009SBarry Smith     ierr = VecDestroy(x);CHKERRQ(ierr);
281064c8009SBarry Smith     ierr = VecDestroy(y);CHKERRQ(ierr);
282064c8009SBarry Smith   }
283064c8009SBarry Smith #endif
284064c8009SBarry Smith 
285064c8009SBarry Smith   ierr = MatDestroy(Aii);CHKERRQ(ierr);
286064c8009SBarry Smith   ierr = MatDestroy(Ais);CHKERRQ(ierr);
287064c8009SBarry Smith   ierr = MatDestroy(Asi);CHKERRQ(ierr);
288064c8009SBarry Smith   ierr = MatDestroy(A);CHKERRQ(ierr);
289064c8009SBarry Smith   ierr = ISDestroy(is);CHKERRQ(ierr);
290064c8009SBarry Smith   ierr = ISDestroy(isint);CHKERRQ(ierr);
291064c8009SBarry Smith   ierr = ISDestroy(issurf);CHKERRQ(ierr);
292064c8009SBarry Smith   ierr = MatDestroy(Xint);CHKERRQ(ierr);
293064c8009SBarry Smith   ierr = MatDestroy(Xsurf);CHKERRQ(ierr);
294064c8009SBarry Smith   PetscFunctionReturn(0);
295064c8009SBarry Smith }
296064c8009SBarry Smith 
297064c8009SBarry Smith #undef __FUNCT__
298064c8009SBarry Smith #define __FUNCT__ "DAGetFaceInterpolation"
299064c8009SBarry Smith /*
300064c8009SBarry Smith       DAGetFaceInterpolation - Gets the interpolation for a face based coarse space
301064c8009SBarry Smith 
302064c8009SBarry Smith */
3036c699258SBarry Smith PetscErrorCode DAGetFaceInterpolation(DA da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P)
304064c8009SBarry Smith {
305064c8009SBarry Smith   PetscErrorCode         ierr;
306064c8009SBarry 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;
307064c8009SBarry Smith   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[6],*globals,Ng,*IIint,*IIsurf;
308064c8009SBarry Smith   Mat                    Xint, Xsurf,Xint_tmp;
309064c8009SBarry Smith   IS                     isint,issurf,is,row,col;
310064c8009SBarry Smith   ISLocalToGlobalMapping ltg;
311064c8009SBarry Smith   MPI_Comm               comm;
312064c8009SBarry Smith   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
313064c8009SBarry Smith   MatFactorInfo          info;
314064c8009SBarry Smith   PetscScalar            *xsurf,*xint;
315064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
316064c8009SBarry Smith   PetscScalar            tmp;
317064c8009SBarry Smith #endif
318064c8009SBarry Smith   PetscTable             ht;
319064c8009SBarry Smith 
320064c8009SBarry Smith   PetscFunctionBegin;
321064c8009SBarry Smith   ierr = DAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0);CHKERRQ(ierr);
322e7e72b3dSBarry Smith   if (dof != 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only for single field problems");
323e7e72b3dSBarry Smith   if (dim != 3) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only coded for 3d problems");
324064c8009SBarry Smith   ierr = DAGetCorners(da,0,0,0,&m,&n,&p);CHKERRQ(ierr);
325064c8009SBarry Smith   ierr = DAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr);
326064c8009SBarry Smith   istart = istart ? -1 : 0;
327064c8009SBarry Smith   jstart = jstart ? -1 : 0;
328064c8009SBarry Smith   kstart = kstart ? -1 : 0;
329064c8009SBarry Smith 
330064c8009SBarry Smith   /*
331064c8009SBarry Smith     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge)
332064c8009SBarry Smith     to all the local degrees of freedom (this includes the vertices, edges and faces).
333064c8009SBarry Smith 
334064c8009SBarry Smith     Xint are the subset of the interpolation into the interior
335064c8009SBarry Smith 
336064c8009SBarry Smith     Xface are the interpolation onto faces but not into the interior
337064c8009SBarry Smith 
338064c8009SBarry Smith     Xsurf are the interpolation onto the vertices and edges (the surfbasket)
339064c8009SBarry Smith                                         Xint
340064c8009SBarry Smith     Symbolically one could write P = (  Xface  ) after interchanging the rows to match the natural ordering on the domain
341064c8009SBarry Smith                                         Xsurf
342064c8009SBarry Smith   */
343064c8009SBarry Smith   N     = (m - istart)*(n - jstart)*(p - kstart);
344064c8009SBarry Smith   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
345064c8009SBarry Smith   Nface = 2*( (m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart) );
346064c8009SBarry Smith   Nwire = 4*( (m-2-istart) + (n-2-jstart) + (p-2-kstart) ) + 8;
347064c8009SBarry Smith   Nsurf = Nface + Nwire;
348064c8009SBarry Smith   ierr = MatCreateSeqDense(MPI_COMM_SELF,Nint,6,PETSC_NULL,&Xint);CHKERRQ(ierr);
349064c8009SBarry Smith   ierr = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,6,PETSC_NULL,&Xsurf);CHKERRQ(ierr);
350064c8009SBarry Smith   ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
351064c8009SBarry Smith 
352064c8009SBarry Smith   /*
353064c8009SBarry Smith      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
354064c8009SBarry Smith      Xsurf will be all zero (thus making the coarse matrix singular).
355064c8009SBarry Smith   */
356e32f2f54SBarry 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");
357e32f2f54SBarry 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");
358e32f2f54SBarry 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");
359064c8009SBarry Smith 
360064c8009SBarry Smith   cnt = 0;
361064c8009SBarry Smith   for (j=1;j<n-1-jstart;j++) {  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 0*Nsurf] = 1;} }
362064c8009SBarry Smith    for (k=1;k<p-1-kstart;k++) {
363064c8009SBarry Smith     for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 1*Nsurf] = 1;}
364064c8009SBarry Smith     for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 2*Nsurf] = 1; /* these are the interior nodes */ xsurf[cnt++ + 3*Nsurf] = 1;}
365064c8009SBarry Smith     for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 4*Nsurf] = 1;}
366064c8009SBarry Smith   }
367064c8009SBarry Smith   for (j=1;j<n-1-jstart;j++) {for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 5*Nsurf] = 1;} }
368064c8009SBarry Smith 
369064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
370064c8009SBarry Smith   for (i=0; i<Nsurf; i++) {
371064c8009SBarry Smith     tmp = 0.0;
372064c8009SBarry Smith     for (j=0; j<6; j++) {
373064c8009SBarry Smith       tmp += xsurf[i+j*Nsurf];
374064c8009SBarry Smith     }
375e32f2f54SBarry 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));
376064c8009SBarry Smith   }
377064c8009SBarry Smith #endif
378064c8009SBarry Smith   ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
379064c8009SBarry Smith   /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/
380064c8009SBarry Smith 
381064c8009SBarry Smith 
382064c8009SBarry Smith   /*
383064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
384064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
385064c8009SBarry Smith             (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it
386064c8009SBarry Smith              is NOT the local DA ordering.)
387064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
388064c8009SBarry Smith   */
389064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
390064c8009SBarry Smith   ierr = PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);CHKERRQ(ierr);
391064c8009SBarry Smith   ierr = PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);CHKERRQ(ierr);
392064c8009SBarry Smith   for (k=0; k<p-kstart; k++) {
393064c8009SBarry Smith     for (j=0; j<n-jstart; j++) {
394064c8009SBarry Smith       for (i=0; i<m-istart; i++) {
395064c8009SBarry Smith         II[c++] = i + j*mwidth + k*mwidth*nwidth;
396064c8009SBarry Smith 
397064c8009SBarry Smith         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
398064c8009SBarry Smith           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
399064c8009SBarry Smith           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
400064c8009SBarry Smith         } else {
401064c8009SBarry Smith           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
402064c8009SBarry Smith           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
403064c8009SBarry Smith         }
404064c8009SBarry Smith       }
405064c8009SBarry Smith     }
406064c8009SBarry Smith   }
407e32f2f54SBarry Smith   if (c != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N");
408e32f2f54SBarry Smith   if (cint != Nint) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint");
409e32f2f54SBarry Smith   if (csurf != Nsurf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf");
410064c8009SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltg);CHKERRQ(ierr);
411064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr);
412064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr);
413064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr);
414064c8009SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
41570b3c8c7SBarry Smith   ierr = ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
41670b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint);CHKERRQ(ierr);
41770b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf);CHKERRQ(ierr);
418064c8009SBarry Smith   ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr);
419064c8009SBarry Smith 
420064c8009SBarry Smith   ierr = MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr);
421064c8009SBarry Smith   A    = *Aholder;
422064c8009SBarry Smith   ierr = PetscFree(Aholder);CHKERRQ(ierr);
423064c8009SBarry Smith 
424064c8009SBarry Smith   ierr = MatGetSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr);
425064c8009SBarry Smith   ierr = MatGetSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr);
426064c8009SBarry Smith   ierr = MatGetSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr);
427064c8009SBarry Smith 
428064c8009SBarry Smith   /*
429064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
430064c8009SBarry Smith   */
431064c8009SBarry Smith   ierr = MatDuplicate(Xint,MAT_DO_NOT_COPY_VALUES,&Xint_tmp);CHKERRQ(ierr);
432064c8009SBarry Smith   ierr = MatMatMult(Ais,Xsurf,MAT_REUSE_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr);
433064c8009SBarry Smith   ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr);
434064c8009SBarry Smith 
4358e722e36SBarry Smith   if (exotic->directSolve) {
4362692d6eeSBarry Smith     ierr = MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr);
437064c8009SBarry Smith     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
4382692d6eeSBarry Smith     ierr = MatGetOrdering(Aii,MATORDERINGND,&row,&col);CHKERRQ(ierr);
439064c8009SBarry Smith     ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr);
440064c8009SBarry Smith     ierr = ISDestroy(row);CHKERRQ(ierr);
441064c8009SBarry Smith     ierr = ISDestroy(col);CHKERRQ(ierr);
442064c8009SBarry Smith     ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr);
443064c8009SBarry Smith     ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr);
444064c8009SBarry Smith     ierr = MatDestroy(iAii);CHKERRQ(ierr);
445064c8009SBarry Smith   } else {
446064c8009SBarry Smith     Vec         b,x;
447064c8009SBarry Smith     PetscScalar *xint_tmp;
448064c8009SBarry Smith 
449064c8009SBarry Smith     ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr);
450064c8009SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&x);CHKERRQ(ierr);
451064c8009SBarry Smith     ierr = MatGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
452064c8009SBarry Smith     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&b);CHKERRQ(ierr);
4538e722e36SBarry Smith     ierr = KSPSetOperators(exotic->ksp,Aii,Aii,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
454064c8009SBarry Smith     for (i=0; i<6; i++) {
455064c8009SBarry Smith       ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr);
456064c8009SBarry Smith       ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr);
4578e722e36SBarry Smith       ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr);
458064c8009SBarry Smith       ierr = VecResetArray(x);CHKERRQ(ierr);
459064c8009SBarry Smith       ierr = VecResetArray(b);CHKERRQ(ierr);
460064c8009SBarry Smith     }
461064c8009SBarry Smith     ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr);
462064c8009SBarry Smith     ierr = MatRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr);
463064c8009SBarry Smith     ierr = VecDestroy(x);CHKERRQ(ierr);
464064c8009SBarry Smith     ierr = VecDestroy(b);CHKERRQ(ierr);
465064c8009SBarry Smith   }
466064c8009SBarry Smith   ierr = MatDestroy(Xint_tmp);CHKERRQ(ierr);
467064c8009SBarry Smith 
468064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
469064c8009SBarry Smith   ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr);
470064c8009SBarry Smith   for (i=0; i<Nint; i++) {
471064c8009SBarry Smith     tmp = 0.0;
472064c8009SBarry Smith     for (j=0; j<6; j++) {
473064c8009SBarry Smith       tmp += xint[i+j*Nint];
474064c8009SBarry Smith     }
475e32f2f54SBarry 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));
476064c8009SBarry Smith   }
477064c8009SBarry Smith   ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr);
478064c8009SBarry Smith   /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
479064c8009SBarry Smith #endif
480064c8009SBarry Smith 
481064c8009SBarry Smith 
482064c8009SBarry Smith   /*         total faces    */
483064c8009SBarry Smith   Ntotal =  mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1);
484064c8009SBarry Smith 
485064c8009SBarry Smith   /*
486064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
487064c8009SBarry Smith   */
488064c8009SBarry Smith   cnt = 0;
489064c8009SBarry Smith   { gl[cnt++] = mwidth+1;}
490064c8009SBarry Smith   {
491064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth+1;}
492064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
493064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;}
494064c8009SBarry Smith   }
495064c8009SBarry Smith   { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;}
496064c8009SBarry Smith 
497064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
498064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
499064c8009SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltg,6,gl,gl);CHKERRQ(ierr);
500064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
501064c8009SBarry Smith   ierr = PetscMalloc(6*mp*np*pp*sizeof(PetscInt),&globals);CHKERRQ(ierr);
502064c8009SBarry Smith   ierr = MPI_Allgather(gl,6,MPIU_INT,globals,6,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
503064c8009SBarry Smith 
504064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
505064c8009SBarry Smith   ierr = PetscTableCreate(Ntotal/3,&ht);CHKERRQ(ierr);
506064c8009SBarry Smith   for (i=0; i<6*mp*np*pp; i++){
507064c8009SBarry Smith     ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr);
508064c8009SBarry Smith   }
509064c8009SBarry Smith   ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr);
510e32f2f54SBarry 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);
511064c8009SBarry Smith   ierr = PetscFree(globals);CHKERRQ(ierr);
512064c8009SBarry Smith   for (i=0; i<6; i++) {
513064c8009SBarry Smith     ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr);
514064c8009SBarry Smith     gl[i]--;
515064c8009SBarry Smith   }
516064c8009SBarry Smith   ierr = PetscTableDestroy(ht);CHKERRQ(ierr);
517064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
518064c8009SBarry Smith 
519064c8009SBarry Smith   /* construct global interpolation matrix */
520064c8009SBarry Smith   ierr = MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);CHKERRQ(ierr);
521064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
522064c8009SBarry Smith     ierr = MatCreateMPIAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint,PETSC_NULL,P);CHKERRQ(ierr);
523064c8009SBarry Smith   } else {
524064c8009SBarry Smith     ierr = MatZeroEntries(*P);CHKERRQ(ierr);
525064c8009SBarry Smith   }
526064c8009SBarry Smith   ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
527064c8009SBarry Smith   ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr);
528064c8009SBarry Smith   ierr = MatSetValues(*P,Nint,IIint,6,gl,xint,INSERT_VALUES);CHKERRQ(ierr);
529064c8009SBarry Smith   ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr);
530064c8009SBarry Smith   ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr);
531064c8009SBarry Smith   ierr = MatSetValues(*P,Nsurf,IIsurf,6,gl,xsurf,INSERT_VALUES);CHKERRQ(ierr);
532064c8009SBarry Smith   ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr);
533064c8009SBarry Smith   ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
534064c8009SBarry Smith   ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
535064c8009SBarry Smith   ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr);
536064c8009SBarry Smith 
537064c8009SBarry Smith 
538064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
539064c8009SBarry Smith   {
540064c8009SBarry Smith     Vec         x,y;
541064c8009SBarry Smith     PetscScalar *yy;
542064c8009SBarry Smith     ierr = VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr);
543064c8009SBarry Smith     ierr = VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr);
544064c8009SBarry Smith     ierr = VecSet(x,1.0);CHKERRQ(ierr);
545064c8009SBarry Smith     ierr = MatMult(*P,x,y);CHKERRQ(ierr);
546064c8009SBarry Smith     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
547064c8009SBarry Smith     for (i=0; i<Ng; i++) {
548e32f2f54SBarry 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]));
549064c8009SBarry Smith     }
550064c8009SBarry Smith     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
551064c8009SBarry Smith     ierr = VecDestroy(x);CHKERRQ(ierr);
552064c8009SBarry Smith     ierr = VecDestroy(y);CHKERRQ(ierr);
553064c8009SBarry Smith   }
554064c8009SBarry Smith #endif
555064c8009SBarry Smith 
556064c8009SBarry Smith   ierr = MatDestroy(Aii);CHKERRQ(ierr);
557064c8009SBarry Smith   ierr = MatDestroy(Ais);CHKERRQ(ierr);
558064c8009SBarry Smith   ierr = MatDestroy(Asi);CHKERRQ(ierr);
559064c8009SBarry Smith   ierr = MatDestroy(A);CHKERRQ(ierr);
560064c8009SBarry Smith   ierr = ISDestroy(is);CHKERRQ(ierr);
561064c8009SBarry Smith   ierr = ISDestroy(isint);CHKERRQ(ierr);
562064c8009SBarry Smith   ierr = ISDestroy(issurf);CHKERRQ(ierr);
563064c8009SBarry Smith   ierr = MatDestroy(Xint);CHKERRQ(ierr);
564064c8009SBarry Smith   ierr = MatDestroy(Xsurf);CHKERRQ(ierr);
565064c8009SBarry Smith   PetscFunctionReturn(0);
566064c8009SBarry Smith }
567174b6946SBarry Smith 
5687233f9f0SBarry Smith 
569174b6946SBarry Smith #undef __FUNCT__
5707233f9f0SBarry Smith #define __FUNCT__ "PCExoticSetType"
5717233f9f0SBarry Smith /*@
5727233f9f0SBarry Smith    PCExoticSetType - Sets the type of coarse grid interpolation to use
5737233f9f0SBarry Smith 
5743f9fe445SBarry Smith    Logically Collective on PC
5757233f9f0SBarry Smith 
5767233f9f0SBarry Smith    Input Parameters:
5777233f9f0SBarry Smith +  pc - the preconditioner context
5787233f9f0SBarry Smith -  type - either PC_EXOTIC_FACE or PC_EXOTIC_WIREBASKET (defaults to face)
5797233f9f0SBarry Smith 
580563e08c6SBarry Smith    Notes: The face based interpolation has 1 degree of freedom per face and ignores the
581563e08c6SBarry Smith      edge and vertex values completely in the coarse problem. For any seven point
582563e08c6SBarry Smith      stencil the interpolation of a constant on all faces into the interior is that constant.
583563e08c6SBarry Smith 
584563e08c6SBarry Smith      The wirebasket interpolation has 1 degree of freedom per vertex, per edge and
585563e08c6SBarry Smith      per face. A constant on the subdomain boundary is interpolated as that constant
586563e08c6SBarry Smith      in the interior of the domain.
587563e08c6SBarry Smith 
588563e08c6SBarry Smith      The coarse grid matrix is obtained via the Galerkin computation A_c = R A R^T, hence
589563e08c6SBarry Smith      if A is nonsingular A_c is also nonsingular.
590563e08c6SBarry Smith 
591563e08c6SBarry Smith      Both interpolations are suitable for only scalar problems.
592563e08c6SBarry Smith 
5937233f9f0SBarry Smith    Level: intermediate
5947233f9f0SBarry Smith 
5957233f9f0SBarry Smith 
5967233f9f0SBarry Smith .seealso: PCEXOTIC, PCExoticType()
5977233f9f0SBarry Smith @*/
5987233f9f0SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCExoticSetType(PC pc,PCExoticType type)
5997233f9f0SBarry Smith {
6004ac538c5SBarry Smith   PetscErrorCode ierr;
6017233f9f0SBarry Smith 
6027233f9f0SBarry Smith   PetscFunctionBegin;
6030700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
604c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
6054ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCExoticSetType_C",(PC,PCExoticType),(pc,type));CHKERRQ(ierr);
6067233f9f0SBarry Smith   PetscFunctionReturn(0);
6077233f9f0SBarry Smith }
6087233f9f0SBarry Smith 
6097233f9f0SBarry Smith #undef __FUNCT__
6107233f9f0SBarry Smith #define __FUNCT__ "PCExoticSetType_Exotic"
6117233f9f0SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCExoticSetType_Exotic(PC pc,PCExoticType type)
6127233f9f0SBarry Smith {
613f3fbd535SBarry Smith   PC_MG     *mg = (PC_MG*)pc->data;
61431567311SBarry Smith   PC_Exotic *ctx = (PC_Exotic*) mg->innerctx;
6157233f9f0SBarry Smith 
6167233f9f0SBarry Smith   PetscFunctionBegin;
6177233f9f0SBarry Smith   ctx->type = type;
6187233f9f0SBarry Smith   PetscFunctionReturn(0);
6197233f9f0SBarry Smith }
6207233f9f0SBarry Smith 
6217233f9f0SBarry Smith #undef __FUNCT__
6227233f9f0SBarry Smith #define __FUNCT__ "PCSetUp_Exotic"
6237233f9f0SBarry Smith PetscErrorCode PCSetUp_Exotic(PC pc)
624174b6946SBarry Smith {
625174b6946SBarry Smith   PetscErrorCode ierr;
62696bdf778SBarry Smith   Mat            A;
627f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
62831567311SBarry Smith   PC_Exotic      *ex = (PC_Exotic*) mg->innerctx;
62996bdf778SBarry Smith   MatReuse       reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX;
630174b6946SBarry Smith 
631174b6946SBarry Smith   PetscFunctionBegin;
632e7e72b3dSBarry Smith   if (!pc->dm) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to call PCSetDM() before using this PC");
633174b6946SBarry Smith   ierr = PCGetOperators(pc,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr);
6347233f9f0SBarry Smith   if (ex->type == PC_EXOTIC_FACE) {
6356c699258SBarry Smith     ierr = DAGetFaceInterpolation((DA)pc->dm,ex,A,reuse,&ex->P);CHKERRQ(ierr);
6367233f9f0SBarry Smith   } else if (ex->type == PC_EXOTIC_WIREBASKET) {
6376c699258SBarry Smith     ierr = DAGetWireBasketInterpolation((DA)pc->dm,ex,A,reuse,&ex->P);CHKERRQ(ierr);
63865e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_PLIB,"Unknown exotic coarse space %d",ex->type);
63996bdf778SBarry Smith   ierr = PCMGSetInterpolation(pc,1,ex->P);CHKERRQ(ierr);
6407233f9f0SBarry Smith   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
641174b6946SBarry Smith   PetscFunctionReturn(0);
642174b6946SBarry Smith }
643174b6946SBarry Smith 
644174b6946SBarry Smith #undef __FUNCT__
6457233f9f0SBarry Smith #define __FUNCT__ "PCDestroy_Exotic"
6467233f9f0SBarry Smith PetscErrorCode PCDestroy_Exotic(PC pc)
647174b6946SBarry Smith {
648174b6946SBarry Smith   PetscErrorCode ierr;
649f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
65031567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
651174b6946SBarry Smith 
652174b6946SBarry Smith   PetscFunctionBegin;
65396bdf778SBarry Smith   if (ctx->P) {ierr = MatDestroy(ctx->P);CHKERRQ(ierr);}
6548e722e36SBarry Smith   if (ctx->ksp) {ierr = KSPDestroy(ctx->ksp);CHKERRQ(ierr);}
6557233f9f0SBarry Smith   ierr = PetscFree(ctx);CHKERRQ(ierr);
6567233f9f0SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
657174b6946SBarry Smith   PetscFunctionReturn(0);
658174b6946SBarry Smith }
659174b6946SBarry Smith 
660174b6946SBarry Smith #undef __FUNCT__
6617233f9f0SBarry Smith #define __FUNCT__ "PCView_Exotic"
6627233f9f0SBarry Smith PetscErrorCode PCView_Exotic(PC pc,PetscViewer viewer)
6637233f9f0SBarry Smith {
664f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
6657233f9f0SBarry Smith   PetscErrorCode ierr;
666ace3abfcSBarry Smith   PetscBool      iascii;
66731567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
6687233f9f0SBarry Smith 
6697233f9f0SBarry Smith   PetscFunctionBegin;
6702692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
6717233f9f0SBarry Smith   if (iascii) {
6727233f9f0SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    Exotic type = %s\n",PCExoticTypes[ctx->type]);CHKERRQ(ierr);
6738e722e36SBarry Smith     if (ctx->directSolve) {
6748e722e36SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"      Using direct solver to construct interpolation\n");CHKERRQ(ierr);
6758e722e36SBarry Smith     } else {
6768e722e36SBarry Smith       PetscViewer sviewer;
6778e722e36SBarry Smith       PetscMPIInt rank;
6788e722e36SBarry Smith 
6798e722e36SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"      Using iterative solver to construct interpolation\n");CHKERRQ(ierr);
6808e722e36SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
6818e722e36SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* should not need to push this twice? */
6828e722e36SBarry Smith       ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
6838e722e36SBarry Smith       ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
6848e722e36SBarry Smith       if (!rank) {
6858e722e36SBarry Smith 	ierr = KSPView(ctx->ksp,sviewer);CHKERRQ(ierr);
6868e722e36SBarry Smith       }
6878e722e36SBarry Smith       ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
6888e722e36SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
6898e722e36SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
6908e722e36SBarry Smith     }
6917233f9f0SBarry Smith   }
6927233f9f0SBarry Smith   ierr = PCView_MG(pc,viewer);CHKERRQ(ierr);
6937233f9f0SBarry Smith   PetscFunctionReturn(0);
6947233f9f0SBarry Smith }
6957233f9f0SBarry Smith 
6967233f9f0SBarry Smith #undef __FUNCT__
6977233f9f0SBarry Smith #define __FUNCT__ "PCSetFromOptions_Exotic"
6987233f9f0SBarry Smith PetscErrorCode PCSetFromOptions_Exotic(PC pc)
6997233f9f0SBarry Smith {
7007233f9f0SBarry Smith   PetscErrorCode ierr;
701ace3abfcSBarry Smith   PetscBool      flg;
702f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
7037233f9f0SBarry Smith   PCExoticType   mgctype;
70431567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
7057233f9f0SBarry Smith 
7067233f9f0SBarry Smith   PetscFunctionBegin;
7077233f9f0SBarry Smith   ierr = PetscOptionsHead("Exotic coarse space options");CHKERRQ(ierr);
7087233f9f0SBarry Smith     ierr = PetscOptionsEnum("-pc_exotic_type","face or wirebasket","PCExoticSetType",PCExoticTypes,(PetscEnum)ctx->type,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
7097233f9f0SBarry Smith     if (flg) {
7107233f9f0SBarry Smith       ierr = PCExoticSetType(pc,mgctype);CHKERRQ(ierr);
7117233f9f0SBarry Smith     }
712*acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_exotic_direct_solver","use direct solver to construct interpolation","None",ctx->directSolve,&ctx->directSolve,PETSC_NULL);CHKERRQ(ierr);
7138e722e36SBarry Smith     if (!ctx->directSolve) {
7148e722e36SBarry Smith       if (!ctx->ksp) {
7158e722e36SBarry Smith         const char *prefix;
7168e722e36SBarry Smith         ierr = KSPCreate(PETSC_COMM_SELF,&ctx->ksp);CHKERRQ(ierr);
7178e722e36SBarry Smith         ierr = PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
7188e722e36SBarry Smith         ierr = PetscLogObjectParent(pc,ctx->ksp);CHKERRQ(ierr);
7198e722e36SBarry Smith         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
7208e722e36SBarry Smith         ierr = KSPSetOptionsPrefix(ctx->ksp,prefix);CHKERRQ(ierr);
7218e722e36SBarry Smith         ierr = KSPAppendOptionsPrefix(ctx->ksp,"exotic_");CHKERRQ(ierr);
7228e722e36SBarry Smith       }
7238e722e36SBarry Smith       ierr = KSPSetFromOptions(ctx->ksp);CHKERRQ(ierr);
7248e722e36SBarry Smith     }
7257233f9f0SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7267233f9f0SBarry Smith   PetscFunctionReturn(0);
7277233f9f0SBarry Smith }
7287233f9f0SBarry Smith 
729174b6946SBarry Smith 
730174b6946SBarry Smith /*MC
7317233f9f0SBarry Smith      PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces
732174b6946SBarry Smith 
7337233f9f0SBarry Smith      This uses the PCMG infrastructure restricted to two levels and the face and wirebasket based coarse
73424c3aa18SBarry Smith    grid spaces.
73524c3aa18SBarry Smith 
73624c3aa18SBarry 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
73724c3aa18SBarry Smith 
73824c3aa18SBarry Smith    References: These coarse grid spaces originate in the work of Bramble, Pasciak  and Schatz, "The Construction
7393b65e785SBarry Smith    of Preconditioners for Elliptic Problems by Substructing IV", Mathematics of Computation, volume 53 pages 1--24, 1989.
7403b65e785SBarry Smith    They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith,
7413b65e785SBarry Smith    New York University, 1990. They were then explored in great detail in Dryja, Smith, Widlund, "Schwarz Analysis
7423b65e785SBarry Smith    of Iterative Substructuring Methods for Elliptic Problems in Three Dimensions, SIAM Journal on Numerical
7433b65e785SBarry Smith    Analysis, volume 31. pages 1662-1694, 1994. These were developed in the context of iterative substructuring preconditioners.
7443b65e785SBarry Smith    They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund.
7453b65e785SBarry Smith    They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example,
7463b65e785SBarry Smith    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Extending theory for domain decomposition algorithms to irregular subdomains. In Ulrich Langer, Marco
7473b65e785SBarry Smith    Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
7483b65e785SBarry Smith    of the 17th International Conference on Domain Decomposition Methods in
7493b65e785SBarry Smith    Science and Engineering, held in Strobl, Austria, July 3-7, 2006, number 60 in
7503b65e785SBarry Smith    Springer-Verlag, Lecture Notes in Computational Science and Engineering, pages 255-261, 2007.
7513b65e785SBarry Smith    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. A family of energy min-
7523b65e785SBarry Smith    imizing coarse spaces for overlapping Schwarz preconditioners. In Ulrich Langer,
7533b65e785SBarry Smith    Marco Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
7543b65e785SBarry Smith    of the 17th International Conference on Domain Decomposition Methods
7553b65e785SBarry Smith    in Science and Engineering, held in Strobl, Austria, July 3-7, 2006, number 60 in
7563b65e785SBarry Smith    Springer-Verlag, Lecture Notes in Computational Science and Engineering, pages 247-254, 2007
7573b65e785SBarry Smith    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Domain decomposition
7583b65e785SBarry Smith    for less regular subdomains: Overlapping Schwarz in two dimensions. SIAM J.
7593b65e785SBarry Smith    Numer. Anal., 46(4):2153-2168, 2008.
7603b65e785SBarry Smith    Clark R. Dohrmann and Olof B. Widlund. An overlapping Schwarz
7613b65e785SBarry Smith    algorithm for almost incompressible elasticity. Technical Report
7623b65e785SBarry Smith    TR2008-912, Department of Computer Science, Courant Institute
7633b65e785SBarry Smith    of Mathematical Sciences, New York University, May 2008. URL:
7643b65e785SBarry Smith    http://cs.nyu.edu/csweb/Research/TechReports/TR2008-912/TR2008-912.pdf
7657233f9f0SBarry Smith 
7667233f9f0SBarry Smith    Options Database: The usual PCMG options are supported, such as -mg_levels_pc_type <type> -mg_coarse_pc_type <type>
7677233f9f0SBarry Smith       -pc_mg_type <type>
7687233f9f0SBarry Smith 
76925a35f6fSSatish Balay    Level: advanced
770174b6946SBarry Smith 
7716c699258SBarry Smith .seealso:  PCMG, PCSetDM(), PCExoticType, PCExoticSetType()
772174b6946SBarry Smith M*/
773174b6946SBarry Smith 
774174b6946SBarry Smith EXTERN_C_BEGIN
775174b6946SBarry Smith #undef __FUNCT__
7767233f9f0SBarry Smith #define __FUNCT__ "PCCreate_Exotic"
7777233f9f0SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Exotic(PC pc)
778174b6946SBarry Smith {
779174b6946SBarry Smith   PetscErrorCode ierr;
7807233f9f0SBarry Smith   PC_Exotic      *ex;
781f3fbd535SBarry Smith   PC_MG          *mg;
782174b6946SBarry Smith 
783174b6946SBarry Smith   PetscFunctionBegin;
784f91d8e95SBarry Smith   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
785f91d8e95SBarry Smith   if (pc->ops->destroy) { ierr =  (*pc->ops->destroy)(pc);CHKERRQ(ierr); pc->data = 0;}
786503cfb0cSBarry Smith   ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr);
787f91d8e95SBarry Smith   ((PetscObject)pc)->type_name = 0;
788f91d8e95SBarry Smith 
789174b6946SBarry Smith   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr);
790174b6946SBarry Smith   ierr = PCMGSetLevels(pc,2,PETSC_NULL);CHKERRQ(ierr);
791174b6946SBarry Smith   ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
79296bdf778SBarry Smith   ierr = PetscNew(PC_Exotic,&ex);CHKERRQ(ierr);\
7937233f9f0SBarry Smith   ex->type = PC_EXOTIC_FACE;
794f3fbd535SBarry Smith   mg = (PC_MG*) pc->data;
79531567311SBarry Smith   mg->innerctx = ex;
7967233f9f0SBarry Smith 
7977233f9f0SBarry Smith 
7987233f9f0SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_Exotic;
7997233f9f0SBarry Smith   pc->ops->view           = PCView_Exotic;
8007233f9f0SBarry Smith   pc->ops->destroy        = PCDestroy_Exotic;
8016c699258SBarry Smith   pc->ops->setup          = PCSetUp_Exotic;
8027233f9f0SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCExoticSetType_C","PCExoticSetType_Exotic",PCExoticSetType_Exotic);CHKERRQ(ierr);
803174b6946SBarry Smith   PetscFunctionReturn(0);
804174b6946SBarry Smith }
805174b6946SBarry Smith EXTERN_C_END
806