xref: /petsc/src/ksp/pc/impls/wb/wb.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
1174b6946SBarry Smith 
2c6db04a5SJed Brown #include <petscdmda.h>   /*I "petscdmda.h" I*/
3af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h>   /*I "petscksp.h" I*/
482c86c8fSBarry Smith #include <petscctable.h>
57233f9f0SBarry Smith 
68e722e36SBarry Smith typedef struct {
78e722e36SBarry Smith   PCExoticType type;
88e722e36SBarry Smith   Mat          P;            /* the constructed interpolation matrix */
9ace3abfcSBarry Smith   PetscBool    directSolve;  /* use direct LU factorization to construct interpolation */
108e722e36SBarry Smith   KSP          ksp;
118e722e36SBarry Smith } PC_Exotic;
12174b6946SBarry Smith 
130a545947SLisandro Dalcin const char *const PCExoticTypes[] = {"face","wirebasket","PCExoticType","PC_Exotic",NULL};
14064c8009SBarry Smith 
15064c8009SBarry Smith /*
16aa219208SBarry Smith       DMDAGetWireBasketInterpolation - Gets the interpolation for a wirebasket based coarse space
17064c8009SBarry Smith 
18064c8009SBarry Smith */
19c0decd05SBarry Smith PetscErrorCode DMDAGetWireBasketInterpolation(PC pc,DM da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P)
20064c8009SBarry Smith {
21064c8009SBarry 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;
2228d20b34SBarry Smith   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[26],*globals,Ng,*IIint,*IIsurf,Nt;
23064c8009SBarry Smith   Mat                    Xint, Xsurf,Xint_tmp;
24064c8009SBarry Smith   IS                     isint,issurf,is,row,col;
25064c8009SBarry Smith   ISLocalToGlobalMapping ltg;
26064c8009SBarry Smith   MPI_Comm               comm;
27064c8009SBarry Smith   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
28064c8009SBarry Smith   MatFactorInfo          info;
29064c8009SBarry Smith   PetscScalar            *xsurf,*xint;
301683a169SBarry Smith   const PetscScalar      *rxint;
318e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
32064c8009SBarry Smith   PetscScalar            tmp;
33064c8009SBarry Smith #endif
34064c8009SBarry Smith   PetscTable             ht;
35064c8009SBarry Smith 
36064c8009SBarry Smith   PetscFunctionBegin;
379566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,&mp,&np,&pp,&dof,NULL,NULL,NULL,NULL,NULL));
3808401ef6SPierre Jolivet   PetscCheck(dof == 1,PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only for single field problems");
3908401ef6SPierre Jolivet   PetscCheck(dim == 3,PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only coded for 3d problems");
409566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,NULL,NULL,NULL,&m,&n,&p));
419566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth));
42064c8009SBarry Smith   istart = istart ? -1 : 0;
43064c8009SBarry Smith   jstart = jstart ? -1 : 0;
44064c8009SBarry Smith   kstart = kstart ? -1 : 0;
45064c8009SBarry Smith 
46064c8009SBarry Smith   /*
47064c8009SBarry Smith     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge)
48064c8009SBarry Smith     to all the local degrees of freedom (this includes the vertices, edges and faces).
49064c8009SBarry Smith 
50064c8009SBarry Smith     Xint are the subset of the interpolation into the interior
51064c8009SBarry Smith 
52064c8009SBarry Smith     Xface are the interpolation onto faces but not into the interior
53064c8009SBarry Smith 
54064c8009SBarry Smith     Xsurf are the interpolation onto the vertices and edges (the surfbasket)
55064c8009SBarry Smith                                       Xint
56064c8009SBarry Smith     Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain
57064c8009SBarry Smith                                       Xsurf
58064c8009SBarry Smith   */
59064c8009SBarry Smith   N     = (m - istart)*(n - jstart)*(p - kstart);
60064c8009SBarry Smith   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
61064c8009SBarry Smith   Nface = 2*((m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart));
62064c8009SBarry Smith   Nwire = 4*((m-2-istart) + (n-2-jstart) + (p-2-kstart)) + 8;
63064c8009SBarry Smith   Nsurf = Nface + Nwire;
649566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(MPI_COMM_SELF,Nint,26,NULL,&Xint));
659566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(MPI_COMM_SELF,Nsurf,26,NULL,&Xsurf));
669566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(Xsurf,&xsurf));
67064c8009SBarry Smith 
68064c8009SBarry Smith   /*
69064c8009SBarry Smith      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
70064c8009SBarry Smith      Xsurf will be all zero (thus making the coarse matrix singular).
71064c8009SBarry Smith   */
7208401ef6SPierre Jolivet   PetscCheck(m-istart >= 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in X direction must be at least 3");
7308401ef6SPierre Jolivet   PetscCheck(n-jstart >= 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in Y direction must be at least 3");
7408401ef6SPierre Jolivet   PetscCheck(p-kstart >= 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in Z direction must be at least 3");
75064c8009SBarry Smith 
76064c8009SBarry Smith   cnt = 0;
772fa5cd67SKarl Rupp 
782fa5cd67SKarl Rupp   xsurf[cnt++] = 1;
792fa5cd67SKarl Rupp   for (i=1; i<m-istart-1; i++) xsurf[cnt++ + Nsurf] = 1;
802fa5cd67SKarl Rupp   xsurf[cnt++ + 2*Nsurf] = 1;
812fa5cd67SKarl Rupp 
822fa5cd67SKarl Rupp   for (j=1; j<n-1-jstart; j++) {
832fa5cd67SKarl Rupp     xsurf[cnt++ + 3*Nsurf] = 1;
842fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 4*Nsurf] = 1;
852fa5cd67SKarl Rupp     xsurf[cnt++ + 5*Nsurf] = 1;
86064c8009SBarry Smith   }
872fa5cd67SKarl Rupp 
882fa5cd67SKarl Rupp   xsurf[cnt++ + 6*Nsurf] = 1;
892fa5cd67SKarl Rupp   for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 7*Nsurf] = 1;
902fa5cd67SKarl Rupp   xsurf[cnt++ + 8*Nsurf] = 1;
912fa5cd67SKarl Rupp 
922fa5cd67SKarl Rupp   for (k=1; k<p-1-kstart; k++) {
932fa5cd67SKarl Rupp     xsurf[cnt++ + 9*Nsurf] = 1;
942fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 10*Nsurf] = 1;
952fa5cd67SKarl Rupp     xsurf[cnt++ + 11*Nsurf] = 1;
962fa5cd67SKarl Rupp 
972fa5cd67SKarl Rupp     for (j=1; j<n-1-jstart; j++) {
982fa5cd67SKarl Rupp       xsurf[cnt++ + 12*Nsurf] = 1;
992fa5cd67SKarl Rupp       /* these are the interior nodes */
1002fa5cd67SKarl Rupp       xsurf[cnt++ + 13*Nsurf] = 1;
1012fa5cd67SKarl Rupp     }
1022fa5cd67SKarl Rupp 
1032fa5cd67SKarl Rupp     xsurf[cnt++ + 14*Nsurf] = 1;
1042fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 15*Nsurf] = 1;
1052fa5cd67SKarl Rupp     xsurf[cnt++ + 16*Nsurf] = 1;
1062fa5cd67SKarl Rupp   }
1072fa5cd67SKarl Rupp 
1082fa5cd67SKarl Rupp   xsurf[cnt++ + 17*Nsurf] = 1;
1092fa5cd67SKarl Rupp   for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 18*Nsurf] = 1;
1102fa5cd67SKarl Rupp   xsurf[cnt++ + 19*Nsurf] = 1;
1112fa5cd67SKarl Rupp 
1122fa5cd67SKarl Rupp   for (j=1;j<n-1-jstart;j++) {
1132fa5cd67SKarl Rupp     xsurf[cnt++ + 20*Nsurf] = 1;
1142fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 21*Nsurf] = 1;
1152fa5cd67SKarl Rupp     xsurf[cnt++ + 22*Nsurf] = 1;
1162fa5cd67SKarl Rupp   }
1172fa5cd67SKarl Rupp 
1182fa5cd67SKarl Rupp   xsurf[cnt++ + 23*Nsurf] = 1;
1192fa5cd67SKarl Rupp   for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 24*Nsurf] = 1;
1202fa5cd67SKarl Rupp   xsurf[cnt++ + 25*Nsurf] = 1;
1212fa5cd67SKarl Rupp 
1228e722e36SBarry Smith   /* interpolations only sum to 1 when using direct solver */
1238e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
124064c8009SBarry Smith   for (i=0; i<Nsurf; i++) {
125064c8009SBarry Smith     tmp = 0.0;
1262fa5cd67SKarl Rupp     for (j=0; j<26; j++) tmp += xsurf[i+j*Nsurf];
12763a3b9bcSJacob Faibussowitsch     PetscCheck(PetscAbsScalar(tmp-1.0) <= 1.e-10,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %" PetscInt_FMT " value %g",i,(double)PetscAbsScalar(tmp));
128064c8009SBarry Smith   }
129064c8009SBarry Smith #endif
1309566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(Xsurf,&xsurf));
1319566063dSJacob Faibussowitsch   /* PetscCall(MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD));*/
132064c8009SBarry Smith 
133064c8009SBarry Smith   /*
134064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
135064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
1367dae84e0SHong Zhang             (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it
137aa219208SBarry Smith              is NOT the local DMDA ordering.)
138064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
139064c8009SBarry Smith   */
140064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
1419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(N,&II,Nint,&Iint,Nsurf,&Isurf));
1429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nint,&IIint,Nsurf,&IIsurf));
143064c8009SBarry Smith   for (k=0; k<p-kstart; k++) {
144064c8009SBarry Smith     for (j=0; j<n-jstart; j++) {
145064c8009SBarry Smith       for (i=0; i<m-istart; i++) {
146064c8009SBarry Smith         II[c++] = i + j*mwidth + k*mwidth*nwidth;
147064c8009SBarry Smith 
148064c8009SBarry Smith         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
149064c8009SBarry Smith           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
150064c8009SBarry Smith           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
151064c8009SBarry Smith         } else {
152064c8009SBarry Smith           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
153064c8009SBarry Smith           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
154064c8009SBarry Smith         }
155064c8009SBarry Smith       }
156064c8009SBarry Smith     }
157064c8009SBarry Smith   }
15808401ef6SPierre Jolivet   PetscCheck(c == N,PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N");
15908401ef6SPierre Jolivet   PetscCheck(cint == Nint,PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint");
16008401ef6SPierre Jolivet   PetscCheck(csurf == Nsurf,PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf");
1619566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(da,&ltg));
1629566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg,N,II,II));
1639566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint));
1649566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf));
1659566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da,&comm));
1669566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is));
1679566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint));
1689566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf));
1699566063dSJacob Faibussowitsch   PetscCall(PetscFree3(II,Iint,Isurf));
170064c8009SBarry Smith 
1719566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder));
172064c8009SBarry Smith   A    = *Aholder;
1739566063dSJacob Faibussowitsch   PetscCall(PetscFree(Aholder));
174064c8009SBarry Smith 
1759566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii));
1769566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais));
1779566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi));
178064c8009SBarry Smith 
179064c8009SBarry Smith   /*
180064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
181064c8009SBarry Smith   */
1829566063dSJacob Faibussowitsch   PetscCall(MatMatMult(Ais,Xsurf,MAT_INITIAL_MATRIX,PETSC_DETERMINE,&Xint_tmp));
1839566063dSJacob Faibussowitsch   PetscCall(MatScale(Xint_tmp,-1.0));
1848e722e36SBarry Smith   if (exotic->directSolve) {
1859566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii));
1869566063dSJacob Faibussowitsch     PetscCall(MatFactorInfoInitialize(&info));
1879566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(Aii,MATORDERINGND,&row,&col));
1889566063dSJacob Faibussowitsch     PetscCall(MatLUFactorSymbolic(iAii,Aii,row,col,&info));
1899566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&row));
1909566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&col));
1919566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(iAii,Aii,&info));
1929566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(iAii,Xint_tmp,Xint));
1939566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&iAii));
1948e722e36SBarry Smith   } else {
1958e722e36SBarry Smith     Vec         b,x;
1968e722e36SBarry Smith     PetscScalar *xint_tmp;
197064c8009SBarry Smith 
1989566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(Xint,&xint));
1999566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,NULL,&x));
2009566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(Xint_tmp,&xint_tmp));
2019566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,NULL,&b));
2029566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(exotic->ksp,Aii,Aii));
2038e722e36SBarry Smith     for (i=0; i<26; i++) {
2049566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(x,xint+i*Nint));
2059566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(b,xint_tmp+i*Nint));
2069566063dSJacob Faibussowitsch       PetscCall(KSPSolve(exotic->ksp,b,x));
2079566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(exotic->ksp,pc,x));
2089566063dSJacob Faibussowitsch       PetscCall(VecResetArray(x));
2099566063dSJacob Faibussowitsch       PetscCall(VecResetArray(b));
2108e722e36SBarry Smith     }
2119566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(Xint,&xint));
2129566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(Xint_tmp,&xint_tmp));
2139566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
2149566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&b));
2158e722e36SBarry Smith   }
2169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xint_tmp));
2178e722e36SBarry Smith 
2188e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
2199566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xint,&rxint));
220064c8009SBarry Smith   for (i=0; i<Nint; i++) {
221064c8009SBarry Smith     tmp = 0.0;
2221683a169SBarry Smith     for (j=0; j<26; j++) tmp += rxint[i+j*Nint];
2232fa5cd67SKarl Rupp 
22463a3b9bcSJacob Faibussowitsch     PetscCheck(PetscAbsScalar(tmp-1.0) <= 1.e-10,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xint interpolation at i %" PetscInt_FMT " value %g",i,(double)PetscAbsScalar(tmp));
225064c8009SBarry Smith   }
2269566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xint,&rxint));
2279566063dSJacob Faibussowitsch   /* PetscCall(MatView(Xint,PETSC_VIEWER_STDOUT_WORLD)); */
228064c8009SBarry Smith #endif
229064c8009SBarry Smith 
230064c8009SBarry Smith   /*         total vertices             total faces                                  total edges */
231064c8009SBarry 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);
232064c8009SBarry Smith 
233064c8009SBarry Smith   /*
234064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
235064c8009SBarry Smith   */
236064c8009SBarry Smith   cnt = 0;
2372fa5cd67SKarl Rupp 
238064c8009SBarry Smith   gl[cnt++] = 0;  { gl[cnt++] = 1;} gl[cnt++] = m-istart-1;
239064c8009SBarry Smith   { gl[cnt++] = mwidth;  { gl[cnt++] = mwidth+1;} gl[cnt++] = mwidth + m-istart-1;}
240064c8009SBarry Smith   gl[cnt++] = mwidth*(n-jstart-1);  { gl[cnt++] = mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*(n-jstart-1) + m-istart-1;
241064c8009SBarry Smith   {
242064c8009SBarry Smith     gl[cnt++] = mwidth*nwidth;  { gl[cnt++] = mwidth*nwidth+1;}  gl[cnt++] = mwidth*nwidth+ m-istart-1;
243064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
244064c8009SBarry 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;
245064c8009SBarry Smith   }
246064c8009SBarry 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;
247064c8009SBarry 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;}
248064c8009SBarry 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;
249064c8009SBarry Smith 
250064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
251064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
2529566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg,26,gl,gl));
253064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
2549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(26*mp*np*pp,&globals));
2559566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(gl,26,MPIU_INT,globals,26,MPIU_INT,PetscObjectComm((PetscObject)da)));
256064c8009SBarry Smith 
257064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
2589566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Aglobal,&Nt,NULL));
2599566063dSJacob Faibussowitsch   PetscCall(PetscTableCreate(Ntotal/3,Nt+1,&ht));
260064c8009SBarry Smith   for (i=0; i<26*mp*np*pp; i++) {
2619566063dSJacob Faibussowitsch     PetscCall(PetscTableAddCount(ht,globals[i]+1));
262064c8009SBarry Smith   }
2639566063dSJacob Faibussowitsch   PetscCall(PetscTableGetCount(ht,&cnt));
26463a3b9bcSJacob Faibussowitsch   PetscCheck(cnt == Ntotal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hash table size %" PetscInt_FMT " not equal to total number coarse grid points %" PetscInt_FMT,cnt,Ntotal);
2659566063dSJacob Faibussowitsch   PetscCall(PetscFree(globals));
266064c8009SBarry Smith   for (i=0; i<26; i++) {
2679566063dSJacob Faibussowitsch     PetscCall(PetscTableFind(ht,gl[i]+1,&gl[i]));
268064c8009SBarry Smith     gl[i]--;
269064c8009SBarry Smith   }
2709566063dSJacob Faibussowitsch   PetscCall(PetscTableDestroy(&ht));
271064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
272064c8009SBarry Smith 
273064c8009SBarry Smith   /* construct global interpolation matrix */
2749566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Aglobal,&Ng,NULL));
275064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
2769566063dSJacob Faibussowitsch     PetscCall(MatCreateAIJ(PetscObjectComm((PetscObject)da),Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,NULL,Nint+Nsurf,NULL,P));
277064c8009SBarry Smith   } else {
2789566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(*P));
279064c8009SBarry Smith   }
2809566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE));
2819566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xint,&rxint));
2829566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*P,Nint,IIint,26,gl,rxint,INSERT_VALUES));
2839566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xint,&rxint));
2849566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xsurf,&rxint));
2859566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*P,Nsurf,IIsurf,26,gl,rxint,INSERT_VALUES));
2869566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xsurf,&rxint));
2879566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY));
2889566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY));
2899566063dSJacob Faibussowitsch   PetscCall(PetscFree2(IIint,IIsurf));
290064c8009SBarry Smith 
2918e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
292064c8009SBarry Smith   {
293064c8009SBarry Smith     Vec         x,y;
294064c8009SBarry Smith     PetscScalar *yy;
2959566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da),Ng,PETSC_DETERMINE,&y));
2969566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da),PETSC_DETERMINE,Ntotal,&x));
2979566063dSJacob Faibussowitsch     PetscCall(VecSet(x,1.0));
2989566063dSJacob Faibussowitsch     PetscCall(MatMult(*P,x,y));
2999566063dSJacob Faibussowitsch     PetscCall(VecGetArray(y,&yy));
300064c8009SBarry Smith     for (i=0; i<Ng; i++) {
30163a3b9bcSJacob Faibussowitsch       PetscCheck(PetscAbsScalar(yy[i]-1.0) <= 1.e-10,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong p interpolation at i %" PetscInt_FMT " value %g",i,(double)PetscAbsScalar(yy[i]));
302064c8009SBarry Smith     }
3039566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(y,&yy));
3049566063dSJacob Faibussowitsch     PetscCall(VecDestroy(x));
3059566063dSJacob Faibussowitsch     PetscCall(VecDestroy(y));
306064c8009SBarry Smith   }
307064c8009SBarry Smith #endif
308064c8009SBarry Smith 
3099566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Aii));
3109566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ais));
3119566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Asi));
3129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3139566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
3149566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isint));
3159566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&issurf));
3169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xint));
3179566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xsurf));
318064c8009SBarry Smith   PetscFunctionReturn(0);
319064c8009SBarry Smith }
320064c8009SBarry Smith 
321064c8009SBarry Smith /*
322aa219208SBarry Smith       DMDAGetFaceInterpolation - Gets the interpolation for a face based coarse space
323064c8009SBarry Smith 
324064c8009SBarry Smith */
325c0decd05SBarry Smith PetscErrorCode DMDAGetFaceInterpolation(PC pc,DM da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P)
326064c8009SBarry Smith {
327064c8009SBarry 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;
32828d20b34SBarry Smith   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[6],*globals,Ng,*IIint,*IIsurf,Nt;
329064c8009SBarry Smith   Mat                    Xint, Xsurf,Xint_tmp;
330064c8009SBarry Smith   IS                     isint,issurf,is,row,col;
331064c8009SBarry Smith   ISLocalToGlobalMapping ltg;
332064c8009SBarry Smith   MPI_Comm               comm;
333064c8009SBarry Smith   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
334064c8009SBarry Smith   MatFactorInfo          info;
335064c8009SBarry Smith   PetscScalar            *xsurf,*xint;
3361683a169SBarry Smith   const PetscScalar      *rxint;
337064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
338064c8009SBarry Smith   PetscScalar            tmp;
339064c8009SBarry Smith #endif
340064c8009SBarry Smith   PetscTable             ht;
341064c8009SBarry Smith 
342064c8009SBarry Smith   PetscFunctionBegin;
3439566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,&mp,&np,&pp,&dof,NULL,NULL,NULL,NULL,NULL));
34408401ef6SPierre Jolivet   PetscCheck(dof == 1,PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only for single field problems");
34508401ef6SPierre Jolivet   PetscCheck(dim == 3,PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only coded for 3d problems");
3469566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,NULL,NULL,NULL,&m,&n,&p));
3479566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth));
348064c8009SBarry Smith   istart = istart ? -1 : 0;
349064c8009SBarry Smith   jstart = jstart ? -1 : 0;
350064c8009SBarry Smith   kstart = kstart ? -1 : 0;
351064c8009SBarry Smith 
352064c8009SBarry Smith   /*
353064c8009SBarry Smith     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge)
354064c8009SBarry Smith     to all the local degrees of freedom (this includes the vertices, edges and faces).
355064c8009SBarry Smith 
356064c8009SBarry Smith     Xint are the subset of the interpolation into the interior
357064c8009SBarry Smith 
358064c8009SBarry Smith     Xface are the interpolation onto faces but not into the interior
359064c8009SBarry Smith 
360064c8009SBarry Smith     Xsurf are the interpolation onto the vertices and edges (the surfbasket)
361064c8009SBarry Smith                                       Xint
362064c8009SBarry Smith     Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain
363064c8009SBarry Smith                                       Xsurf
364064c8009SBarry Smith   */
365064c8009SBarry Smith   N     = (m - istart)*(n - jstart)*(p - kstart);
366064c8009SBarry Smith   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
367064c8009SBarry Smith   Nface = 2*((m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart));
368064c8009SBarry Smith   Nwire = 4*((m-2-istart) + (n-2-jstart) + (p-2-kstart)) + 8;
369064c8009SBarry Smith   Nsurf = Nface + Nwire;
3709566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(MPI_COMM_SELF,Nint,6,NULL,&Xint));
3719566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(MPI_COMM_SELF,Nsurf,6,NULL,&Xsurf));
3729566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(Xsurf,&xsurf));
373064c8009SBarry Smith 
374064c8009SBarry Smith   /*
375064c8009SBarry Smith      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
376064c8009SBarry Smith      Xsurf will be all zero (thus making the coarse matrix singular).
377064c8009SBarry Smith   */
37808401ef6SPierre Jolivet   PetscCheck(m-istart >= 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in X direction must be at least 3");
37908401ef6SPierre Jolivet   PetscCheck(n-jstart >= 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in Y direction must be at least 3");
38008401ef6SPierre Jolivet   PetscCheck(p-kstart >= 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in Z direction must be at least 3");
381064c8009SBarry Smith 
382064c8009SBarry Smith   cnt = 0;
3832fa5cd67SKarl Rupp   for (j=1; j<n-1-jstart; j++) {
3842fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 0*Nsurf] = 1;
385064c8009SBarry Smith   }
3862fa5cd67SKarl Rupp 
3872fa5cd67SKarl Rupp   for (k=1; k<p-1-kstart; k++) {
3882fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 1*Nsurf] = 1;
3892fa5cd67SKarl Rupp     for (j=1; j<n-1-jstart; j++) {
3902fa5cd67SKarl Rupp       xsurf[cnt++ + 2*Nsurf] = 1;
3912fa5cd67SKarl Rupp       /* these are the interior nodes */
3922fa5cd67SKarl Rupp       xsurf[cnt++ + 3*Nsurf] = 1;
3932fa5cd67SKarl Rupp     }
3942fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 4*Nsurf] = 1;
3952fa5cd67SKarl Rupp   }
3962fa5cd67SKarl Rupp   for (j=1;j<n-1-jstart;j++) {
3972fa5cd67SKarl Rupp     for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 5*Nsurf] = 1;
3982fa5cd67SKarl Rupp   }
399064c8009SBarry Smith 
400064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
401064c8009SBarry Smith   for (i=0; i<Nsurf; i++) {
402064c8009SBarry Smith     tmp = 0.0;
4032fa5cd67SKarl Rupp     for (j=0; j<6; j++) tmp += xsurf[i+j*Nsurf];
4042fa5cd67SKarl Rupp 
40563a3b9bcSJacob Faibussowitsch     PetscCheck(PetscAbsScalar(tmp-1.0) <= 1.e-10,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %" PetscInt_FMT " value %g",i,(double)PetscAbsScalar(tmp));
406064c8009SBarry Smith   }
407064c8009SBarry Smith #endif
4089566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(Xsurf,&xsurf));
4099566063dSJacob Faibussowitsch   /* PetscCall(MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD));*/
410064c8009SBarry Smith 
411064c8009SBarry Smith   /*
412064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
413064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
4147dae84e0SHong Zhang             (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it
415aa219208SBarry Smith              is NOT the local DMDA ordering.)
416064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
417064c8009SBarry Smith   */
418064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
4199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(N,&II,Nint,&Iint,Nsurf,&Isurf));
4209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nint,&IIint,Nsurf,&IIsurf));
421064c8009SBarry Smith   for (k=0; k<p-kstart; k++) {
422064c8009SBarry Smith     for (j=0; j<n-jstart; j++) {
423064c8009SBarry Smith       for (i=0; i<m-istart; i++) {
424064c8009SBarry Smith         II[c++] = i + j*mwidth + k*mwidth*nwidth;
425064c8009SBarry Smith 
426064c8009SBarry Smith         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
427064c8009SBarry Smith           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
428064c8009SBarry Smith           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
429064c8009SBarry Smith         } else {
430064c8009SBarry Smith           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
431064c8009SBarry Smith           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
432064c8009SBarry Smith         }
433064c8009SBarry Smith       }
434064c8009SBarry Smith     }
435064c8009SBarry Smith   }
43608401ef6SPierre Jolivet   PetscCheck(c == N,PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N");
43708401ef6SPierre Jolivet   PetscCheck(cint == Nint,PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint");
43808401ef6SPierre Jolivet   PetscCheck(csurf == Nsurf,PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf");
4399566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(da,&ltg));
4409566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg,N,II,II));
4419566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint));
4429566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf));
4439566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da,&comm));
4449566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is));
4459566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint));
4469566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf));
4479566063dSJacob Faibussowitsch   PetscCall(PetscFree3(II,Iint,Isurf));
448064c8009SBarry Smith 
4499566063dSJacob Faibussowitsch   PetscCall(ISSort(is));
4509566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder));
451064c8009SBarry Smith   A    = *Aholder;
4529566063dSJacob Faibussowitsch   PetscCall(PetscFree(Aholder));
453064c8009SBarry Smith 
4549566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii));
4559566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais));
4569566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi));
457064c8009SBarry Smith 
458064c8009SBarry Smith   /*
459064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
460064c8009SBarry Smith   */
4619566063dSJacob Faibussowitsch   PetscCall(MatMatMult(Ais,Xsurf,MAT_INITIAL_MATRIX,PETSC_DETERMINE,&Xint_tmp));
4629566063dSJacob Faibussowitsch   PetscCall(MatScale(Xint_tmp,-1.0));
463064c8009SBarry Smith 
4648e722e36SBarry Smith   if (exotic->directSolve) {
4659566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii));
4669566063dSJacob Faibussowitsch     PetscCall(MatFactorInfoInitialize(&info));
4679566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(Aii,MATORDERINGND,&row,&col));
4689566063dSJacob Faibussowitsch     PetscCall(MatLUFactorSymbolic(iAii,Aii,row,col,&info));
4699566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&row));
4709566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&col));
4719566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(iAii,Aii,&info));
4729566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(iAii,Xint_tmp,Xint));
4739566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&iAii));
474064c8009SBarry Smith   } else {
475064c8009SBarry Smith     Vec         b,x;
476064c8009SBarry Smith     PetscScalar *xint_tmp;
477064c8009SBarry Smith 
4789566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(Xint,&xint));
4799566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,NULL,&x));
4809566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(Xint_tmp,&xint_tmp));
4819566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,NULL,&b));
4829566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(exotic->ksp,Aii,Aii));
483064c8009SBarry Smith     for (i=0; i<6; i++) {
4849566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(x,xint+i*Nint));
4859566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(b,xint_tmp+i*Nint));
4869566063dSJacob Faibussowitsch       PetscCall(KSPSolve(exotic->ksp,b,x));
4879566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(exotic->ksp,pc,x));
4889566063dSJacob Faibussowitsch       PetscCall(VecResetArray(x));
4899566063dSJacob Faibussowitsch       PetscCall(VecResetArray(b));
490064c8009SBarry Smith     }
4919566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(Xint,&xint));
4929566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(Xint_tmp,&xint_tmp));
4939566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
4949566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&b));
495064c8009SBarry Smith   }
4969566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xint_tmp));
497064c8009SBarry Smith 
498064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
4999566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xint,&rxint));
500064c8009SBarry Smith   for (i=0; i<Nint; i++) {
501064c8009SBarry Smith     tmp = 0.0;
5021683a169SBarry Smith     for (j=0; j<6; j++) tmp += rxint[i+j*Nint];
5032fa5cd67SKarl Rupp 
50463a3b9bcSJacob Faibussowitsch     PetscCheck(PetscAbsScalar(tmp-1.0) <= 1.e-10,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xint interpolation at i %" PetscInt_FMT " value %g",i,(double)PetscAbsScalar(tmp));
505064c8009SBarry Smith   }
5069566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xint,&rxint));
5079566063dSJacob Faibussowitsch   /* PetscCall(MatView(Xint,PETSC_VIEWER_STDOUT_WORLD)); */
508064c8009SBarry Smith #endif
509064c8009SBarry Smith 
510064c8009SBarry Smith   /*         total faces    */
511064c8009SBarry Smith   Ntotal =  mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1);
512064c8009SBarry Smith 
513064c8009SBarry Smith   /*
514064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
515064c8009SBarry Smith   */
516064c8009SBarry Smith   cnt = 0;
517064c8009SBarry Smith   { gl[cnt++] = mwidth+1;}
518064c8009SBarry Smith   {
519064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth+1;}
520064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
521064c8009SBarry Smith     { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;}
522064c8009SBarry Smith   }
523064c8009SBarry Smith   { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;}
524064c8009SBarry Smith 
525064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
526064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
5279566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg,6,gl,gl));
528064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
5299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(6*mp*np*pp,&globals));
5309566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(gl,6,MPIU_INT,globals,6,MPIU_INT,PetscObjectComm((PetscObject)da)));
531064c8009SBarry Smith 
532064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
5339566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Aglobal,&Nt,NULL));
5349566063dSJacob Faibussowitsch   PetscCall(PetscTableCreate(Ntotal/3,Nt+1,&ht));
535064c8009SBarry Smith   for (i=0; i<6*mp*np*pp; i++) {
5369566063dSJacob Faibussowitsch     PetscCall(PetscTableAddCount(ht,globals[i]+1));
537064c8009SBarry Smith   }
5389566063dSJacob Faibussowitsch   PetscCall(PetscTableGetCount(ht,&cnt));
53963a3b9bcSJacob Faibussowitsch   PetscCheck(cnt == Ntotal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hash table size %" PetscInt_FMT " not equal to total number coarse grid points %" PetscInt_FMT,cnt,Ntotal);
5409566063dSJacob Faibussowitsch   PetscCall(PetscFree(globals));
541064c8009SBarry Smith   for (i=0; i<6; i++) {
5429566063dSJacob Faibussowitsch     PetscCall(PetscTableFind(ht,gl[i]+1,&gl[i]));
543064c8009SBarry Smith     gl[i]--;
544064c8009SBarry Smith   }
5459566063dSJacob Faibussowitsch   PetscCall(PetscTableDestroy(&ht));
546064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
547064c8009SBarry Smith 
548064c8009SBarry Smith   /* construct global interpolation matrix */
5499566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Aglobal,&Ng,NULL));
550064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
5519566063dSJacob Faibussowitsch     PetscCall(MatCreateAIJ(PetscObjectComm((PetscObject)da),Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,NULL,Nint,NULL,P));
552064c8009SBarry Smith   } else {
5539566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(*P));
554064c8009SBarry Smith   }
5559566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE));
5569566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xint,&rxint));
5579566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*P,Nint,IIint,6,gl,rxint,INSERT_VALUES));
5589566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xint,&rxint));
5599566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xsurf,&rxint));
5609566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*P,Nsurf,IIsurf,6,gl,rxint,INSERT_VALUES));
5619566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xsurf,&rxint));
5629566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY));
5639566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY));
5649566063dSJacob Faibussowitsch   PetscCall(PetscFree2(IIint,IIsurf));
565064c8009SBarry Smith 
566064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
567064c8009SBarry Smith   {
568064c8009SBarry Smith     Vec         x,y;
569064c8009SBarry Smith     PetscScalar *yy;
5709566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da),Ng,PETSC_DETERMINE,&y));
5719566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da),PETSC_DETERMINE,Ntotal,&x));
5729566063dSJacob Faibussowitsch     PetscCall(VecSet(x,1.0));
5739566063dSJacob Faibussowitsch     PetscCall(MatMult(*P,x,y));
5749566063dSJacob Faibussowitsch     PetscCall(VecGetArray(y,&yy));
575064c8009SBarry Smith     for (i=0; i<Ng; i++) {
57663a3b9bcSJacob Faibussowitsch       PetscCheck(PetscAbsScalar(yy[i]-1.0) <= 1.e-10,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong p interpolation at i %" PetscInt_FMT " value %g",i,(double)PetscAbsScalar(yy[i]));
577064c8009SBarry Smith     }
5789566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(y,&yy));
5799566063dSJacob Faibussowitsch     PetscCall(VecDestroy(x));
5809566063dSJacob Faibussowitsch     PetscCall(VecDestroy(y));
581064c8009SBarry Smith   }
582064c8009SBarry Smith #endif
583064c8009SBarry Smith 
5849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Aii));
5859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ais));
5869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Asi));
5879566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
5889566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
5899566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isint));
5909566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&issurf));
5919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xint));
5929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xsurf));
593064c8009SBarry Smith   PetscFunctionReturn(0);
594064c8009SBarry Smith }
595174b6946SBarry Smith 
5967233f9f0SBarry Smith /*@
5977233f9f0SBarry Smith    PCExoticSetType - Sets the type of coarse grid interpolation to use
5987233f9f0SBarry Smith 
5993f9fe445SBarry Smith    Logically Collective on PC
6007233f9f0SBarry Smith 
6017233f9f0SBarry Smith    Input Parameters:
6027233f9f0SBarry Smith +  pc - the preconditioner context
6037233f9f0SBarry Smith -  type - either PC_EXOTIC_FACE or PC_EXOTIC_WIREBASKET (defaults to face)
6047233f9f0SBarry Smith 
60595452b02SPatrick Sanan    Notes:
60695452b02SPatrick Sanan     The face based interpolation has 1 degree of freedom per face and ignores the
607563e08c6SBarry Smith      edge and vertex values completely in the coarse problem. For any seven point
608563e08c6SBarry Smith      stencil the interpolation of a constant on all faces into the interior is that constant.
609563e08c6SBarry Smith 
610563e08c6SBarry Smith      The wirebasket interpolation has 1 degree of freedom per vertex, per edge and
611563e08c6SBarry Smith      per face. A constant on the subdomain boundary is interpolated as that constant
612563e08c6SBarry Smith      in the interior of the domain.
613563e08c6SBarry Smith 
614563e08c6SBarry Smith      The coarse grid matrix is obtained via the Galerkin computation A_c = R A R^T, hence
615563e08c6SBarry Smith      if A is nonsingular A_c is also nonsingular.
616563e08c6SBarry Smith 
617563e08c6SBarry Smith      Both interpolations are suitable for only scalar problems.
618563e08c6SBarry Smith 
6197233f9f0SBarry Smith    Level: intermediate
6207233f9f0SBarry Smith 
621db781477SPatrick Sanan .seealso: `PCEXOTIC`, `PCExoticType()`
6227233f9f0SBarry Smith @*/
6237087cfbeSBarry Smith PetscErrorCode  PCExoticSetType(PC pc,PCExoticType type)
6247233f9f0SBarry Smith {
6257233f9f0SBarry Smith   PetscFunctionBegin;
6260700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
627c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
628cac4c232SBarry Smith   PetscTryMethod(pc,"PCExoticSetType_C",(PC,PCExoticType),(pc,type));
6297233f9f0SBarry Smith   PetscFunctionReturn(0);
6307233f9f0SBarry Smith }
6317233f9f0SBarry Smith 
632f7a08781SBarry Smith static PetscErrorCode  PCExoticSetType_Exotic(PC pc,PCExoticType type)
6337233f9f0SBarry Smith {
634f3fbd535SBarry Smith   PC_MG     *mg  = (PC_MG*)pc->data;
63531567311SBarry Smith   PC_Exotic *ctx = (PC_Exotic*) mg->innerctx;
6367233f9f0SBarry Smith 
6377233f9f0SBarry Smith   PetscFunctionBegin;
6387233f9f0SBarry Smith   ctx->type = type;
6397233f9f0SBarry Smith   PetscFunctionReturn(0);
6407233f9f0SBarry Smith }
6417233f9f0SBarry Smith 
6427233f9f0SBarry Smith PetscErrorCode PCSetUp_Exotic(PC pc)
643174b6946SBarry Smith {
64496bdf778SBarry Smith   Mat            A;
645f3fbd535SBarry Smith   PC_MG          *mg   = (PC_MG*)pc->data;
64631567311SBarry Smith   PC_Exotic      *ex   = (PC_Exotic*) mg->innerctx;
64796bdf778SBarry Smith   MatReuse       reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX;
648174b6946SBarry Smith 
649174b6946SBarry Smith   PetscFunctionBegin;
65028b400f6SJacob Faibussowitsch   PetscCheck(pc->dm,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Need to call PCSetDM() before using this PC");
6519566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc,NULL,&A));
6527233f9f0SBarry Smith   if (ex->type == PC_EXOTIC_FACE) {
6539566063dSJacob Faibussowitsch     PetscCall(DMDAGetFaceInterpolation(pc,pc->dm,ex,A,reuse,&ex->P));
6547233f9f0SBarry Smith   } else if (ex->type == PC_EXOTIC_WIREBASKET) {
6559566063dSJacob Faibussowitsch     PetscCall(DMDAGetWireBasketInterpolation(pc,pc->dm,ex,A,reuse,&ex->P));
65698921bdaSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unknown exotic coarse space %d",ex->type);
6579566063dSJacob Faibussowitsch   PetscCall(PCMGSetInterpolation(pc,1,ex->P));
658d0660788SBarry Smith   /* if PC has attached DM we must remove it or the PCMG will use it to compute incorrect sized vectors and interpolations */
6599566063dSJacob Faibussowitsch   PetscCall(PCSetDM(pc,NULL));
6609566063dSJacob Faibussowitsch   PetscCall(PCSetUp_MG(pc));
661174b6946SBarry Smith   PetscFunctionReturn(0);
662174b6946SBarry Smith }
663174b6946SBarry Smith 
6647233f9f0SBarry Smith PetscErrorCode PCDestroy_Exotic(PC pc)
665174b6946SBarry Smith {
666f3fbd535SBarry Smith   PC_MG          *mg  = (PC_MG*)pc->data;
66731567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
668174b6946SBarry Smith 
669174b6946SBarry Smith   PetscFunctionBegin;
6709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->P));
6719566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ctx->ksp));
6729566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
673*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCExoticSetType_C",NULL));
6749566063dSJacob Faibussowitsch   PetscCall(PCDestroy_MG(pc));
675174b6946SBarry Smith   PetscFunctionReturn(0);
676174b6946SBarry Smith }
677174b6946SBarry Smith 
6787233f9f0SBarry Smith PetscErrorCode PCView_Exotic(PC pc,PetscViewer viewer)
6797233f9f0SBarry Smith {
680f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
681ace3abfcSBarry Smith   PetscBool      iascii;
68231567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
6837233f9f0SBarry Smith 
6847233f9f0SBarry Smith   PetscFunctionBegin;
6859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
6867233f9f0SBarry Smith   if (iascii) {
6879566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"    Exotic type = %s\n",PCExoticTypes[ctx->type]));
6888e722e36SBarry Smith     if (ctx->directSolve) {
6899566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"      Using direct solver to construct interpolation\n"));
6908e722e36SBarry Smith     } else {
6918e722e36SBarry Smith       PetscViewer sviewer;
6928e722e36SBarry Smith       PetscMPIInt rank;
6938e722e36SBarry Smith 
6949566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"      Using iterative solver to construct interpolation\n"));
6959566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
6969566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));  /* should not need to push this twice? */
6979566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
6989566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank));
699dd400576SPatrick Sanan       if (rank == 0) {
7009566063dSJacob Faibussowitsch         PetscCall(KSPView(ctx->ksp,sviewer));
7018e722e36SBarry Smith       }
7029566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
7039566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
7049566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
7058e722e36SBarry Smith     }
7067233f9f0SBarry Smith   }
7079566063dSJacob Faibussowitsch   PetscCall(PCView_MG(pc,viewer));
7087233f9f0SBarry Smith   PetscFunctionReturn(0);
7097233f9f0SBarry Smith }
7107233f9f0SBarry Smith 
7114416b707SBarry Smith PetscErrorCode PCSetFromOptions_Exotic(PetscOptionItems *PetscOptionsObject,PC pc)
7127233f9f0SBarry Smith {
713ace3abfcSBarry Smith   PetscBool      flg;
714f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
7157233f9f0SBarry Smith   PCExoticType   mgctype;
71631567311SBarry Smith   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;
7177233f9f0SBarry Smith 
7187233f9f0SBarry Smith   PetscFunctionBegin;
719d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"Exotic coarse space options");
7209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_exotic_type","face or wirebasket","PCExoticSetType",PCExoticTypes,(PetscEnum)ctx->type,(PetscEnum*)&mgctype,&flg));
7217233f9f0SBarry Smith   if (flg) {
7229566063dSJacob Faibussowitsch     PetscCall(PCExoticSetType(pc,mgctype));
7237233f9f0SBarry Smith   }
7249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_exotic_direct_solver","use direct solver to construct interpolation","None",ctx->directSolve,&ctx->directSolve,NULL));
7258e722e36SBarry Smith   if (!ctx->directSolve) {
7268e722e36SBarry Smith     if (!ctx->ksp) {
7278e722e36SBarry Smith       const char *prefix;
7289566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PETSC_COMM_SELF,&ctx->ksp));
7299566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(ctx->ksp,pc->erroriffailure));
7309566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)pc,1));
7319566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ctx->ksp));
7329566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc,&prefix));
7339566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(ctx->ksp,prefix));
7349566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(ctx->ksp,"exotic_"));
7358e722e36SBarry Smith     }
7369566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ctx->ksp));
7378e722e36SBarry Smith   }
738d0609cedSBarry Smith   PetscOptionsHeadEnd();
7397233f9f0SBarry Smith   PetscFunctionReturn(0);
7407233f9f0SBarry Smith }
7417233f9f0SBarry Smith 
742174b6946SBarry Smith /*MC
7437233f9f0SBarry Smith      PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces
744174b6946SBarry Smith 
7457233f9f0SBarry Smith      This uses the PCMG infrastructure restricted to two levels and the face and wirebasket based coarse
74624c3aa18SBarry Smith    grid spaces.
74724c3aa18SBarry Smith 
74895452b02SPatrick Sanan    Notes:
74995452b02SPatrick Sanan     By default this uses GMRES on the fine grid smoother so this should be used with KSPFGMRES or the smoother changed to not use GMRES
75024c3aa18SBarry Smith 
75196a0c994SBarry Smith    References:
752606c0280SSatish Balay +  * - These coarse grid spaces originate in the work of Bramble, Pasciak  and Schatz, "The Construction
75396a0c994SBarry Smith    of Preconditioners for Elliptic Problems by Substructing IV", Mathematics of Computation, volume 53, 1989.
754606c0280SSatish Balay .  * - They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith,
7554f02bc6aSBarry Smith    New York University, 1990.
756606c0280SSatish Balay .  * - They were then explored in great detail in Dryja, Smith, Widlund, "Schwarz Analysis
7573b65e785SBarry Smith    of Iterative Substructuring Methods for Elliptic Problems in Three Dimensions, SIAM Journal on Numerical
75896a0c994SBarry Smith    Analysis, volume 31. 1994. These were developed in the context of iterative substructuring preconditioners.
759606c0280SSatish Balay .  * - They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund.
7603b65e785SBarry Smith    They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example,
7613b65e785SBarry Smith    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Extending theory for domain decomposition algorithms to irregular subdomains. In Ulrich Langer, Marco
7623b65e785SBarry Smith    Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
7633b65e785SBarry Smith    of the 17th International Conference on Domain Decomposition Methods in
76496a0c994SBarry Smith    Science and Engineering, held in Strobl, Austria, 2006, number 60 in
76596a0c994SBarry Smith    Springer Verlag, Lecture Notes in Computational Science and Engineering, 2007.
766606c0280SSatish Balay .  * -  Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. A family of energy minimizing coarse spaces for overlapping Schwarz preconditioners. In Ulrich Langer,
7673b65e785SBarry Smith    Marco Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
7683b65e785SBarry Smith    of the 17th International Conference on Domain Decomposition Methods
76996a0c994SBarry Smith    in Science and Engineering, held in Strobl, Austria, 2006, number 60 in
77096a0c994SBarry Smith    Springer Verlag, Lecture Notes in Computational Science and Engineering, 2007
771606c0280SSatish Balay .  * - Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Domain decomposition
7723b65e785SBarry Smith    for less regular subdomains: Overlapping Schwarz in two dimensions. SIAM J.
77396a0c994SBarry Smith    Numer. Anal., 46(4), 2008.
774606c0280SSatish Balay -  * - Clark R. Dohrmann and Olof B. Widlund. An overlapping Schwarz
7753b65e785SBarry Smith    algorithm for almost incompressible elasticity. Technical Report
77696a0c994SBarry Smith    TR2008 912, Department of Computer Science, Courant Institute
7773b65e785SBarry Smith    of Mathematical Sciences, New York University, May 2008. URL:
7787233f9f0SBarry Smith 
7797233f9f0SBarry Smith    Options Database: The usual PCMG options are supported, such as -mg_levels_pc_type <type> -mg_coarse_pc_type <type>
7807233f9f0SBarry Smith       -pc_mg_type <type>
7817233f9f0SBarry Smith 
78225a35f6fSSatish Balay    Level: advanced
783174b6946SBarry Smith 
784db781477SPatrick Sanan .seealso: `PCMG`, `PCSetDM()`, `PCExoticType`, `PCExoticSetType()`
785174b6946SBarry Smith M*/
786174b6946SBarry Smith 
7878cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Exotic(PC pc)
788174b6946SBarry Smith {
7897233f9f0SBarry Smith   PC_Exotic      *ex;
790f3fbd535SBarry Smith   PC_MG          *mg;
791174b6946SBarry Smith 
792174b6946SBarry Smith   PetscFunctionBegin;
793f91d8e95SBarry Smith   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
7942fa5cd67SKarl Rupp   if (pc->ops->destroy) {
7959566063dSJacob Faibussowitsch     PetscCall((*pc->ops->destroy)(pc));
7960a545947SLisandro Dalcin     pc->data = NULL;
7972fa5cd67SKarl Rupp   }
7989566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PetscObject)pc)->type_name));
7990a545947SLisandro Dalcin   ((PetscObject)pc)->type_name = NULL;
800f91d8e95SBarry Smith 
8019566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc,PCMG));
8029566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pc,2,NULL));
8039566063dSJacob Faibussowitsch   PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_PMAT));
8049566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ex)); \
8057233f9f0SBarry Smith   ex->type     = PC_EXOTIC_FACE;
806f3fbd535SBarry Smith   mg           = (PC_MG*) pc->data;
80731567311SBarry Smith   mg->innerctx = ex;
8087233f9f0SBarry Smith 
8097233f9f0SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_Exotic;
8107233f9f0SBarry Smith   pc->ops->view           = PCView_Exotic;
8117233f9f0SBarry Smith   pc->ops->destroy        = PCDestroy_Exotic;
8126c699258SBarry Smith   pc->ops->setup          = PCSetUp_Exotic;
8132fa5cd67SKarl Rupp 
8149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCExoticSetType_C",PCExoticSetType_Exotic));
815174b6946SBarry Smith   PetscFunctionReturn(0);
816174b6946SBarry Smith }
817