xref: /petsc/src/ksp/pc/impls/wb/wb.c (revision eec179cf895b6fbcd6a0b58694319392b06e361b)
1174b6946SBarry Smith 
2c6db04a5SJed Brown #include <petscdmda.h>              /*I "petscdmda.h" I*/
3af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/
4*eec179cfSJacob Faibussowitsch #include <petsc/private/hashmapi.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 */
19d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetWireBasketInterpolation(PC pc, DM da, PC_Exotic *exotic, Mat Aglobal, MatReuse reuse, Mat *P)
20d71ae5a4SJacob Faibussowitsch {
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;
22*eec179cfSJacob Faibussowitsch   PetscInt               mwidth, nwidth, pwidth, cnt, mp, np, pp, Ntotal, gl[26], *globals, Ng, *IIint, *IIsurf;
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
34*eec179cfSJacob Faibussowitsch   PetscHMapI ht = NULL;
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   }
158*eec179cfSJacob Faibussowitsch #undef Endpoint
15908401ef6SPierre Jolivet   PetscCheck(c == N, PETSC_COMM_SELF, PETSC_ERR_PLIB, "c != N");
16008401ef6SPierre Jolivet   PetscCheck(cint == Nint, PETSC_COMM_SELF, PETSC_ERR_PLIB, "cint != Nint");
16108401ef6SPierre Jolivet   PetscCheck(csurf == Nsurf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "csurf != Nsurf");
1629566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(da, &ltg));
1639566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg, N, II, II));
1649566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg, Nint, IIint, IIint));
1659566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg, Nsurf, IIsurf, IIsurf));
1669566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
1679566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, N, II, PETSC_COPY_VALUES, &is));
1689566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nint, Iint, PETSC_COPY_VALUES, &isint));
1699566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nsurf, Isurf, PETSC_COPY_VALUES, &issurf));
1709566063dSJacob Faibussowitsch   PetscCall(PetscFree3(II, Iint, Isurf));
171064c8009SBarry Smith 
1729566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Aglobal, 1, &is, &is, MAT_INITIAL_MATRIX, &Aholder));
173064c8009SBarry Smith   A = *Aholder;
1749566063dSJacob Faibussowitsch   PetscCall(PetscFree(Aholder));
175064c8009SBarry Smith 
1769566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, isint, isint, MAT_INITIAL_MATRIX, &Aii));
1779566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, isint, issurf, MAT_INITIAL_MATRIX, &Ais));
1789566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, issurf, isint, MAT_INITIAL_MATRIX, &Asi));
179064c8009SBarry Smith 
180064c8009SBarry Smith   /*
181064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
182064c8009SBarry Smith   */
1839566063dSJacob Faibussowitsch   PetscCall(MatMatMult(Ais, Xsurf, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &Xint_tmp));
1849566063dSJacob Faibussowitsch   PetscCall(MatScale(Xint_tmp, -1.0));
1858e722e36SBarry Smith   if (exotic->directSolve) {
1869566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(Aii, MATSOLVERPETSC, MAT_FACTOR_LU, &iAii));
1879566063dSJacob Faibussowitsch     PetscCall(MatFactorInfoInitialize(&info));
1889566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(Aii, MATORDERINGND, &row, &col));
1899566063dSJacob Faibussowitsch     PetscCall(MatLUFactorSymbolic(iAii, Aii, row, col, &info));
1909566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&row));
1919566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&col));
1929566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(iAii, Aii, &info));
1939566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(iAii, Xint_tmp, Xint));
1949566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&iAii));
1958e722e36SBarry Smith   } else {
1968e722e36SBarry Smith     Vec          b, x;
1978e722e36SBarry Smith     PetscScalar *xint_tmp;
198064c8009SBarry Smith 
1999566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(Xint, &xint));
2009566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &x));
2019566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(Xint_tmp, &xint_tmp));
2029566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &b));
2039566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(exotic->ksp, Aii, Aii));
2048e722e36SBarry Smith     for (i = 0; i < 26; i++) {
2059566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(x, xint + i * Nint));
2069566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(b, xint_tmp + i * Nint));
2079566063dSJacob Faibussowitsch       PetscCall(KSPSolve(exotic->ksp, b, x));
2089566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(exotic->ksp, pc, x));
2099566063dSJacob Faibussowitsch       PetscCall(VecResetArray(x));
2109566063dSJacob Faibussowitsch       PetscCall(VecResetArray(b));
2118e722e36SBarry Smith     }
2129566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(Xint, &xint));
2139566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(Xint_tmp, &xint_tmp));
2149566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
2159566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&b));
2168e722e36SBarry Smith   }
2179566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xint_tmp));
2188e722e36SBarry Smith 
2198e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
2209566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xint, &rxint));
221064c8009SBarry Smith   for (i = 0; i < Nint; i++) {
222064c8009SBarry Smith     tmp = 0.0;
2231683a169SBarry Smith     for (j = 0; j < 26; j++) tmp += rxint[i + j * Nint];
2242fa5cd67SKarl Rupp 
22563a3b9bcSJacob 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));
226064c8009SBarry Smith   }
2279566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xint, &rxint));
2289566063dSJacob Faibussowitsch   /* PetscCall(MatView(Xint,PETSC_VIEWER_STDOUT_WORLD)); */
229064c8009SBarry Smith #endif
230064c8009SBarry Smith 
231064c8009SBarry Smith   /*         total vertices             total faces                                  total edges */
232064c8009SBarry 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);
233064c8009SBarry Smith 
234064c8009SBarry Smith   /*
235064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
236064c8009SBarry Smith   */
237064c8009SBarry Smith   cnt = 0;
2382fa5cd67SKarl Rupp 
2399371c9d4SSatish Balay   gl[cnt++] = 0;
240d71ae5a4SJacob Faibussowitsch   {
241d71ae5a4SJacob Faibussowitsch     gl[cnt++] = 1;
242d71ae5a4SJacob Faibussowitsch   }
2439371c9d4SSatish Balay   gl[cnt++] = m - istart - 1;
244064c8009SBarry Smith   {
2459371c9d4SSatish Balay     gl[cnt++] = mwidth;
246d71ae5a4SJacob Faibussowitsch     {
247d71ae5a4SJacob Faibussowitsch       gl[cnt++] = mwidth + 1;
248d71ae5a4SJacob Faibussowitsch     }
2499371c9d4SSatish Balay     gl[cnt++] = mwidth + m - istart - 1;
250064c8009SBarry Smith   }
2519371c9d4SSatish Balay   gl[cnt++] = mwidth * (n - jstart - 1);
252d71ae5a4SJacob Faibussowitsch   {
253d71ae5a4SJacob Faibussowitsch     gl[cnt++] = mwidth * (n - jstart - 1) + 1;
254d71ae5a4SJacob Faibussowitsch   }
2559371c9d4SSatish Balay   gl[cnt++] = mwidth * (n - jstart - 1) + m - istart - 1;
2569371c9d4SSatish Balay   {
2579371c9d4SSatish Balay     gl[cnt++] = mwidth * nwidth;
258d71ae5a4SJacob Faibussowitsch     {
259d71ae5a4SJacob Faibussowitsch       gl[cnt++] = mwidth * nwidth + 1;
260d71ae5a4SJacob Faibussowitsch     }
2619371c9d4SSatish Balay     gl[cnt++] = mwidth * nwidth + m - istart - 1;
2629371c9d4SSatish Balay     {
2639371c9d4SSatish Balay       gl[cnt++] = mwidth * nwidth + mwidth; /* these are the interior nodes */
2649371c9d4SSatish Balay       gl[cnt++] = mwidth * nwidth + mwidth + m - istart - 1;
2659371c9d4SSatish Balay     }
2669371c9d4SSatish Balay     gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1);
267d71ae5a4SJacob Faibussowitsch     {
268d71ae5a4SJacob Faibussowitsch       gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1) + 1;
269d71ae5a4SJacob Faibussowitsch     }
2709371c9d4SSatish Balay     gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1) + m - istart - 1;
2719371c9d4SSatish Balay   }
2729371c9d4SSatish Balay   gl[cnt++] = mwidth * nwidth * (p - kstart - 1);
273d71ae5a4SJacob Faibussowitsch   {
274d71ae5a4SJacob Faibussowitsch     gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + 1;
275d71ae5a4SJacob Faibussowitsch   }
2769371c9d4SSatish Balay   gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + m - istart - 1;
2779371c9d4SSatish Balay   {
2789371c9d4SSatish Balay     gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth;
279d71ae5a4SJacob Faibussowitsch     {
280d71ae5a4SJacob Faibussowitsch       gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth + 1;
281d71ae5a4SJacob Faibussowitsch     }
2829371c9d4SSatish Balay     gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth + m - istart - 1;
2839371c9d4SSatish Balay   }
2849371c9d4SSatish Balay   gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth * (n - jstart - 1);
285d71ae5a4SJacob Faibussowitsch   {
286d71ae5a4SJacob Faibussowitsch     gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth * (n - jstart - 1) + 1;
287d71ae5a4SJacob Faibussowitsch   }
2889371c9d4SSatish Balay   gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth * (n - jstart - 1) + m - istart - 1;
289064c8009SBarry Smith 
290064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
291064c8009SBarry Smith   /* convert that to global numbering and get them on all processes */
2929566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg, 26, gl, gl));
293064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
2949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(26 * mp * np * pp, &globals));
2959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(gl, 26, MPIU_INT, globals, 26, MPIU_INT, PetscObjectComm((PetscObject)da)));
296064c8009SBarry Smith 
297064c8009SBarry Smith   /* Number the coarse grid points from 0 to Ntotal */
298*eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapICreateWithSize(Ntotal / 3, &ht));
299*eec179cfSJacob Faibussowitsch   for (i = 0, cnt = 0; i < 26 * mp * np * pp; i++) {
300*eec179cfSJacob Faibussowitsch     PetscHashIter it      = 0;
301*eec179cfSJacob Faibussowitsch     PetscBool     missing = PETSC_TRUE;
302*eec179cfSJacob Faibussowitsch 
303*eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIPut(ht, globals[i] + 1, &it, &missing));
304*eec179cfSJacob Faibussowitsch     if (missing) {
305*eec179cfSJacob Faibussowitsch       ++cnt;
306*eec179cfSJacob Faibussowitsch       PetscCall(PetscHMapIIterSet(ht, it, cnt));
307*eec179cfSJacob Faibussowitsch     }
308*eec179cfSJacob Faibussowitsch   }
30963a3b9bcSJacob 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);
3109566063dSJacob Faibussowitsch   PetscCall(PetscFree(globals));
311064c8009SBarry Smith   for (i = 0; i < 26; i++) {
312*eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIGetWithDefault(ht, gl[i] + 1, 0, gl + i));
313*eec179cfSJacob Faibussowitsch     --(gl[i]);
314064c8009SBarry Smith   }
315*eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&ht));
316064c8009SBarry Smith   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
317064c8009SBarry Smith 
318064c8009SBarry Smith   /* construct global interpolation matrix */
3199566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Aglobal, &Ng, NULL));
320064c8009SBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
3219566063dSJacob Faibussowitsch     PetscCall(MatCreateAIJ(PetscObjectComm((PetscObject)da), Ng, PETSC_DECIDE, PETSC_DECIDE, Ntotal, Nint + Nsurf, NULL, Nint + Nsurf, NULL, P));
322064c8009SBarry Smith   } else {
3239566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(*P));
324064c8009SBarry Smith   }
3259566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*P, MAT_ROW_ORIENTED, PETSC_FALSE));
3269566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xint, &rxint));
3279566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*P, Nint, IIint, 26, gl, rxint, INSERT_VALUES));
3289566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xint, &rxint));
3299566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xsurf, &rxint));
3309566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*P, Nsurf, IIsurf, 26, gl, rxint, INSERT_VALUES));
3319566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xsurf, &rxint));
3329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
3339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
3349566063dSJacob Faibussowitsch   PetscCall(PetscFree2(IIint, IIsurf));
335064c8009SBarry Smith 
3368e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
337064c8009SBarry Smith   {
338064c8009SBarry Smith     Vec          x, y;
339064c8009SBarry Smith     PetscScalar *yy;
3409566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), Ng, PETSC_DETERMINE, &y));
3419566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), PETSC_DETERMINE, Ntotal, &x));
3429566063dSJacob Faibussowitsch     PetscCall(VecSet(x, 1.0));
3439566063dSJacob Faibussowitsch     PetscCall(MatMult(*P, x, y));
3449566063dSJacob Faibussowitsch     PetscCall(VecGetArray(y, &yy));
345ad540459SPierre 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]));
3469566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(y, &yy));
3479566063dSJacob Faibussowitsch     PetscCall(VecDestroy(x));
3489566063dSJacob Faibussowitsch     PetscCall(VecDestroy(y));
349064c8009SBarry Smith   }
350064c8009SBarry Smith #endif
351064c8009SBarry Smith 
3529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Aii));
3539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ais));
3549566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Asi));
3559566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3569566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
3579566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isint));
3589566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&issurf));
3599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xint));
3609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xsurf));
361064c8009SBarry Smith   PetscFunctionReturn(0);
362064c8009SBarry Smith }
363064c8009SBarry Smith 
364064c8009SBarry Smith /*
365aa219208SBarry Smith       DMDAGetFaceInterpolation - Gets the interpolation for a face based coarse space
366064c8009SBarry Smith 
367064c8009SBarry Smith */
368d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetFaceInterpolation(PC pc, DM da, PC_Exotic *exotic, Mat Aglobal, MatReuse reuse, Mat *P)
369d71ae5a4SJacob Faibussowitsch {
370064c8009SBarry 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;
371*eec179cfSJacob Faibussowitsch   PetscInt               mwidth, nwidth, pwidth, cnt, mp, np, pp, Ntotal, gl[6], *globals, Ng, *IIint, *IIsurf;
372064c8009SBarry Smith   Mat                    Xint, Xsurf, Xint_tmp;
373064c8009SBarry Smith   IS                     isint, issurf, is, row, col;
374064c8009SBarry Smith   ISLocalToGlobalMapping ltg;
375064c8009SBarry Smith   MPI_Comm               comm;
376064c8009SBarry Smith   Mat                    A, Aii, Ais, Asi, *Aholder, iAii;
377064c8009SBarry Smith   MatFactorInfo          info;
378064c8009SBarry Smith   PetscScalar           *xsurf, *xint;
3791683a169SBarry Smith   const PetscScalar     *rxint;
380064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
381064c8009SBarry Smith   PetscScalar tmp;
382064c8009SBarry Smith #endif
383*eec179cfSJacob Faibussowitsch   PetscHMapI ht;
384064c8009SBarry Smith 
385064c8009SBarry Smith   PetscFunctionBegin;
3869566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, &mp, &np, &pp, &dof, NULL, NULL, NULL, NULL, NULL));
38708401ef6SPierre Jolivet   PetscCheck(dof == 1, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only for single field problems");
38808401ef6SPierre Jolivet   PetscCheck(dim == 3, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only coded for 3d problems");
3899566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, NULL, NULL, NULL, &m, &n, &p));
3909566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &istart, &jstart, &kstart, &mwidth, &nwidth, &pwidth));
391064c8009SBarry Smith   istart = istart ? -1 : 0;
392064c8009SBarry Smith   jstart = jstart ? -1 : 0;
393064c8009SBarry Smith   kstart = kstart ? -1 : 0;
394064c8009SBarry Smith 
395064c8009SBarry Smith   /*
396064c8009SBarry Smith     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge)
397064c8009SBarry Smith     to all the local degrees of freedom (this includes the vertices, edges and faces).
398064c8009SBarry Smith 
399064c8009SBarry Smith     Xint are the subset of the interpolation into the interior
400064c8009SBarry Smith 
401064c8009SBarry Smith     Xface are the interpolation onto faces but not into the interior
402064c8009SBarry Smith 
403064c8009SBarry Smith     Xsurf are the interpolation onto the vertices and edges (the surfbasket)
404064c8009SBarry Smith                                       Xint
405064c8009SBarry Smith     Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain
406064c8009SBarry Smith                                       Xsurf
407064c8009SBarry Smith   */
408064c8009SBarry Smith   N     = (m - istart) * (n - jstart) * (p - kstart);
409064c8009SBarry Smith   Nint  = (m - 2 - istart) * (n - 2 - jstart) * (p - 2 - kstart);
410064c8009SBarry Smith   Nface = 2 * ((m - 2 - istart) * (n - 2 - jstart) + (m - 2 - istart) * (p - 2 - kstart) + (n - 2 - jstart) * (p - 2 - kstart));
411064c8009SBarry Smith   Nwire = 4 * ((m - 2 - istart) + (n - 2 - jstart) + (p - 2 - kstart)) + 8;
412064c8009SBarry Smith   Nsurf = Nface + Nwire;
4139566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nint, 6, NULL, &Xint));
4149566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nsurf, 6, NULL, &Xsurf));
4159566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(Xsurf, &xsurf));
416064c8009SBarry Smith 
417064c8009SBarry Smith   /*
418064c8009SBarry Smith      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
419064c8009SBarry Smith      Xsurf will be all zero (thus making the coarse matrix singular).
420064c8009SBarry Smith   */
42108401ef6SPierre 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");
42208401ef6SPierre 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");
42308401ef6SPierre 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");
424064c8009SBarry Smith 
425064c8009SBarry Smith   cnt = 0;
4262fa5cd67SKarl Rupp   for (j = 1; j < n - 1 - jstart; j++) {
4272fa5cd67SKarl Rupp     for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 0 * Nsurf] = 1;
428064c8009SBarry Smith   }
4292fa5cd67SKarl Rupp 
4302fa5cd67SKarl Rupp   for (k = 1; k < p - 1 - kstart; k++) {
4312fa5cd67SKarl Rupp     for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 1 * Nsurf] = 1;
4322fa5cd67SKarl Rupp     for (j = 1; j < n - 1 - jstart; j++) {
4332fa5cd67SKarl Rupp       xsurf[cnt++ + 2 * Nsurf] = 1;
4342fa5cd67SKarl Rupp       /* these are the interior nodes */
4352fa5cd67SKarl Rupp       xsurf[cnt++ + 3 * Nsurf] = 1;
4362fa5cd67SKarl Rupp     }
4372fa5cd67SKarl Rupp     for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 4 * Nsurf] = 1;
4382fa5cd67SKarl Rupp   }
4392fa5cd67SKarl Rupp   for (j = 1; j < n - 1 - jstart; j++) {
4402fa5cd67SKarl Rupp     for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 5 * Nsurf] = 1;
4412fa5cd67SKarl Rupp   }
442064c8009SBarry Smith 
443064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
444064c8009SBarry Smith   for (i = 0; i < Nsurf; i++) {
445064c8009SBarry Smith     tmp = 0.0;
4462fa5cd67SKarl Rupp     for (j = 0; j < 6; j++) tmp += xsurf[i + j * Nsurf];
4472fa5cd67SKarl Rupp 
44863a3b9bcSJacob 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));
449064c8009SBarry Smith   }
450064c8009SBarry Smith #endif
4519566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(Xsurf, &xsurf));
4529566063dSJacob Faibussowitsch   /* PetscCall(MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD));*/
453064c8009SBarry Smith 
454064c8009SBarry Smith   /*
455064c8009SBarry Smith        I are the indices for all the needed vertices (in global numbering)
456064c8009SBarry Smith        Iint are the indices for the interior values, I surf for the surface values
4577dae84e0SHong Zhang             (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it
458aa219208SBarry Smith              is NOT the local DMDA ordering.)
459064c8009SBarry Smith        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
460064c8009SBarry Smith   */
461064c8009SBarry Smith #define Endpoint(a, start, b) (a == 0 || a == (b - 1 - start))
4629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(N, &II, Nint, &Iint, Nsurf, &Isurf));
4639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nint, &IIint, Nsurf, &IIsurf));
464064c8009SBarry Smith   for (k = 0; k < p - kstart; k++) {
465064c8009SBarry Smith     for (j = 0; j < n - jstart; j++) {
466064c8009SBarry Smith       for (i = 0; i < m - istart; i++) {
467064c8009SBarry Smith         II[c++] = i + j * mwidth + k * mwidth * nwidth;
468064c8009SBarry Smith 
469064c8009SBarry Smith         if (!Endpoint(i, istart, m) && !Endpoint(j, jstart, n) && !Endpoint(k, kstart, p)) {
470064c8009SBarry Smith           IIint[cint]  = i + j * mwidth + k * mwidth * nwidth;
471064c8009SBarry Smith           Iint[cint++] = i + j * (m - istart) + k * (m - istart) * (n - jstart);
472064c8009SBarry Smith         } else {
473064c8009SBarry Smith           IIsurf[csurf]  = i + j * mwidth + k * mwidth * nwidth;
474064c8009SBarry Smith           Isurf[csurf++] = i + j * (m - istart) + k * (m - istart) * (n - jstart);
475064c8009SBarry Smith         }
476064c8009SBarry Smith       }
477064c8009SBarry Smith     }
478064c8009SBarry Smith   }
479*eec179cfSJacob Faibussowitsch #undef Endpoint
48008401ef6SPierre Jolivet   PetscCheck(c == N, PETSC_COMM_SELF, PETSC_ERR_PLIB, "c != N");
48108401ef6SPierre Jolivet   PetscCheck(cint == Nint, PETSC_COMM_SELF, PETSC_ERR_PLIB, "cint != Nint");
48208401ef6SPierre Jolivet   PetscCheck(csurf == Nsurf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "csurf != Nsurf");
4839566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(da, &ltg));
4849566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg, N, II, II));
4859566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg, Nint, IIint, IIint));
4869566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(ltg, Nsurf, IIsurf, IIsurf));
4879566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
4889566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, N, II, PETSC_COPY_VALUES, &is));
4899566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nint, Iint, PETSC_COPY_VALUES, &isint));
4909566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nsurf, Isurf, PETSC_COPY_VALUES, &issurf));
4919566063dSJacob Faibussowitsch   PetscCall(PetscFree3(II, Iint, Isurf));
492064c8009SBarry Smith 
4939566063dSJacob Faibussowitsch   PetscCall(ISSort(is));
4949566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Aglobal, 1, &is, &is, MAT_INITIAL_MATRIX, &Aholder));
495064c8009SBarry Smith   A = *Aholder;
4969566063dSJacob Faibussowitsch   PetscCall(PetscFree(Aholder));
497064c8009SBarry Smith 
4989566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, isint, isint, MAT_INITIAL_MATRIX, &Aii));
4999566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, isint, issurf, MAT_INITIAL_MATRIX, &Ais));
5009566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, issurf, isint, MAT_INITIAL_MATRIX, &Asi));
501064c8009SBarry Smith 
502064c8009SBarry Smith   /*
503064c8009SBarry Smith      Solve for the interpolation onto the interior Xint
504064c8009SBarry Smith   */
5059566063dSJacob Faibussowitsch   PetscCall(MatMatMult(Ais, Xsurf, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &Xint_tmp));
5069566063dSJacob Faibussowitsch   PetscCall(MatScale(Xint_tmp, -1.0));
507064c8009SBarry Smith 
5088e722e36SBarry Smith   if (exotic->directSolve) {
5099566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(Aii, MATSOLVERPETSC, MAT_FACTOR_LU, &iAii));
5109566063dSJacob Faibussowitsch     PetscCall(MatFactorInfoInitialize(&info));
5119566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(Aii, MATORDERINGND, &row, &col));
5129566063dSJacob Faibussowitsch     PetscCall(MatLUFactorSymbolic(iAii, Aii, row, col, &info));
5139566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&row));
5149566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&col));
5159566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(iAii, Aii, &info));
5169566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(iAii, Xint_tmp, Xint));
5179566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&iAii));
518064c8009SBarry Smith   } else {
519064c8009SBarry Smith     Vec          b, x;
520064c8009SBarry Smith     PetscScalar *xint_tmp;
521064c8009SBarry Smith 
5229566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(Xint, &xint));
5239566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &x));
5249566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(Xint_tmp, &xint_tmp));
5259566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &b));
5269566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(exotic->ksp, Aii, Aii));
527064c8009SBarry Smith     for (i = 0; i < 6; i++) {
5289566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(x, xint + i * Nint));
5299566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(b, xint_tmp + i * Nint));
5309566063dSJacob Faibussowitsch       PetscCall(KSPSolve(exotic->ksp, b, x));
5319566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(exotic->ksp, pc, x));
5329566063dSJacob Faibussowitsch       PetscCall(VecResetArray(x));
5339566063dSJacob Faibussowitsch       PetscCall(VecResetArray(b));
534064c8009SBarry Smith     }
5359566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(Xint, &xint));
5369566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(Xint_tmp, &xint_tmp));
5379566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
5389566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&b));
539064c8009SBarry Smith   }
5409566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Xint_tmp));
541064c8009SBarry Smith 
542064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
5439566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Xint, &rxint));
544064c8009SBarry Smith   for (i = 0; i < Nint; i++) {
545064c8009SBarry Smith     tmp = 0.0;
5461683a169SBarry Smith     for (j = 0; j < 6; j++) tmp += rxint[i + j * Nint];
5472fa5cd67SKarl Rupp 
54863a3b9bcSJacob 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));
549064c8009SBarry Smith   }
5509566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Xint, &rxint));
5519566063dSJacob Faibussowitsch   /* PetscCall(MatView(Xint,PETSC_VIEWER_STDOUT_WORLD)); */
552064c8009SBarry Smith #endif
553064c8009SBarry Smith 
554064c8009SBarry Smith   /*         total faces    */
555064c8009SBarry Smith   Ntotal = mp * np * (pp + 1) + mp * pp * (np + 1) + np * pp * (mp + 1);
556064c8009SBarry Smith 
557064c8009SBarry Smith   /*
558064c8009SBarry Smith       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
559064c8009SBarry Smith   */
560064c8009SBarry Smith   cnt = 0;
561d71ae5a4SJacob Faibussowitsch   {
562d71ae5a4SJacob Faibussowitsch     gl[cnt++] = mwidth + 1;
563d71ae5a4SJacob Faibussowitsch   }
5649371c9d4SSatish Balay   {{gl[cnt++] = mwidth * nwidth + 1;
5659371c9d4SSatish Balay }
566064c8009SBarry Smith {
5679371c9d4SSatish Balay   gl[cnt++] = mwidth * nwidth + mwidth; /* these are the interior nodes */
5689371c9d4SSatish Balay   gl[cnt++] = mwidth * nwidth + mwidth + m - istart - 1;
5699371c9d4SSatish Balay }
570d71ae5a4SJacob Faibussowitsch {
571d71ae5a4SJacob Faibussowitsch   gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1) + 1;
572064c8009SBarry Smith }
573d71ae5a4SJacob Faibussowitsch }
574d71ae5a4SJacob Faibussowitsch {
575d71ae5a4SJacob Faibussowitsch   gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth + 1;
576d71ae5a4SJacob Faibussowitsch }
577064c8009SBarry Smith 
578064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
579064c8009SBarry Smith /* convert that to global numbering and get them on all processes */
5809566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, 6, gl, gl));
581064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
5829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(6 * mp * np * pp, &globals));
5839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(gl, 6, MPIU_INT, globals, 6, MPIU_INT, PetscObjectComm((PetscObject)da)));
584064c8009SBarry Smith 
585064c8009SBarry Smith /* Number the coarse grid points from 0 to Ntotal */
586*eec179cfSJacob Faibussowitsch PetscCall(PetscHMapICreateWithSize(Ntotal / 3, &ht));
587*eec179cfSJacob Faibussowitsch for (i = 0, cnt = 0; i < 6 * mp * np * pp; i++) {
588*eec179cfSJacob Faibussowitsch   PetscHashIter it      = 0;
589*eec179cfSJacob Faibussowitsch   PetscBool     missing = PETSC_TRUE;
590*eec179cfSJacob Faibussowitsch 
591*eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapIPut(ht, globals[i] + 1, &it, &missing));
592*eec179cfSJacob Faibussowitsch   if (missing) {
593*eec179cfSJacob Faibussowitsch     ++cnt;
594*eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIIterSet(ht, it, cnt));
595*eec179cfSJacob Faibussowitsch   }
596*eec179cfSJacob Faibussowitsch }
59763a3b9bcSJacob 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);
5989566063dSJacob Faibussowitsch PetscCall(PetscFree(globals));
599064c8009SBarry Smith for (i = 0; i < 6; i++) {
600*eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapIGetWithDefault(ht, gl[i] + 1, 0, gl + i));
601*eec179cfSJacob Faibussowitsch   --(gl[i]);
602064c8009SBarry Smith }
603*eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&ht));
604064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
605064c8009SBarry Smith 
606064c8009SBarry Smith /* construct global interpolation matrix */
6079566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Aglobal, &Ng, NULL));
608064c8009SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {
6099566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(PetscObjectComm((PetscObject)da), Ng, PETSC_DECIDE, PETSC_DECIDE, Ntotal, Nint + Nsurf, NULL, Nint, NULL, P));
610064c8009SBarry Smith } else {
6119566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(*P));
612064c8009SBarry Smith }
6139566063dSJacob Faibussowitsch PetscCall(MatSetOption(*P, MAT_ROW_ORIENTED, PETSC_FALSE));
6149566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xint, &rxint));
6159566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, Nint, IIint, 6, gl, rxint, INSERT_VALUES));
6169566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xint, &rxint));
6179566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xsurf, &rxint));
6189566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, Nsurf, IIsurf, 6, gl, rxint, INSERT_VALUES));
6199566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xsurf, &rxint));
6209566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
6219566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
6229566063dSJacob Faibussowitsch PetscCall(PetscFree2(IIint, IIsurf));
623064c8009SBarry Smith 
624064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
625064c8009SBarry Smith {
626064c8009SBarry Smith   Vec          x, y;
627064c8009SBarry Smith   PetscScalar *yy;
6289566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), Ng, PETSC_DETERMINE, &y));
6299566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), PETSC_DETERMINE, Ntotal, &x));
6309566063dSJacob Faibussowitsch   PetscCall(VecSet(x, 1.0));
6319566063dSJacob Faibussowitsch   PetscCall(MatMult(*P, x, y));
6329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
633ad540459SPierre 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]));
6349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
6359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(x));
6369566063dSJacob Faibussowitsch   PetscCall(VecDestroy(y));
637064c8009SBarry Smith }
638064c8009SBarry Smith #endif
639064c8009SBarry Smith 
6409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aii));
6419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ais));
6429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Asi));
6439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
6449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is));
6459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isint));
6469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&issurf));
6479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xint));
6489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xsurf));
649064c8009SBarry Smith PetscFunctionReturn(0);
650064c8009SBarry Smith }
651174b6946SBarry Smith 
6527233f9f0SBarry Smith /*@
6537233f9f0SBarry Smith    PCExoticSetType - Sets the type of coarse grid interpolation to use
6547233f9f0SBarry Smith 
655f1580f4eSBarry Smith    Logically Collective on pc
6567233f9f0SBarry Smith 
6577233f9f0SBarry Smith    Input Parameters:
6587233f9f0SBarry Smith +  pc - the preconditioner context
6597233f9f0SBarry Smith -  type - either PC_EXOTIC_FACE or PC_EXOTIC_WIREBASKET (defaults to face)
6607233f9f0SBarry Smith 
66195452b02SPatrick Sanan    Notes:
66295452b02SPatrick Sanan     The face based interpolation has 1 degree of freedom per face and ignores the
663563e08c6SBarry Smith      edge and vertex values completely in the coarse problem. For any seven point
664563e08c6SBarry Smith      stencil the interpolation of a constant on all faces into the interior is that constant.
665563e08c6SBarry Smith 
666563e08c6SBarry Smith      The wirebasket interpolation has 1 degree of freedom per vertex, per edge and
667563e08c6SBarry Smith      per face. A constant on the subdomain boundary is interpolated as that constant
668563e08c6SBarry Smith      in the interior of the domain.
669563e08c6SBarry Smith 
670563e08c6SBarry Smith      The coarse grid matrix is obtained via the Galerkin computation A_c = R A R^T, hence
671563e08c6SBarry Smith      if A is nonsingular A_c is also nonsingular.
672563e08c6SBarry Smith 
673563e08c6SBarry Smith      Both interpolations are suitable for only scalar problems.
674563e08c6SBarry Smith 
6757233f9f0SBarry Smith    Level: intermediate
6767233f9f0SBarry Smith 
677db781477SPatrick Sanan .seealso: `PCEXOTIC`, `PCExoticType()`
6787233f9f0SBarry Smith @*/
679d71ae5a4SJacob Faibussowitsch PetscErrorCode PCExoticSetType(PC pc, PCExoticType type)
680d71ae5a4SJacob Faibussowitsch {
6817233f9f0SBarry Smith   PetscFunctionBegin;
6820700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
683c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc, type, 2);
684cac4c232SBarry Smith   PetscTryMethod(pc, "PCExoticSetType_C", (PC, PCExoticType), (pc, type));
6857233f9f0SBarry Smith   PetscFunctionReturn(0);
6867233f9f0SBarry Smith }
6877233f9f0SBarry Smith 
688d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCExoticSetType_Exotic(PC pc, PCExoticType type)
689d71ae5a4SJacob Faibussowitsch {
690f3fbd535SBarry Smith   PC_MG     *mg  = (PC_MG *)pc->data;
69131567311SBarry Smith   PC_Exotic *ctx = (PC_Exotic *)mg->innerctx;
6927233f9f0SBarry Smith 
6937233f9f0SBarry Smith   PetscFunctionBegin;
6947233f9f0SBarry Smith   ctx->type = type;
6957233f9f0SBarry Smith   PetscFunctionReturn(0);
6967233f9f0SBarry Smith }
6977233f9f0SBarry Smith 
698d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp_Exotic(PC pc)
699d71ae5a4SJacob Faibussowitsch {
70096bdf778SBarry Smith   Mat        A;
701f3fbd535SBarry Smith   PC_MG     *mg    = (PC_MG *)pc->data;
70231567311SBarry Smith   PC_Exotic *ex    = (PC_Exotic *)mg->innerctx;
70396bdf778SBarry Smith   MatReuse   reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX;
704174b6946SBarry Smith 
705174b6946SBarry Smith   PetscFunctionBegin;
70628b400f6SJacob Faibussowitsch   PetscCheck(pc->dm, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Need to call PCSetDM() before using this PC");
7079566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &A));
7080fdf79fbSJacob 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);
7097233f9f0SBarry Smith   if (ex->type == PC_EXOTIC_FACE) {
7109566063dSJacob Faibussowitsch     PetscCall(DMDAGetFaceInterpolation(pc, pc->dm, ex, A, reuse, &ex->P));
7110fdf79fbSJacob Faibussowitsch   } else /* if (ex->type == PC_EXOTIC_WIREBASKET) */ {
7129566063dSJacob Faibussowitsch     PetscCall(DMDAGetWireBasketInterpolation(pc, pc->dm, ex, A, reuse, &ex->P));
7130fdf79fbSJacob Faibussowitsch   }
7149566063dSJacob Faibussowitsch   PetscCall(PCMGSetInterpolation(pc, 1, ex->P));
715d0660788SBarry Smith   /* if PC has attached DM we must remove it or the PCMG will use it to compute incorrect sized vectors and interpolations */
7169566063dSJacob Faibussowitsch   PetscCall(PCSetDM(pc, NULL));
7179566063dSJacob Faibussowitsch   PetscCall(PCSetUp_MG(pc));
718174b6946SBarry Smith   PetscFunctionReturn(0);
719174b6946SBarry Smith }
720174b6946SBarry Smith 
721d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy_Exotic(PC pc)
722d71ae5a4SJacob Faibussowitsch {
723f3fbd535SBarry Smith   PC_MG     *mg  = (PC_MG *)pc->data;
72431567311SBarry Smith   PC_Exotic *ctx = (PC_Exotic *)mg->innerctx;
725174b6946SBarry Smith 
726174b6946SBarry Smith   PetscFunctionBegin;
7279566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->P));
7289566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ctx->ksp));
7299566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
7302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCExoticSetType_C", NULL));
7319566063dSJacob Faibussowitsch   PetscCall(PCDestroy_MG(pc));
732174b6946SBarry Smith   PetscFunctionReturn(0);
733174b6946SBarry Smith }
734174b6946SBarry Smith 
735d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView_Exotic(PC pc, PetscViewer viewer)
736d71ae5a4SJacob Faibussowitsch {
737f3fbd535SBarry Smith   PC_MG     *mg = (PC_MG *)pc->data;
738ace3abfcSBarry Smith   PetscBool  iascii;
73931567311SBarry Smith   PC_Exotic *ctx = (PC_Exotic *)mg->innerctx;
7407233f9f0SBarry Smith 
7417233f9f0SBarry Smith   PetscFunctionBegin;
7429566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
7437233f9f0SBarry Smith   if (iascii) {
7449566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Exotic type = %s\n", PCExoticTypes[ctx->type]));
7458e722e36SBarry Smith     if (ctx->directSolve) {
7469566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "      Using direct solver to construct interpolation\n"));
7478e722e36SBarry Smith     } else {
7488e722e36SBarry Smith       PetscViewer sviewer;
7498e722e36SBarry Smith       PetscMPIInt rank;
7508e722e36SBarry Smith 
7519566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "      Using iterative solver to construct interpolation\n"));
7529566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
7539566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer)); /* should not need to push this twice? */
7549566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
7559566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
75648a46eb9SPierre Jolivet       if (rank == 0) PetscCall(KSPView(ctx->ksp, sviewer));
7579566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
7589566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
7599566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
7608e722e36SBarry Smith     }
7617233f9f0SBarry Smith   }
7629566063dSJacob Faibussowitsch   PetscCall(PCView_MG(pc, viewer));
7637233f9f0SBarry Smith   PetscFunctionReturn(0);
7647233f9f0SBarry Smith }
7657233f9f0SBarry Smith 
766d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFromOptions_Exotic(PC pc, PetscOptionItems *PetscOptionsObject)
767d71ae5a4SJacob Faibussowitsch {
768ace3abfcSBarry Smith   PetscBool    flg;
769f3fbd535SBarry Smith   PC_MG       *mg = (PC_MG *)pc->data;
7707233f9f0SBarry Smith   PCExoticType mgctype;
77131567311SBarry Smith   PC_Exotic   *ctx = (PC_Exotic *)mg->innerctx;
7727233f9f0SBarry Smith 
7737233f9f0SBarry Smith   PetscFunctionBegin;
774d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Exotic coarse space options");
7759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_exotic_type", "face or wirebasket", "PCExoticSetType", PCExoticTypes, (PetscEnum)ctx->type, (PetscEnum *)&mgctype, &flg));
7761baa6e33SBarry Smith   if (flg) PetscCall(PCExoticSetType(pc, mgctype));
7779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_exotic_direct_solver", "use direct solver to construct interpolation", "None", ctx->directSolve, &ctx->directSolve, NULL));
7788e722e36SBarry Smith   if (!ctx->directSolve) {
7798e722e36SBarry Smith     if (!ctx->ksp) {
7808e722e36SBarry Smith       const char *prefix;
7819566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PETSC_COMM_SELF, &ctx->ksp));
7829566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(ctx->ksp, pc->erroriffailure));
7839566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->ksp, (PetscObject)pc, 1));
7849566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
7859566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(ctx->ksp, prefix));
7869566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(ctx->ksp, "exotic_"));
7878e722e36SBarry Smith     }
7889566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ctx->ksp));
7898e722e36SBarry Smith   }
790d0609cedSBarry Smith   PetscOptionsHeadEnd();
7917233f9f0SBarry Smith   PetscFunctionReturn(0);
7927233f9f0SBarry Smith }
7937233f9f0SBarry Smith 
794174b6946SBarry Smith /*MC
7957233f9f0SBarry Smith      PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces
796174b6946SBarry Smith 
797f1580f4eSBarry Smith      This uses the `PCMG` infrastructure restricted to two levels and the face and wirebasket based coarse
79824c3aa18SBarry Smith    grid spaces.
79924c3aa18SBarry Smith 
800f1580f4eSBarry Smith    Level: advanced
801f1580f4eSBarry Smith 
80295452b02SPatrick Sanan    Notes:
803f1580f4eSBarry Smith    Must be used with `DMDA`
804f1580f4eSBarry Smith 
805f1580f4eSBarry 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`
80624c3aa18SBarry Smith 
80796a0c994SBarry Smith    References:
808606c0280SSatish Balay +  * - These coarse grid spaces originate in the work of Bramble, Pasciak  and Schatz, "The Construction
80996a0c994SBarry Smith    of Preconditioners for Elliptic Problems by Substructing IV", Mathematics of Computation, volume 53, 1989.
810606c0280SSatish Balay .  * - They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith,
8114f02bc6aSBarry Smith    New York University, 1990.
812606c0280SSatish Balay .  * - They were then explored in great detail in Dryja, Smith, Widlund, "Schwarz Analysis
8133b65e785SBarry Smith    of Iterative Substructuring Methods for Elliptic Problems in Three Dimensions, SIAM Journal on Numerical
81496a0c994SBarry Smith    Analysis, volume 31. 1994. These were developed in the context of iterative substructuring preconditioners.
815606c0280SSatish Balay .  * - They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund.
8163b65e785SBarry Smith    They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example,
8173b65e785SBarry Smith    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Extending theory for domain decomposition algorithms to irregular subdomains. In Ulrich Langer, Marco
8183b65e785SBarry Smith    Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
8193b65e785SBarry Smith    of the 17th International Conference on Domain Decomposition Methods in
82096a0c994SBarry Smith    Science and Engineering, held in Strobl, Austria, 2006, number 60 in
82196a0c994SBarry Smith    Springer Verlag, Lecture Notes in Computational Science and Engineering, 2007.
822606c0280SSatish Balay .  * -  Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. A family of energy minimizing coarse spaces for overlapping Schwarz preconditioners. In Ulrich Langer,
8233b65e785SBarry Smith    Marco Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
8243b65e785SBarry Smith    of the 17th International Conference on Domain Decomposition Methods
82596a0c994SBarry Smith    in Science and Engineering, held in Strobl, Austria, 2006, number 60 in
82696a0c994SBarry Smith    Springer Verlag, Lecture Notes in Computational Science and Engineering, 2007
827606c0280SSatish Balay .  * - Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Domain decomposition
8283b65e785SBarry Smith    for less regular subdomains: Overlapping Schwarz in two dimensions. SIAM J.
82996a0c994SBarry Smith    Numer. Anal., 46(4), 2008.
830606c0280SSatish Balay -  * - Clark R. Dohrmann and Olof B. Widlund. An overlapping Schwarz
8313b65e785SBarry Smith    algorithm for almost incompressible elasticity. Technical Report
83296a0c994SBarry Smith    TR2008 912, Department of Computer Science, Courant Institute
8333b65e785SBarry Smith    of Mathematical Sciences, New York University, May 2008. URL:
8347233f9f0SBarry Smith 
835f1580f4eSBarry Smith    The usual `PCMG` options are supported, such as -mg_levels_pc_type <type> -mg_coarse_pc_type <type> and  -pc_mg_type <type>
836174b6946SBarry Smith 
837db781477SPatrick Sanan .seealso: `PCMG`, `PCSetDM()`, `PCExoticType`, `PCExoticSetType()`
838174b6946SBarry Smith M*/
839174b6946SBarry Smith 
840d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Exotic(PC pc)
841d71ae5a4SJacob Faibussowitsch {
8427233f9f0SBarry Smith   PC_Exotic *ex;
843f3fbd535SBarry Smith   PC_MG     *mg;
844174b6946SBarry Smith 
845174b6946SBarry Smith   PetscFunctionBegin;
846f91d8e95SBarry Smith   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
847dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, destroy);
8480a545947SLisandro Dalcin   pc->data = NULL;
849dbbe0bcdSBarry Smith 
8509566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PetscObject)pc)->type_name));
8510a545947SLisandro Dalcin   ((PetscObject)pc)->type_name = NULL;
852f91d8e95SBarry Smith 
8539566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCMG));
8549566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pc, 2, NULL));
8559566063dSJacob Faibussowitsch   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_PMAT));
8569371c9d4SSatish Balay   PetscCall(PetscNew(&ex));
8577233f9f0SBarry Smith   ex->type     = PC_EXOTIC_FACE;
858f3fbd535SBarry Smith   mg           = (PC_MG *)pc->data;
85931567311SBarry Smith   mg->innerctx = ex;
8607233f9f0SBarry Smith 
8617233f9f0SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_Exotic;
8627233f9f0SBarry Smith   pc->ops->view           = PCView_Exotic;
8637233f9f0SBarry Smith   pc->ops->destroy        = PCDestroy_Exotic;
8646c699258SBarry Smith   pc->ops->setup          = PCSetUp_Exotic;
8652fa5cd67SKarl Rupp 
8669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCExoticSetType_C", PCExoticSetType_Exotic));
867174b6946SBarry Smith   PetscFunctionReturn(0);
868174b6946SBarry Smith }
869