147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 3aa219208SBarry Smith Code for interpolating between grids represented by DMDAs 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 6d52bd9f3SBarry Smith /* 7d52bd9f3SBarry Smith For linear elements there are two branches of code to compute the interpolation. They should compute the same results but may not. The "new version" does 8d52bd9f3SBarry Smith not work for periodic domains, the old does. Change NEWVERSION to 1 to compile in the new version. Eventually when we are sure the two produce identical results 92adb9dcfSBarry Smith we will remove/merge the new version. Based on current tests, these both produce the same results. We are leaving NEWVERSION for now in the code since some 102adb9dcfSBarry Smith consider it cleaner, but old version is turned on since it handles periodic case. 11d52bd9f3SBarry Smith */ 129314d695SBarry Smith #define NEWVERSION 0 1385afcc9aSBarry Smith 14af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 1547c6ae99SBarry Smith 16fdc842d1SBarry Smith /* 17fdc842d1SBarry Smith Since the interpolation uses MATMAIJ for dof > 0 we convert request for non-MATAIJ baseded matrices to MATAIJ. 18fdc842d1SBarry Smith This is a bit of a hack, the reason for it is partially because -dm_mat_type defines the 19fdc842d1SBarry Smith matrix type for both the operator matrices and the interpolation matrices so that users 20fdc842d1SBarry Smith can select matrix types of base MATAIJ for accelerators 21fdc842d1SBarry Smith */ 22d71ae5a4SJacob Faibussowitsch static PetscErrorCode ConvertToAIJ(MatType intype, MatType *outtype) 23d71ae5a4SJacob Faibussowitsch { 24fdc842d1SBarry Smith PetscInt i; 25fdc842d1SBarry Smith char const *types[3] = {MATAIJ, MATSEQAIJ, MATMPIAIJ}; 26fdc842d1SBarry Smith PetscBool flg; 27fdc842d1SBarry Smith 28fdc842d1SBarry Smith PetscFunctionBegin; 2933a4d587SStefano Zampini *outtype = MATAIJ; 30fdc842d1SBarry Smith for (i = 0; i < 3; i++) { 319566063dSJacob Faibussowitsch PetscCall(PetscStrbeginswith(intype, types[i], &flg)); 32fdc842d1SBarry Smith if (flg) { 33fdc842d1SBarry Smith *outtype = intype; 343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35fdc842d1SBarry Smith } 36fdc842d1SBarry Smith } 373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38fdc842d1SBarry Smith } 39fdc842d1SBarry Smith 40d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInterpolation_DA_1D_Q1(DM dac, DM daf, Mat *A) 41d71ae5a4SJacob Faibussowitsch { 428ea3bf28SBarry Smith PetscInt i, i_start, m_f, Mx; 438ea3bf28SBarry Smith const PetscInt *idx_f, *idx_c; 448ea3bf28SBarry Smith PetscInt m_ghost, m_ghost_c; 4547c6ae99SBarry Smith PetscInt row, col, i_start_ghost, mx, m_c, nc, ratio; 4647c6ae99SBarry Smith PetscInt i_c, i_start_c, i_start_ghost_c, cols[2], dof; 4785afcc9aSBarry Smith PetscScalar v[2], x; 4847c6ae99SBarry Smith Mat mat; 49bff4a2f0SMatthew G. Knepley DMBoundaryType bx; 50e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c; 51fdc842d1SBarry Smith MatType mattype; 5247c6ae99SBarry Smith 5347c6ae99SBarry Smith PetscFunctionBegin; 549566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bx, NULL, NULL, NULL)); 559566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 56bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) { 5747c6ae99SBarry Smith ratio = mx / Mx; 5863a3b9bcSJacob Faibussowitsch PetscCheck(ratio * Mx == mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: mx/Mx must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 5947c6ae99SBarry Smith } else { 6047c6ae99SBarry Smith ratio = (mx - 1) / (Mx - 1); 611dca8a05SBarry Smith PetscCheck(ratio * (Mx - 1) == mx - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 6247c6ae99SBarry Smith } 6347c6ae99SBarry Smith 649566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, NULL, NULL, &m_f, NULL, NULL)); 659566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, NULL, NULL, &m_ghost, NULL, NULL)); 669566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f)); 679566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f)); 6847c6ae99SBarry Smith 699566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, NULL, NULL, &m_c, NULL, NULL)); 709566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, NULL, NULL, &m_ghost_c, NULL, NULL)); 719566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c)); 729566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c)); 7347c6ae99SBarry Smith 7447c6ae99SBarry Smith /* create interpolation matrix */ 759566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dac), &mat)); 76fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA) 77fdc842d1SBarry Smith /* 78fdc842d1SBarry Smith Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU 79fdc842d1SBarry Smith we don't want the original unconverted matrix copied to the GPU 80fdc842d1SBarry Smith */ 8148a46eb9SPierre Jolivet if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE)); 82fdc842d1SBarry Smith #endif 839566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, m_f, m_c, mx, Mx)); 849566063dSJacob Faibussowitsch PetscCall(ConvertToAIJ(dac->mattype, &mattype)); 859566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, mattype)); 869566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 2, NULL)); 879566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 2, NULL, 1, NULL)); 8847c6ae99SBarry Smith 8947c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 9085afcc9aSBarry Smith if (!NEWVERSION) { 9147c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) { 9247c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 93e3fbd1f4SBarry Smith row = idx_f[i - i_start_ghost]; 9447c6ae99SBarry Smith 9547c6ae99SBarry Smith i_c = (i / ratio); /* coarse grid node to left of fine grid node */ 9608401ef6SPierre Jolivet PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 979371c9d4SSatish Balay i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT, 989371c9d4SSatish Balay i_start, i_c, i_start_ghost_c); 9947c6ae99SBarry Smith 10047c6ae99SBarry Smith /* 10147c6ae99SBarry Smith Only include those interpolation points that are truly 10247c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 10347c6ae99SBarry Smith in x direction; since they have no right neighbor 10447c6ae99SBarry Smith */ 1056712e2f1SBarry Smith x = ((PetscReal)(i - i_c * ratio)) / ((PetscReal)ratio); 10647c6ae99SBarry Smith nc = 0; 10747c6ae99SBarry Smith /* one left and below; or we are right on it */ 108e3fbd1f4SBarry Smith col = (i_c - i_start_ghost_c); 109e3fbd1f4SBarry Smith cols[nc] = idx_c[col]; 11047c6ae99SBarry Smith v[nc++] = -x + 1.0; 11147c6ae99SBarry Smith /* one right? */ 11247c6ae99SBarry Smith if (i_c * ratio != i) { 113e3fbd1f4SBarry Smith cols[nc] = idx_c[col + 1]; 11447c6ae99SBarry Smith v[nc++] = x; 11547c6ae99SBarry Smith } 1169566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES)); 11747c6ae99SBarry Smith } 118b3bd94feSDave May 119b3bd94feSDave May } else { 120b3bd94feSDave May PetscScalar *xi; 121b3bd94feSDave May PetscInt li, nxi, n; 122b3bd94feSDave May PetscScalar Ni[2]; 123b3bd94feSDave May 124b3bd94feSDave May /* compute local coordinate arrays */ 125b3bd94feSDave May nxi = ratio + 1; 1269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nxi, &xi)); 127ad540459SPierre Jolivet for (li = 0; li < nxi; li++) xi[li] = -1.0 + (PetscScalar)li * (2.0 / (PetscScalar)(nxi - 1)); 128b3bd94feSDave May 129b3bd94feSDave May for (i = i_start; i < i_start + m_f; i++) { 130b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */ 131e3fbd1f4SBarry Smith row = idx_f[(i - i_start_ghost)]; 132b3bd94feSDave May 133b3bd94feSDave May i_c = (i / ratio); /* coarse grid node to left of fine grid node */ 13408401ef6SPierre Jolivet PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 1359371c9d4SSatish Balay i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT, 1369371c9d4SSatish Balay i_start, i_c, i_start_ghost_c); 137b3bd94feSDave May 138b3bd94feSDave May /* remainders */ 139b3bd94feSDave May li = i - ratio * (i / ratio); 1408865f1eaSKarl Rupp if (i == mx - 1) li = nxi - 1; 141b3bd94feSDave May 142b3bd94feSDave May /* corners */ 143e3fbd1f4SBarry Smith col = (i_c - i_start_ghost_c); 144e3fbd1f4SBarry Smith cols[0] = idx_c[col]; 145b3bd94feSDave May Ni[0] = 1.0; 146b3bd94feSDave May if ((li == 0) || (li == nxi - 1)) { 1479566063dSJacob Faibussowitsch PetscCall(MatSetValue(mat, row, cols[0], Ni[0], INSERT_VALUES)); 148b3bd94feSDave May continue; 149b3bd94feSDave May } 150b3bd94feSDave May 151b3bd94feSDave May /* edges + interior */ 152b3bd94feSDave May /* remainders */ 1538865f1eaSKarl Rupp if (i == mx - 1) i_c--; 154b3bd94feSDave May 155e3fbd1f4SBarry Smith col = (i_c - i_start_ghost_c); 156e3fbd1f4SBarry Smith cols[0] = idx_c[col]; /* one left and below; or we are right on it */ 157e3fbd1f4SBarry Smith cols[1] = idx_c[col + 1]; 158b3bd94feSDave May 159b3bd94feSDave May Ni[0] = 0.5 * (1.0 - xi[li]); 160b3bd94feSDave May Ni[1] = 0.5 * (1.0 + xi[li]); 161b3bd94feSDave May for (n = 0; n < 2; n++) { 1628865f1eaSKarl Rupp if (PetscAbsScalar(Ni[n]) < 1.0e-32) cols[n] = -1; 163b3bd94feSDave May } 1649566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, 2, cols, Ni, INSERT_VALUES)); 165b3bd94feSDave May } 1669566063dSJacob Faibussowitsch PetscCall(PetscFree(xi)); 167b3bd94feSDave May } 1689566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f)); 1699566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c)); 1709566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 1719566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 1729566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(mat, dof, A)); 1739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17547c6ae99SBarry Smith } 17647c6ae99SBarry Smith 177d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInterpolation_DA_1D_Q0(DM dac, DM daf, Mat *A) 178d71ae5a4SJacob Faibussowitsch { 1798ea3bf28SBarry Smith PetscInt i, i_start, m_f, Mx; 1808ea3bf28SBarry Smith const PetscInt *idx_f, *idx_c; 181e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c; 1828ea3bf28SBarry Smith PetscInt m_ghost, m_ghost_c; 18347c6ae99SBarry Smith PetscInt row, col, i_start_ghost, mx, m_c, nc, ratio; 18447c6ae99SBarry Smith PetscInt i_c, i_start_c, i_start_ghost_c, cols[2], dof; 18547c6ae99SBarry Smith PetscScalar v[2], x; 18647c6ae99SBarry Smith Mat mat; 187bff4a2f0SMatthew G. Knepley DMBoundaryType bx; 188fdc842d1SBarry Smith MatType mattype; 18947c6ae99SBarry Smith 19047c6ae99SBarry Smith PetscFunctionBegin; 1919566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bx, NULL, NULL, NULL)); 1929566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 193bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) { 19463a3b9bcSJacob Faibussowitsch PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx); 19547c6ae99SBarry Smith ratio = mx / Mx; 19663a3b9bcSJacob Faibussowitsch PetscCheck(ratio * Mx == mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: mx/Mx must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 19747c6ae99SBarry Smith } else { 19863a3b9bcSJacob Faibussowitsch PetscCheck(Mx >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be greater than 1", Mx); 19947c6ae99SBarry Smith ratio = (mx - 1) / (Mx - 1); 2001dca8a05SBarry Smith PetscCheck(ratio * (Mx - 1) == mx - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 20147c6ae99SBarry Smith } 20247c6ae99SBarry Smith 2039566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, NULL, NULL, &m_f, NULL, NULL)); 2049566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, NULL, NULL, &m_ghost, NULL, NULL)); 2059566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f)); 2069566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f)); 20747c6ae99SBarry Smith 2089566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, NULL, NULL, &m_c, NULL, NULL)); 2099566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, NULL, NULL, &m_ghost_c, NULL, NULL)); 2109566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c)); 2119566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c)); 21247c6ae99SBarry Smith 21347c6ae99SBarry Smith /* create interpolation matrix */ 2149566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dac), &mat)); 215fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA) 216fdc842d1SBarry Smith /* 217fdc842d1SBarry Smith Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU 218fdc842d1SBarry Smith we don't want the original unconverted matrix copied to the GPU 219fdc842d1SBarry Smith */ 22048a46eb9SPierre Jolivet if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE)); 221fdc842d1SBarry Smith #endif 2229566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, m_f, m_c, mx, Mx)); 2239566063dSJacob Faibussowitsch PetscCall(ConvertToAIJ(dac->mattype, &mattype)); 2249566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, mattype)); 2259566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 2, NULL)); 2269566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 2, NULL, 0, NULL)); 22747c6ae99SBarry Smith 22847c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 22947c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) { 23047c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 231e3fbd1f4SBarry Smith row = idx_f[(i - i_start_ghost)]; 23247c6ae99SBarry Smith 23347c6ae99SBarry Smith i_c = (i / ratio); /* coarse grid node to left of fine grid node */ 23447c6ae99SBarry Smith 23547c6ae99SBarry Smith /* 23647c6ae99SBarry Smith Only include those interpolation points that are truly 23747c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 23847c6ae99SBarry Smith in x direction; since they have no right neighbor 23947c6ae99SBarry Smith */ 2406712e2f1SBarry Smith x = ((PetscReal)(i - i_c * ratio)) / ((PetscReal)ratio); 24147c6ae99SBarry Smith nc = 0; 24247c6ae99SBarry Smith /* one left and below; or we are right on it */ 243e3fbd1f4SBarry Smith col = (i_c - i_start_ghost_c); 244e3fbd1f4SBarry Smith cols[nc] = idx_c[col]; 24547c6ae99SBarry Smith v[nc++] = -x + 1.0; 24647c6ae99SBarry Smith /* one right? */ 24747c6ae99SBarry Smith if (i_c * ratio != i) { 248e3fbd1f4SBarry Smith cols[nc] = idx_c[col + 1]; 24947c6ae99SBarry Smith v[nc++] = x; 25047c6ae99SBarry Smith } 2519566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES)); 25247c6ae99SBarry Smith } 2539566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f)); 2549566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c)); 2559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 2569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 2579566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(mat, dof, A)); 2589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 2599566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(5.0 * m_f)); 2603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26147c6ae99SBarry Smith } 26247c6ae99SBarry Smith 263d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInterpolation_DA_2D_Q1(DM dac, DM daf, Mat *A) 264d71ae5a4SJacob Faibussowitsch { 2658ea3bf28SBarry Smith PetscInt i, j, i_start, j_start, m_f, n_f, Mx, My, dof; 2668ea3bf28SBarry Smith const PetscInt *idx_c, *idx_f; 267e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c; 2688ea3bf28SBarry Smith PetscInt m_ghost, n_ghost, m_ghost_c, n_ghost_c, *dnz, *onz; 26947c6ae99SBarry Smith PetscInt row, col, i_start_ghost, j_start_ghost, cols[4], mx, m_c, my, nc, ratioi, ratioj; 27047c6ae99SBarry Smith PetscInt i_c, j_c, i_start_c, j_start_c, n_c, i_start_ghost_c, j_start_ghost_c, col_shift, col_scale; 27147c6ae99SBarry Smith PetscMPIInt size_c, size_f, rank_f; 27247c6ae99SBarry Smith PetscScalar v[4], x, y; 27347c6ae99SBarry Smith Mat mat; 274bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by; 275fdc842d1SBarry Smith MatType mattype; 27647c6ae99SBarry Smith 27747c6ae99SBarry Smith PetscFunctionBegin; 2789566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, NULL, NULL, NULL, NULL, NULL, NULL, &bx, &by, NULL, NULL)); 2799566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 280bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) { 28163a3b9bcSJacob Faibussowitsch PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx); 28247c6ae99SBarry Smith ratioi = mx / Mx; 28363a3b9bcSJacob Faibussowitsch PetscCheck(ratioi * Mx == mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: mx/Mx must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 28447c6ae99SBarry Smith } else { 28563a3b9bcSJacob Faibussowitsch PetscCheck(Mx >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be greater than 1", Mx); 28647c6ae99SBarry Smith ratioi = (mx - 1) / (Mx - 1); 2871dca8a05SBarry Smith PetscCheck(ratioi * (Mx - 1) == mx - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 28847c6ae99SBarry Smith } 289bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) { 29063a3b9bcSJacob Faibussowitsch PetscCheck(My, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be positive", My); 29147c6ae99SBarry Smith ratioj = my / My; 29263a3b9bcSJacob Faibussowitsch PetscCheck(ratioj * My == my, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: my/My must be integer: my %" PetscInt_FMT " My %" PetscInt_FMT, my, My); 29347c6ae99SBarry Smith } else { 29463a3b9bcSJacob Faibussowitsch PetscCheck(My >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be greater than 1", My); 29547c6ae99SBarry Smith ratioj = (my - 1) / (My - 1); 2961dca8a05SBarry Smith PetscCheck(ratioj * (My - 1) == my - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (my - 1)/(My - 1) must be integer: my %" PetscInt_FMT " My %" PetscInt_FMT, my, My); 29747c6ae99SBarry Smith } 29847c6ae99SBarry Smith 2999566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, NULL, &m_f, &n_f, NULL)); 3009566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, NULL, &m_ghost, &n_ghost, NULL)); 3019566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f)); 3029566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f)); 30347c6ae99SBarry Smith 3049566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, NULL, &m_c, &n_c, NULL)); 3059566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, NULL, &m_ghost_c, &n_ghost_c, NULL)); 3069566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c)); 3079566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c)); 30847c6ae99SBarry Smith 30947c6ae99SBarry Smith /* 310aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 31147c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 31247c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 31347c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 31447c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 31547c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 31647c6ae99SBarry Smith 31747c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 31847c6ae99SBarry Smith */ 3199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dac), &size_c)); 3209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)daf), &size_f)); 3219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)daf), &rank_f)); 32247c6ae99SBarry Smith col_scale = size_f / size_c; 32347c6ae99SBarry Smith col_shift = Mx * My * (rank_f / size_c); 32447c6ae99SBarry Smith 325d0609cedSBarry Smith MatPreallocateBegin(PetscObjectComm((PetscObject)daf), m_f * n_f, col_scale * m_c * n_c, dnz, onz); 32647c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) { 32747c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) { 32847c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 329e3fbd1f4SBarry Smith row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))]; 33047c6ae99SBarry Smith 33147c6ae99SBarry Smith i_c = (i / ratioi); /* coarse grid node to left of fine grid node */ 33247c6ae99SBarry Smith j_c = (j / ratioj); /* coarse grid node below fine grid node */ 33347c6ae99SBarry Smith 33408401ef6SPierre Jolivet PetscCheck(j_c >= j_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 3359371c9d4SSatish Balay j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT, 3369371c9d4SSatish Balay j_start, j_c, j_start_ghost_c); 33708401ef6SPierre Jolivet PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 3389371c9d4SSatish Balay i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT, 3399371c9d4SSatish Balay i_start, i_c, i_start_ghost_c); 34047c6ae99SBarry Smith 34147c6ae99SBarry Smith /* 34247c6ae99SBarry Smith Only include those interpolation points that are truly 34347c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 34447c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 34547c6ae99SBarry Smith */ 34647c6ae99SBarry Smith nc = 0; 34747c6ae99SBarry Smith /* one left and below; or we are right on it */ 348e3fbd1f4SBarry Smith col = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c)); 349e3fbd1f4SBarry Smith cols[nc++] = col_shift + idx_c[col]; 35047c6ae99SBarry Smith /* one right and below */ 351e3fbd1f4SBarry Smith if (i_c * ratioi != i) cols[nc++] = col_shift + idx_c[col + 1]; 35247c6ae99SBarry Smith /* one left and above */ 353e3fbd1f4SBarry Smith if (j_c * ratioj != j) cols[nc++] = col_shift + idx_c[col + m_ghost_c]; 35447c6ae99SBarry Smith /* one right and above */ 355e3fbd1f4SBarry Smith if (i_c * ratioi != i && j_c * ratioj != j) cols[nc++] = col_shift + idx_c[col + (m_ghost_c + 1)]; 3569566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, nc, cols, dnz, onz)); 35747c6ae99SBarry Smith } 35847c6ae99SBarry Smith } 3599566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)daf), &mat)); 360fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA) 361fdc842d1SBarry Smith /* 362fdc842d1SBarry Smith Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU 363fdc842d1SBarry Smith we don't want the original unconverted matrix copied to the GPU 364fdc842d1SBarry Smith */ 36548a46eb9SPierre Jolivet if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE)); 366fdc842d1SBarry Smith #endif 3679566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, m_f * n_f, col_scale * m_c * n_c, mx * my, col_scale * Mx * My)); 3689566063dSJacob Faibussowitsch PetscCall(ConvertToAIJ(dac->mattype, &mattype)); 3699566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, mattype)); 3709566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 0, dnz)); 3719566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 0, dnz, 0, onz)); 372d0609cedSBarry Smith MatPreallocateEnd(dnz, onz); 37347c6ae99SBarry Smith 37447c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 37585afcc9aSBarry Smith if (!NEWVERSION) { 37647c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) { 37747c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) { 37847c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 379e3fbd1f4SBarry Smith row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))]; 38047c6ae99SBarry Smith 38147c6ae99SBarry Smith i_c = (i / ratioi); /* coarse grid node to left of fine grid node */ 38247c6ae99SBarry Smith j_c = (j / ratioj); /* coarse grid node below fine grid node */ 38347c6ae99SBarry Smith 38447c6ae99SBarry Smith /* 38547c6ae99SBarry Smith Only include those interpolation points that are truly 38647c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 38747c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 38847c6ae99SBarry Smith */ 3896712e2f1SBarry Smith x = ((PetscReal)(i - i_c * ratioi)) / ((PetscReal)ratioi); 3906712e2f1SBarry Smith y = ((PetscReal)(j - j_c * ratioj)) / ((PetscReal)ratioj); 391b3bd94feSDave May 39247c6ae99SBarry Smith nc = 0; 39347c6ae99SBarry Smith /* one left and below; or we are right on it */ 394e3fbd1f4SBarry Smith col = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c)); 395e3fbd1f4SBarry Smith cols[nc] = col_shift + idx_c[col]; 39647c6ae99SBarry Smith v[nc++] = x * y - x - y + 1.0; 39747c6ae99SBarry Smith /* one right and below */ 39847c6ae99SBarry Smith if (i_c * ratioi != i) { 399e3fbd1f4SBarry Smith cols[nc] = col_shift + idx_c[col + 1]; 40047c6ae99SBarry Smith v[nc++] = -x * y + x; 40147c6ae99SBarry Smith } 40247c6ae99SBarry Smith /* one left and above */ 40347c6ae99SBarry Smith if (j_c * ratioj != j) { 404e3fbd1f4SBarry Smith cols[nc] = col_shift + idx_c[col + m_ghost_c]; 40547c6ae99SBarry Smith v[nc++] = -x * y + y; 40647c6ae99SBarry Smith } 40747c6ae99SBarry Smith /* one right and above */ 40847c6ae99SBarry Smith if (j_c * ratioj != j && i_c * ratioi != i) { 409e3fbd1f4SBarry Smith cols[nc] = col_shift + idx_c[col + (m_ghost_c + 1)]; 41047c6ae99SBarry Smith v[nc++] = x * y; 41147c6ae99SBarry Smith } 4129566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES)); 41347c6ae99SBarry Smith } 41447c6ae99SBarry Smith } 415b3bd94feSDave May 416b3bd94feSDave May } else { 417b3bd94feSDave May PetscScalar Ni[4]; 418b3bd94feSDave May PetscScalar *xi, *eta; 419b3bd94feSDave May PetscInt li, nxi, lj, neta; 420b3bd94feSDave May 421b3bd94feSDave May /* compute local coordinate arrays */ 422b3bd94feSDave May nxi = ratioi + 1; 423b3bd94feSDave May neta = ratioj + 1; 4249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nxi, &xi)); 4259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(neta, &eta)); 426ad540459SPierre Jolivet for (li = 0; li < nxi; li++) xi[li] = -1.0 + (PetscScalar)li * (2.0 / (PetscScalar)(nxi - 1)); 427ad540459SPierre Jolivet for (lj = 0; lj < neta; lj++) eta[lj] = -1.0 + (PetscScalar)lj * (2.0 / (PetscScalar)(neta - 1)); 428b3bd94feSDave May 429b3bd94feSDave May /* loop over local fine grid nodes setting interpolation for those*/ 430b3bd94feSDave May for (j = j_start; j < j_start + n_f; j++) { 431b3bd94feSDave May for (i = i_start; i < i_start + m_f; i++) { 432b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */ 433e3fbd1f4SBarry Smith row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))]; 434b3bd94feSDave May 435b3bd94feSDave May i_c = (i / ratioi); /* coarse grid node to left of fine grid node */ 436b3bd94feSDave May j_c = (j / ratioj); /* coarse grid node below fine grid node */ 437b3bd94feSDave May 438b3bd94feSDave May /* remainders */ 439b3bd94feSDave May li = i - ratioi * (i / ratioi); 4408865f1eaSKarl Rupp if (i == mx - 1) li = nxi - 1; 441b3bd94feSDave May lj = j - ratioj * (j / ratioj); 4428865f1eaSKarl Rupp if (j == my - 1) lj = neta - 1; 443b3bd94feSDave May 444b3bd94feSDave May /* corners */ 445e3fbd1f4SBarry Smith col = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c)); 446e3fbd1f4SBarry Smith cols[0] = col_shift + idx_c[col]; /* left, below */ 447b3bd94feSDave May Ni[0] = 1.0; 448b3bd94feSDave May if ((li == 0) || (li == nxi - 1)) { 449b3bd94feSDave May if ((lj == 0) || (lj == neta - 1)) { 4509566063dSJacob Faibussowitsch PetscCall(MatSetValue(mat, row, cols[0], Ni[0], INSERT_VALUES)); 451b3bd94feSDave May continue; 452b3bd94feSDave May } 453b3bd94feSDave May } 454b3bd94feSDave May 455b3bd94feSDave May /* edges + interior */ 456b3bd94feSDave May /* remainders */ 4578865f1eaSKarl Rupp if (i == mx - 1) i_c--; 4588865f1eaSKarl Rupp if (j == my - 1) j_c--; 459b3bd94feSDave May 460e3fbd1f4SBarry Smith col = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c)); 461e3fbd1f4SBarry Smith cols[0] = col_shift + idx_c[col]; /* left, below */ 462e3fbd1f4SBarry Smith cols[1] = col_shift + idx_c[col + 1]; /* right, below */ 463e3fbd1f4SBarry Smith cols[2] = col_shift + idx_c[col + m_ghost_c]; /* left, above */ 464e3fbd1f4SBarry Smith cols[3] = col_shift + idx_c[col + (m_ghost_c + 1)]; /* right, above */ 465b3bd94feSDave May 466b3bd94feSDave May Ni[0] = 0.25 * (1.0 - xi[li]) * (1.0 - eta[lj]); 467b3bd94feSDave May Ni[1] = 0.25 * (1.0 + xi[li]) * (1.0 - eta[lj]); 468b3bd94feSDave May Ni[2] = 0.25 * (1.0 - xi[li]) * (1.0 + eta[lj]); 469b3bd94feSDave May Ni[3] = 0.25 * (1.0 + xi[li]) * (1.0 + eta[lj]); 470b3bd94feSDave May 471b3bd94feSDave May nc = 0; 4728865f1eaSKarl Rupp if (PetscAbsScalar(Ni[0]) < 1.0e-32) cols[0] = -1; 4738865f1eaSKarl Rupp if (PetscAbsScalar(Ni[1]) < 1.0e-32) cols[1] = -1; 4748865f1eaSKarl Rupp if (PetscAbsScalar(Ni[2]) < 1.0e-32) cols[2] = -1; 4758865f1eaSKarl Rupp if (PetscAbsScalar(Ni[3]) < 1.0e-32) cols[3] = -1; 476b3bd94feSDave May 4779566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, 4, cols, Ni, INSERT_VALUES)); 478b3bd94feSDave May } 479b3bd94feSDave May } 4809566063dSJacob Faibussowitsch PetscCall(PetscFree(xi)); 4819566063dSJacob Faibussowitsch PetscCall(PetscFree(eta)); 482b3bd94feSDave May } 4839566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f)); 4849566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c)); 4859566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 4869566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 4879566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(mat, dof, A)); 4889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 4893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49047c6ae99SBarry Smith } 49147c6ae99SBarry Smith 49247c6ae99SBarry Smith /* 49347c6ae99SBarry Smith Contributed by Andrei Draganescu <aidraga@sandia.gov> 49447c6ae99SBarry Smith */ 495d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInterpolation_DA_2D_Q0(DM dac, DM daf, Mat *A) 496d71ae5a4SJacob Faibussowitsch { 4978ea3bf28SBarry Smith PetscInt i, j, i_start, j_start, m_f, n_f, Mx, My, dof; 4988ea3bf28SBarry Smith const PetscInt *idx_c, *idx_f; 499e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c; 5008ea3bf28SBarry Smith PetscInt m_ghost, n_ghost, m_ghost_c, n_ghost_c, *dnz, *onz; 50147c6ae99SBarry Smith PetscInt row, col, i_start_ghost, j_start_ghost, cols[4], mx, m_c, my, nc, ratioi, ratioj; 50247c6ae99SBarry Smith PetscInt i_c, j_c, i_start_c, j_start_c, n_c, i_start_ghost_c, j_start_ghost_c, col_shift, col_scale; 50347c6ae99SBarry Smith PetscMPIInt size_c, size_f, rank_f; 50447c6ae99SBarry Smith PetscScalar v[4]; 50547c6ae99SBarry Smith Mat mat; 506bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by; 507fdc842d1SBarry Smith MatType mattype; 50847c6ae99SBarry Smith 50947c6ae99SBarry Smith PetscFunctionBegin; 5109566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, NULL, NULL, NULL, NULL, NULL, NULL, &bx, &by, NULL, NULL)); 5119566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 51263a3b9bcSJacob Faibussowitsch PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx); 51363a3b9bcSJacob Faibussowitsch PetscCheck(My, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be positive", My); 51447c6ae99SBarry Smith ratioi = mx / Mx; 51547c6ae99SBarry Smith ratioj = my / My; 51608401ef6SPierre Jolivet PetscCheck(ratioi * Mx == mx, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in x"); 51708401ef6SPierre Jolivet PetscCheck(ratioj * My == my, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in y"); 51808401ef6SPierre Jolivet PetscCheck(ratioi == 2, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in x must be 2"); 51908401ef6SPierre Jolivet PetscCheck(ratioj == 2, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in y must be 2"); 52047c6ae99SBarry Smith 5219566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, NULL, &m_f, &n_f, NULL)); 5229566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, NULL, &m_ghost, &n_ghost, NULL)); 5239566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f)); 5249566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f)); 52547c6ae99SBarry Smith 5269566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, NULL, &m_c, &n_c, NULL)); 5279566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, NULL, &m_ghost_c, &n_ghost_c, NULL)); 5289566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c)); 5299566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c)); 53047c6ae99SBarry Smith 53147c6ae99SBarry Smith /* 532aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 53347c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 53447c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 53547c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 53647c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 53747c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 53847c6ae99SBarry Smith 53947c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 54047c6ae99SBarry Smith */ 5419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dac), &size_c)); 5429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)daf), &size_f)); 5439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)daf), &rank_f)); 54447c6ae99SBarry Smith col_scale = size_f / size_c; 54547c6ae99SBarry Smith col_shift = Mx * My * (rank_f / size_c); 54647c6ae99SBarry Smith 547d0609cedSBarry Smith MatPreallocateBegin(PetscObjectComm((PetscObject)daf), m_f * n_f, col_scale * m_c * n_c, dnz, onz); 54847c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) { 54947c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) { 55047c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 551e3fbd1f4SBarry Smith row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))]; 55247c6ae99SBarry Smith 55347c6ae99SBarry Smith i_c = (i / ratioi); /* coarse grid node to left of fine grid node */ 55447c6ae99SBarry Smith j_c = (j / ratioj); /* coarse grid node below fine grid node */ 55547c6ae99SBarry Smith 55608401ef6SPierre Jolivet PetscCheck(j_c >= j_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 5579371c9d4SSatish Balay j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT, 5589371c9d4SSatish Balay j_start, j_c, j_start_ghost_c); 55908401ef6SPierre Jolivet PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 5609371c9d4SSatish Balay i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT, 5619371c9d4SSatish Balay i_start, i_c, i_start_ghost_c); 56247c6ae99SBarry Smith 56347c6ae99SBarry Smith /* 56447c6ae99SBarry Smith Only include those interpolation points that are truly 56547c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 56647c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 56747c6ae99SBarry Smith */ 56847c6ae99SBarry Smith nc = 0; 56947c6ae99SBarry Smith /* one left and below; or we are right on it */ 570e3fbd1f4SBarry Smith col = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c)); 571e3fbd1f4SBarry Smith cols[nc++] = col_shift + idx_c[col]; 5729566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, nc, cols, dnz, onz)); 57347c6ae99SBarry Smith } 57447c6ae99SBarry Smith } 5759566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)daf), &mat)); 576fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA) 577fdc842d1SBarry Smith /* 578fdc842d1SBarry Smith Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU 579fdc842d1SBarry Smith we don't want the original unconverted matrix copied to the GPU 580fdc842d1SBarry Smith */ 58148a46eb9SPierre Jolivet if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE)); 582fdc842d1SBarry Smith #endif 5839566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, m_f * n_f, col_scale * m_c * n_c, mx * my, col_scale * Mx * My)); 5849566063dSJacob Faibussowitsch PetscCall(ConvertToAIJ(dac->mattype, &mattype)); 5859566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, mattype)); 5869566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 0, dnz)); 5879566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 0, dnz, 0, onz)); 588d0609cedSBarry Smith MatPreallocateEnd(dnz, onz); 58947c6ae99SBarry Smith 59047c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 59147c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) { 59247c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) { 59347c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 594e3fbd1f4SBarry Smith row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))]; 59547c6ae99SBarry Smith 59647c6ae99SBarry Smith i_c = (i / ratioi); /* coarse grid node to left of fine grid node */ 59747c6ae99SBarry Smith j_c = (j / ratioj); /* coarse grid node below fine grid node */ 59847c6ae99SBarry Smith nc = 0; 59947c6ae99SBarry Smith /* one left and below; or we are right on it */ 600e3fbd1f4SBarry Smith col = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c)); 601e3fbd1f4SBarry Smith cols[nc] = col_shift + idx_c[col]; 60247c6ae99SBarry Smith v[nc++] = 1.0; 60347c6ae99SBarry Smith 6049566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES)); 60547c6ae99SBarry Smith } 60647c6ae99SBarry Smith } 6079566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f)); 6089566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c)); 6099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 6109566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 6119566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(mat, dof, A)); 6129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 6139566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(13.0 * m_f * n_f)); 6143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 61547c6ae99SBarry Smith } 61647c6ae99SBarry Smith 61747c6ae99SBarry Smith /* 61847c6ae99SBarry Smith Contributed by Jianming Yang <jianming-yang@uiowa.edu> 61947c6ae99SBarry Smith */ 620d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInterpolation_DA_3D_Q0(DM dac, DM daf, Mat *A) 621d71ae5a4SJacob Faibussowitsch { 6228ea3bf28SBarry Smith PetscInt i, j, l, i_start, j_start, l_start, m_f, n_f, p_f, Mx, My, Mz, dof; 6238ea3bf28SBarry Smith const PetscInt *idx_c, *idx_f; 624e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c; 6258ea3bf28SBarry Smith PetscInt m_ghost, n_ghost, p_ghost, m_ghost_c, n_ghost_c, p_ghost_c, nc, *dnz, *onz; 62647c6ae99SBarry Smith PetscInt row, col, i_start_ghost, j_start_ghost, l_start_ghost, cols[8], mx, m_c, my, n_c, mz, p_c, ratioi, ratioj, ratiol; 62747c6ae99SBarry Smith PetscInt i_c, j_c, l_c, i_start_c, j_start_c, l_start_c, i_start_ghost_c, j_start_ghost_c, l_start_ghost_c, col_shift, col_scale; 62847c6ae99SBarry Smith PetscMPIInt size_c, size_f, rank_f; 62947c6ae99SBarry Smith PetscScalar v[8]; 63047c6ae99SBarry Smith Mat mat; 631bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by, bz; 632fdc842d1SBarry Smith MatType mattype; 63347c6ae99SBarry Smith 63447c6ae99SBarry Smith PetscFunctionBegin; 6359566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, &Mz, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL)); 63663a3b9bcSJacob Faibussowitsch PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx); 63763a3b9bcSJacob Faibussowitsch PetscCheck(My, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be positive", My); 63863a3b9bcSJacob Faibussowitsch PetscCheck(Mz, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of z coarse grid points %" PetscInt_FMT " must be positive", Mz); 6399566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, &mz, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 64047c6ae99SBarry Smith ratioi = mx / Mx; 64147c6ae99SBarry Smith ratioj = my / My; 64247c6ae99SBarry Smith ratiol = mz / Mz; 64308401ef6SPierre Jolivet PetscCheck(ratioi * Mx == mx, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in x"); 64408401ef6SPierre Jolivet PetscCheck(ratioj * My == my, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in y"); 64508401ef6SPierre Jolivet PetscCheck(ratiol * Mz == mz, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in z"); 6461dca8a05SBarry Smith PetscCheck(ratioi == 2 || ratioi == 1, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in x must be 1 or 2"); 6471dca8a05SBarry Smith PetscCheck(ratioj == 2 || ratioj == 1, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in y must be 1 or 2"); 6481dca8a05SBarry Smith PetscCheck(ratiol == 2 || ratiol == 1, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in z must be 1 or 2"); 64947c6ae99SBarry Smith 6509566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, &l_start, &m_f, &n_f, &p_f)); 6519566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, &l_start_ghost, &m_ghost, &n_ghost, &p_ghost)); 6529566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f)); 6539566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f)); 65447c6ae99SBarry Smith 6559566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, &l_start_c, &m_c, &n_c, &p_c)); 6569566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, &l_start_ghost_c, &m_ghost_c, &n_ghost_c, &p_ghost_c)); 6579566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c)); 6589566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c)); 659e3fbd1f4SBarry Smith 66047c6ae99SBarry Smith /* 661aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 66247c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 66347c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 66447c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 66547c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 66647c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 66747c6ae99SBarry Smith 66847c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 66947c6ae99SBarry Smith */ 6709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dac), &size_c)); 6719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)daf), &size_f)); 6729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)daf), &rank_f)); 67347c6ae99SBarry Smith col_scale = size_f / size_c; 67447c6ae99SBarry Smith col_shift = Mx * My * Mz * (rank_f / size_c); 67547c6ae99SBarry Smith 676d0609cedSBarry Smith MatPreallocateBegin(PetscObjectComm((PetscObject)daf), m_f * n_f * p_f, col_scale * m_c * n_c * p_c, dnz, onz); 67747c6ae99SBarry Smith for (l = l_start; l < l_start + p_f; l++) { 67847c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) { 67947c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) { 68047c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 681e3fbd1f4SBarry Smith row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))]; 68247c6ae99SBarry Smith 68347c6ae99SBarry Smith i_c = (i / ratioi); /* coarse grid node to left of fine grid node */ 68447c6ae99SBarry Smith j_c = (j / ratioj); /* coarse grid node below fine grid node */ 68547c6ae99SBarry Smith l_c = (l / ratiol); 68647c6ae99SBarry Smith 68708401ef6SPierre Jolivet PetscCheck(l_c >= l_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 6889371c9d4SSatish Balay l_start %" PetscInt_FMT " l_c %" PetscInt_FMT " l_start_ghost_c %" PetscInt_FMT, 6899371c9d4SSatish Balay l_start, l_c, l_start_ghost_c); 69008401ef6SPierre Jolivet PetscCheck(j_c >= j_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 6919371c9d4SSatish Balay j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT, 6929371c9d4SSatish Balay j_start, j_c, j_start_ghost_c); 69308401ef6SPierre Jolivet PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 6949371c9d4SSatish Balay i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT, 6959371c9d4SSatish Balay i_start, i_c, i_start_ghost_c); 69647c6ae99SBarry Smith 69747c6ae99SBarry Smith /* 69847c6ae99SBarry Smith Only include those interpolation points that are truly 69947c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 70047c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 70147c6ae99SBarry Smith */ 70247c6ae99SBarry Smith nc = 0; 70347c6ae99SBarry Smith /* one left and below; or we are right on it */ 704e3fbd1f4SBarry Smith col = (m_ghost_c * n_ghost_c * (l_c - l_start_ghost_c) + m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c)); 705e3fbd1f4SBarry Smith cols[nc++] = col_shift + idx_c[col]; 7069566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, nc, cols, dnz, onz)); 70747c6ae99SBarry Smith } 70847c6ae99SBarry Smith } 70947c6ae99SBarry Smith } 7109566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)daf), &mat)); 711fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA) 712fdc842d1SBarry Smith /* 713fdc842d1SBarry Smith Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU 714fdc842d1SBarry Smith we don't want the original unconverted matrix copied to the GPU 715fdc842d1SBarry Smith */ 71648a46eb9SPierre Jolivet if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE)); 717fdc842d1SBarry Smith #endif 7189566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, m_f * n_f * p_f, col_scale * m_c * n_c * p_c, mx * my * mz, col_scale * Mx * My * Mz)); 7199566063dSJacob Faibussowitsch PetscCall(ConvertToAIJ(dac->mattype, &mattype)); 7209566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, mattype)); 7219566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 0, dnz)); 7229566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 0, dnz, 0, onz)); 723d0609cedSBarry Smith MatPreallocateEnd(dnz, onz); 72447c6ae99SBarry Smith 72547c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 72647c6ae99SBarry Smith for (l = l_start; l < l_start + p_f; l++) { 72747c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) { 72847c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) { 72947c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 730e3fbd1f4SBarry Smith row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))]; 73147c6ae99SBarry Smith 73247c6ae99SBarry Smith i_c = (i / ratioi); /* coarse grid node to left of fine grid node */ 73347c6ae99SBarry Smith j_c = (j / ratioj); /* coarse grid node below fine grid node */ 73447c6ae99SBarry Smith l_c = (l / ratiol); 73547c6ae99SBarry Smith nc = 0; 73647c6ae99SBarry Smith /* one left and below; or we are right on it */ 737e3fbd1f4SBarry Smith col = (m_ghost_c * n_ghost_c * (l_c - l_start_ghost_c) + m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c)); 738e3fbd1f4SBarry Smith cols[nc] = col_shift + idx_c[col]; 73947c6ae99SBarry Smith v[nc++] = 1.0; 74047c6ae99SBarry Smith 7419566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES)); 74247c6ae99SBarry Smith } 74347c6ae99SBarry Smith } 74447c6ae99SBarry Smith } 7459566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f)); 7469566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c)); 7479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 7489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 7499566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(mat, dof, A)); 7509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 7519566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(13.0 * m_f * n_f * p_f)); 7523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 75347c6ae99SBarry Smith } 75447c6ae99SBarry Smith 755d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInterpolation_DA_3D_Q1(DM dac, DM daf, Mat *A) 756d71ae5a4SJacob Faibussowitsch { 7578ea3bf28SBarry Smith PetscInt i, j, i_start, j_start, m_f, n_f, Mx, My, dof, l; 7588ea3bf28SBarry Smith const PetscInt *idx_c, *idx_f; 759e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c; 7608ea3bf28SBarry Smith PetscInt m_ghost, n_ghost, m_ghost_c, n_ghost_c, Mz, mz; 76147c6ae99SBarry Smith PetscInt row, col, i_start_ghost, j_start_ghost, cols[8], mx, m_c, my, nc, ratioi, ratioj, ratiok; 76247c6ae99SBarry Smith PetscInt i_c, j_c, i_start_c, j_start_c, n_c, i_start_ghost_c, j_start_ghost_c; 76347c6ae99SBarry Smith PetscInt l_start, p_f, l_start_ghost, p_ghost, l_start_c, p_c; 76447c6ae99SBarry Smith PetscInt l_start_ghost_c, p_ghost_c, l_c, *dnz, *onz; 76547c6ae99SBarry Smith PetscScalar v[8], x, y, z; 76647c6ae99SBarry Smith Mat mat; 767bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by, bz; 768fdc842d1SBarry Smith MatType mattype; 76947c6ae99SBarry Smith 77047c6ae99SBarry Smith PetscFunctionBegin; 7719566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, &Mz, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL)); 7729566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, &mz, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 77347c6ae99SBarry Smith if (mx == Mx) { 77447c6ae99SBarry Smith ratioi = 1; 775bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_PERIODIC) { 77663a3b9bcSJacob Faibussowitsch PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx); 77747c6ae99SBarry Smith ratioi = mx / Mx; 77863a3b9bcSJacob Faibussowitsch PetscCheck(ratioi * Mx == mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: mx/Mx must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 77947c6ae99SBarry Smith } else { 78063a3b9bcSJacob Faibussowitsch PetscCheck(Mx >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be greater than 1", Mx); 78147c6ae99SBarry Smith ratioi = (mx - 1) / (Mx - 1); 7821dca8a05SBarry Smith PetscCheck(ratioi * (Mx - 1) == mx - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 78347c6ae99SBarry Smith } 78447c6ae99SBarry Smith if (my == My) { 78547c6ae99SBarry Smith ratioj = 1; 786bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_PERIODIC) { 78763a3b9bcSJacob Faibussowitsch PetscCheck(My, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be positive", My); 78847c6ae99SBarry Smith ratioj = my / My; 78963a3b9bcSJacob Faibussowitsch PetscCheck(ratioj * My == my, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: my/My must be integer: my %" PetscInt_FMT " My %" PetscInt_FMT, my, My); 79047c6ae99SBarry Smith } else { 79163a3b9bcSJacob Faibussowitsch PetscCheck(My >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be greater than 1", My); 79247c6ae99SBarry Smith ratioj = (my - 1) / (My - 1); 7931dca8a05SBarry Smith PetscCheck(ratioj * (My - 1) == my - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (my - 1)/(My - 1) must be integer: my %" PetscInt_FMT " My %" PetscInt_FMT, my, My); 79447c6ae99SBarry Smith } 79547c6ae99SBarry Smith if (mz == Mz) { 79647c6ae99SBarry Smith ratiok = 1; 797bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_PERIODIC) { 79863a3b9bcSJacob Faibussowitsch PetscCheck(Mz, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of z coarse grid points %" PetscInt_FMT " must be positive", Mz); 79947c6ae99SBarry Smith ratiok = mz / Mz; 80063a3b9bcSJacob Faibussowitsch PetscCheck(ratiok * Mz == mz, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: mz/Mz must be integer: mz %" PetscInt_FMT " Mz %" PetscInt_FMT, mz, Mz); 80147c6ae99SBarry Smith } else { 80263a3b9bcSJacob Faibussowitsch PetscCheck(Mz >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of z coarse grid points %" PetscInt_FMT " must be greater than 1", Mz); 80347c6ae99SBarry Smith ratiok = (mz - 1) / (Mz - 1); 8041dca8a05SBarry Smith PetscCheck(ratiok * (Mz - 1) == mz - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (mz - 1)/(Mz - 1) must be integer: mz %" PetscInt_FMT " Mz %" PetscInt_FMT, mz, Mz); 80547c6ae99SBarry Smith } 80647c6ae99SBarry Smith 8079566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, &l_start, &m_f, &n_f, &p_f)); 8089566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, &l_start_ghost, &m_ghost, &n_ghost, &p_ghost)); 8099566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f)); 8109566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f)); 81147c6ae99SBarry Smith 8129566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, &l_start_c, &m_c, &n_c, &p_c)); 8139566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, &l_start_ghost_c, &m_ghost_c, &n_ghost_c, &p_ghost_c)); 8149566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c)); 8159566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c)); 81647c6ae99SBarry Smith 81747c6ae99SBarry Smith /* create interpolation matrix, determining exact preallocation */ 818d0609cedSBarry Smith MatPreallocateBegin(PetscObjectComm((PetscObject)dac), m_f * n_f * p_f, m_c * n_c * p_c, dnz, onz); 81947c6ae99SBarry Smith /* loop over local fine grid nodes counting interpolating points */ 82047c6ae99SBarry Smith for (l = l_start; l < l_start + p_f; l++) { 82147c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) { 82247c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) { 82347c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 824e3fbd1f4SBarry Smith row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))]; 82547c6ae99SBarry Smith i_c = (i / ratioi); 82647c6ae99SBarry Smith j_c = (j / ratioj); 82747c6ae99SBarry Smith l_c = (l / ratiok); 82808401ef6SPierre Jolivet PetscCheck(l_c >= l_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 8299371c9d4SSatish Balay l_start %" PetscInt_FMT " l_c %" PetscInt_FMT " l_start_ghost_c %" PetscInt_FMT, 8309371c9d4SSatish Balay l_start, l_c, l_start_ghost_c); 83108401ef6SPierre Jolivet PetscCheck(j_c >= j_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 8329371c9d4SSatish Balay j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT, 8339371c9d4SSatish Balay j_start, j_c, j_start_ghost_c); 83408401ef6SPierre Jolivet PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 8359371c9d4SSatish Balay i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT, 8369371c9d4SSatish Balay i_start, i_c, i_start_ghost_c); 83747c6ae99SBarry Smith 83847c6ae99SBarry Smith /* 83947c6ae99SBarry Smith Only include those interpolation points that are truly 84047c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 84147c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 84247c6ae99SBarry Smith */ 84347c6ae99SBarry Smith nc = 0; 844e3fbd1f4SBarry Smith col = (m_ghost_c * n_ghost_c * (l_c - l_start_ghost_c) + m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c)); 845e3fbd1f4SBarry Smith cols[nc++] = idx_c[col]; 846ad540459SPierre Jolivet if (i_c * ratioi != i) cols[nc++] = idx_c[col + 1]; 847ad540459SPierre Jolivet if (j_c * ratioj != j) cols[nc++] = idx_c[col + m_ghost_c]; 848ad540459SPierre Jolivet if (l_c * ratiok != l) cols[nc++] = idx_c[col + m_ghost_c * n_ghost_c]; 849ad540459SPierre Jolivet if (j_c * ratioj != j && i_c * ratioi != i) cols[nc++] = idx_c[col + (m_ghost_c + 1)]; 850ad540459SPierre Jolivet if (j_c * ratioj != j && l_c * ratiok != l) cols[nc++] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c)]; 851ad540459SPierre Jolivet if (i_c * ratioi != i && l_c * ratiok != l) cols[nc++] = idx_c[col + (m_ghost_c * n_ghost_c + 1)]; 852ad540459SPierre Jolivet if (i_c * ratioi != i && l_c * ratiok != l && j_c * ratioj != j) cols[nc++] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c + 1)]; 8539566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, nc, cols, dnz, onz)); 85447c6ae99SBarry Smith } 85547c6ae99SBarry Smith } 85647c6ae99SBarry Smith } 8579566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dac), &mat)); 858fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA) 859fdc842d1SBarry Smith /* 860fdc842d1SBarry Smith Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU 861fdc842d1SBarry Smith we don't want the original unconverted matrix copied to the GPU 862fdc842d1SBarry Smith */ 86348a46eb9SPierre Jolivet if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE)); 864fdc842d1SBarry Smith #endif 8659566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, m_f * n_f * p_f, m_c * n_c * p_c, mx * my * mz, Mx * My * Mz)); 8669566063dSJacob Faibussowitsch PetscCall(ConvertToAIJ(dac->mattype, &mattype)); 8679566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, mattype)); 8689566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 0, dnz)); 8699566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 0, dnz, 0, onz)); 870d0609cedSBarry Smith MatPreallocateEnd(dnz, onz); 87147c6ae99SBarry Smith 87247c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 8732adb9dcfSBarry Smith if (!NEWVERSION) { 87447c6ae99SBarry Smith for (l = l_start; l < l_start + p_f; l++) { 87547c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) { 87647c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) { 87747c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 878e3fbd1f4SBarry Smith row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))]; 87947c6ae99SBarry Smith 88047c6ae99SBarry Smith i_c = (i / ratioi); 88147c6ae99SBarry Smith j_c = (j / ratioj); 88247c6ae99SBarry Smith l_c = (l / ratiok); 88347c6ae99SBarry Smith 88447c6ae99SBarry Smith /* 88547c6ae99SBarry Smith Only include those interpolation points that are truly 88647c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 88747c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 88847c6ae99SBarry Smith */ 8896712e2f1SBarry Smith x = ((PetscReal)(i - i_c * ratioi)) / ((PetscReal)ratioi); 8906712e2f1SBarry Smith y = ((PetscReal)(j - j_c * ratioj)) / ((PetscReal)ratioj); 8916712e2f1SBarry Smith z = ((PetscReal)(l - l_c * ratiok)) / ((PetscReal)ratiok); 892b3bd94feSDave May 89347c6ae99SBarry Smith nc = 0; 89447c6ae99SBarry Smith /* one left and below; or we are right on it */ 895e3fbd1f4SBarry Smith col = (m_ghost_c * n_ghost_c * (l_c - l_start_ghost_c) + m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c)); 89647c6ae99SBarry Smith 897e3fbd1f4SBarry Smith cols[nc] = idx_c[col]; 89847c6ae99SBarry Smith v[nc++] = .125 * (1. - (2.0 * x - 1.)) * (1. - (2.0 * y - 1.)) * (1. - (2.0 * z - 1.)); 89947c6ae99SBarry Smith 90047c6ae99SBarry Smith if (i_c * ratioi != i) { 901e3fbd1f4SBarry Smith cols[nc] = idx_c[col + 1]; 90247c6ae99SBarry Smith v[nc++] = .125 * (1. + (2.0 * x - 1.)) * (1. - (2.0 * y - 1.)) * (1. - (2.0 * z - 1.)); 90347c6ae99SBarry Smith } 90447c6ae99SBarry Smith 90547c6ae99SBarry Smith if (j_c * ratioj != j) { 906e3fbd1f4SBarry Smith cols[nc] = idx_c[col + m_ghost_c]; 90747c6ae99SBarry Smith v[nc++] = .125 * (1. - (2.0 * x - 1.)) * (1. + (2.0 * y - 1.)) * (1. - (2.0 * z - 1.)); 90847c6ae99SBarry Smith } 90947c6ae99SBarry Smith 91047c6ae99SBarry Smith if (l_c * ratiok != l) { 911e3fbd1f4SBarry Smith cols[nc] = idx_c[col + m_ghost_c * n_ghost_c]; 91247c6ae99SBarry Smith v[nc++] = .125 * (1. - (2.0 * x - 1.)) * (1. - (2.0 * y - 1.)) * (1. + (2.0 * z - 1.)); 91347c6ae99SBarry Smith } 91447c6ae99SBarry Smith 91547c6ae99SBarry Smith if (j_c * ratioj != j && i_c * ratioi != i) { 916e3fbd1f4SBarry Smith cols[nc] = idx_c[col + (m_ghost_c + 1)]; 91747c6ae99SBarry Smith v[nc++] = .125 * (1. + (2.0 * x - 1.)) * (1. + (2.0 * y - 1.)) * (1. - (2.0 * z - 1.)); 91847c6ae99SBarry Smith } 91947c6ae99SBarry Smith 92047c6ae99SBarry Smith if (j_c * ratioj != j && l_c * ratiok != l) { 921e3fbd1f4SBarry Smith cols[nc] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c)]; 92247c6ae99SBarry Smith v[nc++] = .125 * (1. - (2.0 * x - 1.)) * (1. + (2.0 * y - 1.)) * (1. + (2.0 * z - 1.)); 92347c6ae99SBarry Smith } 92447c6ae99SBarry Smith 92547c6ae99SBarry Smith if (i_c * ratioi != i && l_c * ratiok != l) { 926e3fbd1f4SBarry Smith cols[nc] = idx_c[col + (m_ghost_c * n_ghost_c + 1)]; 92747c6ae99SBarry Smith v[nc++] = .125 * (1. + (2.0 * x - 1.)) * (1. - (2.0 * y - 1.)) * (1. + (2.0 * z - 1.)); 92847c6ae99SBarry Smith } 92947c6ae99SBarry Smith 93047c6ae99SBarry Smith if (i_c * ratioi != i && l_c * ratiok != l && j_c * ratioj != j) { 931e3fbd1f4SBarry Smith cols[nc] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c + 1)]; 93247c6ae99SBarry Smith v[nc++] = .125 * (1. + (2.0 * x - 1.)) * (1. + (2.0 * y - 1.)) * (1. + (2.0 * z - 1.)); 93347c6ae99SBarry Smith } 9349566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES)); 93547c6ae99SBarry Smith } 93647c6ae99SBarry Smith } 93747c6ae99SBarry Smith } 938b3bd94feSDave May 939b3bd94feSDave May } else { 940b3bd94feSDave May PetscScalar *xi, *eta, *zeta; 941b3bd94feSDave May PetscInt li, nxi, lj, neta, lk, nzeta, n; 942b3bd94feSDave May PetscScalar Ni[8]; 943b3bd94feSDave May 944b3bd94feSDave May /* compute local coordinate arrays */ 945b3bd94feSDave May nxi = ratioi + 1; 946b3bd94feSDave May neta = ratioj + 1; 947b3bd94feSDave May nzeta = ratiok + 1; 9489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nxi, &xi)); 9499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(neta, &eta)); 9509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeta, &zeta)); 9518865f1eaSKarl Rupp for (li = 0; li < nxi; li++) xi[li] = -1.0 + (PetscScalar)li * (2.0 / (PetscScalar)(nxi - 1)); 9528865f1eaSKarl Rupp for (lj = 0; lj < neta; lj++) eta[lj] = -1.0 + (PetscScalar)lj * (2.0 / (PetscScalar)(neta - 1)); 9538865f1eaSKarl Rupp for (lk = 0; lk < nzeta; lk++) zeta[lk] = -1.0 + (PetscScalar)lk * (2.0 / (PetscScalar)(nzeta - 1)); 954b3bd94feSDave May 955b3bd94feSDave May for (l = l_start; l < l_start + p_f; l++) { 956b3bd94feSDave May for (j = j_start; j < j_start + n_f; j++) { 957b3bd94feSDave May for (i = i_start; i < i_start + m_f; i++) { 958b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */ 959e3fbd1f4SBarry Smith row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))]; 960b3bd94feSDave May 961b3bd94feSDave May i_c = (i / ratioi); 962b3bd94feSDave May j_c = (j / ratioj); 963b3bd94feSDave May l_c = (l / ratiok); 964b3bd94feSDave May 965b3bd94feSDave May /* remainders */ 966b3bd94feSDave May li = i - ratioi * (i / ratioi); 9678865f1eaSKarl Rupp if (i == mx - 1) li = nxi - 1; 968b3bd94feSDave May lj = j - ratioj * (j / ratioj); 9698865f1eaSKarl Rupp if (j == my - 1) lj = neta - 1; 970b3bd94feSDave May lk = l - ratiok * (l / ratiok); 9718865f1eaSKarl Rupp if (l == mz - 1) lk = nzeta - 1; 972b3bd94feSDave May 973b3bd94feSDave May /* corners */ 974e3fbd1f4SBarry Smith col = (m_ghost_c * n_ghost_c * (l_c - l_start_ghost_c) + m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c)); 975e3fbd1f4SBarry Smith cols[0] = idx_c[col]; 976b3bd94feSDave May Ni[0] = 1.0; 977b3bd94feSDave May if ((li == 0) || (li == nxi - 1)) { 978b3bd94feSDave May if ((lj == 0) || (lj == neta - 1)) { 979b3bd94feSDave May if ((lk == 0) || (lk == nzeta - 1)) { 9809566063dSJacob Faibussowitsch PetscCall(MatSetValue(mat, row, cols[0], Ni[0], INSERT_VALUES)); 981b3bd94feSDave May continue; 982b3bd94feSDave May } 983b3bd94feSDave May } 984b3bd94feSDave May } 985b3bd94feSDave May 986b3bd94feSDave May /* edges + interior */ 987b3bd94feSDave May /* remainders */ 9888865f1eaSKarl Rupp if (i == mx - 1) i_c--; 9898865f1eaSKarl Rupp if (j == my - 1) j_c--; 9908865f1eaSKarl Rupp if (l == mz - 1) l_c--; 991b3bd94feSDave May 992e3fbd1f4SBarry Smith col = (m_ghost_c * n_ghost_c * (l_c - l_start_ghost_c) + m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c)); 993e3fbd1f4SBarry Smith cols[0] = idx_c[col]; /* one left and below; or we are right on it */ 994e3fbd1f4SBarry Smith cols[1] = idx_c[col + 1]; /* one right and below */ 995e3fbd1f4SBarry Smith cols[2] = idx_c[col + m_ghost_c]; /* one left and above */ 996e3fbd1f4SBarry Smith cols[3] = idx_c[col + (m_ghost_c + 1)]; /* one right and above */ 997b3bd94feSDave May 998e3fbd1f4SBarry Smith cols[4] = idx_c[col + m_ghost_c * n_ghost_c]; /* one left and below and front; or we are right on it */ 999e3fbd1f4SBarry Smith cols[5] = idx_c[col + (m_ghost_c * n_ghost_c + 1)]; /* one right and below, and front */ 1000e3fbd1f4SBarry Smith cols[6] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c)]; /* one left and above and front*/ 1001e3fbd1f4SBarry Smith cols[7] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c + 1)]; /* one right and above and front */ 1002b3bd94feSDave May 1003b3bd94feSDave May Ni[0] = 0.125 * (1.0 - xi[li]) * (1.0 - eta[lj]) * (1.0 - zeta[lk]); 1004b3bd94feSDave May Ni[1] = 0.125 * (1.0 + xi[li]) * (1.0 - eta[lj]) * (1.0 - zeta[lk]); 1005b3bd94feSDave May Ni[2] = 0.125 * (1.0 - xi[li]) * (1.0 + eta[lj]) * (1.0 - zeta[lk]); 1006b3bd94feSDave May Ni[3] = 0.125 * (1.0 + xi[li]) * (1.0 + eta[lj]) * (1.0 - zeta[lk]); 1007b3bd94feSDave May 1008b3bd94feSDave May Ni[4] = 0.125 * (1.0 - xi[li]) * (1.0 - eta[lj]) * (1.0 + zeta[lk]); 1009b3bd94feSDave May Ni[5] = 0.125 * (1.0 + xi[li]) * (1.0 - eta[lj]) * (1.0 + zeta[lk]); 1010b3bd94feSDave May Ni[6] = 0.125 * (1.0 - xi[li]) * (1.0 + eta[lj]) * (1.0 + zeta[lk]); 1011b3bd94feSDave May Ni[7] = 0.125 * (1.0 + xi[li]) * (1.0 + eta[lj]) * (1.0 + zeta[lk]); 1012b3bd94feSDave May 1013b3bd94feSDave May for (n = 0; n < 8; n++) { 10148865f1eaSKarl Rupp if (PetscAbsScalar(Ni[n]) < 1.0e-32) cols[n] = -1; 1015b3bd94feSDave May } 10169566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, 8, cols, Ni, INSERT_VALUES)); 1017b3bd94feSDave May } 1018b3bd94feSDave May } 1019b3bd94feSDave May } 10209566063dSJacob Faibussowitsch PetscCall(PetscFree(xi)); 10219566063dSJacob Faibussowitsch PetscCall(PetscFree(eta)); 10229566063dSJacob Faibussowitsch PetscCall(PetscFree(zeta)); 1023b3bd94feSDave May } 10249566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f)); 10259566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c)); 10269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 10279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 102847c6ae99SBarry Smith 10299566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(mat, dof, A)); 10309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 10313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 103247c6ae99SBarry Smith } 103347c6ae99SBarry Smith 1034d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInterpolation_DA(DM dac, DM daf, Mat *A, Vec *scale) 1035d71ae5a4SJacob Faibussowitsch { 103647c6ae99SBarry Smith PetscInt dimc, Mc, Nc, Pc, mc, nc, pc, dofc, sc, dimf, Mf, Nf, Pf, mf, nf, pf, doff, sf; 1037bff4a2f0SMatthew G. Knepley DMBoundaryType bxc, byc, bzc, bxf, byf, bzf; 1038aa219208SBarry Smith DMDAStencilType stc, stf; 103947c6ae99SBarry Smith DM_DA *ddc = (DM_DA *)dac->data; 104047c6ae99SBarry Smith 104147c6ae99SBarry Smith PetscFunctionBegin; 104247c6ae99SBarry Smith PetscValidHeaderSpecific(dac, DM_CLASSID, 1); 104347c6ae99SBarry Smith PetscValidHeaderSpecific(daf, DM_CLASSID, 2); 10444f572ea9SToby Isaac PetscAssertPointer(A, 3); 10454f572ea9SToby Isaac if (scale) PetscAssertPointer(scale, 4); 104647c6ae99SBarry Smith 10479566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, &dimc, &Mc, &Nc, &Pc, &mc, &nc, &pc, &dofc, &sc, &bxc, &byc, &bzc, &stc)); 10489566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, &dimf, &Mf, &Nf, &Pf, &mf, &nf, &pf, &doff, &sf, &bxf, &byf, &bzf, &stf)); 104963a3b9bcSJacob Faibussowitsch PetscCheck(dimc == dimf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Dimensions of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dimc, dimf); 105063a3b9bcSJacob Faibussowitsch PetscCheck(dofc == doff, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "DOF of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dofc, doff); 105163a3b9bcSJacob Faibussowitsch PetscCheck(sc == sf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil width of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, sc, sf); 10521dca8a05SBarry Smith PetscCheck(bxc == bxf && byc == byf && bzc == bzf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Boundary type different in two DMDAs"); 105308401ef6SPierre Jolivet PetscCheck(stc == stf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil type different in two DMDAs"); 10541dca8a05SBarry Smith PetscCheck(Mc >= 2 || Mf <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in x direction"); 10551dca8a05SBarry Smith PetscCheck(dimc <= 1 || Nc >= 2 || Nf <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in y direction"); 10561dca8a05SBarry Smith PetscCheck(dimc <= 2 || Pc >= 2 || Pf <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in z direction"); 105747c6ae99SBarry Smith 1058aa219208SBarry Smith if (ddc->interptype == DMDA_Q1) { 105947c6ae99SBarry Smith if (dimc == 1) { 10609566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation_DA_1D_Q1(dac, daf, A)); 106147c6ae99SBarry Smith } else if (dimc == 2) { 10629566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation_DA_2D_Q1(dac, daf, A)); 106347c6ae99SBarry Smith } else if (dimc == 3) { 10649566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation_DA_3D_Q1(dac, daf, A)); 106563a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)daf), PETSC_ERR_SUP, "No support for this DMDA dimension %" PetscInt_FMT " for interpolation type %d", dimc, (int)ddc->interptype); 1066aa219208SBarry Smith } else if (ddc->interptype == DMDA_Q0) { 106747c6ae99SBarry Smith if (dimc == 1) { 10689566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation_DA_1D_Q0(dac, daf, A)); 106947c6ae99SBarry Smith } else if (dimc == 2) { 10709566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation_DA_2D_Q0(dac, daf, A)); 107147c6ae99SBarry Smith } else if (dimc == 3) { 10729566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation_DA_3D_Q0(dac, daf, A)); 107363a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)daf), PETSC_ERR_SUP, "No support for this DMDA dimension %" PetscInt_FMT " for interpolation type %d", dimc, (int)ddc->interptype); 107447c6ae99SBarry Smith } 10751baa6e33SBarry Smith if (scale) PetscCall(DMCreateInterpolationScale((DM)dac, (DM)daf, *A, scale)); 10763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 107747c6ae99SBarry Smith } 107847c6ae99SBarry Smith 1079d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInjection_DA_1D(DM dac, DM daf, VecScatter *inject) 1080d71ae5a4SJacob Faibussowitsch { 10818ea3bf28SBarry Smith PetscInt i, i_start, m_f, Mx, dof; 10828ea3bf28SBarry Smith const PetscInt *idx_f; 1083e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f; 10848ea3bf28SBarry Smith PetscInt m_ghost, m_ghost_c; 10850bb2b966SJungho Lee PetscInt row, i_start_ghost, mx, m_c, nc, ratioi; 10860bb2b966SJungho Lee PetscInt i_start_c, i_start_ghost_c; 10870bb2b966SJungho Lee PetscInt *cols; 1088bff4a2f0SMatthew G. Knepley DMBoundaryType bx; 10890bb2b966SJungho Lee Vec vecf, vecc; 10900bb2b966SJungho Lee IS isf; 10910bb2b966SJungho Lee 10920bb2b966SJungho Lee PetscFunctionBegin; 10939566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bx, NULL, NULL, NULL)); 10949566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 1095bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) { 10960bb2b966SJungho Lee ratioi = mx / Mx; 109763a3b9bcSJacob Faibussowitsch PetscCheck(ratioi * Mx == mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: mx/Mx must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 10980bb2b966SJungho Lee } else { 10990bb2b966SJungho Lee ratioi = (mx - 1) / (Mx - 1); 11001dca8a05SBarry Smith PetscCheck(ratioi * (Mx - 1) == mx - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 11010bb2b966SJungho Lee } 11020bb2b966SJungho Lee 11039566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, NULL, NULL, &m_f, NULL, NULL)); 11049566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, NULL, NULL, &m_ghost, NULL, NULL)); 11059566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f)); 11069566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f)); 11070bb2b966SJungho Lee 11089566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, NULL, NULL, &m_c, NULL, NULL)); 11099566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, NULL, NULL, &m_ghost_c, NULL, NULL)); 11100bb2b966SJungho Lee 11110bb2b966SJungho Lee /* loop over local fine grid nodes setting interpolation for those*/ 11120bb2b966SJungho Lee nc = 0; 11139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m_f, &cols)); 11140bb2b966SJungho Lee 11150bb2b966SJungho Lee for (i = i_start_c; i < i_start_c + m_c; i++) { 11160bb2b966SJungho Lee PetscInt i_f = i * ratioi; 11170bb2b966SJungho Lee 11181dca8a05SBarry Smith PetscCheck(i_f >= i_start_ghost && i_f < i_start_ghost + m_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\ni_c %" PetscInt_FMT " i_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]", i, i_f, i_start_ghost, i_start_ghost + m_ghost); 11198865f1eaSKarl Rupp 1120e3fbd1f4SBarry Smith row = idx_f[(i_f - i_start_ghost)]; 1121e3fbd1f4SBarry Smith cols[nc++] = row; 11220bb2b966SJungho Lee } 11230bb2b966SJungho Lee 11249566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f)); 11259566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)daf), dof, nc, cols, PETSC_OWN_POINTER, &isf)); 11269566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dac, &vecc)); 11279566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(daf, &vecf)); 11289566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vecf, isf, vecc, NULL, inject)); 11299566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dac, &vecc)); 11309566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(daf, &vecf)); 11319566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isf)); 11323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11330bb2b966SJungho Lee } 11340bb2b966SJungho Lee 1135d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInjection_DA_2D(DM dac, DM daf, VecScatter *inject) 1136d71ae5a4SJacob Faibussowitsch { 11378ea3bf28SBarry Smith PetscInt i, j, i_start, j_start, m_f, n_f, Mx, My, dof; 11388ea3bf28SBarry Smith const PetscInt *idx_c, *idx_f; 1139e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c; 11408ea3bf28SBarry Smith PetscInt m_ghost, n_ghost, m_ghost_c, n_ghost_c; 114147c6ae99SBarry Smith PetscInt row, i_start_ghost, j_start_ghost, mx, m_c, my, nc, ratioi, ratioj; 1142076202ddSJed Brown PetscInt i_start_c, j_start_c, n_c, i_start_ghost_c, j_start_ghost_c; 114347c6ae99SBarry Smith PetscInt *cols; 1144bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by; 114547c6ae99SBarry Smith Vec vecf, vecc; 114647c6ae99SBarry Smith IS isf; 114747c6ae99SBarry Smith 114847c6ae99SBarry Smith PetscFunctionBegin; 11499566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, NULL, NULL, NULL, NULL, NULL, NULL, &bx, &by, NULL, NULL)); 11509566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 1151bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) { 115247c6ae99SBarry Smith ratioi = mx / Mx; 115363a3b9bcSJacob Faibussowitsch PetscCheck(ratioi * Mx == mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: mx/Mx must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 115447c6ae99SBarry Smith } else { 115547c6ae99SBarry Smith ratioi = (mx - 1) / (Mx - 1); 11561dca8a05SBarry Smith PetscCheck(ratioi * (Mx - 1) == mx - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 115747c6ae99SBarry Smith } 1158bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) { 115947c6ae99SBarry Smith ratioj = my / My; 116063a3b9bcSJacob Faibussowitsch PetscCheck(ratioj * My == my, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: my/My must be integer: my %" PetscInt_FMT " My %" PetscInt_FMT, my, My); 116147c6ae99SBarry Smith } else { 116247c6ae99SBarry Smith ratioj = (my - 1) / (My - 1); 11631dca8a05SBarry Smith PetscCheck(ratioj * (My - 1) == my - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (my - 1)/(My - 1) must be integer: my %" PetscInt_FMT " My %" PetscInt_FMT, my, My); 116447c6ae99SBarry Smith } 116547c6ae99SBarry Smith 11669566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, NULL, &m_f, &n_f, NULL)); 11679566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, NULL, &m_ghost, &n_ghost, NULL)); 11689566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f)); 11699566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f)); 117047c6ae99SBarry Smith 11719566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, NULL, &m_c, &n_c, NULL)); 11729566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, NULL, &m_ghost_c, &n_ghost_c, NULL)); 11739566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c)); 11749566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c)); 117547c6ae99SBarry Smith 117647c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 117747c6ae99SBarry Smith nc = 0; 11789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_f * m_f, &cols)); 1179076202ddSJed Brown for (j = j_start_c; j < j_start_c + n_c; j++) { 1180076202ddSJed Brown for (i = i_start_c; i < i_start_c + m_c; i++) { 1181076202ddSJed Brown PetscInt i_f = i * ratioi, j_f = j * ratioj; 11821dca8a05SBarry Smith PetscCheck(j_f >= j_start_ghost && j_f < j_start_ghost + n_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 11839371c9d4SSatish Balay j_c %" PetscInt_FMT " j_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]", 11849371c9d4SSatish Balay j, j_f, j_start_ghost, j_start_ghost + n_ghost); 11851dca8a05SBarry Smith PetscCheck(i_f >= i_start_ghost && i_f < i_start_ghost + m_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\ 11869371c9d4SSatish Balay i_c %" PetscInt_FMT " i_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]", 11879371c9d4SSatish Balay i, i_f, i_start_ghost, i_start_ghost + m_ghost); 1188e3fbd1f4SBarry Smith row = idx_f[(m_ghost * (j_f - j_start_ghost) + (i_f - i_start_ghost))]; 1189e3fbd1f4SBarry Smith cols[nc++] = row; 119047c6ae99SBarry Smith } 119147c6ae99SBarry Smith } 11929566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f)); 11939566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c)); 119447c6ae99SBarry Smith 11959566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)daf), dof, nc, cols, PETSC_OWN_POINTER, &isf)); 11969566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dac, &vecc)); 11979566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(daf, &vecf)); 11989566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vecf, isf, vecc, NULL, inject)); 11999566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dac, &vecc)); 12009566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(daf, &vecf)); 12019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isf)); 12023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 120347c6ae99SBarry Smith } 120447c6ae99SBarry Smith 1205d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInjection_DA_3D(DM dac, DM daf, VecScatter *inject) 1206d71ae5a4SJacob Faibussowitsch { 1207b1757049SJed Brown PetscInt i, j, k, i_start, j_start, k_start, m_f, n_f, p_f, Mx, My, Mz; 1208b1757049SJed Brown PetscInt m_ghost, n_ghost, p_ghost, m_ghost_c, n_ghost_c, p_ghost_c; 1209b1757049SJed Brown PetscInt i_start_ghost, j_start_ghost, k_start_ghost; 1210b1757049SJed Brown PetscInt mx, my, mz, ratioi, ratioj, ratiok; 1211b1757049SJed Brown PetscInt i_start_c, j_start_c, k_start_c; 1212b1757049SJed Brown PetscInt m_c, n_c, p_c; 1213b1757049SJed Brown PetscInt i_start_ghost_c, j_start_ghost_c, k_start_ghost_c; 1214b1757049SJed Brown PetscInt row, nc, dof; 12158ea3bf28SBarry Smith const PetscInt *idx_c, *idx_f; 1216e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c; 1217b1757049SJed Brown PetscInt *cols; 1218bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by, bz; 1219b1757049SJed Brown Vec vecf, vecc; 1220b1757049SJed Brown IS isf; 1221b1757049SJed Brown 1222b1757049SJed Brown PetscFunctionBegin; 12239566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, &Mz, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL)); 12249566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, &mz, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 1225b1757049SJed Brown 1226bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) { 1227b1757049SJed Brown ratioi = mx / Mx; 122863a3b9bcSJacob Faibussowitsch PetscCheck(ratioi * Mx == mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: mx/Mx must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 1229b1757049SJed Brown } else { 1230b1757049SJed Brown ratioi = (mx - 1) / (Mx - 1); 12311dca8a05SBarry Smith PetscCheck(ratioi * (Mx - 1) == mx - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT, mx, Mx); 1232b1757049SJed Brown } 1233bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) { 1234b1757049SJed Brown ratioj = my / My; 123563a3b9bcSJacob Faibussowitsch PetscCheck(ratioj * My == my, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: my/My must be integer: my %" PetscInt_FMT " My %" PetscInt_FMT, my, My); 1236b1757049SJed Brown } else { 1237b1757049SJed Brown ratioj = (my - 1) / (My - 1); 12381dca8a05SBarry Smith PetscCheck(ratioj * (My - 1) == my - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (my - 1)/(My - 1) must be integer: my %" PetscInt_FMT " My %" PetscInt_FMT, my, My); 1239b1757049SJed Brown } 1240bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC) { 1241b1757049SJed Brown ratiok = mz / Mz; 124263a3b9bcSJacob Faibussowitsch PetscCheck(ratiok * Mz == mz, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: mz/Mz must be integer: mz %" PetscInt_FMT " My %" PetscInt_FMT, mz, Mz); 1243b1757049SJed Brown } else { 1244b1757049SJed Brown ratiok = (mz - 1) / (Mz - 1); 12451dca8a05SBarry Smith PetscCheck(ratiok * (Mz - 1) == mz - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Ratio between levels: (mz - 1)/(Mz - 1) must be integer: mz %" PetscInt_FMT " Mz %" PetscInt_FMT, mz, Mz); 1246b1757049SJed Brown } 1247b1757049SJed Brown 12489566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, &k_start, &m_f, &n_f, &p_f)); 12499566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, &k_start_ghost, &m_ghost, &n_ghost, &p_ghost)); 12509566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f)); 12519566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f)); 1252b1757049SJed Brown 12539566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, &k_start_c, &m_c, &n_c, &p_c)); 12549566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, &k_start_ghost_c, &m_ghost_c, &n_ghost_c, &p_ghost_c)); 12559566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c)); 12569566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c)); 1257b1757049SJed Brown 1258b1757049SJed Brown /* loop over local fine grid nodes setting interpolation for those*/ 1259b1757049SJed Brown nc = 0; 12609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_f * m_f * p_f, &cols)); 1261b1757049SJed Brown for (k = k_start_c; k < k_start_c + p_c; k++) { 1262b1757049SJed Brown for (j = j_start_c; j < j_start_c + n_c; j++) { 1263b1757049SJed Brown for (i = i_start_c; i < i_start_c + m_c; i++) { 1264b1757049SJed Brown PetscInt i_f = i * ratioi, j_f = j * ratioj, k_f = k * ratiok; 12659371c9d4SSatish Balay PetscCheck(k_f >= k_start_ghost && k_f < k_start_ghost + p_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 12669371c9d4SSatish Balay "Processor's coarse DMDA must lie over fine DMDA " 12679371c9d4SSatish Balay "k_c %" PetscInt_FMT " k_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]", 12689371c9d4SSatish Balay k, k_f, k_start_ghost, k_start_ghost + p_ghost); 12699371c9d4SSatish Balay PetscCheck(j_f >= j_start_ghost && j_f < j_start_ghost + n_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 12709371c9d4SSatish Balay "Processor's coarse DMDA must lie over fine DMDA " 12719371c9d4SSatish Balay "j_c %" PetscInt_FMT " j_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]", 12729371c9d4SSatish Balay j, j_f, j_start_ghost, j_start_ghost + n_ghost); 12739371c9d4SSatish Balay PetscCheck(i_f >= i_start_ghost && i_f < i_start_ghost + m_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 12749371c9d4SSatish Balay "Processor's coarse DMDA must lie over fine DMDA " 12759371c9d4SSatish Balay "i_c %" PetscInt_FMT " i_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]", 12769371c9d4SSatish Balay i, i_f, i_start_ghost, i_start_ghost + m_ghost); 1277e3fbd1f4SBarry Smith row = idx_f[(m_ghost * n_ghost * (k_f - k_start_ghost) + m_ghost * (j_f - j_start_ghost) + (i_f - i_start_ghost))]; 1278e3fbd1f4SBarry Smith cols[nc++] = row; 1279b1757049SJed Brown } 1280b1757049SJed Brown } 1281b1757049SJed Brown } 12829566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f)); 12839566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c)); 1284b1757049SJed Brown 12859566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)daf), dof, nc, cols, PETSC_OWN_POINTER, &isf)); 12869566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dac, &vecc)); 12879566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(daf, &vecf)); 12889566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vecf, isf, vecc, NULL, inject)); 12899566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dac, &vecc)); 12909566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(daf, &vecf)); 12919566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isf)); 12923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1293b1757049SJed Brown } 1294b1757049SJed Brown 1295d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInjection_DA(DM dac, DM daf, Mat *mat) 1296d71ae5a4SJacob Faibussowitsch { 129747c6ae99SBarry Smith PetscInt dimc, Mc, Nc, Pc, mc, nc, pc, dofc, sc, dimf, Mf, Nf, Pf, mf, nf, pf, doff, sf; 1298bff4a2f0SMatthew G. Knepley DMBoundaryType bxc, byc, bzc, bxf, byf, bzf; 1299aa219208SBarry Smith DMDAStencilType stc, stf; 130060c16ac1SBarry Smith VecScatter inject = NULL; 130147c6ae99SBarry Smith 130247c6ae99SBarry Smith PetscFunctionBegin; 130347c6ae99SBarry Smith PetscValidHeaderSpecific(dac, DM_CLASSID, 1); 130447c6ae99SBarry Smith PetscValidHeaderSpecific(daf, DM_CLASSID, 2); 13054f572ea9SToby Isaac PetscAssertPointer(mat, 3); 130647c6ae99SBarry Smith 13079566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, &dimc, &Mc, &Nc, &Pc, &mc, &nc, &pc, &dofc, &sc, &bxc, &byc, &bzc, &stc)); 13089566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, &dimf, &Mf, &Nf, &Pf, &mf, &nf, &pf, &doff, &sf, &bxf, &byf, &bzf, &stf)); 130963a3b9bcSJacob Faibussowitsch PetscCheck(dimc == dimf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Dimensions of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dimc, dimf); 131063a3b9bcSJacob Faibussowitsch PetscCheck(dofc == doff, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "DOF of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dofc, doff); 131163a3b9bcSJacob Faibussowitsch PetscCheck(sc == sf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil width of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, sc, sf); 13121dca8a05SBarry Smith PetscCheck(bxc == bxf && byc == byf && bzc == bzf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Boundary type different in two DMDAs"); 131308401ef6SPierre Jolivet PetscCheck(stc == stf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil type different in two DMDAs"); 131408401ef6SPierre Jolivet PetscCheck(Mc >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in x direction"); 13151dca8a05SBarry Smith PetscCheck(dimc <= 1 || Nc >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in y direction"); 13161dca8a05SBarry Smith PetscCheck(dimc <= 2 || Pc >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in z direction"); 131747c6ae99SBarry Smith 13180bb2b966SJungho Lee if (dimc == 1) { 13199566063dSJacob Faibussowitsch PetscCall(DMCreateInjection_DA_1D(dac, daf, &inject)); 13200bb2b966SJungho Lee } else if (dimc == 2) { 13219566063dSJacob Faibussowitsch PetscCall(DMCreateInjection_DA_2D(dac, daf, &inject)); 1322b1757049SJed Brown } else if (dimc == 3) { 13239566063dSJacob Faibussowitsch PetscCall(DMCreateInjection_DA_3D(dac, daf, &inject)); 132447c6ae99SBarry Smith } 13259566063dSJacob Faibussowitsch PetscCall(MatCreateScatter(PetscObjectComm((PetscObject)inject), inject, mat)); 13269566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&inject)); 13273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 132847c6ae99SBarry Smith } 132947c6ae99SBarry Smith 1330*a4e35b19SJacob Faibussowitsch // PetscClangLinter pragma ignore: -fdoc-* 133197779f9aSLisandro Dalcin /*@ 133297779f9aSLisandro Dalcin DMCreateAggregates - Deprecated, see DMDACreateAggregates. 13336c877ef6SSatish Balay 13346c877ef6SSatish Balay Level: intermediate 133597779f9aSLisandro Dalcin @*/ 1336d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateAggregates(DM dac, DM daf, Mat *mat) 1337d71ae5a4SJacob Faibussowitsch { 1338a5bc1bf3SBarry Smith return DMDACreateAggregates(dac, daf, mat); 133997779f9aSLisandro Dalcin } 134097779f9aSLisandro Dalcin 134197779f9aSLisandro Dalcin /*@ 134297779f9aSLisandro Dalcin DMDACreateAggregates - Gets the aggregates that map between 1343dce8aebaSBarry Smith grids associated with two `DMDA` 134497779f9aSLisandro Dalcin 134520f4b53cSBarry Smith Collective 134697779f9aSLisandro Dalcin 134797779f9aSLisandro Dalcin Input Parameters: 134860225df5SJacob Faibussowitsch + dac - the coarse grid `DMDA` 134960225df5SJacob Faibussowitsch - daf - the fine grid `DMDA` 135097779f9aSLisandro Dalcin 13512fe279fdSBarry Smith Output Parameter: 135297779f9aSLisandro Dalcin . rest - the restriction matrix (transpose of the projection matrix) 135397779f9aSLisandro Dalcin 135497779f9aSLisandro Dalcin Level: intermediate 135597779f9aSLisandro Dalcin 1356dce8aebaSBarry Smith Note: 1357dce8aebaSBarry Smith This routine is not used by PETSc. 135897779f9aSLisandro Dalcin It is not clear what its use case is and it may be removed in a future release. 135997779f9aSLisandro Dalcin Users should contact petsc-maint@mcs.anl.gov if they plan to use it. 136097779f9aSLisandro Dalcin 1361db781477SPatrick Sanan .seealso: `DMRefine()`, `DMCreateInjection()`, `DMCreateInterpolation()` 136297779f9aSLisandro Dalcin @*/ 1363d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDACreateAggregates(DM dac, DM daf, Mat *rest) 1364d71ae5a4SJacob Faibussowitsch { 136547c6ae99SBarry Smith PetscInt dimc, Mc, Nc, Pc, mc, nc, pc, dofc, sc; 136647c6ae99SBarry Smith PetscInt dimf, Mf, Nf, Pf, mf, nf, pf, doff, sf; 1367bff4a2f0SMatthew G. Knepley DMBoundaryType bxc, byc, bzc, bxf, byf, bzf; 1368aa219208SBarry Smith DMDAStencilType stc, stf; 136947c6ae99SBarry Smith PetscInt i, j, l; 137047c6ae99SBarry Smith PetscInt i_start, j_start, l_start, m_f, n_f, p_f; 137147c6ae99SBarry Smith PetscInt i_start_ghost, j_start_ghost, l_start_ghost, m_ghost, n_ghost, p_ghost; 13728ea3bf28SBarry Smith const PetscInt *idx_f; 137347c6ae99SBarry Smith PetscInt i_c, j_c, l_c; 137447c6ae99SBarry Smith PetscInt i_start_c, j_start_c, l_start_c, m_c, n_c, p_c; 137547c6ae99SBarry Smith PetscInt i_start_ghost_c, j_start_ghost_c, l_start_ghost_c, m_ghost_c, n_ghost_c, p_ghost_c; 13768ea3bf28SBarry Smith const PetscInt *idx_c; 137747c6ae99SBarry Smith PetscInt d; 137847c6ae99SBarry Smith PetscInt a; 137947c6ae99SBarry Smith PetscInt max_agg_size; 138047c6ae99SBarry Smith PetscInt *fine_nodes; 138147c6ae99SBarry Smith PetscScalar *one_vec; 138247c6ae99SBarry Smith PetscInt fn_idx; 1383565245c5SBarry Smith ISLocalToGlobalMapping ltogmf, ltogmc; 138447c6ae99SBarry Smith 138547c6ae99SBarry Smith PetscFunctionBegin; 138697779f9aSLisandro Dalcin PetscValidHeaderSpecificType(dac, DM_CLASSID, 1, DMDA); 138797779f9aSLisandro Dalcin PetscValidHeaderSpecificType(daf, DM_CLASSID, 2, DMDA); 13884f572ea9SToby Isaac PetscAssertPointer(rest, 3); 138947c6ae99SBarry Smith 13909566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, &dimc, &Mc, &Nc, &Pc, &mc, &nc, &pc, &dofc, &sc, &bxc, &byc, &bzc, &stc)); 13919566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, &dimf, &Mf, &Nf, &Pf, &mf, &nf, &pf, &doff, &sf, &bxf, &byf, &bzf, &stf)); 139263a3b9bcSJacob Faibussowitsch PetscCheck(dimc == dimf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Dimensions of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dimc, dimf); 139363a3b9bcSJacob Faibussowitsch PetscCheck(dofc == doff, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "DOF of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dofc, doff); 139463a3b9bcSJacob Faibussowitsch PetscCheck(sc == sf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil width of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, sc, sf); 13951dca8a05SBarry Smith PetscCheck(bxc == bxf && byc == byf && bzc == bzf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Boundary type different in two DMDAs"); 139608401ef6SPierre Jolivet PetscCheck(stc == stf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil type different in two DMDAs"); 139747c6ae99SBarry Smith 139863a3b9bcSJacob Faibussowitsch PetscCheck(Mf >= Mc, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Coarse grid has more points than fine grid, Mc %" PetscInt_FMT ", Mf %" PetscInt_FMT, Mc, Mf); 139963a3b9bcSJacob Faibussowitsch PetscCheck(Nf >= Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Coarse grid has more points than fine grid, Nc %" PetscInt_FMT ", Nf %" PetscInt_FMT, Nc, Nf); 140063a3b9bcSJacob Faibussowitsch PetscCheck(Pf >= Pc, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Coarse grid has more points than fine grid, Pc %" PetscInt_FMT ", Pf %" PetscInt_FMT, Pc, Pf); 140147c6ae99SBarry Smith 140247c6ae99SBarry Smith if (Pc < 0) Pc = 1; 140347c6ae99SBarry Smith if (Pf < 0) Pf = 1; 140447c6ae99SBarry Smith if (Nc < 0) Nc = 1; 140547c6ae99SBarry Smith if (Nf < 0) Nf = 1; 140647c6ae99SBarry Smith 14079566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, &l_start, &m_f, &n_f, &p_f)); 14089566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, &l_start_ghost, &m_ghost, &n_ghost, &p_ghost)); 1409565245c5SBarry Smith 14109566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <ogmf)); 14119566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltogmf, &idx_f)); 141247c6ae99SBarry Smith 14139566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, &l_start_c, &m_c, &n_c, &p_c)); 14149566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, &l_start_ghost_c, &m_ghost_c, &n_ghost_c, &p_ghost_c)); 1415565245c5SBarry Smith 14169566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <ogmc)); 14179566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltogmc, &idx_c)); 141847c6ae99SBarry Smith 141947c6ae99SBarry Smith /* 142047c6ae99SBarry Smith Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios 142147c6ae99SBarry Smith for dimension 1 and 2 respectively. 142247c6ae99SBarry Smith Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1) 1423da81f932SPierre Jolivet and r_y*j and r_y*(j+1) will be grouped into the same coarse grid aggregate. 142447c6ae99SBarry Smith Each specific dof on the fine grid is mapped to one dof on the coarse grid. 142547c6ae99SBarry Smith */ 142647c6ae99SBarry Smith 142747c6ae99SBarry Smith max_agg_size = (Mf / Mc + 1) * (Nf / Nc + 1) * (Pf / Pc + 1); 142847c6ae99SBarry Smith 142947c6ae99SBarry Smith /* create the matrix that will contain the restriction operator */ 14309371c9d4SSatish Balay PetscCall(MatCreateAIJ(PetscObjectComm((PetscObject)daf), m_c * n_c * p_c * dofc, m_f * n_f * p_f * doff, Mc * Nc * Pc * dofc, Mf * Nf * Pf * doff, max_agg_size, NULL, max_agg_size, NULL, rest)); 143147c6ae99SBarry Smith 143247c6ae99SBarry Smith /* store nodes in the fine grid here */ 14339566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(max_agg_size, &one_vec, max_agg_size, &fine_nodes)); 143447c6ae99SBarry Smith for (i = 0; i < max_agg_size; i++) one_vec[i] = 1.0; 143547c6ae99SBarry Smith 143647c6ae99SBarry Smith /* loop over all coarse nodes */ 143747c6ae99SBarry Smith for (l_c = l_start_c; l_c < l_start_c + p_c; l_c++) { 143847c6ae99SBarry Smith for (j_c = j_start_c; j_c < j_start_c + n_c; j_c++) { 143947c6ae99SBarry Smith for (i_c = i_start_c; i_c < i_start_c + m_c; i_c++) { 144047c6ae99SBarry Smith for (d = 0; d < dofc; d++) { 144147c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 144247c6ae99SBarry Smith a = idx_c[dofc * (m_ghost_c * n_ghost_c * (l_c - l_start_ghost_c) + m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c))] + d; 144347c6ae99SBarry Smith 144447c6ae99SBarry Smith fn_idx = 0; 144547c6ae99SBarry Smith /* Corresponding fine points are all points (i_f, j_f, l_f) such that 144647c6ae99SBarry Smith i_c*Mf/Mc <= i_f < (i_c+1)*Mf/Mc 144747c6ae99SBarry Smith (same for other dimensions) 144847c6ae99SBarry Smith */ 144947c6ae99SBarry Smith for (l = l_c * Pf / Pc; l < PetscMin((l_c + 1) * Pf / Pc, Pf); l++) { 145047c6ae99SBarry Smith for (j = j_c * Nf / Nc; j < PetscMin((j_c + 1) * Nf / Nc, Nf); j++) { 145147c6ae99SBarry Smith for (i = i_c * Mf / Mc; i < PetscMin((i_c + 1) * Mf / Mc, Mf); i++) { 145247c6ae99SBarry Smith fine_nodes[fn_idx] = idx_f[doff * (m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))] + d; 145347c6ae99SBarry Smith fn_idx++; 145447c6ae99SBarry Smith } 145547c6ae99SBarry Smith } 145647c6ae99SBarry Smith } 145747c6ae99SBarry Smith /* add all these points to one aggregate */ 14589566063dSJacob Faibussowitsch PetscCall(MatSetValues(*rest, 1, &a, fn_idx, fine_nodes, one_vec, INSERT_VALUES)); 145947c6ae99SBarry Smith } 146047c6ae99SBarry Smith } 146147c6ae99SBarry Smith } 146247c6ae99SBarry Smith } 14639566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogmf, &idx_f)); 14649566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogmc, &idx_c)); 14659566063dSJacob Faibussowitsch PetscCall(PetscFree2(one_vec, fine_nodes)); 14669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*rest, MAT_FINAL_ASSEMBLY)); 14679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*rest, MAT_FINAL_ASSEMBLY)); 14683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 146947c6ae99SBarry Smith } 1470