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 */
DMDAGetWireBasketInterpolation(PC pc,DM da,PC_Exotic * exotic,Mat Aglobal,MatReuse reuse,Mat * P)1766976f2fSJacob Faibussowitsch static PetscErrorCode DMDAGetWireBasketInterpolation(PC pc, DM da, PC_Exotic *exotic, Mat Aglobal, MatReuse reuse, Mat *P)
18d71ae5a4SJacob Faibussowitsch {
19064c8009SBarry 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;
20eec179cfSJacob Faibussowitsch PetscInt mwidth, nwidth, pwidth, cnt, mp, np, pp, Ntotal, gl[26], *globals, Ng, *IIint, *IIsurf;
21064c8009SBarry Smith Mat Xint, Xsurf, Xint_tmp;
22064c8009SBarry Smith IS isint, issurf, is, row, col;
23064c8009SBarry Smith ISLocalToGlobalMapping ltg;
24064c8009SBarry Smith MPI_Comm comm;
25064c8009SBarry Smith Mat A, Aii, Ais, Asi, *Aholder, iAii;
26064c8009SBarry Smith MatFactorInfo info;
27064c8009SBarry Smith PetscScalar *xsurf, *xint;
281683a169SBarry Smith const PetscScalar *rxint;
298e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
30064c8009SBarry Smith PetscScalar tmp;
31064c8009SBarry Smith #endif
32eec179cfSJacob Faibussowitsch PetscHMapI ht = NULL;
33064c8009SBarry Smith
34064c8009SBarry Smith PetscFunctionBegin;
359566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, &mp, &np, &pp, &dof, NULL, NULL, NULL, NULL, NULL));
3608401ef6SPierre Jolivet PetscCheck(dof == 1, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only for single field problems");
3708401ef6SPierre Jolivet PetscCheck(dim == 3, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only coded for 3d problems");
389566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, NULL, NULL, NULL, &m, &n, &p));
399566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &istart, &jstart, &kstart, &mwidth, &nwidth, &pwidth));
40064c8009SBarry Smith istart = istart ? -1 : 0;
41064c8009SBarry Smith jstart = jstart ? -1 : 0;
42064c8009SBarry Smith kstart = kstart ? -1 : 0;
43064c8009SBarry Smith
44064c8009SBarry Smith /*
454bc3789fSBarry Smith the columns of P are the interpolation of each coarse (wirebasket) grid point (one for each face, vertex and edge)
46064c8009SBarry Smith to all the local degrees of freedom (this includes the vertices, edges and faces).
47064c8009SBarry Smith
48064c8009SBarry Smith Xint are the subset of the interpolation into the interior
49064c8009SBarry Smith
50064c8009SBarry Smith Xface are the interpolation onto faces but not into the interior
51064c8009SBarry Smith
524bc3789fSBarry Smith Xsurf are the interpolation onto the vertices and edges (the wirebasket)
53064c8009SBarry Smith Xint
54064c8009SBarry Smith Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain
55064c8009SBarry Smith Xsurf
56064c8009SBarry Smith */
57064c8009SBarry Smith N = (m - istart) * (n - jstart) * (p - kstart);
58064c8009SBarry Smith Nint = (m - 2 - istart) * (n - 2 - jstart) * (p - 2 - kstart);
59064c8009SBarry Smith Nface = 2 * ((m - 2 - istart) * (n - 2 - jstart) + (m - 2 - istart) * (p - 2 - kstart) + (n - 2 - jstart) * (p - 2 - kstart));
60064c8009SBarry Smith Nwire = 4 * ((m - 2 - istart) + (n - 2 - jstart) + (p - 2 - kstart)) + 8;
61064c8009SBarry Smith Nsurf = Nface + Nwire;
629566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nint, 26, NULL, &Xint));
639566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nsurf, 26, NULL, &Xsurf));
649566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(Xsurf, &xsurf));
65064c8009SBarry Smith
66064c8009SBarry Smith /*
67064c8009SBarry Smith Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
68064c8009SBarry Smith Xsurf will be all zero (thus making the coarse matrix singular).
69064c8009SBarry Smith */
7008401ef6SPierre 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");
7108401ef6SPierre 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");
7208401ef6SPierre 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");
73064c8009SBarry Smith
74064c8009SBarry Smith cnt = 0;
752fa5cd67SKarl Rupp
762fa5cd67SKarl Rupp xsurf[cnt++] = 1;
772fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + Nsurf] = 1;
782fa5cd67SKarl Rupp xsurf[cnt++ + 2 * Nsurf] = 1;
792fa5cd67SKarl Rupp
802fa5cd67SKarl Rupp for (j = 1; j < n - 1 - jstart; j++) {
812fa5cd67SKarl Rupp xsurf[cnt++ + 3 * Nsurf] = 1;
822fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 4 * Nsurf] = 1;
832fa5cd67SKarl Rupp xsurf[cnt++ + 5 * Nsurf] = 1;
84064c8009SBarry Smith }
852fa5cd67SKarl Rupp
862fa5cd67SKarl Rupp xsurf[cnt++ + 6 * Nsurf] = 1;
872fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 7 * Nsurf] = 1;
882fa5cd67SKarl Rupp xsurf[cnt++ + 8 * Nsurf] = 1;
892fa5cd67SKarl Rupp
902fa5cd67SKarl Rupp for (k = 1; k < p - 1 - kstart; k++) {
912fa5cd67SKarl Rupp xsurf[cnt++ + 9 * Nsurf] = 1;
922fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 10 * Nsurf] = 1;
932fa5cd67SKarl Rupp xsurf[cnt++ + 11 * Nsurf] = 1;
942fa5cd67SKarl Rupp
952fa5cd67SKarl Rupp for (j = 1; j < n - 1 - jstart; j++) {
962fa5cd67SKarl Rupp xsurf[cnt++ + 12 * Nsurf] = 1;
972fa5cd67SKarl Rupp /* these are the interior nodes */
982fa5cd67SKarl Rupp xsurf[cnt++ + 13 * Nsurf] = 1;
992fa5cd67SKarl Rupp }
1002fa5cd67SKarl Rupp
1012fa5cd67SKarl Rupp xsurf[cnt++ + 14 * Nsurf] = 1;
1022fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 15 * Nsurf] = 1;
1032fa5cd67SKarl Rupp xsurf[cnt++ + 16 * Nsurf] = 1;
1042fa5cd67SKarl Rupp }
1052fa5cd67SKarl Rupp
1062fa5cd67SKarl Rupp xsurf[cnt++ + 17 * Nsurf] = 1;
1072fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 18 * Nsurf] = 1;
1082fa5cd67SKarl Rupp xsurf[cnt++ + 19 * Nsurf] = 1;
1092fa5cd67SKarl Rupp
1102fa5cd67SKarl Rupp for (j = 1; j < n - 1 - jstart; j++) {
1112fa5cd67SKarl Rupp xsurf[cnt++ + 20 * Nsurf] = 1;
1122fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 21 * Nsurf] = 1;
1132fa5cd67SKarl Rupp xsurf[cnt++ + 22 * Nsurf] = 1;
1142fa5cd67SKarl Rupp }
1152fa5cd67SKarl Rupp
1162fa5cd67SKarl Rupp xsurf[cnt++ + 23 * Nsurf] = 1;
1172fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 24 * Nsurf] = 1;
1182fa5cd67SKarl Rupp xsurf[cnt++ + 25 * Nsurf] = 1;
1192fa5cd67SKarl Rupp
1208e722e36SBarry Smith /* interpolations only sum to 1 when using direct solver */
1218e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
122064c8009SBarry Smith for (i = 0; i < Nsurf; i++) {
123064c8009SBarry Smith tmp = 0.0;
1242fa5cd67SKarl Rupp for (j = 0; j < 26; j++) tmp += xsurf[i + j * Nsurf];
12563a3b9bcSJacob 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));
126064c8009SBarry Smith }
127064c8009SBarry Smith #endif
1289566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(Xsurf, &xsurf));
1299566063dSJacob Faibussowitsch /* PetscCall(MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD));*/
130064c8009SBarry Smith
131064c8009SBarry Smith /*
132064c8009SBarry Smith I are the indices for all the needed vertices (in global numbering)
133064c8009SBarry Smith Iint are the indices for the interior values, I surf for the surface values
1347dae84e0SHong Zhang (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it
135aa219208SBarry Smith is NOT the local DMDA ordering.)
136064c8009SBarry Smith IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
137064c8009SBarry Smith */
138064c8009SBarry Smith #define Endpoint(a, start, b) (a == 0 || a == (b - 1 - start))
1399566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(N, &II, Nint, &Iint, Nsurf, &Isurf));
1409566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nint, &IIint, Nsurf, &IIsurf));
141064c8009SBarry Smith for (k = 0; k < p - kstart; k++) {
142064c8009SBarry Smith for (j = 0; j < n - jstart; j++) {
143064c8009SBarry Smith for (i = 0; i < m - istart; i++) {
144064c8009SBarry Smith II[c++] = i + j * mwidth + k * mwidth * nwidth;
145064c8009SBarry Smith
146064c8009SBarry Smith if (!Endpoint(i, istart, m) && !Endpoint(j, jstart, n) && !Endpoint(k, kstart, p)) {
147064c8009SBarry Smith IIint[cint] = i + j * mwidth + k * mwidth * nwidth;
148064c8009SBarry Smith Iint[cint++] = i + j * (m - istart) + k * (m - istart) * (n - jstart);
149064c8009SBarry Smith } else {
150064c8009SBarry Smith IIsurf[csurf] = i + j * mwidth + k * mwidth * nwidth;
151064c8009SBarry Smith Isurf[csurf++] = i + j * (m - istart) + k * (m - istart) * (n - jstart);
152064c8009SBarry Smith }
153064c8009SBarry Smith }
154064c8009SBarry Smith }
155064c8009SBarry Smith }
156eec179cfSJacob Faibussowitsch #undef Endpoint
15708401ef6SPierre Jolivet PetscCheck(c == N, PETSC_COMM_SELF, PETSC_ERR_PLIB, "c != N");
15808401ef6SPierre Jolivet PetscCheck(cint == Nint, PETSC_COMM_SELF, PETSC_ERR_PLIB, "cint != Nint");
15908401ef6SPierre Jolivet PetscCheck(csurf == Nsurf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "csurf != Nsurf");
1609566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(da, <g));
1619566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, N, II, II));
1629566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, Nint, IIint, IIint));
1639566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, Nsurf, IIsurf, IIsurf));
1649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
1659566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, N, II, PETSC_COPY_VALUES, &is));
1669566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nint, Iint, PETSC_COPY_VALUES, &isint));
1679566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nsurf, Isurf, PETSC_COPY_VALUES, &issurf));
1689566063dSJacob Faibussowitsch PetscCall(PetscFree3(II, Iint, Isurf));
169064c8009SBarry Smith
1709566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(Aglobal, 1, &is, &is, MAT_INITIAL_MATRIX, &Aholder));
171064c8009SBarry Smith A = *Aholder;
1729566063dSJacob Faibussowitsch PetscCall(PetscFree(Aholder));
173064c8009SBarry Smith
1749566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, isint, isint, MAT_INITIAL_MATRIX, &Aii));
1759566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, isint, issurf, MAT_INITIAL_MATRIX, &Ais));
1769566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, issurf, isint, MAT_INITIAL_MATRIX, &Asi));
177064c8009SBarry Smith
178064c8009SBarry Smith /*
179064c8009SBarry Smith Solve for the interpolation onto the interior Xint
180064c8009SBarry Smith */
1819566063dSJacob Faibussowitsch PetscCall(MatMatMult(Ais, Xsurf, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &Xint_tmp));
1829566063dSJacob Faibussowitsch PetscCall(MatScale(Xint_tmp, -1.0));
1838e722e36SBarry Smith if (exotic->directSolve) {
1849566063dSJacob Faibussowitsch PetscCall(MatGetFactor(Aii, MATSOLVERPETSC, MAT_FACTOR_LU, &iAii));
1859566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info));
1869566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(Aii, MATORDERINGND, &row, &col));
1879566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(iAii, Aii, row, col, &info));
1889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row));
1899566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col));
1909566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(iAii, Aii, &info));
1919566063dSJacob Faibussowitsch PetscCall(MatMatSolve(iAii, Xint_tmp, Xint));
1929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&iAii));
1938e722e36SBarry Smith } else {
1948e722e36SBarry Smith Vec b, x;
1958e722e36SBarry Smith PetscScalar *xint_tmp;
196064c8009SBarry Smith
1979566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(Xint, &xint));
1989566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &x));
1999566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(Xint_tmp, &xint_tmp));
2009566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &b));
2019566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(exotic->ksp, Aii, Aii));
2028e722e36SBarry Smith for (i = 0; i < 26; i++) {
2039566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(x, xint + i * Nint));
2049566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(b, xint_tmp + i * Nint));
2059566063dSJacob Faibussowitsch PetscCall(KSPSolve(exotic->ksp, b, x));
2069566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(exotic->ksp, pc, x));
2079566063dSJacob Faibussowitsch PetscCall(VecResetArray(x));
2089566063dSJacob Faibussowitsch PetscCall(VecResetArray(b));
2098e722e36SBarry Smith }
2109566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(Xint, &xint));
2119566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(Xint_tmp, &xint_tmp));
2129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
2139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b));
2148e722e36SBarry Smith }
2159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xint_tmp));
2168e722e36SBarry Smith
2178e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
2189566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xint, &rxint));
219064c8009SBarry Smith for (i = 0; i < Nint; i++) {
220064c8009SBarry Smith tmp = 0.0;
2211683a169SBarry Smith for (j = 0; j < 26; j++) tmp += rxint[i + j * Nint];
2222fa5cd67SKarl Rupp
22363a3b9bcSJacob 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));
224064c8009SBarry Smith }
2259566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xint, &rxint));
2269566063dSJacob Faibussowitsch /* PetscCall(MatView(Xint,PETSC_VIEWER_STDOUT_WORLD)); */
227064c8009SBarry Smith #endif
228064c8009SBarry Smith
229064c8009SBarry Smith /* total vertices total faces total edges */
230064c8009SBarry 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);
231064c8009SBarry Smith
232064c8009SBarry Smith /*
233064c8009SBarry Smith For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
234064c8009SBarry Smith */
235064c8009SBarry Smith cnt = 0;
2362fa5cd67SKarl Rupp
2379371c9d4SSatish Balay gl[cnt++] = 0;
238d71ae5a4SJacob Faibussowitsch {
239d71ae5a4SJacob Faibussowitsch gl[cnt++] = 1;
240d71ae5a4SJacob Faibussowitsch }
2419371c9d4SSatish Balay gl[cnt++] = m - istart - 1;
242064c8009SBarry Smith {
2439371c9d4SSatish Balay gl[cnt++] = mwidth;
244d71ae5a4SJacob Faibussowitsch {
245d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth + 1;
246d71ae5a4SJacob Faibussowitsch }
2479371c9d4SSatish Balay gl[cnt++] = mwidth + m - istart - 1;
248064c8009SBarry Smith }
2499371c9d4SSatish Balay gl[cnt++] = mwidth * (n - jstart - 1);
250d71ae5a4SJacob Faibussowitsch {
251d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * (n - jstart - 1) + 1;
252d71ae5a4SJacob Faibussowitsch }
2539371c9d4SSatish Balay gl[cnt++] = mwidth * (n - jstart - 1) + m - istart - 1;
2549371c9d4SSatish Balay {
2559371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth;
256d71ae5a4SJacob Faibussowitsch {
257d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth + 1;
258d71ae5a4SJacob Faibussowitsch }
2599371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + m - istart - 1;
2609371c9d4SSatish Balay {
2619371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + mwidth; /* these are the interior nodes */
2629371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + mwidth + m - istart - 1;
2639371c9d4SSatish Balay }
2649371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1);
265d71ae5a4SJacob Faibussowitsch {
266d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1) + 1;
267d71ae5a4SJacob Faibussowitsch }
2689371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1) + m - istart - 1;
2699371c9d4SSatish Balay }
2709371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth * (p - kstart - 1);
271d71ae5a4SJacob Faibussowitsch {
272d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + 1;
273d71ae5a4SJacob Faibussowitsch }
2749371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + m - istart - 1;
2759371c9d4SSatish Balay {
2769371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth;
277d71ae5a4SJacob Faibussowitsch {
278d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth + 1;
279d71ae5a4SJacob Faibussowitsch }
2809371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth + m - istart - 1;
2819371c9d4SSatish Balay }
2829371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth * (n - jstart - 1);
283d71ae5a4SJacob Faibussowitsch {
284d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth * (n - jstart - 1) + 1;
285d71ae5a4SJacob Faibussowitsch }
2869371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth * (n - jstart - 1) + m - istart - 1;
287064c8009SBarry Smith
288064c8009SBarry Smith /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
289064c8009SBarry Smith /* convert that to global numbering and get them on all processes */
2909566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, 26, gl, gl));
291064c8009SBarry Smith /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
2929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(26 * mp * np * pp, &globals));
2939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(gl, 26, MPIU_INT, globals, 26, MPIU_INT, PetscObjectComm((PetscObject)da)));
294064c8009SBarry Smith
295064c8009SBarry Smith /* Number the coarse grid points from 0 to Ntotal */
296eec179cfSJacob Faibussowitsch PetscCall(PetscHMapICreateWithSize(Ntotal / 3, &ht));
297eec179cfSJacob Faibussowitsch for (i = 0, cnt = 0; i < 26 * mp * np * pp; i++) {
298eec179cfSJacob Faibussowitsch PetscHashIter it = 0;
299eec179cfSJacob Faibussowitsch PetscBool missing = PETSC_TRUE;
300eec179cfSJacob Faibussowitsch
301eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIPut(ht, globals[i] + 1, &it, &missing));
302eec179cfSJacob Faibussowitsch if (missing) {
303eec179cfSJacob Faibussowitsch ++cnt;
304eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIIterSet(ht, it, cnt));
305eec179cfSJacob Faibussowitsch }
306eec179cfSJacob Faibussowitsch }
30763a3b9bcSJacob 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);
3089566063dSJacob Faibussowitsch PetscCall(PetscFree(globals));
309064c8009SBarry Smith for (i = 0; i < 26; i++) {
310eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(ht, gl[i] + 1, 0, gl + i));
311*663e1fa8SPierre Jolivet --gl[i];
312064c8009SBarry Smith }
313eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&ht));
314064c8009SBarry Smith /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
315064c8009SBarry Smith
316064c8009SBarry Smith /* construct global interpolation matrix */
3179566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Aglobal, &Ng, NULL));
318064c8009SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {
3199566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PetscObjectComm((PetscObject)da), Ng, PETSC_DECIDE, PETSC_DECIDE, Ntotal, Nint + Nsurf, NULL, Nint + Nsurf, NULL, P));
320064c8009SBarry Smith } else {
3219566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(*P));
322064c8009SBarry Smith }
3239566063dSJacob Faibussowitsch PetscCall(MatSetOption(*P, MAT_ROW_ORIENTED, PETSC_FALSE));
3249566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xint, &rxint));
3259566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, Nint, IIint, 26, gl, rxint, INSERT_VALUES));
3269566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xint, &rxint));
3279566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xsurf, &rxint));
3289566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, Nsurf, IIsurf, 26, gl, rxint, INSERT_VALUES));
3299566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xsurf, &rxint));
3309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
3319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
3329566063dSJacob Faibussowitsch PetscCall(PetscFree2(IIint, IIsurf));
333064c8009SBarry Smith
3348e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
335064c8009SBarry Smith {
336064c8009SBarry Smith Vec x, y;
337064c8009SBarry Smith PetscScalar *yy;
3389566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), Ng, PETSC_DETERMINE, &y));
3399566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), PETSC_DETERMINE, Ntotal, &x));
3409566063dSJacob Faibussowitsch PetscCall(VecSet(x, 1.0));
3419566063dSJacob Faibussowitsch PetscCall(MatMult(*P, x, y));
3429566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy));
343ad540459SPierre 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]));
3449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy));
3459566063dSJacob Faibussowitsch PetscCall(VecDestroy(x));
3469566063dSJacob Faibussowitsch PetscCall(VecDestroy(y));
347064c8009SBarry Smith }
348064c8009SBarry Smith #endif
349064c8009SBarry Smith
3509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aii));
3519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ais));
3529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Asi));
3539566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
3549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is));
3559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isint));
3569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&issurf));
3579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xint));
3589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xsurf));
3593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
360064c8009SBarry Smith }
361064c8009SBarry Smith
362064c8009SBarry Smith /*
363aa219208SBarry Smith DMDAGetFaceInterpolation - Gets the interpolation for a face based coarse space
364064c8009SBarry Smith */
DMDAGetFaceInterpolation(PC pc,DM da,PC_Exotic * exotic,Mat Aglobal,MatReuse reuse,Mat * P)36566976f2fSJacob Faibussowitsch static PetscErrorCode DMDAGetFaceInterpolation(PC pc, DM da, PC_Exotic *exotic, Mat Aglobal, MatReuse reuse, Mat *P)
366d71ae5a4SJacob Faibussowitsch {
367064c8009SBarry 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;
368eec179cfSJacob Faibussowitsch PetscInt mwidth, nwidth, pwidth, cnt, mp, np, pp, Ntotal, gl[6], *globals, Ng, *IIint, *IIsurf;
369064c8009SBarry Smith Mat Xint, Xsurf, Xint_tmp;
370064c8009SBarry Smith IS isint, issurf, is, row, col;
371064c8009SBarry Smith ISLocalToGlobalMapping ltg;
372064c8009SBarry Smith MPI_Comm comm;
373064c8009SBarry Smith Mat A, Aii, Ais, Asi, *Aholder, iAii;
374064c8009SBarry Smith MatFactorInfo info;
375064c8009SBarry Smith PetscScalar *xsurf, *xint;
3761683a169SBarry Smith const PetscScalar *rxint;
377064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
378064c8009SBarry Smith PetscScalar tmp;
379064c8009SBarry Smith #endif
380eec179cfSJacob Faibussowitsch PetscHMapI ht;
381064c8009SBarry Smith
382064c8009SBarry Smith PetscFunctionBegin;
3839566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, &mp, &np, &pp, &dof, NULL, NULL, NULL, NULL, NULL));
38408401ef6SPierre Jolivet PetscCheck(dof == 1, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only for single field problems");
38508401ef6SPierre Jolivet PetscCheck(dim == 3, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only coded for 3d problems");
3869566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, NULL, NULL, NULL, &m, &n, &p));
3879566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &istart, &jstart, &kstart, &mwidth, &nwidth, &pwidth));
388064c8009SBarry Smith istart = istart ? -1 : 0;
389064c8009SBarry Smith jstart = jstart ? -1 : 0;
390064c8009SBarry Smith kstart = kstart ? -1 : 0;
391064c8009SBarry Smith
392064c8009SBarry Smith /*
3934bc3789fSBarry Smith the columns of P are the interpolation of each coarse (face) grid point (one for each face)
394064c8009SBarry Smith to all the local degrees of freedom (this includes the vertices, edges and faces).
395064c8009SBarry Smith
396064c8009SBarry Smith Xint are the subset of the interpolation into the interior
397064c8009SBarry Smith
398064c8009SBarry Smith Xface are the interpolation onto faces but not into the interior
399064c8009SBarry Smith
4004bc3789fSBarry Smith Xsurf are the interpolation onto the vertices and edges (the wirebasket)
401064c8009SBarry Smith Xint
402064c8009SBarry Smith Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain
403064c8009SBarry Smith Xsurf
404064c8009SBarry Smith */
405064c8009SBarry Smith N = (m - istart) * (n - jstart) * (p - kstart);
406064c8009SBarry Smith Nint = (m - 2 - istart) * (n - 2 - jstart) * (p - 2 - kstart);
407064c8009SBarry Smith Nface = 2 * ((m - 2 - istart) * (n - 2 - jstart) + (m - 2 - istart) * (p - 2 - kstart) + (n - 2 - jstart) * (p - 2 - kstart));
408064c8009SBarry Smith Nwire = 4 * ((m - 2 - istart) + (n - 2 - jstart) + (p - 2 - kstart)) + 8;
409064c8009SBarry Smith Nsurf = Nface + Nwire;
4109566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nint, 6, NULL, &Xint));
4119566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nsurf, 6, NULL, &Xsurf));
4129566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(Xsurf, &xsurf));
413064c8009SBarry Smith
414064c8009SBarry Smith /*
415064c8009SBarry Smith Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of
416064c8009SBarry Smith Xsurf will be all zero (thus making the coarse matrix singular).
417064c8009SBarry Smith */
41808401ef6SPierre 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");
41908401ef6SPierre 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");
42008401ef6SPierre 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");
421064c8009SBarry Smith
422064c8009SBarry Smith cnt = 0;
4232fa5cd67SKarl Rupp for (j = 1; j < n - 1 - jstart; j++) {
4242fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 0 * Nsurf] = 1;
425064c8009SBarry Smith }
4262fa5cd67SKarl Rupp
4272fa5cd67SKarl Rupp for (k = 1; k < p - 1 - kstart; k++) {
4282fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 1 * Nsurf] = 1;
4292fa5cd67SKarl Rupp for (j = 1; j < n - 1 - jstart; j++) {
4302fa5cd67SKarl Rupp xsurf[cnt++ + 2 * Nsurf] = 1;
4312fa5cd67SKarl Rupp /* these are the interior nodes */
4322fa5cd67SKarl Rupp xsurf[cnt++ + 3 * Nsurf] = 1;
4332fa5cd67SKarl Rupp }
4342fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 4 * Nsurf] = 1;
4352fa5cd67SKarl Rupp }
4362fa5cd67SKarl Rupp for (j = 1; j < n - 1 - jstart; j++) {
4372fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 5 * Nsurf] = 1;
4382fa5cd67SKarl Rupp }
439064c8009SBarry Smith
440064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
441064c8009SBarry Smith for (i = 0; i < Nsurf; i++) {
442064c8009SBarry Smith tmp = 0.0;
4432fa5cd67SKarl Rupp for (j = 0; j < 6; j++) tmp += xsurf[i + j * Nsurf];
4442fa5cd67SKarl Rupp
44563a3b9bcSJacob 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));
446064c8009SBarry Smith }
447064c8009SBarry Smith #endif
4489566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(Xsurf, &xsurf));
4499566063dSJacob Faibussowitsch /* PetscCall(MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD));*/
450064c8009SBarry Smith
451064c8009SBarry Smith /*
452064c8009SBarry Smith I are the indices for all the needed vertices (in global numbering)
453064c8009SBarry Smith Iint are the indices for the interior values, I surf for the surface values
4547dae84e0SHong Zhang (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it
455aa219208SBarry Smith is NOT the local DMDA ordering.)
456064c8009SBarry Smith IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
457064c8009SBarry Smith */
458064c8009SBarry Smith #define Endpoint(a, start, b) (a == 0 || a == (b - 1 - start))
4599566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(N, &II, Nint, &Iint, Nsurf, &Isurf));
4609566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nint, &IIint, Nsurf, &IIsurf));
461064c8009SBarry Smith for (k = 0; k < p - kstart; k++) {
462064c8009SBarry Smith for (j = 0; j < n - jstart; j++) {
463064c8009SBarry Smith for (i = 0; i < m - istart; i++) {
464064c8009SBarry Smith II[c++] = i + j * mwidth + k * mwidth * nwidth;
465064c8009SBarry Smith
466064c8009SBarry Smith if (!Endpoint(i, istart, m) && !Endpoint(j, jstart, n) && !Endpoint(k, kstart, p)) {
467064c8009SBarry Smith IIint[cint] = i + j * mwidth + k * mwidth * nwidth;
468064c8009SBarry Smith Iint[cint++] = i + j * (m - istart) + k * (m - istart) * (n - jstart);
469064c8009SBarry Smith } else {
470064c8009SBarry Smith IIsurf[csurf] = i + j * mwidth + k * mwidth * nwidth;
471064c8009SBarry Smith Isurf[csurf++] = i + j * (m - istart) + k * (m - istart) * (n - jstart);
472064c8009SBarry Smith }
473064c8009SBarry Smith }
474064c8009SBarry Smith }
475064c8009SBarry Smith }
476eec179cfSJacob Faibussowitsch #undef Endpoint
47708401ef6SPierre Jolivet PetscCheck(c == N, PETSC_COMM_SELF, PETSC_ERR_PLIB, "c != N");
47808401ef6SPierre Jolivet PetscCheck(cint == Nint, PETSC_COMM_SELF, PETSC_ERR_PLIB, "cint != Nint");
47908401ef6SPierre Jolivet PetscCheck(csurf == Nsurf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "csurf != Nsurf");
4809566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(da, <g));
4819566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, N, II, II));
4829566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, Nint, IIint, IIint));
4839566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, Nsurf, IIsurf, IIsurf));
4849566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
4859566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, N, II, PETSC_COPY_VALUES, &is));
4869566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nint, Iint, PETSC_COPY_VALUES, &isint));
4879566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nsurf, Isurf, PETSC_COPY_VALUES, &issurf));
4889566063dSJacob Faibussowitsch PetscCall(PetscFree3(II, Iint, Isurf));
489064c8009SBarry Smith
4909566063dSJacob Faibussowitsch PetscCall(ISSort(is));
4919566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(Aglobal, 1, &is, &is, MAT_INITIAL_MATRIX, &Aholder));
492064c8009SBarry Smith A = *Aholder;
4939566063dSJacob Faibussowitsch PetscCall(PetscFree(Aholder));
494064c8009SBarry Smith
4959566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, isint, isint, MAT_INITIAL_MATRIX, &Aii));
4969566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, isint, issurf, MAT_INITIAL_MATRIX, &Ais));
4979566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, issurf, isint, MAT_INITIAL_MATRIX, &Asi));
498064c8009SBarry Smith
499064c8009SBarry Smith /*
500064c8009SBarry Smith Solve for the interpolation onto the interior Xint
501064c8009SBarry Smith */
5029566063dSJacob Faibussowitsch PetscCall(MatMatMult(Ais, Xsurf, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &Xint_tmp));
5039566063dSJacob Faibussowitsch PetscCall(MatScale(Xint_tmp, -1.0));
504064c8009SBarry Smith
5058e722e36SBarry Smith if (exotic->directSolve) {
5069566063dSJacob Faibussowitsch PetscCall(MatGetFactor(Aii, MATSOLVERPETSC, MAT_FACTOR_LU, &iAii));
5079566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info));
5089566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(Aii, MATORDERINGND, &row, &col));
5099566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(iAii, Aii, row, col, &info));
5109566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row));
5119566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col));
5129566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(iAii, Aii, &info));
5139566063dSJacob Faibussowitsch PetscCall(MatMatSolve(iAii, Xint_tmp, Xint));
5149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&iAii));
515064c8009SBarry Smith } else {
516064c8009SBarry Smith Vec b, x;
517064c8009SBarry Smith PetscScalar *xint_tmp;
518064c8009SBarry Smith
5199566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(Xint, &xint));
5209566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &x));
5219566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(Xint_tmp, &xint_tmp));
5229566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &b));
5239566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(exotic->ksp, Aii, Aii));
524064c8009SBarry Smith for (i = 0; i < 6; i++) {
5259566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(x, xint + i * Nint));
5269566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(b, xint_tmp + i * Nint));
5279566063dSJacob Faibussowitsch PetscCall(KSPSolve(exotic->ksp, b, x));
5289566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(exotic->ksp, pc, x));
5299566063dSJacob Faibussowitsch PetscCall(VecResetArray(x));
5309566063dSJacob Faibussowitsch PetscCall(VecResetArray(b));
531064c8009SBarry Smith }
5329566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(Xint, &xint));
5339566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(Xint_tmp, &xint_tmp));
5349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
5359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b));
536064c8009SBarry Smith }
5379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xint_tmp));
538064c8009SBarry Smith
539064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
5409566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xint, &rxint));
541064c8009SBarry Smith for (i = 0; i < Nint; i++) {
542064c8009SBarry Smith tmp = 0.0;
5431683a169SBarry Smith for (j = 0; j < 6; j++) tmp += rxint[i + j * Nint];
5442fa5cd67SKarl Rupp
54563a3b9bcSJacob 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));
546064c8009SBarry Smith }
5479566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xint, &rxint));
5489566063dSJacob Faibussowitsch /* PetscCall(MatView(Xint,PETSC_VIEWER_STDOUT_WORLD)); */
549064c8009SBarry Smith #endif
550064c8009SBarry Smith
551064c8009SBarry Smith /* total faces */
552064c8009SBarry Smith Ntotal = mp * np * (pp + 1) + mp * pp * (np + 1) + np * pp * (mp + 1);
553064c8009SBarry Smith
554064c8009SBarry Smith /*
555064c8009SBarry Smith For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points
556064c8009SBarry Smith */
557064c8009SBarry Smith cnt = 0;
558d71ae5a4SJacob Faibussowitsch {
559d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth + 1;
560d71ae5a4SJacob Faibussowitsch }
561fbccb6d4SPierre Jolivet {
562fbccb6d4SPierre Jolivet {
563fbccb6d4SPierre Jolivet gl[cnt++] = mwidth * nwidth + 1;
5649371c9d4SSatish Balay }
565064c8009SBarry Smith {
5669371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + mwidth; /* these are the interior nodes */
5679371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + mwidth + m - istart - 1;
5689371c9d4SSatish Balay }
569d71ae5a4SJacob Faibussowitsch {
570d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1) + 1;
571064c8009SBarry Smith }
572d71ae5a4SJacob Faibussowitsch }
573d71ae5a4SJacob Faibussowitsch {
574d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth + 1;
575d71ae5a4SJacob Faibussowitsch }
576064c8009SBarry Smith
577064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
578064c8009SBarry Smith /* convert that to global numbering and get them on all processes */
5799566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, 6, gl, gl));
580064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
5819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(6 * mp * np * pp, &globals));
5829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(gl, 6, MPIU_INT, globals, 6, MPIU_INT, PetscObjectComm((PetscObject)da)));
583064c8009SBarry Smith
584064c8009SBarry Smith /* Number the coarse grid points from 0 to Ntotal */
585eec179cfSJacob Faibussowitsch PetscCall(PetscHMapICreateWithSize(Ntotal / 3, &ht));
586eec179cfSJacob Faibussowitsch for (i = 0, cnt = 0; i < 6 * mp * np * pp; i++) {
587eec179cfSJacob Faibussowitsch PetscHashIter it = 0;
588eec179cfSJacob Faibussowitsch PetscBool missing = PETSC_TRUE;
589eec179cfSJacob Faibussowitsch
590eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIPut(ht, globals[i] + 1, &it, &missing));
591eec179cfSJacob Faibussowitsch if (missing) {
592eec179cfSJacob Faibussowitsch ++cnt;
593eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIIterSet(ht, it, cnt));
594eec179cfSJacob Faibussowitsch }
595eec179cfSJacob Faibussowitsch }
59663a3b9bcSJacob 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);
5979566063dSJacob Faibussowitsch PetscCall(PetscFree(globals));
598064c8009SBarry Smith for (i = 0; i < 6; i++) {
599eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(ht, gl[i] + 1, 0, gl + i));
600*663e1fa8SPierre Jolivet --gl[i];
601064c8009SBarry Smith }
602eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&ht));
603064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
604064c8009SBarry Smith
605064c8009SBarry Smith /* construct global interpolation matrix */
6069566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Aglobal, &Ng, NULL));
607064c8009SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {
6089566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PetscObjectComm((PetscObject)da), Ng, PETSC_DECIDE, PETSC_DECIDE, Ntotal, Nint + Nsurf, NULL, Nint, NULL, P));
609064c8009SBarry Smith } else {
6109566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(*P));
611064c8009SBarry Smith }
6129566063dSJacob Faibussowitsch PetscCall(MatSetOption(*P, MAT_ROW_ORIENTED, PETSC_FALSE));
6139566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xint, &rxint));
6149566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, Nint, IIint, 6, gl, rxint, INSERT_VALUES));
6159566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xint, &rxint));
6169566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xsurf, &rxint));
6179566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, Nsurf, IIsurf, 6, gl, rxint, INSERT_VALUES));
6189566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xsurf, &rxint));
6199566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
6209566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
6219566063dSJacob Faibussowitsch PetscCall(PetscFree2(IIint, IIsurf));
622064c8009SBarry Smith
623064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo)
624064c8009SBarry Smith {
625064c8009SBarry Smith Vec x, y;
626064c8009SBarry Smith PetscScalar *yy;
6279566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), Ng, PETSC_DETERMINE, &y));
6289566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), PETSC_DETERMINE, Ntotal, &x));
6299566063dSJacob Faibussowitsch PetscCall(VecSet(x, 1.0));
6309566063dSJacob Faibussowitsch PetscCall(MatMult(*P, x, y));
6319566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy));
632ad540459SPierre 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]));
6339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy));
6349566063dSJacob Faibussowitsch PetscCall(VecDestroy(x));
6359566063dSJacob Faibussowitsch PetscCall(VecDestroy(y));
636064c8009SBarry Smith }
637064c8009SBarry Smith #endif
638064c8009SBarry Smith
6399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aii));
6409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ais));
6419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Asi));
6429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
6439566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is));
6449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isint));
6459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&issurf));
6469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xint));
6479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xsurf));
6483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
649064c8009SBarry Smith }
650174b6946SBarry Smith
6517233f9f0SBarry Smith /*@
6527233f9f0SBarry Smith PCExoticSetType - Sets the type of coarse grid interpolation to use
6537233f9f0SBarry Smith
654c3339decSBarry Smith Logically Collective
6557233f9f0SBarry Smith
6567233f9f0SBarry Smith Input Parameters:
6577233f9f0SBarry Smith + pc - the preconditioner context
6584bc3789fSBarry Smith - type - either `PC_EXOTIC_FACE` or `PC_EXOTIC_WIREBASKET` (defaults to face)
6594bc3789fSBarry Smith
6604bc3789fSBarry Smith Options Database Keys:
6614bc3789fSBarry Smith . -pc_exotic_type <face,wirebasket> - use a coarse grid point for each face, or edge and vertex
6627233f9f0SBarry Smith
66395452b02SPatrick Sanan Notes:
66495452b02SPatrick Sanan The face based interpolation has 1 degree of freedom per face and ignores the
665563e08c6SBarry Smith edge and vertex values completely in the coarse problem. For any seven point
666563e08c6SBarry Smith stencil the interpolation of a constant on all faces into the interior is that constant.
667563e08c6SBarry Smith
668563e08c6SBarry Smith The wirebasket interpolation has 1 degree of freedom per vertex, per edge and
669563e08c6SBarry Smith per face. A constant on the subdomain boundary is interpolated as that constant
670563e08c6SBarry Smith in the interior of the domain.
671563e08c6SBarry Smith
6724bc3789fSBarry Smith The coarse grid matrix is obtained via the Galerkin computation $A_c = R A R^T$, hence
6734bc3789fSBarry Smith if $A$ is nonsingular $A_c$ is also nonsingular.
674563e08c6SBarry Smith
675563e08c6SBarry Smith Both interpolations are suitable for only scalar problems.
676563e08c6SBarry Smith
6777233f9f0SBarry Smith Level: intermediate
6787233f9f0SBarry Smith
679562efe2eSBarry Smith .seealso: [](ch_ksp), `PCEXOTIC`, `PCExoticType()`
6807233f9f0SBarry Smith @*/
PCExoticSetType(PC pc,PCExoticType type)681d71ae5a4SJacob Faibussowitsch PetscErrorCode PCExoticSetType(PC pc, PCExoticType type)
682d71ae5a4SJacob Faibussowitsch {
6837233f9f0SBarry Smith PetscFunctionBegin;
6840700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
685c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, type, 2);
686cac4c232SBarry Smith PetscTryMethod(pc, "PCExoticSetType_C", (PC, PCExoticType), (pc, type));
6873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6887233f9f0SBarry Smith }
6897233f9f0SBarry Smith
PCExoticSetType_Exotic(PC pc,PCExoticType type)690d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCExoticSetType_Exotic(PC pc, PCExoticType type)
691d71ae5a4SJacob Faibussowitsch {
692f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data;
69331567311SBarry Smith PC_Exotic *ctx = (PC_Exotic *)mg->innerctx;
6947233f9f0SBarry Smith
6957233f9f0SBarry Smith PetscFunctionBegin;
6967233f9f0SBarry Smith ctx->type = type;
6973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6987233f9f0SBarry Smith }
6997233f9f0SBarry Smith
PCSetUp_Exotic(PC pc)70066976f2fSJacob Faibussowitsch static PetscErrorCode PCSetUp_Exotic(PC pc)
701d71ae5a4SJacob Faibussowitsch {
70296bdf778SBarry Smith Mat A;
703f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data;
70431567311SBarry Smith PC_Exotic *ex = (PC_Exotic *)mg->innerctx;
70596bdf778SBarry Smith MatReuse reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX;
706174b6946SBarry Smith
707174b6946SBarry Smith PetscFunctionBegin;
70828b400f6SJacob Faibussowitsch PetscCheck(pc->dm, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Need to call PCSetDM() before using this PC");
7099566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, NULL, &A));
7100fdf79fbSJacob 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);
7117233f9f0SBarry Smith if (ex->type == PC_EXOTIC_FACE) {
7129566063dSJacob Faibussowitsch PetscCall(DMDAGetFaceInterpolation(pc, pc->dm, ex, A, reuse, &ex->P));
7130fdf79fbSJacob Faibussowitsch } else /* if (ex->type == PC_EXOTIC_WIREBASKET) */ {
7149566063dSJacob Faibussowitsch PetscCall(DMDAGetWireBasketInterpolation(pc, pc->dm, ex, A, reuse, &ex->P));
7150fdf79fbSJacob Faibussowitsch }
7169566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc, 1, ex->P));
717d0660788SBarry Smith /* if PC has attached DM we must remove it or the PCMG will use it to compute incorrect sized vectors and interpolations */
7189566063dSJacob Faibussowitsch PetscCall(PCSetDM(pc, NULL));
7199566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc));
7203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
721174b6946SBarry Smith }
722174b6946SBarry Smith
PCDestroy_Exotic(PC pc)72366976f2fSJacob Faibussowitsch static PetscErrorCode PCDestroy_Exotic(PC pc)
724d71ae5a4SJacob Faibussowitsch {
725f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data;
72631567311SBarry Smith PC_Exotic *ctx = (PC_Exotic *)mg->innerctx;
727174b6946SBarry Smith
728174b6946SBarry Smith PetscFunctionBegin;
7299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->P));
7309566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ctx->ksp));
7319566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx));
7322e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCExoticSetType_C", NULL));
7339566063dSJacob Faibussowitsch PetscCall(PCDestroy_MG(pc));
7343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
735174b6946SBarry Smith }
736174b6946SBarry Smith
PCView_Exotic(PC pc,PetscViewer viewer)73766976f2fSJacob Faibussowitsch static PetscErrorCode PCView_Exotic(PC pc, PetscViewer viewer)
738d71ae5a4SJacob Faibussowitsch {
739f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data;
7409f196a02SMartin Diehl PetscBool isascii;
74131567311SBarry Smith PC_Exotic *ctx = (PC_Exotic *)mg->innerctx;
7427233f9f0SBarry Smith
7437233f9f0SBarry Smith PetscFunctionBegin;
7449f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
7459f196a02SMartin Diehl if (isascii) {
7469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Exotic type = %s\n", PCExoticTypes[ctx->type]));
7478e722e36SBarry Smith if (ctx->directSolve) {
7489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using direct solver to construct interpolation\n"));
7498e722e36SBarry Smith } else {
7508e722e36SBarry Smith PetscViewer sviewer;
7518e722e36SBarry Smith PetscMPIInt rank;
7528e722e36SBarry Smith
7539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using iterative solver to construct interpolation\n"));
7549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer));
7559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); /* should not need to push this twice? */
7569566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
7579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
75848a46eb9SPierre Jolivet if (rank == 0) PetscCall(KSPView(ctx->ksp, sviewer));
7599566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
7609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
7619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
7628e722e36SBarry Smith }
7637233f9f0SBarry Smith }
7649566063dSJacob Faibussowitsch PetscCall(PCView_MG(pc, viewer));
7653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7667233f9f0SBarry Smith }
7677233f9f0SBarry Smith
PCSetFromOptions_Exotic(PC pc,PetscOptionItems PetscOptionsObject)768ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_Exotic(PC pc, PetscOptionItems PetscOptionsObject)
769d71ae5a4SJacob Faibussowitsch {
770ace3abfcSBarry Smith PetscBool flg;
771f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data;
7727233f9f0SBarry Smith PCExoticType mgctype;
77331567311SBarry Smith PC_Exotic *ctx = (PC_Exotic *)mg->innerctx;
7747233f9f0SBarry Smith
7757233f9f0SBarry Smith PetscFunctionBegin;
776d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Exotic coarse space options");
7779566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_exotic_type", "face or wirebasket", "PCExoticSetType", PCExoticTypes, (PetscEnum)ctx->type, (PetscEnum *)&mgctype, &flg));
7781baa6e33SBarry Smith if (flg) PetscCall(PCExoticSetType(pc, mgctype));
7799566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_exotic_direct_solver", "use direct solver to construct interpolation", "None", ctx->directSolve, &ctx->directSolve, NULL));
7808e722e36SBarry Smith if (!ctx->directSolve) {
7818e722e36SBarry Smith if (!ctx->ksp) {
7828e722e36SBarry Smith const char *prefix;
7839566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF, &ctx->ksp));
7843821be0aSBarry Smith PetscCall(KSPSetNestLevel(ctx->ksp, pc->kspnestlevel));
7859566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ctx->ksp, pc->erroriffailure));
7869566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->ksp, (PetscObject)pc, 1));
7879566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix));
7889566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ctx->ksp, prefix));
7899566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ctx->ksp, "exotic_"));
7908e722e36SBarry Smith }
7919566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ctx->ksp));
7928e722e36SBarry Smith }
793d0609cedSBarry Smith PetscOptionsHeadEnd();
7943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7957233f9f0SBarry Smith }
7967233f9f0SBarry Smith
797174b6946SBarry Smith /*MC
7987233f9f0SBarry Smith PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces
799174b6946SBarry Smith
800f1580f4eSBarry Smith This uses the `PCMG` infrastructure restricted to two levels and the face and wirebasket based coarse
80124c3aa18SBarry Smith grid spaces.
80224c3aa18SBarry Smith
8034bc3789fSBarry Smith Options Database Keys:
8044bc3789fSBarry Smith + -pc_exotic_type <face,wirebasket> - use a coarse grid point for each face, or edge and vertex
8054bc3789fSBarry Smith - -pc_exotic_direct_solver - use a direct solver to construct interpolation instead of an iterative solver
8064bc3789fSBarry Smith
807f1580f4eSBarry Smith Level: advanced
808f1580f4eSBarry Smith
80995452b02SPatrick Sanan Notes:
8104bc3789fSBarry Smith Must be used with `DMDA` in three dimensions
811f1580f4eSBarry Smith
812f1580f4eSBarry 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`
81324c3aa18SBarry Smith
8141d27aa22SBarry Smith These coarse grid spaces originate in the work of Bramble, Pasciak and Schatz {cite}`bramble1989construction`.
8151d27aa22SBarry Smith They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith, {cite}`smith1990domain`.
8161d27aa22SBarry Smith They were then explored in great detail in Dryja, Smith, Widlund {cite}`dryja1994schwarz`. These were developed in the context of iterative substructuring preconditioners.
8171d27aa22SBarry Smith
8181d27aa22SBarry Smith They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund.
8191d27aa22SBarry Smith They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example, {cite}`dohrmann2008extending`, {cite}`dohrmann2008family`,
8201d27aa22SBarry Smith {cite}`dohrmann2008domain`, {cite}`dohrmann2009overlapping`.
8217233f9f0SBarry Smith
8224bc3789fSBarry Smith In this code the wirebasket includes a constant for each face, as well as the true "wirebasket". Other wirebasket algorithms exist that
8234bc3789fSBarry Smith only have constants for edges and vertices.
8244bc3789fSBarry Smith
8254bc3789fSBarry Smith 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>
826174b6946SBarry Smith
827562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCSetDM()`, `PCExoticType`, `PCExoticSetType()`
828174b6946SBarry Smith M*/
829174b6946SBarry Smith
PCCreate_Exotic(PC pc)830d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Exotic(PC pc)
831d71ae5a4SJacob Faibussowitsch {
8327233f9f0SBarry Smith PC_Exotic *ex;
833f3fbd535SBarry Smith PC_MG *mg;
834174b6946SBarry Smith
835174b6946SBarry Smith PetscFunctionBegin;
836f91d8e95SBarry Smith /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
837dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, destroy);
8380a545947SLisandro Dalcin pc->data = NULL;
839dbbe0bcdSBarry Smith
8409566063dSJacob Faibussowitsch PetscCall(PetscFree(((PetscObject)pc)->type_name));
8410a545947SLisandro Dalcin ((PetscObject)pc)->type_name = NULL;
842f91d8e95SBarry Smith
8439566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCMG));
8449566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, 2, NULL));
8459566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_PMAT));
8469371c9d4SSatish Balay PetscCall(PetscNew(&ex));
8477233f9f0SBarry Smith ex->type = PC_EXOTIC_FACE;
848f3fbd535SBarry Smith mg = (PC_MG *)pc->data;
84931567311SBarry Smith mg->innerctx = ex;
8507233f9f0SBarry Smith
8517233f9f0SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Exotic;
8527233f9f0SBarry Smith pc->ops->view = PCView_Exotic;
8537233f9f0SBarry Smith pc->ops->destroy = PCDestroy_Exotic;
8546c699258SBarry Smith pc->ops->setup = PCSetUp_Exotic;
8552fa5cd67SKarl Rupp
8569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCExoticSetType_C", PCExoticSetType_Exotic));
8573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
858174b6946SBarry Smith }
859