147c6ae99SBarry Smith /*
2aa219208SBarry Smith Code for interpolating between grids represented by DMDAs
347c6ae99SBarry Smith */
447c6ae99SBarry Smith
5d52bd9f3SBarry Smith /*
6d52bd9f3SBarry 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
7d52bd9f3SBarry 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
82adb9dcfSBarry 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
92adb9dcfSBarry Smith consider it cleaner, but old version is turned on since it handles periodic case.
10d52bd9f3SBarry Smith */
119314d695SBarry Smith #define NEWVERSION 0
1285afcc9aSBarry Smith
13af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/
1447c6ae99SBarry Smith
15fdc842d1SBarry Smith /*
16fdc842d1SBarry Smith Since the interpolation uses MATMAIJ for dof > 0 we convert request for non-MATAIJ baseded matrices to MATAIJ.
17fdc842d1SBarry Smith This is a bit of a hack, the reason for it is partially because -dm_mat_type defines the
18fdc842d1SBarry Smith matrix type for both the operator matrices and the interpolation matrices so that users
19fdc842d1SBarry Smith can select matrix types of base MATAIJ for accelerators
20fdc842d1SBarry Smith */
ConvertToAIJ(MatType intype,MatType * outtype)21d71ae5a4SJacob Faibussowitsch static PetscErrorCode ConvertToAIJ(MatType intype, MatType *outtype)
22d71ae5a4SJacob Faibussowitsch {
23fdc842d1SBarry Smith PetscInt i;
24fdc842d1SBarry Smith char const *types[3] = {MATAIJ, MATSEQAIJ, MATMPIAIJ};
25fdc842d1SBarry Smith PetscBool flg;
26fdc842d1SBarry Smith
27fdc842d1SBarry Smith PetscFunctionBegin;
2833a4d587SStefano Zampini *outtype = MATAIJ;
29fdc842d1SBarry Smith for (i = 0; i < 3; i++) {
309566063dSJacob Faibussowitsch PetscCall(PetscStrbeginswith(intype, types[i], &flg));
31fdc842d1SBarry Smith if (flg) {
32fdc842d1SBarry Smith *outtype = intype;
333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
34fdc842d1SBarry Smith }
35fdc842d1SBarry Smith }
363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
37fdc842d1SBarry Smith }
38fdc842d1SBarry Smith
DMCreateInterpolation_DA_1D_Q1(DM dac,DM daf,Mat * A)3966976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_DA_1D_Q1(DM dac, DM daf, Mat *A)
40d71ae5a4SJacob Faibussowitsch {
418ea3bf28SBarry Smith PetscInt i, i_start, m_f, Mx;
428ea3bf28SBarry Smith const PetscInt *idx_f, *idx_c;
438ea3bf28SBarry Smith PetscInt m_ghost, m_ghost_c;
4447c6ae99SBarry Smith PetscInt row, col, i_start_ghost, mx, m_c, nc, ratio;
4547c6ae99SBarry Smith PetscInt i_c, i_start_c, i_start_ghost_c, cols[2], dof;
4685afcc9aSBarry Smith PetscScalar v[2], x;
4747c6ae99SBarry Smith Mat mat;
48bff4a2f0SMatthew G. Knepley DMBoundaryType bx;
49e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c;
50fdc842d1SBarry Smith MatType mattype;
5147c6ae99SBarry Smith
5247c6ae99SBarry Smith PetscFunctionBegin;
539566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bx, NULL, NULL, NULL));
549566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
55bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) {
5647c6ae99SBarry Smith ratio = mx / Mx;
5763a3b9bcSJacob 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);
5847c6ae99SBarry Smith } else {
5947c6ae99SBarry Smith ratio = (mx - 1) / (Mx - 1);
601dca8a05SBarry 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);
6147c6ae99SBarry Smith }
6247c6ae99SBarry Smith
639566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, NULL, NULL, &m_f, NULL, NULL));
649566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, NULL, NULL, &m_ghost, NULL, NULL));
659566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f));
669566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
6747c6ae99SBarry Smith
689566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, NULL, NULL, &m_c, NULL, NULL));
699566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, NULL, NULL, &m_ghost_c, NULL, NULL));
709566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c));
719566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
7247c6ae99SBarry Smith
7347c6ae99SBarry Smith /* create interpolation matrix */
749566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dac), &mat));
75fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
76fdc842d1SBarry Smith /*
77fdc842d1SBarry Smith Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
78fdc842d1SBarry Smith we don't want the original unconverted matrix copied to the GPU
79fdc842d1SBarry Smith */
8048a46eb9SPierre Jolivet if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE));
81fdc842d1SBarry Smith #endif
829566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, m_f, m_c, mx, Mx));
839566063dSJacob Faibussowitsch PetscCall(ConvertToAIJ(dac->mattype, &mattype));
849566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, mattype));
859566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 2, NULL));
869566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 2, NULL, 1, NULL));
8747c6ae99SBarry Smith
8847c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/
8985afcc9aSBarry Smith if (!NEWVERSION) {
9047c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) {
9147c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */
92e3fbd1f4SBarry Smith row = idx_f[i - i_start_ghost];
9347c6ae99SBarry Smith
9447c6ae99SBarry Smith i_c = (i / ratio); /* coarse grid node to left of fine grid node */
9500045ab3SPierre Jolivet PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA, i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT, i_start, i_c, i_start_ghost_c);
9647c6ae99SBarry Smith
9747c6ae99SBarry Smith /*
9847c6ae99SBarry Smith Only include those interpolation points that are truly
9947c6ae99SBarry Smith nonzero. Note this is very important for final grid lines
10047c6ae99SBarry Smith in x direction; since they have no right neighbor
10147c6ae99SBarry Smith */
1026712e2f1SBarry Smith x = ((PetscReal)(i - i_c * ratio)) / ((PetscReal)ratio);
10347c6ae99SBarry Smith nc = 0;
10447c6ae99SBarry Smith /* one left and below; or we are right on it */
105e3fbd1f4SBarry Smith col = (i_c - i_start_ghost_c);
106e3fbd1f4SBarry Smith cols[nc] = idx_c[col];
10747c6ae99SBarry Smith v[nc++] = -x + 1.0;
10847c6ae99SBarry Smith /* one right? */
10947c6ae99SBarry Smith if (i_c * ratio != i) {
110e3fbd1f4SBarry Smith cols[nc] = idx_c[col + 1];
11147c6ae99SBarry Smith v[nc++] = x;
11247c6ae99SBarry Smith }
1139566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES));
11447c6ae99SBarry Smith }
115b3bd94feSDave May
116b3bd94feSDave May } else {
117b3bd94feSDave May PetscScalar *xi;
118b3bd94feSDave May PetscInt li, nxi, n;
119b3bd94feSDave May PetscScalar Ni[2];
120b3bd94feSDave May
121b3bd94feSDave May /* compute local coordinate arrays */
122b3bd94feSDave May nxi = ratio + 1;
1239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nxi, &xi));
124ad540459SPierre Jolivet for (li = 0; li < nxi; li++) xi[li] = -1.0 + (PetscScalar)li * (2.0 / (PetscScalar)(nxi - 1));
125b3bd94feSDave May
126b3bd94feSDave May for (i = i_start; i < i_start + m_f; i++) {
127b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */
1284ad8454bSPierre Jolivet row = idx_f[i - i_start_ghost];
129b3bd94feSDave May
130b3bd94feSDave May i_c = (i / ratio); /* coarse grid node to left of fine grid node */
13100045ab3SPierre Jolivet PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA, i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT, i_start, i_c, i_start_ghost_c); /* remainders */
132b3bd94feSDave May li = i - ratio * (i / ratio);
1338865f1eaSKarl Rupp if (i == mx - 1) li = nxi - 1;
134b3bd94feSDave May
135b3bd94feSDave May /* corners */
136e3fbd1f4SBarry Smith col = (i_c - i_start_ghost_c);
137e3fbd1f4SBarry Smith cols[0] = idx_c[col];
138b3bd94feSDave May Ni[0] = 1.0;
139b3bd94feSDave May if ((li == 0) || (li == nxi - 1)) {
1409566063dSJacob Faibussowitsch PetscCall(MatSetValue(mat, row, cols[0], Ni[0], INSERT_VALUES));
141b3bd94feSDave May continue;
142b3bd94feSDave May }
143b3bd94feSDave May
144b3bd94feSDave May /* edges + interior */
145b3bd94feSDave May /* remainders */
1468865f1eaSKarl Rupp if (i == mx - 1) i_c--;
147b3bd94feSDave May
148e3fbd1f4SBarry Smith col = (i_c - i_start_ghost_c);
149e3fbd1f4SBarry Smith cols[0] = idx_c[col]; /* one left and below; or we are right on it */
150e3fbd1f4SBarry Smith cols[1] = idx_c[col + 1];
151b3bd94feSDave May
152b3bd94feSDave May Ni[0] = 0.5 * (1.0 - xi[li]);
153b3bd94feSDave May Ni[1] = 0.5 * (1.0 + xi[li]);
154b3bd94feSDave May for (n = 0; n < 2; n++) {
1558865f1eaSKarl Rupp if (PetscAbsScalar(Ni[n]) < 1.0e-32) cols[n] = -1;
156b3bd94feSDave May }
1579566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, 2, cols, Ni, INSERT_VALUES));
158b3bd94feSDave May }
1599566063dSJacob Faibussowitsch PetscCall(PetscFree(xi));
160b3bd94feSDave May }
1619566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
1629566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
1639566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
1649566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
1659566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(mat, dof, A));
1669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat));
1673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
16847c6ae99SBarry Smith }
16947c6ae99SBarry Smith
DMCreateInterpolation_DA_1D_Q0(DM dac,DM daf,Mat * A)17066976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_DA_1D_Q0(DM dac, DM daf, Mat *A)
171d71ae5a4SJacob Faibussowitsch {
1728ea3bf28SBarry Smith PetscInt i, i_start, m_f, Mx;
1738ea3bf28SBarry Smith const PetscInt *idx_f, *idx_c;
174e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c;
1758ea3bf28SBarry Smith PetscInt m_ghost, m_ghost_c;
17647c6ae99SBarry Smith PetscInt row, col, i_start_ghost, mx, m_c, nc, ratio;
17747c6ae99SBarry Smith PetscInt i_c, i_start_c, i_start_ghost_c, cols[2], dof;
17847c6ae99SBarry Smith PetscScalar v[2], x;
17947c6ae99SBarry Smith Mat mat;
180bff4a2f0SMatthew G. Knepley DMBoundaryType bx;
181fdc842d1SBarry Smith MatType mattype;
18247c6ae99SBarry Smith
18347c6ae99SBarry Smith PetscFunctionBegin;
1849566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bx, NULL, NULL, NULL));
1859566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
186bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) {
18763a3b9bcSJacob Faibussowitsch PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx);
18847c6ae99SBarry Smith ratio = mx / Mx;
18963a3b9bcSJacob 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);
19047c6ae99SBarry Smith } else {
19163a3b9bcSJacob 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);
19247c6ae99SBarry Smith ratio = (mx - 1) / (Mx - 1);
1931dca8a05SBarry 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);
19447c6ae99SBarry Smith }
19547c6ae99SBarry Smith
1969566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, NULL, NULL, &m_f, NULL, NULL));
1979566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, NULL, NULL, &m_ghost, NULL, NULL));
1989566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f));
1999566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
20047c6ae99SBarry Smith
2019566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, NULL, NULL, &m_c, NULL, NULL));
2029566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, NULL, NULL, &m_ghost_c, NULL, NULL));
2039566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c));
2049566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
20547c6ae99SBarry Smith
20647c6ae99SBarry Smith /* create interpolation matrix */
2079566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dac), &mat));
208fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
209fdc842d1SBarry Smith /*
210fdc842d1SBarry Smith Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
211fdc842d1SBarry Smith we don't want the original unconverted matrix copied to the GPU
212fdc842d1SBarry Smith */
21348a46eb9SPierre Jolivet if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE));
214fdc842d1SBarry Smith #endif
2159566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, m_f, m_c, mx, Mx));
2169566063dSJacob Faibussowitsch PetscCall(ConvertToAIJ(dac->mattype, &mattype));
2179566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, mattype));
2189566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 2, NULL));
2199566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 2, NULL, 0, NULL));
22047c6ae99SBarry Smith
22147c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/
22247c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) {
22347c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */
2244ad8454bSPierre Jolivet row = idx_f[i - i_start_ghost];
22547c6ae99SBarry Smith
22647c6ae99SBarry Smith i_c = (i / ratio); /* coarse grid node to left of fine grid node */
22747c6ae99SBarry Smith
22847c6ae99SBarry Smith /*
22947c6ae99SBarry Smith Only include those interpolation points that are truly
23047c6ae99SBarry Smith nonzero. Note this is very important for final grid lines
23147c6ae99SBarry Smith in x direction; since they have no right neighbor
23247c6ae99SBarry Smith */
2336712e2f1SBarry Smith x = ((PetscReal)(i - i_c * ratio)) / ((PetscReal)ratio);
23447c6ae99SBarry Smith nc = 0;
23547c6ae99SBarry Smith /* one left and below; or we are right on it */
236e3fbd1f4SBarry Smith col = (i_c - i_start_ghost_c);
237e3fbd1f4SBarry Smith cols[nc] = idx_c[col];
23847c6ae99SBarry Smith v[nc++] = -x + 1.0;
23947c6ae99SBarry Smith /* one right? */
24047c6ae99SBarry Smith if (i_c * ratio != i) {
241e3fbd1f4SBarry Smith cols[nc] = idx_c[col + 1];
24247c6ae99SBarry Smith v[nc++] = x;
24347c6ae99SBarry Smith }
2449566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES));
24547c6ae99SBarry Smith }
2469566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
2479566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
2489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
2499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
2509566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(mat, dof, A));
2519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat));
2529566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(5.0 * m_f));
2533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
25447c6ae99SBarry Smith }
25547c6ae99SBarry Smith
DMCreateInterpolation_DA_2D_Q1(DM dac,DM daf,Mat * A)25666976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_DA_2D_Q1(DM dac, DM daf, Mat *A)
257d71ae5a4SJacob Faibussowitsch {
2588ea3bf28SBarry Smith PetscInt i, j, i_start, j_start, m_f, n_f, Mx, My, dof;
2598ea3bf28SBarry Smith const PetscInt *idx_c, *idx_f;
260e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c;
2618ea3bf28SBarry Smith PetscInt m_ghost, n_ghost, m_ghost_c, n_ghost_c, *dnz, *onz;
26247c6ae99SBarry Smith PetscInt row, col, i_start_ghost, j_start_ghost, cols[4], mx, m_c, my, nc, ratioi, ratioj;
26347c6ae99SBarry 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;
26447c6ae99SBarry Smith PetscMPIInt size_c, size_f, rank_f;
26547c6ae99SBarry Smith PetscScalar v[4], x, y;
26647c6ae99SBarry Smith Mat mat;
267bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by;
268fdc842d1SBarry Smith MatType mattype;
26947c6ae99SBarry Smith
27047c6ae99SBarry Smith PetscFunctionBegin;
2719566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, NULL, NULL, NULL, NULL, NULL, NULL, &bx, &by, NULL, NULL));
2729566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
273bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) {
27463a3b9bcSJacob Faibussowitsch PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx);
27547c6ae99SBarry Smith ratioi = mx / Mx;
27663a3b9bcSJacob 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);
27747c6ae99SBarry Smith } else {
27863a3b9bcSJacob 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);
27947c6ae99SBarry Smith ratioi = (mx - 1) / (Mx - 1);
2801dca8a05SBarry 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);
28147c6ae99SBarry Smith }
282bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) {
28363a3b9bcSJacob Faibussowitsch PetscCheck(My, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be positive", My);
28447c6ae99SBarry Smith ratioj = my / My;
28563a3b9bcSJacob 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);
28647c6ae99SBarry Smith } else {
28763a3b9bcSJacob 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);
28847c6ae99SBarry Smith ratioj = (my - 1) / (My - 1);
2891dca8a05SBarry 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);
29047c6ae99SBarry Smith }
29147c6ae99SBarry Smith
2929566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, NULL, &m_f, &n_f, NULL));
2939566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, NULL, &m_ghost, &n_ghost, NULL));
2949566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f));
2959566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
29647c6ae99SBarry Smith
2979566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, NULL, &m_c, &n_c, NULL));
2989566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, NULL, &m_ghost_c, &n_ghost_c, NULL));
2999566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c));
3009566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
30147c6ae99SBarry Smith
30247c6ae99SBarry Smith /*
303aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
30447c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
30547c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is
30647c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes.
30747c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS
30847c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better.
30947c6ae99SBarry Smith
31047c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
31147c6ae99SBarry Smith */
3129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dac), &size_c));
3139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)daf), &size_f));
3149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)daf), &rank_f));
31547c6ae99SBarry Smith col_scale = size_f / size_c;
31647c6ae99SBarry Smith col_shift = Mx * My * (rank_f / size_c);
31747c6ae99SBarry Smith
318d0609cedSBarry Smith MatPreallocateBegin(PetscObjectComm((PetscObject)daf), m_f * n_f, col_scale * m_c * n_c, dnz, onz);
31947c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) {
32047c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) {
32147c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */
322e3fbd1f4SBarry Smith row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
32347c6ae99SBarry Smith
32447c6ae99SBarry Smith i_c = (i / ratioi); /* coarse grid node to left of fine grid node */
32547c6ae99SBarry Smith j_c = (j / ratioj); /* coarse grid node below fine grid node */
32647c6ae99SBarry Smith
32700045ab3SPierre Jolivet PetscCheck(j_c >= j_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA, j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT, j_start, j_c, j_start_ghost_c);
32800045ab3SPierre Jolivet PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA, i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT, i_start, i_c, i_start_ghost_c);
32947c6ae99SBarry Smith
33047c6ae99SBarry Smith /*
33147c6ae99SBarry Smith Only include those interpolation points that are truly
33247c6ae99SBarry Smith nonzero. Note this is very important for final grid lines
33347c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors
33447c6ae99SBarry Smith */
33547c6ae99SBarry Smith nc = 0;
33647c6ae99SBarry Smith /* one left and below; or we are right on it */
337e3fbd1f4SBarry Smith col = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c));
338e3fbd1f4SBarry Smith cols[nc++] = col_shift + idx_c[col];
33947c6ae99SBarry Smith /* one right and below */
340e3fbd1f4SBarry Smith if (i_c * ratioi != i) cols[nc++] = col_shift + idx_c[col + 1];
34147c6ae99SBarry Smith /* one left and above */
342e3fbd1f4SBarry Smith if (j_c * ratioj != j) cols[nc++] = col_shift + idx_c[col + m_ghost_c];
34347c6ae99SBarry Smith /* one right and above */
344e3fbd1f4SBarry Smith if (i_c * ratioi != i && j_c * ratioj != j) cols[nc++] = col_shift + idx_c[col + (m_ghost_c + 1)];
3459566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, nc, cols, dnz, onz));
34647c6ae99SBarry Smith }
34747c6ae99SBarry Smith }
3489566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)daf), &mat));
349fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
350fdc842d1SBarry Smith /*
351fdc842d1SBarry Smith Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
352fdc842d1SBarry Smith we don't want the original unconverted matrix copied to the GPU
353fdc842d1SBarry Smith */
35448a46eb9SPierre Jolivet if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE));
355fdc842d1SBarry Smith #endif
3569566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, m_f * n_f, col_scale * m_c * n_c, mx * my, col_scale * Mx * My));
3579566063dSJacob Faibussowitsch PetscCall(ConvertToAIJ(dac->mattype, &mattype));
3589566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, mattype));
3599566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 0, dnz));
3609566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 0, dnz, 0, onz));
361d0609cedSBarry Smith MatPreallocateEnd(dnz, onz);
36247c6ae99SBarry Smith
36347c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/
36485afcc9aSBarry Smith if (!NEWVERSION) {
36547c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) {
36647c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) {
36747c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */
368e3fbd1f4SBarry Smith row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
36947c6ae99SBarry Smith
37047c6ae99SBarry Smith i_c = (i / ratioi); /* coarse grid node to left of fine grid node */
37147c6ae99SBarry Smith j_c = (j / ratioj); /* coarse grid node below fine grid node */
37247c6ae99SBarry Smith
37347c6ae99SBarry Smith /*
37447c6ae99SBarry Smith Only include those interpolation points that are truly
37547c6ae99SBarry Smith nonzero. Note this is very important for final grid lines
37647c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors
37747c6ae99SBarry Smith */
3786712e2f1SBarry Smith x = ((PetscReal)(i - i_c * ratioi)) / ((PetscReal)ratioi);
3796712e2f1SBarry Smith y = ((PetscReal)(j - j_c * ratioj)) / ((PetscReal)ratioj);
380b3bd94feSDave May
38147c6ae99SBarry Smith nc = 0;
38247c6ae99SBarry Smith /* one left and below; or we are right on it */
383e3fbd1f4SBarry Smith col = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c));
384e3fbd1f4SBarry Smith cols[nc] = col_shift + idx_c[col];
38547c6ae99SBarry Smith v[nc++] = x * y - x - y + 1.0;
38647c6ae99SBarry Smith /* one right and below */
38747c6ae99SBarry Smith if (i_c * ratioi != i) {
388e3fbd1f4SBarry Smith cols[nc] = col_shift + idx_c[col + 1];
38947c6ae99SBarry Smith v[nc++] = -x * y + x;
39047c6ae99SBarry Smith }
39147c6ae99SBarry Smith /* one left and above */
39247c6ae99SBarry Smith if (j_c * ratioj != j) {
393e3fbd1f4SBarry Smith cols[nc] = col_shift + idx_c[col + m_ghost_c];
39447c6ae99SBarry Smith v[nc++] = -x * y + y;
39547c6ae99SBarry Smith }
39647c6ae99SBarry Smith /* one right and above */
39747c6ae99SBarry Smith if (j_c * ratioj != j && i_c * ratioi != i) {
398e3fbd1f4SBarry Smith cols[nc] = col_shift + idx_c[col + (m_ghost_c + 1)];
39947c6ae99SBarry Smith v[nc++] = x * y;
40047c6ae99SBarry Smith }
4019566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES));
40247c6ae99SBarry Smith }
40347c6ae99SBarry Smith }
404b3bd94feSDave May
405b3bd94feSDave May } else {
406b3bd94feSDave May PetscScalar Ni[4];
407b3bd94feSDave May PetscScalar *xi, *eta;
408b3bd94feSDave May PetscInt li, nxi, lj, neta;
409b3bd94feSDave May
410b3bd94feSDave May /* compute local coordinate arrays */
411b3bd94feSDave May nxi = ratioi + 1;
412b3bd94feSDave May neta = ratioj + 1;
4139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nxi, &xi));
4149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(neta, &eta));
415ad540459SPierre Jolivet for (li = 0; li < nxi; li++) xi[li] = -1.0 + (PetscScalar)li * (2.0 / (PetscScalar)(nxi - 1));
416ad540459SPierre Jolivet for (lj = 0; lj < neta; lj++) eta[lj] = -1.0 + (PetscScalar)lj * (2.0 / (PetscScalar)(neta - 1));
417b3bd94feSDave May
418b3bd94feSDave May /* loop over local fine grid nodes setting interpolation for those*/
419b3bd94feSDave May for (j = j_start; j < j_start + n_f; j++) {
420b3bd94feSDave May for (i = i_start; i < i_start + m_f; i++) {
421b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */
422e3fbd1f4SBarry Smith row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
423b3bd94feSDave May
424b3bd94feSDave May i_c = (i / ratioi); /* coarse grid node to left of fine grid node */
425b3bd94feSDave May j_c = (j / ratioj); /* coarse grid node below fine grid node */
426b3bd94feSDave May
427b3bd94feSDave May /* remainders */
428b3bd94feSDave May li = i - ratioi * (i / ratioi);
4298865f1eaSKarl Rupp if (i == mx - 1) li = nxi - 1;
430b3bd94feSDave May lj = j - ratioj * (j / ratioj);
4318865f1eaSKarl Rupp if (j == my - 1) lj = neta - 1;
432b3bd94feSDave May
433b3bd94feSDave May /* corners */
434e3fbd1f4SBarry Smith col = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c));
435e3fbd1f4SBarry Smith cols[0] = col_shift + idx_c[col]; /* left, below */
436b3bd94feSDave May Ni[0] = 1.0;
437b3bd94feSDave May if ((li == 0) || (li == nxi - 1)) {
438b3bd94feSDave May if ((lj == 0) || (lj == neta - 1)) {
4399566063dSJacob Faibussowitsch PetscCall(MatSetValue(mat, row, cols[0], Ni[0], INSERT_VALUES));
440b3bd94feSDave May continue;
441b3bd94feSDave May }
442b3bd94feSDave May }
443b3bd94feSDave May
444b3bd94feSDave May /* edges + interior */
445b3bd94feSDave May /* remainders */
4468865f1eaSKarl Rupp if (i == mx - 1) i_c--;
4478865f1eaSKarl Rupp if (j == my - 1) j_c--;
448b3bd94feSDave May
449e3fbd1f4SBarry Smith col = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c));
450e3fbd1f4SBarry Smith cols[0] = col_shift + idx_c[col]; /* left, below */
451e3fbd1f4SBarry Smith cols[1] = col_shift + idx_c[col + 1]; /* right, below */
452e3fbd1f4SBarry Smith cols[2] = col_shift + idx_c[col + m_ghost_c]; /* left, above */
453e3fbd1f4SBarry Smith cols[3] = col_shift + idx_c[col + (m_ghost_c + 1)]; /* right, above */
454b3bd94feSDave May
455b3bd94feSDave May Ni[0] = 0.25 * (1.0 - xi[li]) * (1.0 - eta[lj]);
456b3bd94feSDave May Ni[1] = 0.25 * (1.0 + xi[li]) * (1.0 - eta[lj]);
457b3bd94feSDave May Ni[2] = 0.25 * (1.0 - xi[li]) * (1.0 + eta[lj]);
458b3bd94feSDave May Ni[3] = 0.25 * (1.0 + xi[li]) * (1.0 + eta[lj]);
459b3bd94feSDave May
460b3bd94feSDave May nc = 0;
4618865f1eaSKarl Rupp if (PetscAbsScalar(Ni[0]) < 1.0e-32) cols[0] = -1;
4628865f1eaSKarl Rupp if (PetscAbsScalar(Ni[1]) < 1.0e-32) cols[1] = -1;
4638865f1eaSKarl Rupp if (PetscAbsScalar(Ni[2]) < 1.0e-32) cols[2] = -1;
4648865f1eaSKarl Rupp if (PetscAbsScalar(Ni[3]) < 1.0e-32) cols[3] = -1;
465b3bd94feSDave May
4669566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, 4, cols, Ni, INSERT_VALUES));
467b3bd94feSDave May }
468b3bd94feSDave May }
4699566063dSJacob Faibussowitsch PetscCall(PetscFree(xi));
4709566063dSJacob Faibussowitsch PetscCall(PetscFree(eta));
471b3bd94feSDave May }
4729566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
4739566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
4749566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
4759566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
4769566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(mat, dof, A));
4779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat));
4783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
47947c6ae99SBarry Smith }
48047c6ae99SBarry Smith
48147c6ae99SBarry Smith /*
48247c6ae99SBarry Smith Contributed by Andrei Draganescu <aidraga@sandia.gov>
48347c6ae99SBarry Smith */
DMCreateInterpolation_DA_2D_Q0(DM dac,DM daf,Mat * A)48466976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_DA_2D_Q0(DM dac, DM daf, Mat *A)
485d71ae5a4SJacob Faibussowitsch {
4868ea3bf28SBarry Smith PetscInt i, j, i_start, j_start, m_f, n_f, Mx, My, dof;
4878ea3bf28SBarry Smith const PetscInt *idx_c, *idx_f;
488e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c;
4898ea3bf28SBarry Smith PetscInt m_ghost, n_ghost, m_ghost_c, n_ghost_c, *dnz, *onz;
49047c6ae99SBarry Smith PetscInt row, col, i_start_ghost, j_start_ghost, cols[4], mx, m_c, my, nc, ratioi, ratioj;
49147c6ae99SBarry 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;
49247c6ae99SBarry Smith PetscMPIInt size_c, size_f, rank_f;
49347c6ae99SBarry Smith PetscScalar v[4];
49447c6ae99SBarry Smith Mat mat;
495bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by;
496fdc842d1SBarry Smith MatType mattype;
49747c6ae99SBarry Smith
49847c6ae99SBarry Smith PetscFunctionBegin;
4999566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, NULL, NULL, NULL, NULL, NULL, NULL, &bx, &by, NULL, NULL));
5009566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
50163a3b9bcSJacob Faibussowitsch PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx);
50263a3b9bcSJacob Faibussowitsch PetscCheck(My, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be positive", My);
50347c6ae99SBarry Smith ratioi = mx / Mx;
50447c6ae99SBarry Smith ratioj = my / My;
50508401ef6SPierre Jolivet PetscCheck(ratioi * Mx == mx, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in x");
50608401ef6SPierre Jolivet PetscCheck(ratioj * My == my, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in y");
50708401ef6SPierre Jolivet PetscCheck(ratioi == 2, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in x must be 2");
50808401ef6SPierre Jolivet PetscCheck(ratioj == 2, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in y must be 2");
50947c6ae99SBarry Smith
5109566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, NULL, &m_f, &n_f, NULL));
5119566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, NULL, &m_ghost, &n_ghost, NULL));
5129566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f));
5139566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
51447c6ae99SBarry Smith
5159566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, NULL, &m_c, &n_c, NULL));
5169566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, NULL, &m_ghost_c, &n_ghost_c, NULL));
5179566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c));
5189566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
51947c6ae99SBarry Smith
52047c6ae99SBarry Smith /*
521aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
52247c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
52347c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is
52447c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes.
52547c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS
52647c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better.
52747c6ae99SBarry Smith
52847c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
52947c6ae99SBarry Smith */
5309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dac), &size_c));
5319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)daf), &size_f));
5329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)daf), &rank_f));
53347c6ae99SBarry Smith col_scale = size_f / size_c;
53447c6ae99SBarry Smith col_shift = Mx * My * (rank_f / size_c);
53547c6ae99SBarry Smith
536d0609cedSBarry Smith MatPreallocateBegin(PetscObjectComm((PetscObject)daf), m_f * n_f, col_scale * m_c * n_c, dnz, onz);
53747c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) {
53847c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) {
53947c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */
540e3fbd1f4SBarry Smith row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
54147c6ae99SBarry Smith
54247c6ae99SBarry Smith i_c = (i / ratioi); /* coarse grid node to left of fine grid node */
54347c6ae99SBarry Smith j_c = (j / ratioj); /* coarse grid node below fine grid node */
54447c6ae99SBarry Smith
54500045ab3SPierre Jolivet PetscCheck(j_c >= j_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA, j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT, j_start, j_c, j_start_ghost_c);
54600045ab3SPierre Jolivet PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA, i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT, i_start, i_c, i_start_ghost_c);
54747c6ae99SBarry Smith
54847c6ae99SBarry Smith /*
54947c6ae99SBarry Smith Only include those interpolation points that are truly
55047c6ae99SBarry Smith nonzero. Note this is very important for final grid lines
55147c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors
55247c6ae99SBarry Smith */
55347c6ae99SBarry Smith nc = 0;
55447c6ae99SBarry Smith /* one left and below; or we are right on it */
555e3fbd1f4SBarry Smith col = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c));
556e3fbd1f4SBarry Smith cols[nc++] = col_shift + idx_c[col];
5579566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, nc, cols, dnz, onz));
55847c6ae99SBarry Smith }
55947c6ae99SBarry Smith }
5609566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)daf), &mat));
561fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
562fdc842d1SBarry Smith /*
563fdc842d1SBarry Smith Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
564fdc842d1SBarry Smith we don't want the original unconverted matrix copied to the GPU
565fdc842d1SBarry Smith */
56648a46eb9SPierre Jolivet if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE));
567fdc842d1SBarry Smith #endif
5689566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, m_f * n_f, col_scale * m_c * n_c, mx * my, col_scale * Mx * My));
5699566063dSJacob Faibussowitsch PetscCall(ConvertToAIJ(dac->mattype, &mattype));
5709566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, mattype));
5719566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 0, dnz));
5729566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 0, dnz, 0, onz));
573d0609cedSBarry Smith MatPreallocateEnd(dnz, onz);
57447c6ae99SBarry Smith
57547c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/
57647c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) {
57747c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) {
57847c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */
579e3fbd1f4SBarry Smith row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
58047c6ae99SBarry Smith
58147c6ae99SBarry Smith i_c = (i / ratioi); /* coarse grid node to left of fine grid node */
58247c6ae99SBarry Smith j_c = (j / ratioj); /* coarse grid node below fine grid node */
58347c6ae99SBarry Smith nc = 0;
58447c6ae99SBarry Smith /* one left and below; or we are right on it */
585e3fbd1f4SBarry Smith col = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c));
586e3fbd1f4SBarry Smith cols[nc] = col_shift + idx_c[col];
58747c6ae99SBarry Smith v[nc++] = 1.0;
58847c6ae99SBarry Smith
5899566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES));
59047c6ae99SBarry Smith }
59147c6ae99SBarry Smith }
5929566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
5939566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
5949566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
5959566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
5969566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(mat, dof, A));
5979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat));
5989566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(13.0 * m_f * n_f));
5993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
60047c6ae99SBarry Smith }
60147c6ae99SBarry Smith
60247c6ae99SBarry Smith /*
60347c6ae99SBarry Smith Contributed by Jianming Yang <jianming-yang@uiowa.edu>
60447c6ae99SBarry Smith */
DMCreateInterpolation_DA_3D_Q0(DM dac,DM daf,Mat * A)60566976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_DA_3D_Q0(DM dac, DM daf, Mat *A)
606d71ae5a4SJacob Faibussowitsch {
6078ea3bf28SBarry Smith PetscInt i, j, l, i_start, j_start, l_start, m_f, n_f, p_f, Mx, My, Mz, dof;
6088ea3bf28SBarry Smith const PetscInt *idx_c, *idx_f;
609e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c;
6108ea3bf28SBarry Smith PetscInt m_ghost, n_ghost, p_ghost, m_ghost_c, n_ghost_c, p_ghost_c, nc, *dnz, *onz;
61147c6ae99SBarry 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;
61247c6ae99SBarry 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;
61347c6ae99SBarry Smith PetscMPIInt size_c, size_f, rank_f;
61447c6ae99SBarry Smith PetscScalar v[8];
61547c6ae99SBarry Smith Mat mat;
616bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by, bz;
617fdc842d1SBarry Smith MatType mattype;
61847c6ae99SBarry Smith
61947c6ae99SBarry Smith PetscFunctionBegin;
6209566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, &Mz, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL));
62163a3b9bcSJacob Faibussowitsch PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx);
62263a3b9bcSJacob Faibussowitsch PetscCheck(My, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be positive", My);
62363a3b9bcSJacob Faibussowitsch PetscCheck(Mz, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of z coarse grid points %" PetscInt_FMT " must be positive", Mz);
6249566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, &mz, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
62547c6ae99SBarry Smith ratioi = mx / Mx;
62647c6ae99SBarry Smith ratioj = my / My;
62747c6ae99SBarry Smith ratiol = mz / Mz;
62808401ef6SPierre Jolivet PetscCheck(ratioi * Mx == mx, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in x");
62908401ef6SPierre Jolivet PetscCheck(ratioj * My == my, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in y");
63008401ef6SPierre Jolivet PetscCheck(ratiol * Mz == mz, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in z");
6311dca8a05SBarry Smith PetscCheck(ratioi == 2 || ratioi == 1, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in x must be 1 or 2");
6321dca8a05SBarry Smith PetscCheck(ratioj == 2 || ratioj == 1, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in y must be 1 or 2");
6331dca8a05SBarry Smith PetscCheck(ratiol == 2 || ratiol == 1, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in z must be 1 or 2");
63447c6ae99SBarry Smith
6359566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, &l_start, &m_f, &n_f, &p_f));
6369566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, &l_start_ghost, &m_ghost, &n_ghost, &p_ghost));
6379566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f));
6389566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
63947c6ae99SBarry Smith
6409566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, &l_start_c, &m_c, &n_c, &p_c));
6419566063dSJacob 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));
6429566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c));
6439566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
644e3fbd1f4SBarry Smith
64547c6ae99SBarry Smith /*
646aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
64747c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
64847c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is
64947c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes.
65047c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS
65147c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better.
65247c6ae99SBarry Smith
65347c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
65447c6ae99SBarry Smith */
6559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dac), &size_c));
6569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)daf), &size_f));
6579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)daf), &rank_f));
65847c6ae99SBarry Smith col_scale = size_f / size_c;
65947c6ae99SBarry Smith col_shift = Mx * My * Mz * (rank_f / size_c);
66047c6ae99SBarry Smith
661d0609cedSBarry Smith MatPreallocateBegin(PetscObjectComm((PetscObject)daf), m_f * n_f * p_f, col_scale * m_c * n_c * p_c, dnz, onz);
66247c6ae99SBarry Smith for (l = l_start; l < l_start + p_f; l++) {
66347c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) {
66447c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) {
66547c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */
666e3fbd1f4SBarry Smith row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
66747c6ae99SBarry Smith
66847c6ae99SBarry Smith i_c = (i / ratioi); /* coarse grid node to left of fine grid node */
66947c6ae99SBarry Smith j_c = (j / ratioj); /* coarse grid node below fine grid node */
67047c6ae99SBarry Smith l_c = (l / ratiol);
67147c6ae99SBarry Smith
67200045ab3SPierre Jolivet PetscCheck(l_c >= l_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA, l_start %" PetscInt_FMT " l_c %" PetscInt_FMT " l_start_ghost_c %" PetscInt_FMT, l_start, l_c, l_start_ghost_c);
67300045ab3SPierre Jolivet PetscCheck(j_c >= j_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA, j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT, j_start, j_c, j_start_ghost_c);
67400045ab3SPierre Jolivet PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA, i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT, i_start, i_c, i_start_ghost_c);
67547c6ae99SBarry Smith
67647c6ae99SBarry Smith /*
67747c6ae99SBarry Smith Only include those interpolation points that are truly
67847c6ae99SBarry Smith nonzero. Note this is very important for final grid lines
67947c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors
68047c6ae99SBarry Smith */
68147c6ae99SBarry Smith nc = 0;
68247c6ae99SBarry Smith /* one left and below; or we are right on it */
683e3fbd1f4SBarry 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));
684e3fbd1f4SBarry Smith cols[nc++] = col_shift + idx_c[col];
6859566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, nc, cols, dnz, onz));
68647c6ae99SBarry Smith }
68747c6ae99SBarry Smith }
68847c6ae99SBarry Smith }
6899566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)daf), &mat));
690fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
691fdc842d1SBarry Smith /*
692fdc842d1SBarry Smith Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
693fdc842d1SBarry Smith we don't want the original unconverted matrix copied to the GPU
694fdc842d1SBarry Smith */
69548a46eb9SPierre Jolivet if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE));
696fdc842d1SBarry Smith #endif
6979566063dSJacob 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));
6989566063dSJacob Faibussowitsch PetscCall(ConvertToAIJ(dac->mattype, &mattype));
6999566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, mattype));
7009566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 0, dnz));
7019566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 0, dnz, 0, onz));
702d0609cedSBarry Smith MatPreallocateEnd(dnz, onz);
70347c6ae99SBarry Smith
70447c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/
70547c6ae99SBarry Smith for (l = l_start; l < l_start + p_f; l++) {
70647c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) {
70747c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) {
70847c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */
709e3fbd1f4SBarry Smith row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
71047c6ae99SBarry Smith
71147c6ae99SBarry Smith i_c = (i / ratioi); /* coarse grid node to left of fine grid node */
71247c6ae99SBarry Smith j_c = (j / ratioj); /* coarse grid node below fine grid node */
71347c6ae99SBarry Smith l_c = (l / ratiol);
71447c6ae99SBarry Smith nc = 0;
71547c6ae99SBarry Smith /* one left and below; or we are right on it */
716e3fbd1f4SBarry 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));
717e3fbd1f4SBarry Smith cols[nc] = col_shift + idx_c[col];
71847c6ae99SBarry Smith v[nc++] = 1.0;
71947c6ae99SBarry Smith
7209566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES));
72147c6ae99SBarry Smith }
72247c6ae99SBarry Smith }
72347c6ae99SBarry Smith }
7249566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
7259566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
7269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
7279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
7289566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(mat, dof, A));
7299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat));
7309566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(13.0 * m_f * n_f * p_f));
7313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
73247c6ae99SBarry Smith }
73347c6ae99SBarry Smith
DMCreateInterpolation_DA_3D_Q1(DM dac,DM daf,Mat * A)73466976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_DA_3D_Q1(DM dac, DM daf, Mat *A)
735d71ae5a4SJacob Faibussowitsch {
7368ea3bf28SBarry Smith PetscInt i, j, i_start, j_start, m_f, n_f, Mx, My, dof, l;
7378ea3bf28SBarry Smith const PetscInt *idx_c, *idx_f;
738e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c;
7398ea3bf28SBarry Smith PetscInt m_ghost, n_ghost, m_ghost_c, n_ghost_c, Mz, mz;
74047c6ae99SBarry Smith PetscInt row, col, i_start_ghost, j_start_ghost, cols[8], mx, m_c, my, nc, ratioi, ratioj, ratiok;
74147c6ae99SBarry Smith PetscInt i_c, j_c, i_start_c, j_start_c, n_c, i_start_ghost_c, j_start_ghost_c;
74247c6ae99SBarry Smith PetscInt l_start, p_f, l_start_ghost, p_ghost, l_start_c, p_c;
74347c6ae99SBarry Smith PetscInt l_start_ghost_c, p_ghost_c, l_c, *dnz, *onz;
74447c6ae99SBarry Smith PetscScalar v[8], x, y, z;
74547c6ae99SBarry Smith Mat mat;
746bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by, bz;
747fdc842d1SBarry Smith MatType mattype;
74847c6ae99SBarry Smith
74947c6ae99SBarry Smith PetscFunctionBegin;
7509566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, &Mz, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL));
7519566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, &mz, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
75247c6ae99SBarry Smith if (mx == Mx) {
75347c6ae99SBarry Smith ratioi = 1;
754bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_PERIODIC) {
75563a3b9bcSJacob Faibussowitsch PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx);
75647c6ae99SBarry Smith ratioi = mx / Mx;
75763a3b9bcSJacob 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);
75847c6ae99SBarry Smith } else {
75963a3b9bcSJacob 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);
76047c6ae99SBarry Smith ratioi = (mx - 1) / (Mx - 1);
7611dca8a05SBarry 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);
76247c6ae99SBarry Smith }
76347c6ae99SBarry Smith if (my == My) {
76447c6ae99SBarry Smith ratioj = 1;
765bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_PERIODIC) {
76663a3b9bcSJacob Faibussowitsch PetscCheck(My, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be positive", My);
76747c6ae99SBarry Smith ratioj = my / My;
76863a3b9bcSJacob 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);
76947c6ae99SBarry Smith } else {
77063a3b9bcSJacob 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);
77147c6ae99SBarry Smith ratioj = (my - 1) / (My - 1);
7721dca8a05SBarry 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);
77347c6ae99SBarry Smith }
77447c6ae99SBarry Smith if (mz == Mz) {
77547c6ae99SBarry Smith ratiok = 1;
776bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_PERIODIC) {
77763a3b9bcSJacob Faibussowitsch PetscCheck(Mz, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of z coarse grid points %" PetscInt_FMT " must be positive", Mz);
77847c6ae99SBarry Smith ratiok = mz / Mz;
77963a3b9bcSJacob 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);
78047c6ae99SBarry Smith } else {
78163a3b9bcSJacob 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);
78247c6ae99SBarry Smith ratiok = (mz - 1) / (Mz - 1);
7831dca8a05SBarry 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);
78447c6ae99SBarry Smith }
78547c6ae99SBarry Smith
7869566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, &l_start, &m_f, &n_f, &p_f));
7879566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, &l_start_ghost, &m_ghost, &n_ghost, &p_ghost));
7889566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f));
7899566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
79047c6ae99SBarry Smith
7919566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, &l_start_c, &m_c, &n_c, &p_c));
7929566063dSJacob 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));
7939566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c));
7949566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
79547c6ae99SBarry Smith
79647c6ae99SBarry Smith /* create interpolation matrix, determining exact preallocation */
797d0609cedSBarry Smith MatPreallocateBegin(PetscObjectComm((PetscObject)dac), m_f * n_f * p_f, m_c * n_c * p_c, dnz, onz);
79847c6ae99SBarry Smith /* loop over local fine grid nodes counting interpolating points */
79947c6ae99SBarry Smith for (l = l_start; l < l_start + p_f; l++) {
80047c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) {
80147c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) {
80247c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */
803e3fbd1f4SBarry Smith row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
80447c6ae99SBarry Smith i_c = (i / ratioi);
80547c6ae99SBarry Smith j_c = (j / ratioj);
80647c6ae99SBarry Smith l_c = (l / ratiok);
80700045ab3SPierre Jolivet PetscCheck(l_c >= l_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA, l_start %" PetscInt_FMT " l_c %" PetscInt_FMT " l_start_ghost_c %" PetscInt_FMT, l_start, l_c, l_start_ghost_c);
80800045ab3SPierre Jolivet PetscCheck(j_c >= j_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA, j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT, j_start, j_c, j_start_ghost_c);
80900045ab3SPierre Jolivet PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA, i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT, i_start, i_c, i_start_ghost_c);
81047c6ae99SBarry Smith
81147c6ae99SBarry Smith /*
81247c6ae99SBarry Smith Only include those interpolation points that are truly
81347c6ae99SBarry Smith nonzero. Note this is very important for final grid lines
81447c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors
81547c6ae99SBarry Smith */
81647c6ae99SBarry Smith nc = 0;
817e3fbd1f4SBarry 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));
818e3fbd1f4SBarry Smith cols[nc++] = idx_c[col];
819ad540459SPierre Jolivet if (i_c * ratioi != i) cols[nc++] = idx_c[col + 1];
820ad540459SPierre Jolivet if (j_c * ratioj != j) cols[nc++] = idx_c[col + m_ghost_c];
821ad540459SPierre Jolivet if (l_c * ratiok != l) cols[nc++] = idx_c[col + m_ghost_c * n_ghost_c];
822ad540459SPierre Jolivet if (j_c * ratioj != j && i_c * ratioi != i) cols[nc++] = idx_c[col + (m_ghost_c + 1)];
823ad540459SPierre Jolivet if (j_c * ratioj != j && l_c * ratiok != l) cols[nc++] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c)];
824ad540459SPierre Jolivet if (i_c * ratioi != i && l_c * ratiok != l) cols[nc++] = idx_c[col + (m_ghost_c * n_ghost_c + 1)];
825ad540459SPierre 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)];
8269566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, nc, cols, dnz, onz));
82747c6ae99SBarry Smith }
82847c6ae99SBarry Smith }
82947c6ae99SBarry Smith }
8309566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dac), &mat));
831fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
832fdc842d1SBarry Smith /*
833fdc842d1SBarry Smith Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
834fdc842d1SBarry Smith we don't want the original unconverted matrix copied to the GPU
835fdc842d1SBarry Smith */
83648a46eb9SPierre Jolivet if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE));
837fdc842d1SBarry Smith #endif
8389566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, m_f * n_f * p_f, m_c * n_c * p_c, mx * my * mz, Mx * My * Mz));
8399566063dSJacob Faibussowitsch PetscCall(ConvertToAIJ(dac->mattype, &mattype));
8409566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, mattype));
8419566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 0, dnz));
8429566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 0, dnz, 0, onz));
843d0609cedSBarry Smith MatPreallocateEnd(dnz, onz);
84447c6ae99SBarry Smith
84547c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/
8462adb9dcfSBarry Smith if (!NEWVERSION) {
84747c6ae99SBarry Smith for (l = l_start; l < l_start + p_f; l++) {
84847c6ae99SBarry Smith for (j = j_start; j < j_start + n_f; j++) {
84947c6ae99SBarry Smith for (i = i_start; i < i_start + m_f; i++) {
85047c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */
851e3fbd1f4SBarry Smith row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
85247c6ae99SBarry Smith
85347c6ae99SBarry Smith i_c = (i / ratioi);
85447c6ae99SBarry Smith j_c = (j / ratioj);
85547c6ae99SBarry Smith l_c = (l / ratiok);
85647c6ae99SBarry Smith
85747c6ae99SBarry Smith /*
85847c6ae99SBarry Smith Only include those interpolation points that are truly
85947c6ae99SBarry Smith nonzero. Note this is very important for final grid lines
86047c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors
86147c6ae99SBarry Smith */
8626712e2f1SBarry Smith x = ((PetscReal)(i - i_c * ratioi)) / ((PetscReal)ratioi);
8636712e2f1SBarry Smith y = ((PetscReal)(j - j_c * ratioj)) / ((PetscReal)ratioj);
8646712e2f1SBarry Smith z = ((PetscReal)(l - l_c * ratiok)) / ((PetscReal)ratiok);
865b3bd94feSDave May
86647c6ae99SBarry Smith nc = 0;
86747c6ae99SBarry Smith /* one left and below; or we are right on it */
868e3fbd1f4SBarry 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));
86947c6ae99SBarry Smith
870e3fbd1f4SBarry Smith cols[nc] = idx_c[col];
87147c6ae99SBarry Smith v[nc++] = .125 * (1. - (2.0 * x - 1.)) * (1. - (2.0 * y - 1.)) * (1. - (2.0 * z - 1.));
87247c6ae99SBarry Smith
87347c6ae99SBarry Smith if (i_c * ratioi != i) {
874e3fbd1f4SBarry Smith cols[nc] = idx_c[col + 1];
87547c6ae99SBarry Smith v[nc++] = .125 * (1. + (2.0 * x - 1.)) * (1. - (2.0 * y - 1.)) * (1. - (2.0 * z - 1.));
87647c6ae99SBarry Smith }
87747c6ae99SBarry Smith
87847c6ae99SBarry Smith if (j_c * ratioj != j) {
879e3fbd1f4SBarry Smith cols[nc] = idx_c[col + m_ghost_c];
88047c6ae99SBarry Smith v[nc++] = .125 * (1. - (2.0 * x - 1.)) * (1. + (2.0 * y - 1.)) * (1. - (2.0 * z - 1.));
88147c6ae99SBarry Smith }
88247c6ae99SBarry Smith
88347c6ae99SBarry Smith if (l_c * ratiok != l) {
884e3fbd1f4SBarry Smith cols[nc] = idx_c[col + m_ghost_c * n_ghost_c];
88547c6ae99SBarry Smith v[nc++] = .125 * (1. - (2.0 * x - 1.)) * (1. - (2.0 * y - 1.)) * (1. + (2.0 * z - 1.));
88647c6ae99SBarry Smith }
88747c6ae99SBarry Smith
88847c6ae99SBarry Smith if (j_c * ratioj != j && i_c * ratioi != i) {
889e3fbd1f4SBarry Smith cols[nc] = idx_c[col + (m_ghost_c + 1)];
89047c6ae99SBarry Smith v[nc++] = .125 * (1. + (2.0 * x - 1.)) * (1. + (2.0 * y - 1.)) * (1. - (2.0 * z - 1.));
89147c6ae99SBarry Smith }
89247c6ae99SBarry Smith
89347c6ae99SBarry Smith if (j_c * ratioj != j && l_c * ratiok != l) {
894e3fbd1f4SBarry Smith cols[nc] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c)];
89547c6ae99SBarry Smith v[nc++] = .125 * (1. - (2.0 * x - 1.)) * (1. + (2.0 * y - 1.)) * (1. + (2.0 * z - 1.));
89647c6ae99SBarry Smith }
89747c6ae99SBarry Smith
89847c6ae99SBarry Smith if (i_c * ratioi != i && l_c * ratiok != l) {
899e3fbd1f4SBarry Smith cols[nc] = idx_c[col + (m_ghost_c * n_ghost_c + 1)];
90047c6ae99SBarry Smith v[nc++] = .125 * (1. + (2.0 * x - 1.)) * (1. - (2.0 * y - 1.)) * (1. + (2.0 * z - 1.));
90147c6ae99SBarry Smith }
90247c6ae99SBarry Smith
90347c6ae99SBarry Smith if (i_c * ratioi != i && l_c * ratiok != l && j_c * ratioj != j) {
904e3fbd1f4SBarry Smith cols[nc] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c + 1)];
90547c6ae99SBarry Smith v[nc++] = .125 * (1. + (2.0 * x - 1.)) * (1. + (2.0 * y - 1.)) * (1. + (2.0 * z - 1.));
90647c6ae99SBarry Smith }
9079566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES));
90847c6ae99SBarry Smith }
90947c6ae99SBarry Smith }
91047c6ae99SBarry Smith }
911b3bd94feSDave May
912b3bd94feSDave May } else {
913b3bd94feSDave May PetscScalar *xi, *eta, *zeta;
914b3bd94feSDave May PetscInt li, nxi, lj, neta, lk, nzeta, n;
915b3bd94feSDave May PetscScalar Ni[8];
916b3bd94feSDave May
917b3bd94feSDave May /* compute local coordinate arrays */
918b3bd94feSDave May nxi = ratioi + 1;
919b3bd94feSDave May neta = ratioj + 1;
920b3bd94feSDave May nzeta = ratiok + 1;
9219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nxi, &xi));
9229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(neta, &eta));
9239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeta, &zeta));
9248865f1eaSKarl Rupp for (li = 0; li < nxi; li++) xi[li] = -1.0 + (PetscScalar)li * (2.0 / (PetscScalar)(nxi - 1));
9258865f1eaSKarl Rupp for (lj = 0; lj < neta; lj++) eta[lj] = -1.0 + (PetscScalar)lj * (2.0 / (PetscScalar)(neta - 1));
9268865f1eaSKarl Rupp for (lk = 0; lk < nzeta; lk++) zeta[lk] = -1.0 + (PetscScalar)lk * (2.0 / (PetscScalar)(nzeta - 1));
927b3bd94feSDave May
928b3bd94feSDave May for (l = l_start; l < l_start + p_f; l++) {
929b3bd94feSDave May for (j = j_start; j < j_start + n_f; j++) {
930b3bd94feSDave May for (i = i_start; i < i_start + m_f; i++) {
931b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */
932e3fbd1f4SBarry Smith row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
933b3bd94feSDave May
934b3bd94feSDave May i_c = (i / ratioi);
935b3bd94feSDave May j_c = (j / ratioj);
936b3bd94feSDave May l_c = (l / ratiok);
937b3bd94feSDave May
938b3bd94feSDave May /* remainders */
939b3bd94feSDave May li = i - ratioi * (i / ratioi);
9408865f1eaSKarl Rupp if (i == mx - 1) li = nxi - 1;
941b3bd94feSDave May lj = j - ratioj * (j / ratioj);
9428865f1eaSKarl Rupp if (j == my - 1) lj = neta - 1;
943b3bd94feSDave May lk = l - ratiok * (l / ratiok);
9448865f1eaSKarl Rupp if (l == mz - 1) lk = nzeta - 1;
945b3bd94feSDave May
946b3bd94feSDave May /* corners */
947e3fbd1f4SBarry 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));
948e3fbd1f4SBarry Smith cols[0] = idx_c[col];
949b3bd94feSDave May Ni[0] = 1.0;
950b3bd94feSDave May if ((li == 0) || (li == nxi - 1)) {
951b3bd94feSDave May if ((lj == 0) || (lj == neta - 1)) {
952b3bd94feSDave May if ((lk == 0) || (lk == nzeta - 1)) {
9539566063dSJacob Faibussowitsch PetscCall(MatSetValue(mat, row, cols[0], Ni[0], INSERT_VALUES));
954b3bd94feSDave May continue;
955b3bd94feSDave May }
956b3bd94feSDave May }
957b3bd94feSDave May }
958b3bd94feSDave May
959b3bd94feSDave May /* edges + interior */
960b3bd94feSDave May /* remainders */
9618865f1eaSKarl Rupp if (i == mx - 1) i_c--;
9628865f1eaSKarl Rupp if (j == my - 1) j_c--;
9638865f1eaSKarl Rupp if (l == mz - 1) l_c--;
964b3bd94feSDave May
965e3fbd1f4SBarry 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));
966e3fbd1f4SBarry Smith cols[0] = idx_c[col]; /* one left and below; or we are right on it */
967e3fbd1f4SBarry Smith cols[1] = idx_c[col + 1]; /* one right and below */
968e3fbd1f4SBarry Smith cols[2] = idx_c[col + m_ghost_c]; /* one left and above */
969e3fbd1f4SBarry Smith cols[3] = idx_c[col + (m_ghost_c + 1)]; /* one right and above */
970b3bd94feSDave May
971e3fbd1f4SBarry Smith cols[4] = idx_c[col + m_ghost_c * n_ghost_c]; /* one left and below and front; or we are right on it */
972e3fbd1f4SBarry Smith cols[5] = idx_c[col + (m_ghost_c * n_ghost_c + 1)]; /* one right and below, and front */
973e3fbd1f4SBarry Smith cols[6] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c)]; /* one left and above and front*/
974e3fbd1f4SBarry Smith cols[7] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c + 1)]; /* one right and above and front */
975b3bd94feSDave May
976b3bd94feSDave May Ni[0] = 0.125 * (1.0 - xi[li]) * (1.0 - eta[lj]) * (1.0 - zeta[lk]);
977b3bd94feSDave May Ni[1] = 0.125 * (1.0 + xi[li]) * (1.0 - eta[lj]) * (1.0 - zeta[lk]);
978b3bd94feSDave May Ni[2] = 0.125 * (1.0 - xi[li]) * (1.0 + eta[lj]) * (1.0 - zeta[lk]);
979b3bd94feSDave May Ni[3] = 0.125 * (1.0 + xi[li]) * (1.0 + eta[lj]) * (1.0 - zeta[lk]);
980b3bd94feSDave May
981b3bd94feSDave May Ni[4] = 0.125 * (1.0 - xi[li]) * (1.0 - eta[lj]) * (1.0 + zeta[lk]);
982b3bd94feSDave May Ni[5] = 0.125 * (1.0 + xi[li]) * (1.0 - eta[lj]) * (1.0 + zeta[lk]);
983b3bd94feSDave May Ni[6] = 0.125 * (1.0 - xi[li]) * (1.0 + eta[lj]) * (1.0 + zeta[lk]);
984b3bd94feSDave May Ni[7] = 0.125 * (1.0 + xi[li]) * (1.0 + eta[lj]) * (1.0 + zeta[lk]);
985b3bd94feSDave May
986b3bd94feSDave May for (n = 0; n < 8; n++) {
9878865f1eaSKarl Rupp if (PetscAbsScalar(Ni[n]) < 1.0e-32) cols[n] = -1;
988b3bd94feSDave May }
9899566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &row, 8, cols, Ni, INSERT_VALUES));
990b3bd94feSDave May }
991b3bd94feSDave May }
992b3bd94feSDave May }
9939566063dSJacob Faibussowitsch PetscCall(PetscFree(xi));
9949566063dSJacob Faibussowitsch PetscCall(PetscFree(eta));
9959566063dSJacob Faibussowitsch PetscCall(PetscFree(zeta));
996b3bd94feSDave May }
9979566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
9989566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
9999566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
10009566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
100147c6ae99SBarry Smith
10029566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(mat, dof, A));
10039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat));
10043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
100547c6ae99SBarry Smith }
100647c6ae99SBarry Smith
DMCreateInterpolation_DA(DM dac,DM daf,Mat * A,Vec * scale)1007d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInterpolation_DA(DM dac, DM daf, Mat *A, Vec *scale)
1008d71ae5a4SJacob Faibussowitsch {
100947c6ae99SBarry Smith PetscInt dimc, Mc, Nc, Pc, mc, nc, pc, dofc, sc, dimf, Mf, Nf, Pf, mf, nf, pf, doff, sf;
1010bff4a2f0SMatthew G. Knepley DMBoundaryType bxc, byc, bzc, bxf, byf, bzf;
1011aa219208SBarry Smith DMDAStencilType stc, stf;
101247c6ae99SBarry Smith DM_DA *ddc = (DM_DA *)dac->data;
101347c6ae99SBarry Smith
101447c6ae99SBarry Smith PetscFunctionBegin;
101547c6ae99SBarry Smith PetscValidHeaderSpecific(dac, DM_CLASSID, 1);
101647c6ae99SBarry Smith PetscValidHeaderSpecific(daf, DM_CLASSID, 2);
10174f572ea9SToby Isaac PetscAssertPointer(A, 3);
10184f572ea9SToby Isaac if (scale) PetscAssertPointer(scale, 4);
101947c6ae99SBarry Smith
10209566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, &dimc, &Mc, &Nc, &Pc, &mc, &nc, &pc, &dofc, &sc, &bxc, &byc, &bzc, &stc));
10219566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, &dimf, &Mf, &Nf, &Pf, &mf, &nf, &pf, &doff, &sf, &bxf, &byf, &bzf, &stf));
102263a3b9bcSJacob Faibussowitsch PetscCheck(dimc == dimf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Dimensions of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dimc, dimf);
102363a3b9bcSJacob Faibussowitsch PetscCheck(dofc == doff, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "DOF of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dofc, doff);
102463a3b9bcSJacob Faibussowitsch PetscCheck(sc == sf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil width of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, sc, sf);
10251dca8a05SBarry Smith PetscCheck(bxc == bxf && byc == byf && bzc == bzf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Boundary type different in two DMDAs");
102608401ef6SPierre Jolivet PetscCheck(stc == stf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil type different in two DMDAs");
10271dca8a05SBarry Smith PetscCheck(Mc >= 2 || Mf <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in x direction");
10281dca8a05SBarry 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");
10291dca8a05SBarry 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");
103047c6ae99SBarry Smith
1031aa219208SBarry Smith if (ddc->interptype == DMDA_Q1) {
103247c6ae99SBarry Smith if (dimc == 1) {
10339566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation_DA_1D_Q1(dac, daf, A));
103447c6ae99SBarry Smith } else if (dimc == 2) {
10359566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation_DA_2D_Q1(dac, daf, A));
103647c6ae99SBarry Smith } else if (dimc == 3) {
10379566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation_DA_3D_Q1(dac, daf, A));
103863a3b9bcSJacob 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);
1039aa219208SBarry Smith } else if (ddc->interptype == DMDA_Q0) {
104047c6ae99SBarry Smith if (dimc == 1) {
10419566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation_DA_1D_Q0(dac, daf, A));
104247c6ae99SBarry Smith } else if (dimc == 2) {
10439566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation_DA_2D_Q0(dac, daf, A));
104447c6ae99SBarry Smith } else if (dimc == 3) {
10459566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation_DA_3D_Q0(dac, daf, A));
104663a3b9bcSJacob 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);
104747c6ae99SBarry Smith }
1048*835f2295SStefano Zampini if (scale) PetscCall(DMCreateInterpolationScale(dac, daf, *A, scale));
10493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
105047c6ae99SBarry Smith }
105147c6ae99SBarry Smith
DMCreateInjection_DA_1D(DM dac,DM daf,VecScatter * inject)105266976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInjection_DA_1D(DM dac, DM daf, VecScatter *inject)
1053d71ae5a4SJacob Faibussowitsch {
10548ea3bf28SBarry Smith PetscInt i, i_start, m_f, Mx, dof;
10558ea3bf28SBarry Smith const PetscInt *idx_f;
1056e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f;
10578ea3bf28SBarry Smith PetscInt m_ghost, m_ghost_c;
10580bb2b966SJungho Lee PetscInt row, i_start_ghost, mx, m_c, nc, ratioi;
10590bb2b966SJungho Lee PetscInt i_start_c, i_start_ghost_c;
10600bb2b966SJungho Lee PetscInt *cols;
1061bff4a2f0SMatthew G. Knepley DMBoundaryType bx;
10620bb2b966SJungho Lee Vec vecf, vecc;
10630bb2b966SJungho Lee IS isf;
10640bb2b966SJungho Lee
10650bb2b966SJungho Lee PetscFunctionBegin;
10669566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bx, NULL, NULL, NULL));
10679566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
1068bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) {
10690bb2b966SJungho Lee ratioi = mx / Mx;
107063a3b9bcSJacob 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);
10710bb2b966SJungho Lee } else {
10720bb2b966SJungho Lee ratioi = (mx - 1) / (Mx - 1);
10731dca8a05SBarry 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);
10740bb2b966SJungho Lee }
10750bb2b966SJungho Lee
10769566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, NULL, NULL, &m_f, NULL, NULL));
10779566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, NULL, NULL, &m_ghost, NULL, NULL));
10789566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f));
10799566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
10800bb2b966SJungho Lee
10819566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, NULL, NULL, &m_c, NULL, NULL));
10829566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, NULL, NULL, &m_ghost_c, NULL, NULL));
10830bb2b966SJungho Lee
10840bb2b966SJungho Lee /* loop over local fine grid nodes setting interpolation for those*/
10850bb2b966SJungho Lee nc = 0;
10869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m_f, &cols));
10870bb2b966SJungho Lee
10880bb2b966SJungho Lee for (i = i_start_c; i < i_start_c + m_c; i++) {
10890bb2b966SJungho Lee PetscInt i_f = i * ratioi;
10900bb2b966SJungho Lee
109100045ab3SPierre Jolivet 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, i_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);
10928865f1eaSKarl Rupp
10934ad8454bSPierre Jolivet row = idx_f[i_f - i_start_ghost];
1094e3fbd1f4SBarry Smith cols[nc++] = row;
10950bb2b966SJungho Lee }
10960bb2b966SJungho Lee
10979566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
10989566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)daf), dof, nc, cols, PETSC_OWN_POINTER, &isf));
10999566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dac, &vecc));
11009566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(daf, &vecf));
11019566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vecf, isf, vecc, NULL, inject));
11029566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dac, &vecc));
11039566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(daf, &vecf));
11049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isf));
11053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11060bb2b966SJungho Lee }
11070bb2b966SJungho Lee
DMCreateInjection_DA_2D(DM dac,DM daf,VecScatter * inject)110866976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInjection_DA_2D(DM dac, DM daf, VecScatter *inject)
1109d71ae5a4SJacob Faibussowitsch {
11108ea3bf28SBarry Smith PetscInt i, j, i_start, j_start, m_f, n_f, Mx, My, dof;
11118ea3bf28SBarry Smith const PetscInt *idx_c, *idx_f;
1112e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c;
11138ea3bf28SBarry Smith PetscInt m_ghost, n_ghost, m_ghost_c, n_ghost_c;
111447c6ae99SBarry Smith PetscInt row, i_start_ghost, j_start_ghost, mx, m_c, my, nc, ratioi, ratioj;
1115076202ddSJed Brown PetscInt i_start_c, j_start_c, n_c, i_start_ghost_c, j_start_ghost_c;
111647c6ae99SBarry Smith PetscInt *cols;
1117bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by;
111847c6ae99SBarry Smith Vec vecf, vecc;
111947c6ae99SBarry Smith IS isf;
112047c6ae99SBarry Smith
112147c6ae99SBarry Smith PetscFunctionBegin;
11229566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, NULL, NULL, NULL, NULL, NULL, NULL, &bx, &by, NULL, NULL));
11239566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
1124bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) {
112547c6ae99SBarry Smith ratioi = mx / Mx;
112663a3b9bcSJacob 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);
112747c6ae99SBarry Smith } else {
112847c6ae99SBarry Smith ratioi = (mx - 1) / (Mx - 1);
11291dca8a05SBarry 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);
113047c6ae99SBarry Smith }
1131bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) {
113247c6ae99SBarry Smith ratioj = my / My;
113363a3b9bcSJacob 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);
113447c6ae99SBarry Smith } else {
113547c6ae99SBarry Smith ratioj = (my - 1) / (My - 1);
11361dca8a05SBarry 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);
113747c6ae99SBarry Smith }
113847c6ae99SBarry Smith
11399566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, NULL, &m_f, &n_f, NULL));
11409566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, NULL, &m_ghost, &n_ghost, NULL));
11419566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f));
11429566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
114347c6ae99SBarry Smith
11449566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, NULL, &m_c, &n_c, NULL));
11459566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, NULL, &m_ghost_c, &n_ghost_c, NULL));
11469566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c));
11479566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
114847c6ae99SBarry Smith
114947c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/
115047c6ae99SBarry Smith nc = 0;
11519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_f * m_f, &cols));
1152076202ddSJed Brown for (j = j_start_c; j < j_start_c + n_c; j++) {
1153076202ddSJed Brown for (i = i_start_c; i < i_start_c + m_c; i++) {
1154076202ddSJed Brown PetscInt i_f = i * ratioi, j_f = j * ratioj;
115500045ab3SPierre Jolivet 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, j_c %" PetscInt_FMT " j_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]", j, j_f, j_start_ghost, j_start_ghost + n_ghost);
115600045ab3SPierre Jolivet 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, i_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);
1157e3fbd1f4SBarry Smith row = idx_f[(m_ghost * (j_f - j_start_ghost) + (i_f - i_start_ghost))];
1158e3fbd1f4SBarry Smith cols[nc++] = row;
115947c6ae99SBarry Smith }
116047c6ae99SBarry Smith }
11619566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
11629566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
116347c6ae99SBarry Smith
11649566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)daf), dof, nc, cols, PETSC_OWN_POINTER, &isf));
11659566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dac, &vecc));
11669566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(daf, &vecf));
11679566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vecf, isf, vecc, NULL, inject));
11689566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dac, &vecc));
11699566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(daf, &vecf));
11709566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isf));
11713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
117247c6ae99SBarry Smith }
117347c6ae99SBarry Smith
DMCreateInjection_DA_3D(DM dac,DM daf,VecScatter * inject)117466976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInjection_DA_3D(DM dac, DM daf, VecScatter *inject)
1175d71ae5a4SJacob Faibussowitsch {
1176b1757049SJed Brown PetscInt i, j, k, i_start, j_start, k_start, m_f, n_f, p_f, Mx, My, Mz;
1177b1757049SJed Brown PetscInt m_ghost, n_ghost, p_ghost, m_ghost_c, n_ghost_c, p_ghost_c;
1178b1757049SJed Brown PetscInt i_start_ghost, j_start_ghost, k_start_ghost;
1179b1757049SJed Brown PetscInt mx, my, mz, ratioi, ratioj, ratiok;
1180b1757049SJed Brown PetscInt i_start_c, j_start_c, k_start_c;
1181b1757049SJed Brown PetscInt m_c, n_c, p_c;
1182b1757049SJed Brown PetscInt i_start_ghost_c, j_start_ghost_c, k_start_ghost_c;
1183b1757049SJed Brown PetscInt row, nc, dof;
11848ea3bf28SBarry Smith const PetscInt *idx_c, *idx_f;
1185e3fbd1f4SBarry Smith ISLocalToGlobalMapping ltog_f, ltog_c;
1186b1757049SJed Brown PetscInt *cols;
1187bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by, bz;
1188b1757049SJed Brown Vec vecf, vecc;
1189b1757049SJed Brown IS isf;
1190b1757049SJed Brown
1191b1757049SJed Brown PetscFunctionBegin;
11929566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, &Mz, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL));
11939566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, &mz, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
1194b1757049SJed Brown
1195bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) {
1196b1757049SJed Brown ratioi = mx / Mx;
119763a3b9bcSJacob 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);
1198b1757049SJed Brown } else {
1199b1757049SJed Brown ratioi = (mx - 1) / (Mx - 1);
12001dca8a05SBarry 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);
1201b1757049SJed Brown }
1202bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) {
1203b1757049SJed Brown ratioj = my / My;
120463a3b9bcSJacob 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);
1205b1757049SJed Brown } else {
1206b1757049SJed Brown ratioj = (my - 1) / (My - 1);
12071dca8a05SBarry 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);
1208b1757049SJed Brown }
1209bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC) {
1210b1757049SJed Brown ratiok = mz / Mz;
121163a3b9bcSJacob 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);
1212b1757049SJed Brown } else {
1213b1757049SJed Brown ratiok = (mz - 1) / (Mz - 1);
12141dca8a05SBarry 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);
1215b1757049SJed Brown }
1216b1757049SJed Brown
12179566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, &k_start, &m_f, &n_f, &p_f));
12189566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, &k_start_ghost, &m_ghost, &n_ghost, &p_ghost));
12199566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <og_f));
12209566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
1221b1757049SJed Brown
12229566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, &k_start_c, &m_c, &n_c, &p_c));
12239566063dSJacob 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));
12249566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <og_c));
12259566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
1226b1757049SJed Brown
1227b1757049SJed Brown /* loop over local fine grid nodes setting interpolation for those*/
1228b1757049SJed Brown nc = 0;
12299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_f * m_f * p_f, &cols));
1230b1757049SJed Brown for (k = k_start_c; k < k_start_c + p_c; k++) {
1231b1757049SJed Brown for (j = j_start_c; j < j_start_c + n_c; j++) {
1232b1757049SJed Brown for (i = i_start_c; i < i_start_c + m_c; i++) {
1233b1757049SJed Brown PetscInt i_f = i * ratioi, j_f = j * ratioj, k_f = k * ratiok;
12349371c9d4SSatish Balay PetscCheck(k_f >= k_start_ghost && k_f < k_start_ghost + p_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
12359371c9d4SSatish Balay "Processor's coarse DMDA must lie over fine DMDA "
12369371c9d4SSatish Balay "k_c %" PetscInt_FMT " k_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",
12379371c9d4SSatish Balay k, k_f, k_start_ghost, k_start_ghost + p_ghost);
12389371c9d4SSatish Balay PetscCheck(j_f >= j_start_ghost && j_f < j_start_ghost + n_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
12399371c9d4SSatish Balay "Processor's coarse DMDA must lie over fine DMDA "
12409371c9d4SSatish Balay "j_c %" PetscInt_FMT " j_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",
12419371c9d4SSatish Balay j, j_f, j_start_ghost, j_start_ghost + n_ghost);
12429371c9d4SSatish Balay PetscCheck(i_f >= i_start_ghost && i_f < i_start_ghost + m_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
12439371c9d4SSatish Balay "Processor's coarse DMDA must lie over fine DMDA "
12449371c9d4SSatish Balay "i_c %" PetscInt_FMT " i_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",
12459371c9d4SSatish Balay i, i_f, i_start_ghost, i_start_ghost + m_ghost);
1246e3fbd1f4SBarry 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))];
1247e3fbd1f4SBarry Smith cols[nc++] = row;
1248b1757049SJed Brown }
1249b1757049SJed Brown }
1250b1757049SJed Brown }
12519566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
12529566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
1253b1757049SJed Brown
12549566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)daf), dof, nc, cols, PETSC_OWN_POINTER, &isf));
12559566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dac, &vecc));
12569566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(daf, &vecf));
12579566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vecf, isf, vecc, NULL, inject));
12589566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dac, &vecc));
12599566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(daf, &vecf));
12609566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isf));
12613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1262b1757049SJed Brown }
1263b1757049SJed Brown
DMCreateInjection_DA(DM dac,DM daf,Mat * mat)1264d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInjection_DA(DM dac, DM daf, Mat *mat)
1265d71ae5a4SJacob Faibussowitsch {
126647c6ae99SBarry Smith PetscInt dimc, Mc, Nc, Pc, mc, nc, pc, dofc, sc, dimf, Mf, Nf, Pf, mf, nf, pf, doff, sf;
1267bff4a2f0SMatthew G. Knepley DMBoundaryType bxc, byc, bzc, bxf, byf, bzf;
1268aa219208SBarry Smith DMDAStencilType stc, stf;
126960c16ac1SBarry Smith VecScatter inject = NULL;
127047c6ae99SBarry Smith
127147c6ae99SBarry Smith PetscFunctionBegin;
127247c6ae99SBarry Smith PetscValidHeaderSpecific(dac, DM_CLASSID, 1);
127347c6ae99SBarry Smith PetscValidHeaderSpecific(daf, DM_CLASSID, 2);
12744f572ea9SToby Isaac PetscAssertPointer(mat, 3);
127547c6ae99SBarry Smith
12769566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, &dimc, &Mc, &Nc, &Pc, &mc, &nc, &pc, &dofc, &sc, &bxc, &byc, &bzc, &stc));
12779566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, &dimf, &Mf, &Nf, &Pf, &mf, &nf, &pf, &doff, &sf, &bxf, &byf, &bzf, &stf));
127863a3b9bcSJacob Faibussowitsch PetscCheck(dimc == dimf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Dimensions of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dimc, dimf);
127963a3b9bcSJacob Faibussowitsch PetscCheck(dofc == doff, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "DOF of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dofc, doff);
128063a3b9bcSJacob Faibussowitsch PetscCheck(sc == sf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil width of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, sc, sf);
12811dca8a05SBarry Smith PetscCheck(bxc == bxf && byc == byf && bzc == bzf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Boundary type different in two DMDAs");
128208401ef6SPierre Jolivet PetscCheck(stc == stf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil type different in two DMDAs");
128308401ef6SPierre Jolivet PetscCheck(Mc >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in x direction");
12841dca8a05SBarry Smith PetscCheck(dimc <= 1 || Nc >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in y direction");
12851dca8a05SBarry Smith PetscCheck(dimc <= 2 || Pc >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in z direction");
128647c6ae99SBarry Smith
12870bb2b966SJungho Lee if (dimc == 1) {
12889566063dSJacob Faibussowitsch PetscCall(DMCreateInjection_DA_1D(dac, daf, &inject));
12890bb2b966SJungho Lee } else if (dimc == 2) {
12909566063dSJacob Faibussowitsch PetscCall(DMCreateInjection_DA_2D(dac, daf, &inject));
1291b1757049SJed Brown } else if (dimc == 3) {
12929566063dSJacob Faibussowitsch PetscCall(DMCreateInjection_DA_3D(dac, daf, &inject));
129347c6ae99SBarry Smith }
12949566063dSJacob Faibussowitsch PetscCall(MatCreateScatter(PetscObjectComm((PetscObject)inject), inject, mat));
12959566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&inject));
12963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
129747c6ae99SBarry Smith }
129847c6ae99SBarry Smith
1299a4e35b19SJacob Faibussowitsch // PetscClangLinter pragma ignore: -fdoc-*
130097779f9aSLisandro Dalcin /*@
130197779f9aSLisandro Dalcin DMCreateAggregates - Deprecated, see DMDACreateAggregates.
13026c877ef6SSatish Balay
13036c877ef6SSatish Balay Level: intermediate
130497779f9aSLisandro Dalcin @*/
DMCreateAggregates(DM dac,DM daf,Mat * mat)1305d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateAggregates(DM dac, DM daf, Mat *mat)
1306d71ae5a4SJacob Faibussowitsch {
1307a5bc1bf3SBarry Smith return DMDACreateAggregates(dac, daf, mat);
130897779f9aSLisandro Dalcin }
130997779f9aSLisandro Dalcin
131097779f9aSLisandro Dalcin /*@
131197779f9aSLisandro Dalcin DMDACreateAggregates - Gets the aggregates that map between
1312dce8aebaSBarry Smith grids associated with two `DMDA`
131397779f9aSLisandro Dalcin
131420f4b53cSBarry Smith Collective
131597779f9aSLisandro Dalcin
131697779f9aSLisandro Dalcin Input Parameters:
131760225df5SJacob Faibussowitsch + dac - the coarse grid `DMDA`
131860225df5SJacob Faibussowitsch - daf - the fine grid `DMDA`
131997779f9aSLisandro Dalcin
13202fe279fdSBarry Smith Output Parameter:
132197779f9aSLisandro Dalcin . rest - the restriction matrix (transpose of the projection matrix)
132297779f9aSLisandro Dalcin
132397779f9aSLisandro Dalcin Level: intermediate
132497779f9aSLisandro Dalcin
1325dce8aebaSBarry Smith Note:
1326dce8aebaSBarry Smith This routine is not used by PETSc.
132797779f9aSLisandro Dalcin It is not clear what its use case is and it may be removed in a future release.
132897779f9aSLisandro Dalcin Users should contact petsc-maint@mcs.anl.gov if they plan to use it.
132997779f9aSLisandro Dalcin
133012b4a537SBarry Smith .seealso: [](sec_struct), `DMRefine()`, `DMCreateInjection()`, `DMCreateInterpolation()`
133197779f9aSLisandro Dalcin @*/
DMDACreateAggregates(DM dac,DM daf,Mat * rest)1332d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDACreateAggregates(DM dac, DM daf, Mat *rest)
1333d71ae5a4SJacob Faibussowitsch {
133447c6ae99SBarry Smith PetscInt dimc, Mc, Nc, Pc, mc, nc, pc, dofc, sc;
133547c6ae99SBarry Smith PetscInt dimf, Mf, Nf, Pf, mf, nf, pf, doff, sf;
1336bff4a2f0SMatthew G. Knepley DMBoundaryType bxc, byc, bzc, bxf, byf, bzf;
1337aa219208SBarry Smith DMDAStencilType stc, stf;
133847c6ae99SBarry Smith PetscInt i, j, l;
133947c6ae99SBarry Smith PetscInt i_start, j_start, l_start, m_f, n_f, p_f;
134047c6ae99SBarry Smith PetscInt i_start_ghost, j_start_ghost, l_start_ghost, m_ghost, n_ghost, p_ghost;
13418ea3bf28SBarry Smith const PetscInt *idx_f;
134247c6ae99SBarry Smith PetscInt i_c, j_c, l_c;
134347c6ae99SBarry Smith PetscInt i_start_c, j_start_c, l_start_c, m_c, n_c, p_c;
134447c6ae99SBarry Smith PetscInt i_start_ghost_c, j_start_ghost_c, l_start_ghost_c, m_ghost_c, n_ghost_c, p_ghost_c;
13458ea3bf28SBarry Smith const PetscInt *idx_c;
134647c6ae99SBarry Smith PetscInt d;
134747c6ae99SBarry Smith PetscInt a;
134847c6ae99SBarry Smith PetscInt max_agg_size;
134947c6ae99SBarry Smith PetscInt *fine_nodes;
135047c6ae99SBarry Smith PetscScalar *one_vec;
135147c6ae99SBarry Smith PetscInt fn_idx;
1352565245c5SBarry Smith ISLocalToGlobalMapping ltogmf, ltogmc;
135347c6ae99SBarry Smith
135447c6ae99SBarry Smith PetscFunctionBegin;
135597779f9aSLisandro Dalcin PetscValidHeaderSpecificType(dac, DM_CLASSID, 1, DMDA);
135697779f9aSLisandro Dalcin PetscValidHeaderSpecificType(daf, DM_CLASSID, 2, DMDA);
13574f572ea9SToby Isaac PetscAssertPointer(rest, 3);
135847c6ae99SBarry Smith
13599566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, &dimc, &Mc, &Nc, &Pc, &mc, &nc, &pc, &dofc, &sc, &bxc, &byc, &bzc, &stc));
13609566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, &dimf, &Mf, &Nf, &Pf, &mf, &nf, &pf, &doff, &sf, &bxf, &byf, &bzf, &stf));
136163a3b9bcSJacob Faibussowitsch PetscCheck(dimc == dimf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Dimensions of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dimc, dimf);
136263a3b9bcSJacob Faibussowitsch PetscCheck(dofc == doff, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "DOF of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dofc, doff);
136363a3b9bcSJacob Faibussowitsch PetscCheck(sc == sf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil width of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, sc, sf);
13641dca8a05SBarry Smith PetscCheck(bxc == bxf && byc == byf && bzc == bzf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Boundary type different in two DMDAs");
136508401ef6SPierre Jolivet PetscCheck(stc == stf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil type different in two DMDAs");
136647c6ae99SBarry Smith
136763a3b9bcSJacob 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);
136863a3b9bcSJacob 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);
136963a3b9bcSJacob 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);
137047c6ae99SBarry Smith
137147c6ae99SBarry Smith if (Pc < 0) Pc = 1;
137247c6ae99SBarry Smith if (Pf < 0) Pf = 1;
137347c6ae99SBarry Smith if (Nc < 0) Nc = 1;
137447c6ae99SBarry Smith if (Nf < 0) Nf = 1;
137547c6ae99SBarry Smith
13769566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(daf, &i_start, &j_start, &l_start, &m_f, &n_f, &p_f));
13779566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, &l_start_ghost, &m_ghost, &n_ghost, &p_ghost));
1378565245c5SBarry Smith
13799566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(daf, <ogmf));
13809566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltogmf, &idx_f));
138147c6ae99SBarry Smith
13829566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, &l_start_c, &m_c, &n_c, &p_c));
13839566063dSJacob 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));
1384565245c5SBarry Smith
13859566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dac, <ogmc));
13869566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltogmc, &idx_c));
138747c6ae99SBarry Smith
138847c6ae99SBarry Smith /*
138947c6ae99SBarry Smith Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios
139047c6ae99SBarry Smith for dimension 1 and 2 respectively.
139147c6ae99SBarry Smith Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1)
1392da81f932SPierre Jolivet and r_y*j and r_y*(j+1) will be grouped into the same coarse grid aggregate.
139347c6ae99SBarry Smith Each specific dof on the fine grid is mapped to one dof on the coarse grid.
139447c6ae99SBarry Smith */
139547c6ae99SBarry Smith
139647c6ae99SBarry Smith max_agg_size = (Mf / Mc + 1) * (Nf / Nc + 1) * (Pf / Pc + 1);
139747c6ae99SBarry Smith
139847c6ae99SBarry Smith /* create the matrix that will contain the restriction operator */
13999371c9d4SSatish 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));
140047c6ae99SBarry Smith
140147c6ae99SBarry Smith /* store nodes in the fine grid here */
14029566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(max_agg_size, &one_vec, max_agg_size, &fine_nodes));
140347c6ae99SBarry Smith for (i = 0; i < max_agg_size; i++) one_vec[i] = 1.0;
140447c6ae99SBarry Smith
140547c6ae99SBarry Smith /* loop over all coarse nodes */
140647c6ae99SBarry Smith for (l_c = l_start_c; l_c < l_start_c + p_c; l_c++) {
140747c6ae99SBarry Smith for (j_c = j_start_c; j_c < j_start_c + n_c; j_c++) {
140847c6ae99SBarry Smith for (i_c = i_start_c; i_c < i_start_c + m_c; i_c++) {
140947c6ae99SBarry Smith for (d = 0; d < dofc; d++) {
141047c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */
141147c6ae99SBarry 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;
141247c6ae99SBarry Smith
141347c6ae99SBarry Smith fn_idx = 0;
141447c6ae99SBarry Smith /* Corresponding fine points are all points (i_f, j_f, l_f) such that
141547c6ae99SBarry Smith i_c*Mf/Mc <= i_f < (i_c+1)*Mf/Mc
141647c6ae99SBarry Smith (same for other dimensions)
141747c6ae99SBarry Smith */
141847c6ae99SBarry Smith for (l = l_c * Pf / Pc; l < PetscMin((l_c + 1) * Pf / Pc, Pf); l++) {
141947c6ae99SBarry Smith for (j = j_c * Nf / Nc; j < PetscMin((j_c + 1) * Nf / Nc, Nf); j++) {
142047c6ae99SBarry Smith for (i = i_c * Mf / Mc; i < PetscMin((i_c + 1) * Mf / Mc, Mf); i++) {
142147c6ae99SBarry 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;
142247c6ae99SBarry Smith fn_idx++;
142347c6ae99SBarry Smith }
142447c6ae99SBarry Smith }
142547c6ae99SBarry Smith }
142647c6ae99SBarry Smith /* add all these points to one aggregate */
14279566063dSJacob Faibussowitsch PetscCall(MatSetValues(*rest, 1, &a, fn_idx, fine_nodes, one_vec, INSERT_VALUES));
142847c6ae99SBarry Smith }
142947c6ae99SBarry Smith }
143047c6ae99SBarry Smith }
143147c6ae99SBarry Smith }
14329566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogmf, &idx_f));
14339566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogmc, &idx_c));
14349566063dSJacob Faibussowitsch PetscCall(PetscFree2(one_vec, fine_nodes));
14359566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*rest, MAT_FINAL_ASSEMBLY));
14369566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*rest, MAT_FINAL_ASSEMBLY));
14373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
143847c6ae99SBarry Smith }
1439