1174b6946SBarry Smith 2c6db04a5SJed Brown #include <petscdmda.h> /*I "petscdmda.h" I*/ 3af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 482c86c8fSBarry Smith #include <petscctable.h> 57233f9f0SBarry Smith 68e722e36SBarry Smith typedef struct { 78e722e36SBarry Smith PCExoticType type; 88e722e36SBarry Smith Mat P; /* the constructed interpolation matrix */ 9ace3abfcSBarry Smith PetscBool directSolve; /* use direct LU factorization to construct interpolation */ 108e722e36SBarry Smith KSP ksp; 118e722e36SBarry Smith } PC_Exotic; 12174b6946SBarry Smith 130a545947SLisandro Dalcin const char *const PCExoticTypes[] = {"face", "wirebasket", "PCExoticType", "PC_Exotic", NULL}; 14064c8009SBarry Smith 15064c8009SBarry Smith /* 16aa219208SBarry Smith DMDAGetWireBasketInterpolation - Gets the interpolation for a wirebasket based coarse space 17064c8009SBarry Smith 18064c8009SBarry Smith */ 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; 2228d20b34SBarry Smith PetscInt mwidth, nwidth, pwidth, cnt, mp, np, pp, Ntotal, gl[26], *globals, Ng, *IIint, *IIsurf, Nt; 23064c8009SBarry Smith Mat Xint, Xsurf, Xint_tmp; 24064c8009SBarry Smith IS isint, issurf, is, row, col; 25064c8009SBarry Smith ISLocalToGlobalMapping ltg; 26064c8009SBarry Smith MPI_Comm comm; 27064c8009SBarry Smith Mat A, Aii, Ais, Asi, *Aholder, iAii; 28064c8009SBarry Smith MatFactorInfo info; 29064c8009SBarry Smith PetscScalar *xsurf, *xint; 301683a169SBarry Smith const PetscScalar *rxint; 318e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 32064c8009SBarry Smith PetscScalar tmp; 33064c8009SBarry Smith #endif 34064c8009SBarry Smith PetscTable ht; 35064c8009SBarry Smith 36064c8009SBarry Smith PetscFunctionBegin; 379566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, &mp, &np, &pp, &dof, NULL, NULL, NULL, NULL, NULL)); 3808401ef6SPierre Jolivet PetscCheck(dof == 1, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only for single field problems"); 3908401ef6SPierre Jolivet PetscCheck(dim == 3, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only coded for 3d problems"); 409566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, NULL, NULL, NULL, &m, &n, &p)); 419566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &istart, &jstart, &kstart, &mwidth, &nwidth, &pwidth)); 42064c8009SBarry Smith istart = istart ? -1 : 0; 43064c8009SBarry Smith jstart = jstart ? -1 : 0; 44064c8009SBarry Smith kstart = kstart ? -1 : 0; 45064c8009SBarry Smith 46064c8009SBarry Smith /* 47064c8009SBarry Smith the columns of P are the interpolation of each coarse grid point (one for each vertex and edge) 48064c8009SBarry Smith to all the local degrees of freedom (this includes the vertices, edges and faces). 49064c8009SBarry Smith 50064c8009SBarry Smith Xint are the subset of the interpolation into the interior 51064c8009SBarry Smith 52064c8009SBarry Smith Xface are the interpolation onto faces but not into the interior 53064c8009SBarry Smith 54064c8009SBarry Smith Xsurf are the interpolation onto the vertices and edges (the surfbasket) 55064c8009SBarry Smith Xint 56064c8009SBarry Smith Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain 57064c8009SBarry Smith Xsurf 58064c8009SBarry Smith */ 59064c8009SBarry Smith N = (m - istart) * (n - jstart) * (p - kstart); 60064c8009SBarry Smith Nint = (m - 2 - istart) * (n - 2 - jstart) * (p - 2 - kstart); 61064c8009SBarry Smith Nface = 2 * ((m - 2 - istart) * (n - 2 - jstart) + (m - 2 - istart) * (p - 2 - kstart) + (n - 2 - jstart) * (p - 2 - kstart)); 62064c8009SBarry Smith Nwire = 4 * ((m - 2 - istart) + (n - 2 - jstart) + (p - 2 - kstart)) + 8; 63064c8009SBarry Smith Nsurf = Nface + Nwire; 649566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nint, 26, NULL, &Xint)); 659566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nsurf, 26, NULL, &Xsurf)); 669566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(Xsurf, &xsurf)); 67064c8009SBarry Smith 68064c8009SBarry Smith /* 69064c8009SBarry Smith Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of 70064c8009SBarry Smith Xsurf will be all zero (thus making the coarse matrix singular). 71064c8009SBarry Smith */ 7208401ef6SPierre Jolivet PetscCheck(m - istart >= 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Number of grid points per process in X direction must be at least 3"); 7308401ef6SPierre Jolivet PetscCheck(n - jstart >= 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Number of grid points per process in Y direction must be at least 3"); 7408401ef6SPierre Jolivet PetscCheck(p - kstart >= 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Number of grid points per process in Z direction must be at least 3"); 75064c8009SBarry Smith 76064c8009SBarry Smith cnt = 0; 772fa5cd67SKarl Rupp 782fa5cd67SKarl Rupp xsurf[cnt++] = 1; 792fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + Nsurf] = 1; 802fa5cd67SKarl Rupp xsurf[cnt++ + 2 * Nsurf] = 1; 812fa5cd67SKarl Rupp 822fa5cd67SKarl Rupp for (j = 1; j < n - 1 - jstart; j++) { 832fa5cd67SKarl Rupp xsurf[cnt++ + 3 * Nsurf] = 1; 842fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 4 * Nsurf] = 1; 852fa5cd67SKarl Rupp xsurf[cnt++ + 5 * Nsurf] = 1; 86064c8009SBarry Smith } 872fa5cd67SKarl Rupp 882fa5cd67SKarl Rupp xsurf[cnt++ + 6 * Nsurf] = 1; 892fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 7 * Nsurf] = 1; 902fa5cd67SKarl Rupp xsurf[cnt++ + 8 * Nsurf] = 1; 912fa5cd67SKarl Rupp 922fa5cd67SKarl Rupp for (k = 1; k < p - 1 - kstart; k++) { 932fa5cd67SKarl Rupp xsurf[cnt++ + 9 * Nsurf] = 1; 942fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 10 * Nsurf] = 1; 952fa5cd67SKarl Rupp xsurf[cnt++ + 11 * Nsurf] = 1; 962fa5cd67SKarl Rupp 972fa5cd67SKarl Rupp for (j = 1; j < n - 1 - jstart; j++) { 982fa5cd67SKarl Rupp xsurf[cnt++ + 12 * Nsurf] = 1; 992fa5cd67SKarl Rupp /* these are the interior nodes */ 1002fa5cd67SKarl Rupp xsurf[cnt++ + 13 * Nsurf] = 1; 1012fa5cd67SKarl Rupp } 1022fa5cd67SKarl Rupp 1032fa5cd67SKarl Rupp xsurf[cnt++ + 14 * Nsurf] = 1; 1042fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 15 * Nsurf] = 1; 1052fa5cd67SKarl Rupp xsurf[cnt++ + 16 * Nsurf] = 1; 1062fa5cd67SKarl Rupp } 1072fa5cd67SKarl Rupp 1082fa5cd67SKarl Rupp xsurf[cnt++ + 17 * Nsurf] = 1; 1092fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 18 * Nsurf] = 1; 1102fa5cd67SKarl Rupp xsurf[cnt++ + 19 * Nsurf] = 1; 1112fa5cd67SKarl Rupp 1122fa5cd67SKarl Rupp for (j = 1; j < n - 1 - jstart; j++) { 1132fa5cd67SKarl Rupp xsurf[cnt++ + 20 * Nsurf] = 1; 1142fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 21 * Nsurf] = 1; 1152fa5cd67SKarl Rupp xsurf[cnt++ + 22 * Nsurf] = 1; 1162fa5cd67SKarl Rupp } 1172fa5cd67SKarl Rupp 1182fa5cd67SKarl Rupp xsurf[cnt++ + 23 * Nsurf] = 1; 1192fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 24 * Nsurf] = 1; 1202fa5cd67SKarl Rupp xsurf[cnt++ + 25 * Nsurf] = 1; 1212fa5cd67SKarl Rupp 1228e722e36SBarry Smith /* interpolations only sum to 1 when using direct solver */ 1238e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 124064c8009SBarry Smith for (i = 0; i < Nsurf; i++) { 125064c8009SBarry Smith tmp = 0.0; 1262fa5cd67SKarl Rupp for (j = 0; j < 26; j++) tmp += xsurf[i + j * Nsurf]; 12763a3b9bcSJacob Faibussowitsch PetscCheck(PetscAbsScalar(tmp - 1.0) <= 1.e-10, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Wrong Xsurf interpolation at i %" PetscInt_FMT " value %g", i, (double)PetscAbsScalar(tmp)); 128064c8009SBarry Smith } 129064c8009SBarry Smith #endif 1309566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(Xsurf, &xsurf)); 1319566063dSJacob Faibussowitsch /* PetscCall(MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD));*/ 132064c8009SBarry Smith 133064c8009SBarry Smith /* 134064c8009SBarry Smith I are the indices for all the needed vertices (in global numbering) 135064c8009SBarry Smith Iint are the indices for the interior values, I surf for the surface values 1367dae84e0SHong Zhang (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it 137aa219208SBarry Smith is NOT the local DMDA ordering.) 138064c8009SBarry Smith IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering 139064c8009SBarry Smith */ 140064c8009SBarry Smith #define Endpoint(a, start, b) (a == 0 || a == (b - 1 - start)) 1419566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(N, &II, Nint, &Iint, Nsurf, &Isurf)); 1429566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nint, &IIint, Nsurf, &IIsurf)); 143064c8009SBarry Smith for (k = 0; k < p - kstart; k++) { 144064c8009SBarry Smith for (j = 0; j < n - jstart; j++) { 145064c8009SBarry Smith for (i = 0; i < m - istart; i++) { 146064c8009SBarry Smith II[c++] = i + j * mwidth + k * mwidth * nwidth; 147064c8009SBarry Smith 148064c8009SBarry Smith if (!Endpoint(i, istart, m) && !Endpoint(j, jstart, n) && !Endpoint(k, kstart, p)) { 149064c8009SBarry Smith IIint[cint] = i + j * mwidth + k * mwidth * nwidth; 150064c8009SBarry Smith Iint[cint++] = i + j * (m - istart) + k * (m - istart) * (n - jstart); 151064c8009SBarry Smith } else { 152064c8009SBarry Smith IIsurf[csurf] = i + j * mwidth + k * mwidth * nwidth; 153064c8009SBarry Smith Isurf[csurf++] = i + j * (m - istart) + k * (m - istart) * (n - jstart); 154064c8009SBarry Smith } 155064c8009SBarry Smith } 156064c8009SBarry Smith } 157064c8009SBarry Smith } 15808401ef6SPierre Jolivet PetscCheck(c == N, PETSC_COMM_SELF, PETSC_ERR_PLIB, "c != N"); 15908401ef6SPierre Jolivet PetscCheck(cint == Nint, PETSC_COMM_SELF, PETSC_ERR_PLIB, "cint != Nint"); 16008401ef6SPierre Jolivet PetscCheck(csurf == Nsurf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "csurf != Nsurf"); 1619566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(da, <g)); 1629566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, N, II, II)); 1639566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, Nint, IIint, IIint)); 1649566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, Nsurf, IIsurf, IIsurf)); 1659566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 1669566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, N, II, PETSC_COPY_VALUES, &is)); 1679566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nint, Iint, PETSC_COPY_VALUES, &isint)); 1689566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nsurf, Isurf, PETSC_COPY_VALUES, &issurf)); 1699566063dSJacob Faibussowitsch PetscCall(PetscFree3(II, Iint, Isurf)); 170064c8009SBarry Smith 1719566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(Aglobal, 1, &is, &is, MAT_INITIAL_MATRIX, &Aholder)); 172064c8009SBarry Smith A = *Aholder; 1739566063dSJacob Faibussowitsch PetscCall(PetscFree(Aholder)); 174064c8009SBarry Smith 1759566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, isint, isint, MAT_INITIAL_MATRIX, &Aii)); 1769566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, isint, issurf, MAT_INITIAL_MATRIX, &Ais)); 1779566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, issurf, isint, MAT_INITIAL_MATRIX, &Asi)); 178064c8009SBarry Smith 179064c8009SBarry Smith /* 180064c8009SBarry Smith Solve for the interpolation onto the interior Xint 181064c8009SBarry Smith */ 1829566063dSJacob Faibussowitsch PetscCall(MatMatMult(Ais, Xsurf, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &Xint_tmp)); 1839566063dSJacob Faibussowitsch PetscCall(MatScale(Xint_tmp, -1.0)); 1848e722e36SBarry Smith if (exotic->directSolve) { 1859566063dSJacob Faibussowitsch PetscCall(MatGetFactor(Aii, MATSOLVERPETSC, MAT_FACTOR_LU, &iAii)); 1869566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 1879566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(Aii, MATORDERINGND, &row, &col)); 1889566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(iAii, Aii, row, col, &info)); 1899566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 1909566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col)); 1919566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(iAii, Aii, &info)); 1929566063dSJacob Faibussowitsch PetscCall(MatMatSolve(iAii, Xint_tmp, Xint)); 1939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&iAii)); 1948e722e36SBarry Smith } else { 1958e722e36SBarry Smith Vec b, x; 1968e722e36SBarry Smith PetscScalar *xint_tmp; 197064c8009SBarry Smith 1989566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(Xint, &xint)); 1999566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &x)); 2009566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(Xint_tmp, &xint_tmp)); 2019566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &b)); 2029566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(exotic->ksp, Aii, Aii)); 2038e722e36SBarry Smith for (i = 0; i < 26; i++) { 2049566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(x, xint + i * Nint)); 2059566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(b, xint_tmp + i * Nint)); 2069566063dSJacob Faibussowitsch PetscCall(KSPSolve(exotic->ksp, b, x)); 2079566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(exotic->ksp, pc, x)); 2089566063dSJacob Faibussowitsch PetscCall(VecResetArray(x)); 2099566063dSJacob Faibussowitsch PetscCall(VecResetArray(b)); 2108e722e36SBarry Smith } 2119566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(Xint, &xint)); 2129566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(Xint_tmp, &xint_tmp)); 2139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 2158e722e36SBarry Smith } 2169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xint_tmp)); 2178e722e36SBarry Smith 2188e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 2199566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xint, &rxint)); 220064c8009SBarry Smith for (i = 0; i < Nint; i++) { 221064c8009SBarry Smith tmp = 0.0; 2221683a169SBarry Smith for (j = 0; j < 26; j++) tmp += rxint[i + j * Nint]; 2232fa5cd67SKarl Rupp 22463a3b9bcSJacob Faibussowitsch PetscCheck(PetscAbsScalar(tmp - 1.0) <= 1.e-10, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Wrong Xint interpolation at i %" PetscInt_FMT " value %g", i, (double)PetscAbsScalar(tmp)); 225064c8009SBarry Smith } 2269566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xint, &rxint)); 2279566063dSJacob Faibussowitsch /* PetscCall(MatView(Xint,PETSC_VIEWER_STDOUT_WORLD)); */ 228064c8009SBarry Smith #endif 229064c8009SBarry Smith 230064c8009SBarry Smith /* total vertices total faces total edges */ 231064c8009SBarry Smith Ntotal = (mp + 1) * (np + 1) * (pp + 1) + mp * np * (pp + 1) + mp * pp * (np + 1) + np * pp * (mp + 1) + mp * (np + 1) * (pp + 1) + np * (mp + 1) * (pp + 1) + pp * (mp + 1) * (np + 1); 232064c8009SBarry Smith 233064c8009SBarry Smith /* 234064c8009SBarry Smith For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points 235064c8009SBarry Smith */ 236064c8009SBarry Smith cnt = 0; 2372fa5cd67SKarl Rupp 2389371c9d4SSatish Balay gl[cnt++] = 0; 239d71ae5a4SJacob Faibussowitsch { 240d71ae5a4SJacob Faibussowitsch gl[cnt++] = 1; 241d71ae5a4SJacob Faibussowitsch } 2429371c9d4SSatish Balay gl[cnt++] = m - istart - 1; 243064c8009SBarry Smith { 2449371c9d4SSatish Balay gl[cnt++] = mwidth; 245d71ae5a4SJacob Faibussowitsch { 246d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth + 1; 247d71ae5a4SJacob Faibussowitsch } 2489371c9d4SSatish Balay gl[cnt++] = mwidth + m - istart - 1; 249064c8009SBarry Smith } 2509371c9d4SSatish Balay gl[cnt++] = mwidth * (n - jstart - 1); 251d71ae5a4SJacob Faibussowitsch { 252d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * (n - jstart - 1) + 1; 253d71ae5a4SJacob Faibussowitsch } 2549371c9d4SSatish Balay gl[cnt++] = mwidth * (n - jstart - 1) + m - istart - 1; 2559371c9d4SSatish Balay { 2569371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth; 257d71ae5a4SJacob Faibussowitsch { 258d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth + 1; 259d71ae5a4SJacob Faibussowitsch } 2609371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + m - istart - 1; 2619371c9d4SSatish Balay { 2629371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + mwidth; /* these are the interior nodes */ 2639371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + mwidth + m - istart - 1; 2649371c9d4SSatish Balay } 2659371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1); 266d71ae5a4SJacob Faibussowitsch { 267d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1) + 1; 268d71ae5a4SJacob Faibussowitsch } 2699371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1) + m - istart - 1; 2709371c9d4SSatish Balay } 2719371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth * (p - kstart - 1); 272d71ae5a4SJacob Faibussowitsch { 273d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + 1; 274d71ae5a4SJacob Faibussowitsch } 2759371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + m - istart - 1; 2769371c9d4SSatish Balay { 2779371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth; 278d71ae5a4SJacob Faibussowitsch { 279d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth + 1; 280d71ae5a4SJacob Faibussowitsch } 2819371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth + m - istart - 1; 2829371c9d4SSatish Balay } 2839371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth * (n - jstart - 1); 284d71ae5a4SJacob Faibussowitsch { 285d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth * (n - jstart - 1) + 1; 286d71ae5a4SJacob Faibussowitsch } 2879371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth * (n - jstart - 1) + m - istart - 1; 288064c8009SBarry Smith 289064c8009SBarry Smith /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */ 290064c8009SBarry Smith /* convert that to global numbering and get them on all processes */ 2919566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, 26, gl, gl)); 292064c8009SBarry Smith /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */ 2939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(26 * mp * np * pp, &globals)); 2949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(gl, 26, MPIU_INT, globals, 26, MPIU_INT, PetscObjectComm((PetscObject)da))); 295064c8009SBarry Smith 296064c8009SBarry Smith /* Number the coarse grid points from 0 to Ntotal */ 2979566063dSJacob Faibussowitsch PetscCall(MatGetSize(Aglobal, &Nt, NULL)); 2989566063dSJacob Faibussowitsch PetscCall(PetscTableCreate(Ntotal / 3, Nt + 1, &ht)); 29948a46eb9SPierre Jolivet for (i = 0; i < 26 * mp * np * pp; i++) PetscCall(PetscTableAddCount(ht, globals[i] + 1)); 3009566063dSJacob Faibussowitsch PetscCall(PetscTableGetCount(ht, &cnt)); 30163a3b9bcSJacob 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); 3029566063dSJacob Faibussowitsch PetscCall(PetscFree(globals)); 303064c8009SBarry Smith for (i = 0; i < 26; i++) { 3049566063dSJacob Faibussowitsch PetscCall(PetscTableFind(ht, gl[i] + 1, &gl[i])); 305064c8009SBarry Smith gl[i]--; 306064c8009SBarry Smith } 3079566063dSJacob Faibussowitsch PetscCall(PetscTableDestroy(&ht)); 308064c8009SBarry Smith /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */ 309064c8009SBarry Smith 310064c8009SBarry Smith /* construct global interpolation matrix */ 3119566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Aglobal, &Ng, NULL)); 312064c8009SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 3139566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PetscObjectComm((PetscObject)da), Ng, PETSC_DECIDE, PETSC_DECIDE, Ntotal, Nint + Nsurf, NULL, Nint + Nsurf, NULL, P)); 314064c8009SBarry Smith } else { 3159566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(*P)); 316064c8009SBarry Smith } 3179566063dSJacob Faibussowitsch PetscCall(MatSetOption(*P, MAT_ROW_ORIENTED, PETSC_FALSE)); 3189566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xint, &rxint)); 3199566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, Nint, IIint, 26, gl, rxint, INSERT_VALUES)); 3209566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xint, &rxint)); 3219566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xsurf, &rxint)); 3229566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, Nsurf, IIsurf, 26, gl, rxint, INSERT_VALUES)); 3239566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xsurf, &rxint)); 3249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY)); 3259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY)); 3269566063dSJacob Faibussowitsch PetscCall(PetscFree2(IIint, IIsurf)); 327064c8009SBarry Smith 3288e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 329064c8009SBarry Smith { 330064c8009SBarry Smith Vec x, y; 331064c8009SBarry Smith PetscScalar *yy; 3329566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), Ng, PETSC_DETERMINE, &y)); 3339566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), PETSC_DETERMINE, Ntotal, &x)); 3349566063dSJacob Faibussowitsch PetscCall(VecSet(x, 1.0)); 3359566063dSJacob Faibussowitsch PetscCall(MatMult(*P, x, y)); 3369566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 337ad540459SPierre 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])); 3389566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 3399566063dSJacob Faibussowitsch PetscCall(VecDestroy(x)); 3409566063dSJacob Faibussowitsch PetscCall(VecDestroy(y)); 341064c8009SBarry Smith } 342064c8009SBarry Smith #endif 343064c8009SBarry Smith 3449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aii)); 3459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ais)); 3469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Asi)); 3479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 3499566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isint)); 3509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&issurf)); 3519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xint)); 3529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xsurf)); 353064c8009SBarry Smith PetscFunctionReturn(0); 354064c8009SBarry Smith } 355064c8009SBarry Smith 356064c8009SBarry Smith /* 357aa219208SBarry Smith DMDAGetFaceInterpolation - Gets the interpolation for a face based coarse space 358064c8009SBarry Smith 359064c8009SBarry Smith */ 360d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetFaceInterpolation(PC pc, DM da, PC_Exotic *exotic, Mat Aglobal, MatReuse reuse, Mat *P) 361d71ae5a4SJacob Faibussowitsch { 362064c8009SBarry 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; 36328d20b34SBarry Smith PetscInt mwidth, nwidth, pwidth, cnt, mp, np, pp, Ntotal, gl[6], *globals, Ng, *IIint, *IIsurf, Nt; 364064c8009SBarry Smith Mat Xint, Xsurf, Xint_tmp; 365064c8009SBarry Smith IS isint, issurf, is, row, col; 366064c8009SBarry Smith ISLocalToGlobalMapping ltg; 367064c8009SBarry Smith MPI_Comm comm; 368064c8009SBarry Smith Mat A, Aii, Ais, Asi, *Aholder, iAii; 369064c8009SBarry Smith MatFactorInfo info; 370064c8009SBarry Smith PetscScalar *xsurf, *xint; 3711683a169SBarry Smith const PetscScalar *rxint; 372064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 373064c8009SBarry Smith PetscScalar tmp; 374064c8009SBarry Smith #endif 375064c8009SBarry Smith PetscTable ht; 376064c8009SBarry Smith 377064c8009SBarry Smith PetscFunctionBegin; 3789566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, &mp, &np, &pp, &dof, NULL, NULL, NULL, NULL, NULL)); 37908401ef6SPierre Jolivet PetscCheck(dof == 1, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only for single field problems"); 38008401ef6SPierre Jolivet PetscCheck(dim == 3, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Only coded for 3d problems"); 3819566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, NULL, NULL, NULL, &m, &n, &p)); 3829566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &istart, &jstart, &kstart, &mwidth, &nwidth, &pwidth)); 383064c8009SBarry Smith istart = istart ? -1 : 0; 384064c8009SBarry Smith jstart = jstart ? -1 : 0; 385064c8009SBarry Smith kstart = kstart ? -1 : 0; 386064c8009SBarry Smith 387064c8009SBarry Smith /* 388064c8009SBarry Smith the columns of P are the interpolation of each coarse grid point (one for each vertex and edge) 389064c8009SBarry Smith to all the local degrees of freedom (this includes the vertices, edges and faces). 390064c8009SBarry Smith 391064c8009SBarry Smith Xint are the subset of the interpolation into the interior 392064c8009SBarry Smith 393064c8009SBarry Smith Xface are the interpolation onto faces but not into the interior 394064c8009SBarry Smith 395064c8009SBarry Smith Xsurf are the interpolation onto the vertices and edges (the surfbasket) 396064c8009SBarry Smith Xint 397064c8009SBarry Smith Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain 398064c8009SBarry Smith Xsurf 399064c8009SBarry Smith */ 400064c8009SBarry Smith N = (m - istart) * (n - jstart) * (p - kstart); 401064c8009SBarry Smith Nint = (m - 2 - istart) * (n - 2 - jstart) * (p - 2 - kstart); 402064c8009SBarry Smith Nface = 2 * ((m - 2 - istart) * (n - 2 - jstart) + (m - 2 - istart) * (p - 2 - kstart) + (n - 2 - jstart) * (p - 2 - kstart)); 403064c8009SBarry Smith Nwire = 4 * ((m - 2 - istart) + (n - 2 - jstart) + (p - 2 - kstart)) + 8; 404064c8009SBarry Smith Nsurf = Nface + Nwire; 4059566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nint, 6, NULL, &Xint)); 4069566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(MPI_COMM_SELF, Nsurf, 6, NULL, &Xsurf)); 4079566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(Xsurf, &xsurf)); 408064c8009SBarry Smith 409064c8009SBarry Smith /* 410064c8009SBarry Smith Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of 411064c8009SBarry Smith Xsurf will be all zero (thus making the coarse matrix singular). 412064c8009SBarry Smith */ 41308401ef6SPierre 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"); 41408401ef6SPierre 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"); 41508401ef6SPierre 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"); 416064c8009SBarry Smith 417064c8009SBarry Smith cnt = 0; 4182fa5cd67SKarl Rupp for (j = 1; j < n - 1 - jstart; j++) { 4192fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 0 * Nsurf] = 1; 420064c8009SBarry Smith } 4212fa5cd67SKarl Rupp 4222fa5cd67SKarl Rupp for (k = 1; k < p - 1 - kstart; k++) { 4232fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 1 * Nsurf] = 1; 4242fa5cd67SKarl Rupp for (j = 1; j < n - 1 - jstart; j++) { 4252fa5cd67SKarl Rupp xsurf[cnt++ + 2 * Nsurf] = 1; 4262fa5cd67SKarl Rupp /* these are the interior nodes */ 4272fa5cd67SKarl Rupp xsurf[cnt++ + 3 * Nsurf] = 1; 4282fa5cd67SKarl Rupp } 4292fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 4 * Nsurf] = 1; 4302fa5cd67SKarl Rupp } 4312fa5cd67SKarl Rupp for (j = 1; j < n - 1 - jstart; j++) { 4322fa5cd67SKarl Rupp for (i = 1; i < m - istart - 1; i++) xsurf[cnt++ + 5 * Nsurf] = 1; 4332fa5cd67SKarl Rupp } 434064c8009SBarry Smith 435064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 436064c8009SBarry Smith for (i = 0; i < Nsurf; i++) { 437064c8009SBarry Smith tmp = 0.0; 4382fa5cd67SKarl Rupp for (j = 0; j < 6; j++) tmp += xsurf[i + j * Nsurf]; 4392fa5cd67SKarl Rupp 44063a3b9bcSJacob 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)); 441064c8009SBarry Smith } 442064c8009SBarry Smith #endif 4439566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(Xsurf, &xsurf)); 4449566063dSJacob Faibussowitsch /* PetscCall(MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD));*/ 445064c8009SBarry Smith 446064c8009SBarry Smith /* 447064c8009SBarry Smith I are the indices for all the needed vertices (in global numbering) 448064c8009SBarry Smith Iint are the indices for the interior values, I surf for the surface values 4497dae84e0SHong Zhang (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it 450aa219208SBarry Smith is NOT the local DMDA ordering.) 451064c8009SBarry Smith IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering 452064c8009SBarry Smith */ 453064c8009SBarry Smith #define Endpoint(a, start, b) (a == 0 || a == (b - 1 - start)) 4549566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(N, &II, Nint, &Iint, Nsurf, &Isurf)); 4559566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nint, &IIint, Nsurf, &IIsurf)); 456064c8009SBarry Smith for (k = 0; k < p - kstart; k++) { 457064c8009SBarry Smith for (j = 0; j < n - jstart; j++) { 458064c8009SBarry Smith for (i = 0; i < m - istart; i++) { 459064c8009SBarry Smith II[c++] = i + j * mwidth + k * mwidth * nwidth; 460064c8009SBarry Smith 461064c8009SBarry Smith if (!Endpoint(i, istart, m) && !Endpoint(j, jstart, n) && !Endpoint(k, kstart, p)) { 462064c8009SBarry Smith IIint[cint] = i + j * mwidth + k * mwidth * nwidth; 463064c8009SBarry Smith Iint[cint++] = i + j * (m - istart) + k * (m - istart) * (n - jstart); 464064c8009SBarry Smith } else { 465064c8009SBarry Smith IIsurf[csurf] = i + j * mwidth + k * mwidth * nwidth; 466064c8009SBarry Smith Isurf[csurf++] = i + j * (m - istart) + k * (m - istart) * (n - jstart); 467064c8009SBarry Smith } 468064c8009SBarry Smith } 469064c8009SBarry Smith } 470064c8009SBarry Smith } 47108401ef6SPierre Jolivet PetscCheck(c == N, PETSC_COMM_SELF, PETSC_ERR_PLIB, "c != N"); 47208401ef6SPierre Jolivet PetscCheck(cint == Nint, PETSC_COMM_SELF, PETSC_ERR_PLIB, "cint != Nint"); 47308401ef6SPierre Jolivet PetscCheck(csurf == Nsurf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "csurf != Nsurf"); 4749566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(da, <g)); 4759566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, N, II, II)); 4769566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, Nint, IIint, IIint)); 4779566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, Nsurf, IIsurf, IIsurf)); 4789566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 4799566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, N, II, PETSC_COPY_VALUES, &is)); 4809566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nint, Iint, PETSC_COPY_VALUES, &isint)); 4819566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Nsurf, Isurf, PETSC_COPY_VALUES, &issurf)); 4829566063dSJacob Faibussowitsch PetscCall(PetscFree3(II, Iint, Isurf)); 483064c8009SBarry Smith 4849566063dSJacob Faibussowitsch PetscCall(ISSort(is)); 4859566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(Aglobal, 1, &is, &is, MAT_INITIAL_MATRIX, &Aholder)); 486064c8009SBarry Smith A = *Aholder; 4879566063dSJacob Faibussowitsch PetscCall(PetscFree(Aholder)); 488064c8009SBarry Smith 4899566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, isint, isint, MAT_INITIAL_MATRIX, &Aii)); 4909566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, isint, issurf, MAT_INITIAL_MATRIX, &Ais)); 4919566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, issurf, isint, MAT_INITIAL_MATRIX, &Asi)); 492064c8009SBarry Smith 493064c8009SBarry Smith /* 494064c8009SBarry Smith Solve for the interpolation onto the interior Xint 495064c8009SBarry Smith */ 4969566063dSJacob Faibussowitsch PetscCall(MatMatMult(Ais, Xsurf, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &Xint_tmp)); 4979566063dSJacob Faibussowitsch PetscCall(MatScale(Xint_tmp, -1.0)); 498064c8009SBarry Smith 4998e722e36SBarry Smith if (exotic->directSolve) { 5009566063dSJacob Faibussowitsch PetscCall(MatGetFactor(Aii, MATSOLVERPETSC, MAT_FACTOR_LU, &iAii)); 5019566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 5029566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(Aii, MATORDERINGND, &row, &col)); 5039566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(iAii, Aii, row, col, &info)); 5049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 5059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col)); 5069566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(iAii, Aii, &info)); 5079566063dSJacob Faibussowitsch PetscCall(MatMatSolve(iAii, Xint_tmp, Xint)); 5089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&iAii)); 509064c8009SBarry Smith } else { 510064c8009SBarry Smith Vec b, x; 511064c8009SBarry Smith PetscScalar *xint_tmp; 512064c8009SBarry Smith 5139566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(Xint, &xint)); 5149566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &x)); 5159566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(Xint_tmp, &xint_tmp)); 5169566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, Nint, NULL, &b)); 5179566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(exotic->ksp, Aii, Aii)); 518064c8009SBarry Smith for (i = 0; i < 6; i++) { 5199566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(x, xint + i * Nint)); 5209566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(b, xint_tmp + i * Nint)); 5219566063dSJacob Faibussowitsch PetscCall(KSPSolve(exotic->ksp, b, x)); 5229566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(exotic->ksp, pc, x)); 5239566063dSJacob Faibussowitsch PetscCall(VecResetArray(x)); 5249566063dSJacob Faibussowitsch PetscCall(VecResetArray(b)); 525064c8009SBarry Smith } 5269566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(Xint, &xint)); 5279566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(Xint_tmp, &xint_tmp)); 5289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 5299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 530064c8009SBarry Smith } 5319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xint_tmp)); 532064c8009SBarry Smith 533064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 5349566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xint, &rxint)); 535064c8009SBarry Smith for (i = 0; i < Nint; i++) { 536064c8009SBarry Smith tmp = 0.0; 5371683a169SBarry Smith for (j = 0; j < 6; j++) tmp += rxint[i + j * Nint]; 5382fa5cd67SKarl Rupp 53963a3b9bcSJacob 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)); 540064c8009SBarry Smith } 5419566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xint, &rxint)); 5429566063dSJacob Faibussowitsch /* PetscCall(MatView(Xint,PETSC_VIEWER_STDOUT_WORLD)); */ 543064c8009SBarry Smith #endif 544064c8009SBarry Smith 545064c8009SBarry Smith /* total faces */ 546064c8009SBarry Smith Ntotal = mp * np * (pp + 1) + mp * pp * (np + 1) + np * pp * (mp + 1); 547064c8009SBarry Smith 548064c8009SBarry Smith /* 549064c8009SBarry Smith For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points 550064c8009SBarry Smith */ 551064c8009SBarry Smith cnt = 0; 552d71ae5a4SJacob Faibussowitsch { 553d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth + 1; 554d71ae5a4SJacob Faibussowitsch } 5559371c9d4SSatish Balay {{gl[cnt++] = mwidth * nwidth + 1; 5569371c9d4SSatish Balay } 557064c8009SBarry Smith { 5589371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + mwidth; /* these are the interior nodes */ 5599371c9d4SSatish Balay gl[cnt++] = mwidth * nwidth + mwidth + m - istart - 1; 5609371c9d4SSatish Balay } 561d71ae5a4SJacob Faibussowitsch { 562d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth + mwidth * (n - jstart - 1) + 1; 563064c8009SBarry Smith } 564d71ae5a4SJacob Faibussowitsch } 565d71ae5a4SJacob Faibussowitsch { 566d71ae5a4SJacob Faibussowitsch gl[cnt++] = mwidth * nwidth * (p - kstart - 1) + mwidth + 1; 567d71ae5a4SJacob Faibussowitsch } 568064c8009SBarry Smith 569064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */ 570064c8009SBarry Smith /* convert that to global numbering and get them on all processes */ 5719566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltg, 6, gl, gl)); 572064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */ 5739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(6 * mp * np * pp, &globals)); 5749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(gl, 6, MPIU_INT, globals, 6, MPIU_INT, PetscObjectComm((PetscObject)da))); 575064c8009SBarry Smith 576064c8009SBarry Smith /* Number the coarse grid points from 0 to Ntotal */ 5779566063dSJacob Faibussowitsch PetscCall(MatGetSize(Aglobal, &Nt, NULL)); 5789566063dSJacob Faibussowitsch PetscCall(PetscTableCreate(Ntotal / 3, Nt + 1, &ht)); 57948a46eb9SPierre Jolivet for (i = 0; i < 6 * mp * np * pp; i++) PetscCall(PetscTableAddCount(ht, globals[i] + 1)); 5809566063dSJacob Faibussowitsch PetscCall(PetscTableGetCount(ht, &cnt)); 58163a3b9bcSJacob 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); 5829566063dSJacob Faibussowitsch PetscCall(PetscFree(globals)); 583064c8009SBarry Smith for (i = 0; i < 6; i++) { 5849566063dSJacob Faibussowitsch PetscCall(PetscTableFind(ht, gl[i] + 1, &gl[i])); 585064c8009SBarry Smith gl[i]--; 586064c8009SBarry Smith } 5879566063dSJacob Faibussowitsch PetscCall(PetscTableDestroy(&ht)); 588064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */ 589064c8009SBarry Smith 590064c8009SBarry Smith /* construct global interpolation matrix */ 5919566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Aglobal, &Ng, NULL)); 592064c8009SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 5939566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PetscObjectComm((PetscObject)da), Ng, PETSC_DECIDE, PETSC_DECIDE, Ntotal, Nint + Nsurf, NULL, Nint, NULL, P)); 594064c8009SBarry Smith } else { 5959566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(*P)); 596064c8009SBarry Smith } 5979566063dSJacob Faibussowitsch PetscCall(MatSetOption(*P, MAT_ROW_ORIENTED, PETSC_FALSE)); 5989566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xint, &rxint)); 5999566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, Nint, IIint, 6, gl, rxint, INSERT_VALUES)); 6009566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xint, &rxint)); 6019566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Xsurf, &rxint)); 6029566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, Nsurf, IIsurf, 6, gl, rxint, INSERT_VALUES)); 6039566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Xsurf, &rxint)); 6049566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY)); 6059566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY)); 6069566063dSJacob Faibussowitsch PetscCall(PetscFree2(IIint, IIsurf)); 607064c8009SBarry Smith 608064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 609064c8009SBarry Smith { 610064c8009SBarry Smith Vec x, y; 611064c8009SBarry Smith PetscScalar *yy; 6129566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), Ng, PETSC_DETERMINE, &y)); 6139566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), PETSC_DETERMINE, Ntotal, &x)); 6149566063dSJacob Faibussowitsch PetscCall(VecSet(x, 1.0)); 6159566063dSJacob Faibussowitsch PetscCall(MatMult(*P, x, y)); 6169566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 617ad540459SPierre 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])); 6189566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 6199566063dSJacob Faibussowitsch PetscCall(VecDestroy(x)); 6209566063dSJacob Faibussowitsch PetscCall(VecDestroy(y)); 621064c8009SBarry Smith } 622064c8009SBarry Smith #endif 623064c8009SBarry Smith 6249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aii)); 6259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ais)); 6269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Asi)); 6279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 6289566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 6299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isint)); 6309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&issurf)); 6319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xint)); 6329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Xsurf)); 633064c8009SBarry Smith PetscFunctionReturn(0); 634064c8009SBarry Smith } 635174b6946SBarry Smith 6367233f9f0SBarry Smith /*@ 6377233f9f0SBarry Smith PCExoticSetType - Sets the type of coarse grid interpolation to use 6387233f9f0SBarry Smith 639f1580f4eSBarry Smith Logically Collective on pc 6407233f9f0SBarry Smith 6417233f9f0SBarry Smith Input Parameters: 6427233f9f0SBarry Smith + pc - the preconditioner context 6437233f9f0SBarry Smith - type - either PC_EXOTIC_FACE or PC_EXOTIC_WIREBASKET (defaults to face) 6447233f9f0SBarry Smith 64595452b02SPatrick Sanan Notes: 64695452b02SPatrick Sanan The face based interpolation has 1 degree of freedom per face and ignores the 647563e08c6SBarry Smith edge and vertex values completely in the coarse problem. For any seven point 648563e08c6SBarry Smith stencil the interpolation of a constant on all faces into the interior is that constant. 649563e08c6SBarry Smith 650563e08c6SBarry Smith The wirebasket interpolation has 1 degree of freedom per vertex, per edge and 651563e08c6SBarry Smith per face. A constant on the subdomain boundary is interpolated as that constant 652563e08c6SBarry Smith in the interior of the domain. 653563e08c6SBarry Smith 654563e08c6SBarry Smith The coarse grid matrix is obtained via the Galerkin computation A_c = R A R^T, hence 655563e08c6SBarry Smith if A is nonsingular A_c is also nonsingular. 656563e08c6SBarry Smith 657563e08c6SBarry Smith Both interpolations are suitable for only scalar problems. 658563e08c6SBarry Smith 6597233f9f0SBarry Smith Level: intermediate 6607233f9f0SBarry Smith 661db781477SPatrick Sanan .seealso: `PCEXOTIC`, `PCExoticType()` 6627233f9f0SBarry Smith @*/ 663d71ae5a4SJacob Faibussowitsch PetscErrorCode PCExoticSetType(PC pc, PCExoticType type) 664d71ae5a4SJacob Faibussowitsch { 6657233f9f0SBarry Smith PetscFunctionBegin; 6660700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 667c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, type, 2); 668cac4c232SBarry Smith PetscTryMethod(pc, "PCExoticSetType_C", (PC, PCExoticType), (pc, type)); 6697233f9f0SBarry Smith PetscFunctionReturn(0); 6707233f9f0SBarry Smith } 6717233f9f0SBarry Smith 672d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCExoticSetType_Exotic(PC pc, PCExoticType type) 673d71ae5a4SJacob Faibussowitsch { 674f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 67531567311SBarry Smith PC_Exotic *ctx = (PC_Exotic *)mg->innerctx; 6767233f9f0SBarry Smith 6777233f9f0SBarry Smith PetscFunctionBegin; 6787233f9f0SBarry Smith ctx->type = type; 6797233f9f0SBarry Smith PetscFunctionReturn(0); 6807233f9f0SBarry Smith } 6817233f9f0SBarry Smith 682d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp_Exotic(PC pc) 683d71ae5a4SJacob Faibussowitsch { 68496bdf778SBarry Smith Mat A; 685f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 68631567311SBarry Smith PC_Exotic *ex = (PC_Exotic *)mg->innerctx; 68796bdf778SBarry Smith MatReuse reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX; 688174b6946SBarry Smith 689174b6946SBarry Smith PetscFunctionBegin; 69028b400f6SJacob Faibussowitsch PetscCheck(pc->dm, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Need to call PCSetDM() before using this PC"); 6919566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, NULL, &A)); 692*0fdf79fbSJacob 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); 6937233f9f0SBarry Smith if (ex->type == PC_EXOTIC_FACE) { 6949566063dSJacob Faibussowitsch PetscCall(DMDAGetFaceInterpolation(pc, pc->dm, ex, A, reuse, &ex->P)); 695*0fdf79fbSJacob Faibussowitsch } else /* if (ex->type == PC_EXOTIC_WIREBASKET) */ { 6969566063dSJacob Faibussowitsch PetscCall(DMDAGetWireBasketInterpolation(pc, pc->dm, ex, A, reuse, &ex->P)); 697*0fdf79fbSJacob Faibussowitsch } 6989566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc, 1, ex->P)); 699d0660788SBarry Smith /* if PC has attached DM we must remove it or the PCMG will use it to compute incorrect sized vectors and interpolations */ 7009566063dSJacob Faibussowitsch PetscCall(PCSetDM(pc, NULL)); 7019566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 702174b6946SBarry Smith PetscFunctionReturn(0); 703174b6946SBarry Smith } 704174b6946SBarry Smith 705d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy_Exotic(PC pc) 706d71ae5a4SJacob Faibussowitsch { 707f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 70831567311SBarry Smith PC_Exotic *ctx = (PC_Exotic *)mg->innerctx; 709174b6946SBarry Smith 710174b6946SBarry Smith PetscFunctionBegin; 7119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->P)); 7129566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ctx->ksp)); 7139566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 7142e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCExoticSetType_C", NULL)); 7159566063dSJacob Faibussowitsch PetscCall(PCDestroy_MG(pc)); 716174b6946SBarry Smith PetscFunctionReturn(0); 717174b6946SBarry Smith } 718174b6946SBarry Smith 719d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView_Exotic(PC pc, PetscViewer viewer) 720d71ae5a4SJacob Faibussowitsch { 721f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 722ace3abfcSBarry Smith PetscBool iascii; 72331567311SBarry Smith PC_Exotic *ctx = (PC_Exotic *)mg->innerctx; 7247233f9f0SBarry Smith 7257233f9f0SBarry Smith PetscFunctionBegin; 7269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 7277233f9f0SBarry Smith if (iascii) { 7289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Exotic type = %s\n", PCExoticTypes[ctx->type])); 7298e722e36SBarry Smith if (ctx->directSolve) { 7309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using direct solver to construct interpolation\n")); 7318e722e36SBarry Smith } else { 7328e722e36SBarry Smith PetscViewer sviewer; 7338e722e36SBarry Smith PetscMPIInt rank; 7348e722e36SBarry Smith 7359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using iterative solver to construct interpolation\n")); 7369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 7379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); /* should not need to push this twice? */ 7389566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 7399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 74048a46eb9SPierre Jolivet if (rank == 0) PetscCall(KSPView(ctx->ksp, sviewer)); 7419566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 7429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 7439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 7448e722e36SBarry Smith } 7457233f9f0SBarry Smith } 7469566063dSJacob Faibussowitsch PetscCall(PCView_MG(pc, viewer)); 7477233f9f0SBarry Smith PetscFunctionReturn(0); 7487233f9f0SBarry Smith } 7497233f9f0SBarry Smith 750d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFromOptions_Exotic(PC pc, PetscOptionItems *PetscOptionsObject) 751d71ae5a4SJacob Faibussowitsch { 752ace3abfcSBarry Smith PetscBool flg; 753f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 7547233f9f0SBarry Smith PCExoticType mgctype; 75531567311SBarry Smith PC_Exotic *ctx = (PC_Exotic *)mg->innerctx; 7567233f9f0SBarry Smith 7577233f9f0SBarry Smith PetscFunctionBegin; 758d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Exotic coarse space options"); 7599566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_exotic_type", "face or wirebasket", "PCExoticSetType", PCExoticTypes, (PetscEnum)ctx->type, (PetscEnum *)&mgctype, &flg)); 7601baa6e33SBarry Smith if (flg) PetscCall(PCExoticSetType(pc, mgctype)); 7619566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_exotic_direct_solver", "use direct solver to construct interpolation", "None", ctx->directSolve, &ctx->directSolve, NULL)); 7628e722e36SBarry Smith if (!ctx->directSolve) { 7638e722e36SBarry Smith if (!ctx->ksp) { 7648e722e36SBarry Smith const char *prefix; 7659566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF, &ctx->ksp)); 7669566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ctx->ksp, pc->erroriffailure)); 7679566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->ksp, (PetscObject)pc, 1)); 7689566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 7699566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ctx->ksp, prefix)); 7709566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ctx->ksp, "exotic_")); 7718e722e36SBarry Smith } 7729566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ctx->ksp)); 7738e722e36SBarry Smith } 774d0609cedSBarry Smith PetscOptionsHeadEnd(); 7757233f9f0SBarry Smith PetscFunctionReturn(0); 7767233f9f0SBarry Smith } 7777233f9f0SBarry Smith 778174b6946SBarry Smith /*MC 7797233f9f0SBarry Smith PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces 780174b6946SBarry Smith 781f1580f4eSBarry Smith This uses the `PCMG` infrastructure restricted to two levels and the face and wirebasket based coarse 78224c3aa18SBarry Smith grid spaces. 78324c3aa18SBarry Smith 784f1580f4eSBarry Smith Level: advanced 785f1580f4eSBarry Smith 78695452b02SPatrick Sanan Notes: 787f1580f4eSBarry Smith Must be used with `DMDA` 788f1580f4eSBarry Smith 789f1580f4eSBarry 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` 79024c3aa18SBarry Smith 79196a0c994SBarry Smith References: 792606c0280SSatish Balay + * - These coarse grid spaces originate in the work of Bramble, Pasciak and Schatz, "The Construction 79396a0c994SBarry Smith of Preconditioners for Elliptic Problems by Substructing IV", Mathematics of Computation, volume 53, 1989. 794606c0280SSatish Balay . * - They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith, 7954f02bc6aSBarry Smith New York University, 1990. 796606c0280SSatish Balay . * - They were then explored in great detail in Dryja, Smith, Widlund, "Schwarz Analysis 7973b65e785SBarry Smith of Iterative Substructuring Methods for Elliptic Problems in Three Dimensions, SIAM Journal on Numerical 79896a0c994SBarry Smith Analysis, volume 31. 1994. These were developed in the context of iterative substructuring preconditioners. 799606c0280SSatish Balay . * - They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund. 8003b65e785SBarry Smith They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example, 8013b65e785SBarry Smith Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Extending theory for domain decomposition algorithms to irregular subdomains. In Ulrich Langer, Marco 8023b65e785SBarry Smith Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings 8033b65e785SBarry Smith of the 17th International Conference on Domain Decomposition Methods in 80496a0c994SBarry Smith Science and Engineering, held in Strobl, Austria, 2006, number 60 in 80596a0c994SBarry Smith Springer Verlag, Lecture Notes in Computational Science and Engineering, 2007. 806606c0280SSatish Balay . * - Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. A family of energy minimizing coarse spaces for overlapping Schwarz preconditioners. In Ulrich Langer, 8073b65e785SBarry Smith Marco Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings 8083b65e785SBarry Smith of the 17th International Conference on Domain Decomposition Methods 80996a0c994SBarry Smith in Science and Engineering, held in Strobl, Austria, 2006, number 60 in 81096a0c994SBarry Smith Springer Verlag, Lecture Notes in Computational Science and Engineering, 2007 811606c0280SSatish Balay . * - Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Domain decomposition 8123b65e785SBarry Smith for less regular subdomains: Overlapping Schwarz in two dimensions. SIAM J. 81396a0c994SBarry Smith Numer. Anal., 46(4), 2008. 814606c0280SSatish Balay - * - Clark R. Dohrmann and Olof B. Widlund. An overlapping Schwarz 8153b65e785SBarry Smith algorithm for almost incompressible elasticity. Technical Report 81696a0c994SBarry Smith TR2008 912, Department of Computer Science, Courant Institute 8173b65e785SBarry Smith of Mathematical Sciences, New York University, May 2008. URL: 8187233f9f0SBarry Smith 819f1580f4eSBarry Smith The usual `PCMG` options are supported, such as -mg_levels_pc_type <type> -mg_coarse_pc_type <type> and -pc_mg_type <type> 820174b6946SBarry Smith 821db781477SPatrick Sanan .seealso: `PCMG`, `PCSetDM()`, `PCExoticType`, `PCExoticSetType()` 822174b6946SBarry Smith M*/ 823174b6946SBarry Smith 824d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Exotic(PC pc) 825d71ae5a4SJacob Faibussowitsch { 8267233f9f0SBarry Smith PC_Exotic *ex; 827f3fbd535SBarry Smith PC_MG *mg; 828174b6946SBarry Smith 829174b6946SBarry Smith PetscFunctionBegin; 830f91d8e95SBarry Smith /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */ 831dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, destroy); 8320a545947SLisandro Dalcin pc->data = NULL; 833dbbe0bcdSBarry Smith 8349566063dSJacob Faibussowitsch PetscCall(PetscFree(((PetscObject)pc)->type_name)); 8350a545947SLisandro Dalcin ((PetscObject)pc)->type_name = NULL; 836f91d8e95SBarry Smith 8379566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCMG)); 8389566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, 2, NULL)); 8399566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_PMAT)); 8409371c9d4SSatish Balay PetscCall(PetscNew(&ex)); 8417233f9f0SBarry Smith ex->type = PC_EXOTIC_FACE; 842f3fbd535SBarry Smith mg = (PC_MG *)pc->data; 84331567311SBarry Smith mg->innerctx = ex; 8447233f9f0SBarry Smith 8457233f9f0SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Exotic; 8467233f9f0SBarry Smith pc->ops->view = PCView_Exotic; 8477233f9f0SBarry Smith pc->ops->destroy = PCDestroy_Exotic; 8486c699258SBarry Smith pc->ops->setup = PCSetUp_Exotic; 8492fa5cd67SKarl Rupp 8509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCExoticSetType_C", PCExoticSetType_Exotic)); 851174b6946SBarry Smith PetscFunctionReturn(0); 852174b6946SBarry Smith } 853