xref: /petsc/src/ksp/pc/impls/wb/wb.c (revision fbccb6d4aee615a408613f07f0a20245ae912c37)
1c6db04a5SJed Brown #include <petscdmda.h>              /*I "petscdmda.h" I*/
2af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/
3eec179cfSJacob Faibussowitsch #include <petsc/private/hashmapi.h>
47233f9f0SBarry Smith 
58e722e36SBarry Smith typedef struct {
68e722e36SBarry Smith   PCExoticType type;
78e722e36SBarry Smith   Mat          P;           /* the constructed interpolation matrix */
8ace3abfcSBarry Smith   PetscBool    directSolve; /* use direct LU factorization to construct interpolation */
98e722e36SBarry Smith   KSP          ksp;
108e722e36SBarry Smith } PC_Exotic;
11174b6946SBarry Smith 
120a545947SLisandro Dalcin const char *const PCExoticTypes[] = {"face", "wirebasket", "PCExoticType", "PC_Exotic", NULL};
13064c8009SBarry Smith 
14064c8009SBarry Smith /*
15aa219208SBarry Smith       DMDAGetWireBasketInterpolation - Gets the interpolation for a wirebasket based coarse space
16064c8009SBarry Smith 
17064c8009SBarry Smith */
1866976f2fSJacob Faibussowitsch static PetscErrorCode DMDAGetWireBasketInterpolation(PC pc, DM da, PC_Exotic *exotic, Mat Aglobal, MatReuse reuse, Mat *P)
19d71ae5a4SJacob Faibussowitsch {
20064c8009SBarry 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;
21eec179cfSJacob Faibussowitsch   PetscInt               mwidth, nwidth, pwidth, cnt, mp, np, pp, Ntotal, gl[26], *globals, Ng, *IIint, *IIsurf;
22064c8009SBarry Smith   Mat                    Xint, Xsurf, Xint_tmp;
23064c8009SBarry Smith   IS                     isint, issurf, is, row, col;
24064c8009SBarry Smith   ISLocalToGlobalMapping ltg;
25064c8009SBarry Smith   MPI_Comm               comm;
26064c8009SBarry Smith   Mat                    A, Aii, Ais, Asi, *Aholder, iAii;
27064c8009SBarry Smith   MatFactorInfo          info;
28064c8009SBarry Smith   PetscScalar           *xsurf, *xint;
291683a169SBarry Smith   const PetscScalar     *rxint;
308e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
31064c8009SBarry Smith   PetscScalar tmp;
32064c8009SBarry Smith #endif
33eec179cfSJacob Faibussowitsch   PetscHMapI ht = NULL;
34064c8009SBarry Smith 
35064c8009SBarry Smith   PetscFunctionBegin;
369566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, &mp, &np, &pp, &dof, NULL, NULL, NULL, NULL, NULL));
3708401ef6SPierre Jolivet   PetscCheck(dof == 1, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only for single field problems");
3808401ef6SPierre Jolivet   PetscCheck(dim == 3, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only coded for 3d problems");
399566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, NULL, NULL, NULL, &m, &n, &p));
409566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &istart, &jstart, &kstart, &mwidth, &nwidth, &pwidth));
41064c8009SBarry Smith   istart = istart ? -1 : 0;
42064c8009SBarry Smith   jstart = jstart ? -1 : 0;
43064c8009SBarry Smith   kstart = kstart ? -1 : 0;
44064c8009SBarry Smith 
45064c8009SBarry Smith   /*
46064c8009SBarry Smith     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge)
47064c8009SBarry Smith     to all the local degrees of freedom (this includes the vertices, edges and faces).
48064c8009SBarry Smith 
49064c8009SBarry Smith     Xint are the subset of the interpolation into the interior
50064c8009SBarry Smith 
51064c8009SBarry Smith     Xface are the interpolation onto faces but not into the interior
52064c8009SBarry Smith 
53064c8009SBarry Smith     Xsurf are the interpolation onto the vertices and edges (the surfbasket)
54064c8009SBarry Smith                                       Xint
55064c8009SBarry Smith     Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain
56064c8009SBarry Smith                                       Xsurf
57064c8009SBarry Smith   */
58064c8009SBarry Smith   N     = (m - istart) * (n - jstart) * (p - kstart);
59064c8009SBarry Smith   Nint  = (m - 2 - istart) * (n - 2 - jstart) * (p - 2 - kstart);
60064c8009SBarry Smith   Nface = 2 * ((m - 2 - istart) * (n - 2 - jstart) + (m - 2 - istart) * (p - 2 - kstart) + (n - 2 - jstart) * (p - 2 - kstart));
61064c8009SBarry Smith   Nwire = 4 * ((m - 2 - istart) + (n - 2 - jstart) + (p - 2 - kstart)) + 8;
62064c8009SBarry Smith   Nsurf = Nface + Nwire;
639566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nint, 26, NULL, &Xint));
649566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nsurf, 26, NULL, &Xsurf));
659566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(Xsurf, &xsurf));
66064c8009SBarry Smith 
67064c8009SBarry Smith   /*
68064c8009SBarry Smith      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
69064c8009SBarry Smith      Xsurf will be all zero (thus making the coarse matrix singular).
70064c8009SBarry Smith   */
7108401ef6SPierre 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");
7208401ef6SPierre 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");
7308401ef6SPierre 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");
74064c8009SBarry Smith 
75064c8009SBarry Smith   cnt = 0;
762fa5cd67SKarl Rupp 
772fa5cd67SKarl Rupp   xsurf[cnt++] = 1;
782fa5cd67SKarl Rupp   for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + Nsurf] = 1;
792fa5cd67SKarl Rupp   xsurf[cnt++ + 2 * Nsurf] = 1;
802fa5cd67SKarl Rupp 
812fa5cd67SKarl Rupp   for (j = 1; j < n - 1 - jstart; j++) {
822fa5cd67SKarl Rupp     xsurf[cnt++ + 3 * Nsurf] = 1;
832fa5cd67SKarl Rupp     for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 4 * Nsurf] = 1;
842fa5cd67SKarl Rupp     xsurf[cnt++ + 5 * Nsurf] = 1;
85064c8009SBarry Smith   }
862fa5cd67SKarl Rupp 
872fa5cd67SKarl Rupp   xsurf[cnt++ + 6 * Nsurf] = 1;
882fa5cd67SKarl Rupp   for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 7 * Nsurf] = 1;
892fa5cd67SKarl Rupp   xsurf[cnt++ + 8 * Nsurf] = 1;
902fa5cd67SKarl Rupp 
912fa5cd67SKarl Rupp   for (k = 1; k < p - 1 - kstart; k++) {
922fa5cd67SKarl Rupp     xsurf[cnt++ + 9 * Nsurf] = 1;
932fa5cd67SKarl Rupp     for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 10 * Nsurf] = 1;
942fa5cd67SKarl Rupp     xsurf[cnt++ + 11 * Nsurf] = 1;
952fa5cd67SKarl Rupp 
962fa5cd67SKarl Rupp     for (j = 1; j < n - 1 - jstart; j++) {
972fa5cd67SKarl Rupp       xsurf[cnt++ + 12 * Nsurf] = 1;
982fa5cd67SKarl Rupp       /* these are the interior nodes */
992fa5cd67SKarl Rupp       xsurf[cnt++ + 13 * Nsurf] = 1;
1002fa5cd67SKarl Rupp     }
1012fa5cd67SKarl Rupp 
1022fa5cd67SKarl Rupp     xsurf[cnt++ + 14 * Nsurf] = 1;
1032fa5cd67SKarl Rupp     for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 15 * Nsurf] = 1;
1042fa5cd67SKarl Rupp     xsurf[cnt++ + 16 * Nsurf] = 1;
1052fa5cd67SKarl Rupp   }
1062fa5cd67SKarl Rupp 
1072fa5cd67SKarl Rupp   xsurf[cnt++ + 17 * Nsurf] = 1;
1082fa5cd67SKarl Rupp   for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 18 * Nsurf] = 1;
1092fa5cd67SKarl Rupp   xsurf[cnt++ + 19 * Nsurf] = 1;
1102fa5cd67SKarl Rupp 
1112fa5cd67SKarl Rupp   for (j = 1; j < n - 1 - jstart; j++) {
1122fa5cd67SKarl Rupp     xsurf[cnt++ + 20 * Nsurf] = 1;
1132fa5cd67SKarl Rupp     for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 21 * Nsurf] = 1;
1142fa5cd67SKarl Rupp     xsurf[cnt++ + 22 * Nsurf] = 1;
1152fa5cd67SKarl Rupp   }
1162fa5cd67SKarl Rupp 
1172fa5cd67SKarl Rupp   xsurf[cnt++ + 23 * Nsurf] = 1;
1182fa5cd67SKarl Rupp   for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 24 * Nsurf] = 1;
1192fa5cd67SKarl Rupp   xsurf[cnt++ + 25 * Nsurf] = 1;
1202fa5cd67SKarl Rupp 
1218e722e36SBarry Smith   /* interpolations only sum to 1 when using direct solver */
1228e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
123064c8009SBarry Smith   for (i = 0; i < Nsurf; i++) {
124064c8009SBarry Smith     tmp = 0.0;
1252fa5cd67SKarl Rupp     for (j = 0; j < 26; j++) tmp += xsurf[i + j * Nsurf];
12663a3b9bcSJacob 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));
127064c8009SBarry Smith   }
128064c8009SBarry Smith #endif
1299566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(Xsurf, &xsurf));
1309566063dSJacob Faibussowitsch   /* PetscCall(MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD));*/
131064c8009SBarry Smith 
132064c8009SBarry Smith   /*
133064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
134064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
1357dae84e0SHong Zhang             (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it
136aa219208SBarry Smith              is NOT the local DMDA ordering.)
137064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
138064c8009SBarry Smith   */
139064c8009SBarry Smith #define Endpoint(a, start, b) (a == 0 || a == (b - 1 - start))
1409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(N, &II, Nint, &Iint, Nsurf, &Isurf));
1419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nint, &IIint, Nsurf, &IIsurf));
142064c8009SBarry Smith   for (k = 0; k < p - kstart; k++) {
143064c8009SBarry Smith     for (j = 0; j < n - jstart; j++) {
144064c8009SBarry Smith       for (i = 0; i < m - istart; i++) {
145064c8009SBarry Smith         II[c++] = i + j * mwidth + k * mwidth * nwidth;
146064c8009SBarry Smith 
147064c8009SBarry Smith         if (!Endpoint(i, istart, m) && !Endpoint(j, jstart, n) && !Endpoint(k, kstart, p)) {
148064c8009SBarry Smith           IIint[cint]  = i + j * mwidth + k * mwidth * nwidth;
149064c8009SBarry Smith           Iint[cint++] = i + j * (m - istart) + k * (m - istart) * (n - jstart);
150064c8009SBarry Smith         } else {
151064c8009SBarry Smith           IIsurf[csurf]  = i + j * mwidth + k * mwidth * nwidth;
152064c8009SBarry Smith           Isurf[csurf++] = i + j * (m - istart) + k * (m - istart) * (n - jstart);
153064c8009SBarry Smith         }
154064c8009SBarry Smith       }
155064c8009SBarry Smith     }
156064c8009SBarry Smith   }
157eec179cfSJacob Faibussowitsch #undef Endpoint
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 
2389371c9d4SSatish Balay   gl[cnt++] = 0;
239d71ae5a4SJacob Faibussowitsch   {
240d71ae5a4SJacob Faibussowitsch     gl[cnt++] = 1;
241d71ae5a4SJacob Faibussowitsch   }
2429371c9d4SSatish Balay   gl[cnt++] = m - istart - 1;
243064c8009SBarry Smith   {
2449371c9d4SSatish Balay     gl[cnt++] = mwidth;
245d71ae5a4SJacob Faibussowitsch     {
246d71ae5a4SJacob Faibussowitsch       gl[cnt++] = mwidth + 1;
247d71ae5a4SJacob Faibussowitsch     }
2489371c9d4SSatish Balay     gl[cnt++] = mwidth + m - istart - 1;
249064c8009SBarry Smith   }
2509371c9d4SSatish Balay   gl[cnt++] = mwidth * (n - jstart - 1);
251d71ae5a4SJacob Faibussowitsch   {
252d71ae5a4SJacob Faibussowitsch     gl[cnt++] = mwidth * (n - jstart - 1) + 1;
253d71ae5a4SJacob Faibussowitsch   }
2549371c9d4SSatish Balay   gl[cnt++] = mwidth * (n - jstart - 1) + m - istart - 1;
2559371c9d4SSatish Balay   {
2569371c9d4SSatish Balay     gl[cnt++] = mwidth * nwidth;
257d71ae5a4SJacob Faibussowitsch     {
258d71ae5a4SJacob Faibussowitsch       gl[cnt++] = mwidth * nwidth + 1;
259d71ae5a4SJacob Faibussowitsch     }
2609371c9d4SSatish Balay     gl[cnt++] = mwidth * nwidth + m - istart - 1;
2619371c9d4SSatish Balay     {
2629371c9d4SSatish Balay       gl[cnt++] = mwidth * nwidth + mwidth; /* these are the interior nodes */
2639371c9d4SSatish Balay       gl[cnt++] = mwidth * nwidth + mwidth + m - istart - 1;
2649371c9d4SSatish Balay     }
2659371c9d4SSatish Balay     gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1);
266d71ae5a4SJacob Faibussowitsch     {
267d71ae5a4SJacob Faibussowitsch       gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1) + 1;
268d71ae5a4SJacob Faibussowitsch     }
2699371c9d4SSatish Balay     gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1) + m - istart - 1;
2709371c9d4SSatish Balay   }
2719371c9d4SSatish Balay   gl[cnt++] = mwidth * nwidth * (p - kstart - 1);
272d71ae5a4SJacob Faibussowitsch   {
273d71ae5a4SJacob Faibussowitsch     gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + 1;
274d71ae5a4SJacob Faibussowitsch   }
2759371c9d4SSatish Balay   gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + m - istart - 1;
2769371c9d4SSatish Balay   {
2779371c9d4SSatish Balay     gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth;
278d71ae5a4SJacob Faibussowitsch     {
279d71ae5a4SJacob Faibussowitsch       gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth + 1;
280d71ae5a4SJacob Faibussowitsch     }
2819371c9d4SSatish Balay     gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth + m - istart - 1;
2829371c9d4SSatish Balay   }
2839371c9d4SSatish Balay   gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth * (n - jstart - 1);
284d71ae5a4SJacob Faibussowitsch   {
285d71ae5a4SJacob Faibussowitsch     gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth * (n - jstart - 1) + 1;
286d71ae5a4SJacob Faibussowitsch   }
2879371c9d4SSatish Balay   gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth * (n - jstart - 1) + m - istart - 1;
288064c8009SBarry Smith 
289064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
290064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
2919566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg, 26, gl, gl));
292064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
2939566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(26 * mp * np * pp, &globals));
2949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(gl, 26, MPIU_INT, globals, 26, MPIU_INT, PetscObjectComm((PetscObject)da)));
295064c8009SBarry Smith 
296064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
297eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapICreateWithSize(Ntotal / 3, &ht));
298eec179cfSJacob Faibussowitsch   for (i = 0, cnt = 0; i < 26 * mp * np * pp; i++) {
299eec179cfSJacob Faibussowitsch     PetscHashIter it      = 0;
300eec179cfSJacob Faibussowitsch     PetscBool     missing = PETSC_TRUE;
301eec179cfSJacob Faibussowitsch 
302eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIPut(ht, globals[i] + 1, &it, &missing));
303eec179cfSJacob Faibussowitsch     if (missing) {
304eec179cfSJacob Faibussowitsch       ++cnt;
305eec179cfSJacob Faibussowitsch       PetscCall(PetscHMapIIterSet(ht, it, cnt));
306eec179cfSJacob Faibussowitsch     }
307eec179cfSJacob Faibussowitsch   }
30863a3b9bcSJacob 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);
3099566063dSJacob Faibussowitsch   PetscCall(PetscFree(globals));
310064c8009SBarry Smith   for (i = 0; i < 26; i++) {
311eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIGetWithDefault(ht, gl[i] + 1, 0, gl + i));
312eec179cfSJacob Faibussowitsch     --(gl[i]);
313064c8009SBarry Smith   }
314eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&ht));
315064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
316064c8009SBarry Smith 
317064c8009SBarry Smith   /* construct global interpolation matrix */
3189566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Aglobal, &Ng, NULL));
319064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
3209566063dSJacob Faibussowitsch     PetscCall(MatCreateAIJ(PetscObjectComm((PetscObject)da), Ng, PETSC_DECIDE, PETSC_DECIDE, Ntotal, Nint + Nsurf, NULL, Nint + Nsurf, NULL, P));
321064c8009SBarry Smith   } else {
3229566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(*P));
323064c8009SBarry Smith   }
3249566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*P, MAT_ROW_ORIENTED, PETSC_FALSE));
3259566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xint, &rxint));
3269566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*P, Nint, IIint, 26, gl, rxint, INSERT_VALUES));
3279566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xint, &rxint));
3289566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xsurf, &rxint));
3299566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*P, Nsurf, IIsurf, 26, gl, rxint, INSERT_VALUES));
3309566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xsurf, &rxint));
3319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
3329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
3339566063dSJacob Faibussowitsch   PetscCall(PetscFree2(IIint, IIsurf));
334064c8009SBarry Smith 
3358e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
336064c8009SBarry Smith   {
337064c8009SBarry Smith     Vec          x, y;
338064c8009SBarry Smith     PetscScalar *yy;
3399566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), Ng, PETSC_DETERMINE, &y));
3409566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), PETSC_DETERMINE, Ntotal, &x));
3419566063dSJacob Faibussowitsch     PetscCall(VecSet(x, 1.0));
3429566063dSJacob Faibussowitsch     PetscCall(MatMult(*P, x, y));
3439566063dSJacob Faibussowitsch     PetscCall(VecGetArray(y, &yy));
344ad540459SPierre Jolivet     for (i = 0; i < Ng; i++) 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]));
3459566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(y, &yy));
3469566063dSJacob Faibussowitsch     PetscCall(VecDestroy(x));
3479566063dSJacob Faibussowitsch     PetscCall(VecDestroy(y));
348064c8009SBarry Smith   }
349064c8009SBarry Smith #endif
350064c8009SBarry Smith 
3519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Aii));
3529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ais));
3539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Asi));
3549566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3559566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
3569566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isint));
3579566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&issurf));
3589566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xint));
3599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xsurf));
3603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
361064c8009SBarry Smith }
362064c8009SBarry Smith 
363064c8009SBarry Smith /*
364aa219208SBarry Smith       DMDAGetFaceInterpolation - Gets the interpolation for a face based coarse space
365064c8009SBarry Smith 
366064c8009SBarry Smith */
36766976f2fSJacob Faibussowitsch static PetscErrorCode DMDAGetFaceInterpolation(PC pc, DM da, PC_Exotic *exotic, Mat Aglobal, MatReuse reuse, Mat *P)
368d71ae5a4SJacob Faibussowitsch {
369064c8009SBarry 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;
370eec179cfSJacob Faibussowitsch   PetscInt               mwidth, nwidth, pwidth, cnt, mp, np, pp, Ntotal, gl[6], *globals, Ng, *IIint, *IIsurf;
371064c8009SBarry Smith   Mat                    Xint, Xsurf, Xint_tmp;
372064c8009SBarry Smith   IS                     isint, issurf, is, row, col;
373064c8009SBarry Smith   ISLocalToGlobalMapping ltg;
374064c8009SBarry Smith   MPI_Comm               comm;
375064c8009SBarry Smith   Mat                    A, Aii, Ais, Asi, *Aholder, iAii;
376064c8009SBarry Smith   MatFactorInfo          info;
377064c8009SBarry Smith   PetscScalar           *xsurf, *xint;
3781683a169SBarry Smith   const PetscScalar     *rxint;
379064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
380064c8009SBarry Smith   PetscScalar tmp;
381064c8009SBarry Smith #endif
382eec179cfSJacob Faibussowitsch   PetscHMapI ht;
383064c8009SBarry Smith 
384064c8009SBarry Smith   PetscFunctionBegin;
3859566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, &mp, &np, &pp, &dof, NULL, NULL, NULL, NULL, NULL));
38608401ef6SPierre Jolivet   PetscCheck(dof == 1, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only for single field problems");
38708401ef6SPierre Jolivet   PetscCheck(dim == 3, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only coded for 3d problems");
3889566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, NULL, NULL, NULL, &m, &n, &p));
3899566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &istart, &jstart, &kstart, &mwidth, &nwidth, &pwidth));
390064c8009SBarry Smith   istart = istart ? -1 : 0;
391064c8009SBarry Smith   jstart = jstart ? -1 : 0;
392064c8009SBarry Smith   kstart = kstart ? -1 : 0;
393064c8009SBarry Smith 
394064c8009SBarry Smith   /*
395064c8009SBarry Smith     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge)
396064c8009SBarry Smith     to all the local degrees of freedom (this includes the vertices, edges and faces).
397064c8009SBarry Smith 
398064c8009SBarry Smith     Xint are the subset of the interpolation into the interior
399064c8009SBarry Smith 
400064c8009SBarry Smith     Xface are the interpolation onto faces but not into the interior
401064c8009SBarry Smith 
402064c8009SBarry Smith     Xsurf are the interpolation onto the vertices and edges (the surfbasket)
403064c8009SBarry Smith                                       Xint
404064c8009SBarry Smith     Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain
405064c8009SBarry Smith                                       Xsurf
406064c8009SBarry Smith   */
407064c8009SBarry Smith   N     = (m - istart) * (n - jstart) * (p - kstart);
408064c8009SBarry Smith   Nint  = (m - 2 - istart) * (n - 2 - jstart) * (p - 2 - kstart);
409064c8009SBarry Smith   Nface = 2 * ((m - 2 - istart) * (n - 2 - jstart) + (m - 2 - istart) * (p - 2 - kstart) + (n - 2 - jstart) * (p - 2 - kstart));
410064c8009SBarry Smith   Nwire = 4 * ((m - 2 - istart) + (n - 2 - jstart) + (p - 2 - kstart)) + 8;
411064c8009SBarry Smith   Nsurf = Nface + Nwire;
4129566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nint, 6, NULL, &Xint));
4139566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nsurf, 6, NULL, &Xsurf));
4149566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(Xsurf, &xsurf));
415064c8009SBarry Smith 
416064c8009SBarry Smith   /*
417064c8009SBarry Smith      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
418064c8009SBarry Smith      Xsurf will be all zero (thus making the coarse matrix singular).
419064c8009SBarry Smith   */
42008401ef6SPierre 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");
42108401ef6SPierre 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");
42208401ef6SPierre 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");
423064c8009SBarry Smith 
424064c8009SBarry Smith   cnt = 0;
4252fa5cd67SKarl Rupp   for (j = 1; j < n - 1 - jstart; j++) {
4262fa5cd67SKarl Rupp     for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 0 * Nsurf] = 1;
427064c8009SBarry Smith   }
4282fa5cd67SKarl Rupp 
4292fa5cd67SKarl Rupp   for (k = 1; k < p - 1 - kstart; k++) {
4302fa5cd67SKarl Rupp     for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 1 * Nsurf] = 1;
4312fa5cd67SKarl Rupp     for (j = 1; j < n - 1 - jstart; j++) {
4322fa5cd67SKarl Rupp       xsurf[cnt++ + 2 * Nsurf] = 1;
4332fa5cd67SKarl Rupp       /* these are the interior nodes */
4342fa5cd67SKarl Rupp       xsurf[cnt++ + 3 * Nsurf] = 1;
4352fa5cd67SKarl Rupp     }
4362fa5cd67SKarl Rupp     for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 4 * Nsurf] = 1;
4372fa5cd67SKarl Rupp   }
4382fa5cd67SKarl Rupp   for (j = 1; j < n - 1 - jstart; j++) {
4392fa5cd67SKarl Rupp     for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 5 * Nsurf] = 1;
4402fa5cd67SKarl Rupp   }
441064c8009SBarry Smith 
442064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
443064c8009SBarry Smith   for (i = 0; i < Nsurf; i++) {
444064c8009SBarry Smith     tmp = 0.0;
4452fa5cd67SKarl Rupp     for (j = 0; j < 6; j++) tmp += xsurf[i + j * Nsurf];
4462fa5cd67SKarl Rupp 
44763a3b9bcSJacob 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));
448064c8009SBarry Smith   }
449064c8009SBarry Smith #endif
4509566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(Xsurf, &xsurf));
4519566063dSJacob Faibussowitsch   /* PetscCall(MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD));*/
452064c8009SBarry Smith 
453064c8009SBarry Smith   /*
454064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
455064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
4567dae84e0SHong Zhang             (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it
457aa219208SBarry Smith              is NOT the local DMDA ordering.)
458064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
459064c8009SBarry Smith   */
460064c8009SBarry Smith #define Endpoint(a, start, b) (a == 0 || a == (b - 1 - start))
4619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(N, &II, Nint, &Iint, Nsurf, &Isurf));
4629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nint, &IIint, Nsurf, &IIsurf));
463064c8009SBarry Smith   for (k = 0; k < p - kstart; k++) {
464064c8009SBarry Smith     for (j = 0; j < n - jstart; j++) {
465064c8009SBarry Smith       for (i = 0; i < m - istart; i++) {
466064c8009SBarry Smith         II[c++] = i + j * mwidth + k * mwidth * nwidth;
467064c8009SBarry Smith 
468064c8009SBarry Smith         if (!Endpoint(i, istart, m) && !Endpoint(j, jstart, n) && !Endpoint(k, kstart, p)) {
469064c8009SBarry Smith           IIint[cint]  = i + j * mwidth + k * mwidth * nwidth;
470064c8009SBarry Smith           Iint[cint++] = i + j * (m - istart) + k * (m - istart) * (n - jstart);
471064c8009SBarry Smith         } else {
472064c8009SBarry Smith           IIsurf[csurf]  = i + j * mwidth + k * mwidth * nwidth;
473064c8009SBarry Smith           Isurf[csurf++] = i + j * (m - istart) + k * (m - istart) * (n - jstart);
474064c8009SBarry Smith         }
475064c8009SBarry Smith       }
476064c8009SBarry Smith     }
477064c8009SBarry Smith   }
478eec179cfSJacob Faibussowitsch #undef Endpoint
47908401ef6SPierre Jolivet   PetscCheck(c == N, PETSC_COMM_SELF, PETSC_ERR_PLIB, "c != N");
48008401ef6SPierre Jolivet   PetscCheck(cint == Nint, PETSC_COMM_SELF, PETSC_ERR_PLIB, "cint != Nint");
48108401ef6SPierre Jolivet   PetscCheck(csurf == Nsurf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "csurf != Nsurf");
4829566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(da, &ltg));
4839566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg, N, II, II));
4849566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg, Nint, IIint, IIint));
4859566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg, Nsurf, IIsurf, IIsurf));
4869566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
4879566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, N, II, PETSC_COPY_VALUES, &is));
4889566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nint, Iint, PETSC_COPY_VALUES, &isint));
4899566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nsurf, Isurf, PETSC_COPY_VALUES, &issurf));
4909566063dSJacob Faibussowitsch   PetscCall(PetscFree3(II, Iint, Isurf));
491064c8009SBarry Smith 
4929566063dSJacob Faibussowitsch   PetscCall(ISSort(is));
4939566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Aglobal, 1, &is, &is, MAT_INITIAL_MATRIX, &Aholder));
494064c8009SBarry Smith   A = *Aholder;
4959566063dSJacob Faibussowitsch   PetscCall(PetscFree(Aholder));
496064c8009SBarry Smith 
4979566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, isint, isint, MAT_INITIAL_MATRIX, &Aii));
4989566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, isint, issurf, MAT_INITIAL_MATRIX, &Ais));
4999566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, issurf, isint, MAT_INITIAL_MATRIX, &Asi));
500064c8009SBarry Smith 
501064c8009SBarry Smith   /*
502064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
503064c8009SBarry Smith   */
5049566063dSJacob Faibussowitsch   PetscCall(MatMatMult(Ais, Xsurf, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &Xint_tmp));
5059566063dSJacob Faibussowitsch   PetscCall(MatScale(Xint_tmp, -1.0));
506064c8009SBarry Smith 
5078e722e36SBarry Smith   if (exotic->directSolve) {
5089566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(Aii, MATSOLVERPETSC, MAT_FACTOR_LU, &iAii));
5099566063dSJacob Faibussowitsch     PetscCall(MatFactorInfoInitialize(&info));
5109566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(Aii, MATORDERINGND, &row, &col));
5119566063dSJacob Faibussowitsch     PetscCall(MatLUFactorSymbolic(iAii, Aii, row, col, &info));
5129566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&row));
5139566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&col));
5149566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(iAii, Aii, &info));
5159566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(iAii, Xint_tmp, Xint));
5169566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&iAii));
517064c8009SBarry Smith   } else {
518064c8009SBarry Smith     Vec          b, x;
519064c8009SBarry Smith     PetscScalar *xint_tmp;
520064c8009SBarry Smith 
5219566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(Xint, &xint));
5229566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &x));
5239566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(Xint_tmp, &xint_tmp));
5249566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &b));
5259566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(exotic->ksp, Aii, Aii));
526064c8009SBarry Smith     for (i = 0; i < 6; i++) {
5279566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(x, xint + i * Nint));
5289566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(b, xint_tmp + i * Nint));
5299566063dSJacob Faibussowitsch       PetscCall(KSPSolve(exotic->ksp, b, x));
5309566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(exotic->ksp, pc, x));
5319566063dSJacob Faibussowitsch       PetscCall(VecResetArray(x));
5329566063dSJacob Faibussowitsch       PetscCall(VecResetArray(b));
533064c8009SBarry Smith     }
5349566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(Xint, &xint));
5359566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(Xint_tmp, &xint_tmp));
5369566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
5379566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&b));
538064c8009SBarry Smith   }
5399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xint_tmp));
540064c8009SBarry Smith 
541064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
5429566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xint, &rxint));
543064c8009SBarry Smith   for (i = 0; i < Nint; i++) {
544064c8009SBarry Smith     tmp = 0.0;
5451683a169SBarry Smith     for (j = 0; j < 6; j++) tmp += rxint[i + j * Nint];
5462fa5cd67SKarl Rupp 
54763a3b9bcSJacob 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));
548064c8009SBarry Smith   }
5499566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xint, &rxint));
5509566063dSJacob Faibussowitsch   /* PetscCall(MatView(Xint,PETSC_VIEWER_STDOUT_WORLD)); */
551064c8009SBarry Smith #endif
552064c8009SBarry Smith 
553064c8009SBarry Smith   /*         total faces    */
554064c8009SBarry Smith   Ntotal = mp * np * (pp + 1) + mp * pp * (np + 1) + np * pp * (mp + 1);
555064c8009SBarry Smith 
556064c8009SBarry Smith   /*
557064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
558064c8009SBarry Smith   */
559064c8009SBarry Smith   cnt = 0;
560d71ae5a4SJacob Faibussowitsch   {
561d71ae5a4SJacob Faibussowitsch     gl[cnt++] = mwidth + 1;
562d71ae5a4SJacob Faibussowitsch   }
563*fbccb6d4SPierre Jolivet   {
564*fbccb6d4SPierre Jolivet     {
565*fbccb6d4SPierre Jolivet       gl[cnt++] = mwidth * nwidth + 1;
5669371c9d4SSatish Balay     }
567064c8009SBarry Smith     {
5689371c9d4SSatish Balay       gl[cnt++] = mwidth * nwidth + mwidth; /* these are the interior nodes */
5699371c9d4SSatish Balay       gl[cnt++] = mwidth * nwidth + mwidth + m - istart - 1;
5709371c9d4SSatish Balay     }
571d71ae5a4SJacob Faibussowitsch     {
572d71ae5a4SJacob Faibussowitsch       gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1) + 1;
573064c8009SBarry Smith     }
574d71ae5a4SJacob Faibussowitsch   }
575d71ae5a4SJacob Faibussowitsch   {
576d71ae5a4SJacob Faibussowitsch     gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth + 1;
577d71ae5a4SJacob Faibussowitsch   }
578064c8009SBarry Smith 
579064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
580064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
5819566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg, 6, gl, gl));
582064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
5839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(6 * mp * np * pp, &globals));
5849566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(gl, 6, MPIU_INT, globals, 6, MPIU_INT, PetscObjectComm((PetscObject)da)));
585064c8009SBarry Smith 
586064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
587eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapICreateWithSize(Ntotal / 3, &ht));
588eec179cfSJacob Faibussowitsch   for (i = 0, cnt = 0; i < 6 * mp * np * pp; i++) {
589eec179cfSJacob Faibussowitsch     PetscHashIter it      = 0;
590eec179cfSJacob Faibussowitsch     PetscBool     missing = PETSC_TRUE;
591eec179cfSJacob Faibussowitsch 
592eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIPut(ht, globals[i] + 1, &it, &missing));
593eec179cfSJacob Faibussowitsch     if (missing) {
594eec179cfSJacob Faibussowitsch       ++cnt;
595eec179cfSJacob Faibussowitsch       PetscCall(PetscHMapIIterSet(ht, it, cnt));
596eec179cfSJacob Faibussowitsch     }
597eec179cfSJacob Faibussowitsch   }
59863a3b9bcSJacob 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);
5999566063dSJacob Faibussowitsch   PetscCall(PetscFree(globals));
600064c8009SBarry Smith   for (i = 0; i < 6; i++) {
601eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIGetWithDefault(ht, gl[i] + 1, 0, gl + i));
602eec179cfSJacob Faibussowitsch     --(gl[i]);
603064c8009SBarry Smith   }
604eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&ht));
605064c8009SBarry Smith   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
606064c8009SBarry Smith 
607064c8009SBarry Smith   /* construct global interpolation matrix */
6089566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Aglobal, &Ng, NULL));
609064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
6109566063dSJacob Faibussowitsch     PetscCall(MatCreateAIJ(PetscObjectComm((PetscObject)da), Ng, PETSC_DECIDE, PETSC_DECIDE, Ntotal, Nint + Nsurf, NULL, Nint, NULL, P));
611064c8009SBarry Smith   } else {
6129566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(*P));
613064c8009SBarry Smith   }
6149566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*P, MAT_ROW_ORIENTED, PETSC_FALSE));
6159566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xint, &rxint));
6169566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*P, Nint, IIint, 6, gl, rxint, INSERT_VALUES));
6179566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xint, &rxint));
6189566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xsurf, &rxint));
6199566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*P, Nsurf, IIsurf, 6, gl, rxint, INSERT_VALUES));
6209566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xsurf, &rxint));
6219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
6229566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
6239566063dSJacob Faibussowitsch   PetscCall(PetscFree2(IIint, IIsurf));
624064c8009SBarry Smith 
625064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
626064c8009SBarry Smith   {
627064c8009SBarry Smith     Vec          x, y;
628064c8009SBarry Smith     PetscScalar *yy;
6299566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), Ng, PETSC_DETERMINE, &y));
6309566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), PETSC_DETERMINE, Ntotal, &x));
6319566063dSJacob Faibussowitsch     PetscCall(VecSet(x, 1.0));
6329566063dSJacob Faibussowitsch     PetscCall(MatMult(*P, x, y));
6339566063dSJacob Faibussowitsch     PetscCall(VecGetArray(y, &yy));
634ad540459SPierre Jolivet     for (i = 0; i < Ng; i++) 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]));
6359566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(y, &yy));
6369566063dSJacob Faibussowitsch     PetscCall(VecDestroy(x));
6379566063dSJacob Faibussowitsch     PetscCall(VecDestroy(y));
638064c8009SBarry Smith   }
639064c8009SBarry Smith #endif
640064c8009SBarry Smith 
6419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Aii));
6429566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ais));
6439566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Asi));
6449566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
6459566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
6469566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isint));
6479566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&issurf));
6489566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xint));
6499566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xsurf));
6503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
651064c8009SBarry Smith }
652174b6946SBarry Smith 
6537233f9f0SBarry Smith /*@
6547233f9f0SBarry Smith   PCExoticSetType - Sets the type of coarse grid interpolation to use
6557233f9f0SBarry Smith 
656c3339decSBarry Smith   Logically Collective
6577233f9f0SBarry Smith 
6587233f9f0SBarry Smith   Input Parameters:
6597233f9f0SBarry Smith + pc   - the preconditioner context
6607233f9f0SBarry Smith - type - either PC_EXOTIC_FACE or PC_EXOTIC_WIREBASKET (defaults to face)
6617233f9f0SBarry Smith 
66295452b02SPatrick Sanan   Notes:
66395452b02SPatrick Sanan   The face based interpolation has 1 degree of freedom per face and ignores the
664563e08c6SBarry Smith   edge and vertex values completely in the coarse problem. For any seven point
665563e08c6SBarry Smith   stencil the interpolation of a constant on all faces into the interior is that constant.
666563e08c6SBarry Smith 
667563e08c6SBarry Smith   The wirebasket interpolation has 1 degree of freedom per vertex, per edge and
668563e08c6SBarry Smith   per face. A constant on the subdomain boundary is interpolated as that constant
669563e08c6SBarry Smith   in the interior of the domain.
670563e08c6SBarry Smith 
671563e08c6SBarry Smith   The coarse grid matrix is obtained via the Galerkin computation A_c = R A R^T, hence
672563e08c6SBarry Smith   if A is nonsingular A_c is also nonsingular.
673563e08c6SBarry Smith 
674563e08c6SBarry Smith   Both interpolations are suitable for only scalar problems.
675563e08c6SBarry Smith 
6767233f9f0SBarry Smith   Level: intermediate
6777233f9f0SBarry Smith 
678562efe2eSBarry Smith .seealso: [](ch_ksp), `PCEXOTIC`, `PCExoticType()`
6797233f9f0SBarry Smith @*/
680d71ae5a4SJacob Faibussowitsch PetscErrorCode PCExoticSetType(PC pc, PCExoticType type)
681d71ae5a4SJacob Faibussowitsch {
6827233f9f0SBarry Smith   PetscFunctionBegin;
6830700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
684c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc, type, 2);
685cac4c232SBarry Smith   PetscTryMethod(pc, "PCExoticSetType_C", (PC, PCExoticType), (pc, type));
6863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6877233f9f0SBarry Smith }
6887233f9f0SBarry Smith 
689d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCExoticSetType_Exotic(PC pc, PCExoticType type)
690d71ae5a4SJacob Faibussowitsch {
691f3fbd535SBarry Smith   PC_MG     *mg  = (PC_MG *)pc->data;
69231567311SBarry Smith   PC_Exotic *ctx = (PC_Exotic *)mg->innerctx;
6937233f9f0SBarry Smith 
6947233f9f0SBarry Smith   PetscFunctionBegin;
6957233f9f0SBarry Smith   ctx->type = type;
6963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6977233f9f0SBarry Smith }
6987233f9f0SBarry Smith 
69966976f2fSJacob Faibussowitsch static PetscErrorCode PCSetUp_Exotic(PC pc)
700d71ae5a4SJacob Faibussowitsch {
70196bdf778SBarry Smith   Mat        A;
702f3fbd535SBarry Smith   PC_MG     *mg    = (PC_MG *)pc->data;
70331567311SBarry Smith   PC_Exotic *ex    = (PC_Exotic *)mg->innerctx;
70496bdf778SBarry Smith   MatReuse   reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX;
705174b6946SBarry Smith 
706174b6946SBarry Smith   PetscFunctionBegin;
70728b400f6SJacob Faibussowitsch   PetscCheck(pc->dm, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Need to call PCSetDM() before using this PC");
7089566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &A));
7090fdf79fbSJacob Faibussowitsch   PetscCheck(ex->type == PC_EXOTIC_FACE || ex->type == PC_EXOTIC_WIREBASKET, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Unknown exotic coarse space %d", ex->type);
7107233f9f0SBarry Smith   if (ex->type == PC_EXOTIC_FACE) {
7119566063dSJacob Faibussowitsch     PetscCall(DMDAGetFaceInterpolation(pc, pc->dm, ex, A, reuse, &ex->P));
7120fdf79fbSJacob Faibussowitsch   } else /* if (ex->type == PC_EXOTIC_WIREBASKET) */ {
7139566063dSJacob Faibussowitsch     PetscCall(DMDAGetWireBasketInterpolation(pc, pc->dm, ex, A, reuse, &ex->P));
7140fdf79fbSJacob Faibussowitsch   }
7159566063dSJacob Faibussowitsch   PetscCall(PCMGSetInterpolation(pc, 1, ex->P));
716d0660788SBarry Smith   /* if PC has attached DM we must remove it or the PCMG will use it to compute incorrect sized vectors and interpolations */
7179566063dSJacob Faibussowitsch   PetscCall(PCSetDM(pc, NULL));
7189566063dSJacob Faibussowitsch   PetscCall(PCSetUp_MG(pc));
7193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
720174b6946SBarry Smith }
721174b6946SBarry Smith 
72266976f2fSJacob Faibussowitsch static PetscErrorCode PCDestroy_Exotic(PC pc)
723d71ae5a4SJacob Faibussowitsch {
724f3fbd535SBarry Smith   PC_MG     *mg  = (PC_MG *)pc->data;
72531567311SBarry Smith   PC_Exotic *ctx = (PC_Exotic *)mg->innerctx;
726174b6946SBarry Smith 
727174b6946SBarry Smith   PetscFunctionBegin;
7289566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->P));
7299566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ctx->ksp));
7309566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
7312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCExoticSetType_C", NULL));
7329566063dSJacob Faibussowitsch   PetscCall(PCDestroy_MG(pc));
7333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
734174b6946SBarry Smith }
735174b6946SBarry Smith 
73666976f2fSJacob Faibussowitsch static PetscErrorCode PCView_Exotic(PC pc, PetscViewer viewer)
737d71ae5a4SJacob Faibussowitsch {
738f3fbd535SBarry Smith   PC_MG     *mg = (PC_MG *)pc->data;
739ace3abfcSBarry Smith   PetscBool  iascii;
74031567311SBarry Smith   PC_Exotic *ctx = (PC_Exotic *)mg->innerctx;
7417233f9f0SBarry Smith 
7427233f9f0SBarry Smith   PetscFunctionBegin;
7439566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
7447233f9f0SBarry Smith   if (iascii) {
7459566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Exotic type = %s\n", PCExoticTypes[ctx->type]));
7468e722e36SBarry Smith     if (ctx->directSolve) {
7479566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "      Using direct solver to construct interpolation\n"));
7488e722e36SBarry Smith     } else {
7498e722e36SBarry Smith       PetscViewer sviewer;
7508e722e36SBarry Smith       PetscMPIInt rank;
7518e722e36SBarry Smith 
7529566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "      Using iterative solver to construct interpolation\n"));
7539566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
7549566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer)); /* should not need to push this twice? */
7559566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
7569566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
75748a46eb9SPierre Jolivet       if (rank == 0) PetscCall(KSPView(ctx->ksp, sviewer));
7589566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
7599566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
7609566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
7618e722e36SBarry Smith     }
7627233f9f0SBarry Smith   }
7639566063dSJacob Faibussowitsch   PetscCall(PCView_MG(pc, viewer));
7643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7657233f9f0SBarry Smith }
7667233f9f0SBarry Smith 
76766976f2fSJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_Exotic(PC pc, PetscOptionItems *PetscOptionsObject)
768d71ae5a4SJacob Faibussowitsch {
769ace3abfcSBarry Smith   PetscBool    flg;
770f3fbd535SBarry Smith   PC_MG       *mg = (PC_MG *)pc->data;
7717233f9f0SBarry Smith   PCExoticType mgctype;
77231567311SBarry Smith   PC_Exotic   *ctx = (PC_Exotic *)mg->innerctx;
7737233f9f0SBarry Smith 
7747233f9f0SBarry Smith   PetscFunctionBegin;
775d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Exotic coarse space options");
7769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_exotic_type", "face or wirebasket", "PCExoticSetType", PCExoticTypes, (PetscEnum)ctx->type, (PetscEnum *)&mgctype, &flg));
7771baa6e33SBarry Smith   if (flg) PetscCall(PCExoticSetType(pc, mgctype));
7789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_exotic_direct_solver", "use direct solver to construct interpolation", "None", ctx->directSolve, &ctx->directSolve, NULL));
7798e722e36SBarry Smith   if (!ctx->directSolve) {
7808e722e36SBarry Smith     if (!ctx->ksp) {
7818e722e36SBarry Smith       const char *prefix;
7829566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PETSC_COMM_SELF, &ctx->ksp));
7833821be0aSBarry Smith       PetscCall(KSPSetNestLevel(ctx->ksp, pc->kspnestlevel));
7849566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(ctx->ksp, pc->erroriffailure));
7859566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->ksp, (PetscObject)pc, 1));
7869566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
7879566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(ctx->ksp, prefix));
7889566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(ctx->ksp, "exotic_"));
7898e722e36SBarry Smith     }
7909566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ctx->ksp));
7918e722e36SBarry Smith   }
792d0609cedSBarry Smith   PetscOptionsHeadEnd();
7933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7947233f9f0SBarry Smith }
7957233f9f0SBarry Smith 
796174b6946SBarry Smith /*MC
7977233f9f0SBarry Smith      PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces
798174b6946SBarry Smith 
799f1580f4eSBarry Smith      This uses the `PCMG` infrastructure restricted to two levels and the face and wirebasket based coarse
80024c3aa18SBarry Smith    grid spaces.
80124c3aa18SBarry Smith 
802f1580f4eSBarry Smith    Level: advanced
803f1580f4eSBarry Smith 
80495452b02SPatrick Sanan    Notes:
805f1580f4eSBarry Smith    Must be used with `DMDA`
806f1580f4eSBarry Smith 
807f1580f4eSBarry Smith    By default this uses `KSPGMRES` on the fine grid smoother so this should be used with `KSPFGMRES` or the smoother changed to not use `KSPGMRES`
80824c3aa18SBarry Smith 
8091d27aa22SBarry Smith    These coarse grid spaces originate in the work of Bramble, Pasciak  and Schatz {cite}`bramble1989construction`.
8101d27aa22SBarry Smith    They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith, {cite}`smith1990domain`.
8111d27aa22SBarry Smith    They were then explored in great detail in Dryja, Smith, Widlund {cite}`dryja1994schwarz`. These were developed in the context of iterative substructuring preconditioners.
8121d27aa22SBarry Smith 
8131d27aa22SBarry Smith    They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund.
8141d27aa22SBarry Smith    They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example, {cite}`dohrmann2008extending`, {cite}`dohrmann2008family`,
8151d27aa22SBarry Smith    {cite}`dohrmann2008domain`, {cite}`dohrmann2009overlapping`.
8167233f9f0SBarry Smith 
8178f4fb22eSMark Adams    The usual `PCMG` options are supported, such as -mg_levels_pc_type <type> -mg_coarse_pc_type <type> -mg_fine_pc_type <type> and -pc_mg_type <type>
818174b6946SBarry Smith 
819562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCSetDM()`, `PCExoticType`, `PCExoticSetType()`
820174b6946SBarry Smith M*/
821174b6946SBarry Smith 
822d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Exotic(PC pc)
823d71ae5a4SJacob Faibussowitsch {
8247233f9f0SBarry Smith   PC_Exotic *ex;
825f3fbd535SBarry Smith   PC_MG     *mg;
826174b6946SBarry Smith 
827174b6946SBarry Smith   PetscFunctionBegin;
828f91d8e95SBarry Smith   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
829dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, destroy);
8300a545947SLisandro Dalcin   pc->data = NULL;
831dbbe0bcdSBarry Smith 
8329566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PetscObject)pc)->type_name));
8330a545947SLisandro Dalcin   ((PetscObject)pc)->type_name = NULL;
834f91d8e95SBarry Smith 
8359566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCMG));
8369566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pc, 2, NULL));
8379566063dSJacob Faibussowitsch   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_PMAT));
8389371c9d4SSatish Balay   PetscCall(PetscNew(&ex));
8397233f9f0SBarry Smith   ex->type     = PC_EXOTIC_FACE;
840f3fbd535SBarry Smith   mg           = (PC_MG *)pc->data;
84131567311SBarry Smith   mg->innerctx = ex;
8427233f9f0SBarry Smith 
8437233f9f0SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_Exotic;
8447233f9f0SBarry Smith   pc->ops->view           = PCView_Exotic;
8457233f9f0SBarry Smith   pc->ops->destroy        = PCDestroy_Exotic;
8466c699258SBarry Smith   pc->ops->setup          = PCSetUp_Exotic;
8472fa5cd67SKarl Rupp 
8489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCExoticSetType_C", PCExoticSetType_Exotic));
8493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
850174b6946SBarry Smith }
851