xref: /petsc/src/dm/impls/da/dainterp.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
3aa219208SBarry Smith   Code for interpolating between grids represented by DMDAs
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
6d52bd9f3SBarry Smith /*
7d52bd9f3SBarry Smith       For linear elements there are two branches of code to compute the interpolation. They should compute the same results but may not. The "new version" does
8d52bd9f3SBarry Smith    not work for periodic domains, the old does. Change NEWVERSION to 1 to compile in the new version. Eventually when we are sure the two produce identical results
92adb9dcfSBarry Smith    we will remove/merge the new version. Based on current tests, these both produce the same results. We are leaving NEWVERSION for now in the code since some
102adb9dcfSBarry Smith    consider it cleaner, but old version is turned on since it handles periodic case.
11d52bd9f3SBarry Smith */
129314d695SBarry Smith #define NEWVERSION 0
1385afcc9aSBarry Smith 
14af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"   I*/
1547c6ae99SBarry Smith 
16fdc842d1SBarry Smith /*
17fdc842d1SBarry Smith    Since the interpolation uses MATMAIJ for dof > 0 we convert request for non-MATAIJ baseded matrices to MATAIJ.
18fdc842d1SBarry Smith    This is a bit of a hack, the reason for it is partially because -dm_mat_type defines the
19fdc842d1SBarry Smith    matrix type for both the operator matrices and the interpolation matrices so that users
20fdc842d1SBarry Smith    can select matrix types of base MATAIJ for accelerators
21fdc842d1SBarry Smith */
229371c9d4SSatish Balay static PetscErrorCode ConvertToAIJ(MatType intype, MatType *outtype) {
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;
33fdc842d1SBarry Smith       PetscFunctionReturn(0);
34fdc842d1SBarry Smith     }
35fdc842d1SBarry Smith   }
36fdc842d1SBarry Smith   PetscFunctionReturn(0);
37fdc842d1SBarry Smith }
38fdc842d1SBarry Smith 
399371c9d4SSatish Balay PetscErrorCode DMCreateInterpolation_DA_1D_Q1(DM dac, DM daf, Mat *A) {
408ea3bf28SBarry Smith   PetscInt               i, i_start, m_f, Mx;
418ea3bf28SBarry Smith   const PetscInt        *idx_f, *idx_c;
428ea3bf28SBarry Smith   PetscInt               m_ghost, m_ghost_c;
4347c6ae99SBarry Smith   PetscInt               row, col, i_start_ghost, mx, m_c, nc, ratio;
4447c6ae99SBarry Smith   PetscInt               i_c, i_start_c, i_start_ghost_c, cols[2], dof;
4585afcc9aSBarry Smith   PetscScalar            v[2], x;
4647c6ae99SBarry Smith   Mat                    mat;
47bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
48e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f, ltog_c;
49fdc842d1SBarry Smith   MatType                mattype;
5047c6ae99SBarry Smith 
5147c6ae99SBarry Smith   PetscFunctionBegin;
529566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, NULL, &Mx, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bx, NULL, NULL, NULL));
539566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, NULL, &mx, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
54bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
5547c6ae99SBarry Smith     ratio = mx / Mx;
5663a3b9bcSJacob 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);
5747c6ae99SBarry Smith   } else {
5847c6ae99SBarry Smith     ratio = (mx - 1) / (Mx - 1);
591dca8a05SBarry 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);
6047c6ae99SBarry Smith   }
6147c6ae99SBarry Smith 
629566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf, &i_start, NULL, NULL, &m_f, NULL, NULL));
639566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, NULL, NULL, &m_ghost, NULL, NULL));
649566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf, &ltog_f));
659566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
6647c6ae99SBarry Smith 
679566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac, &i_start_c, NULL, NULL, &m_c, NULL, NULL));
689566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, NULL, NULL, &m_ghost_c, NULL, NULL));
699566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac, &ltog_c));
709566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
7147c6ae99SBarry Smith 
7247c6ae99SBarry Smith   /* create interpolation matrix */
739566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dac), &mat));
74fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
75fdc842d1SBarry Smith   /*
76fdc842d1SBarry Smith      Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
77fdc842d1SBarry Smith      we don't want the original unconverted matrix copied to the GPU
78fdc842d1SBarry Smith    */
79*48a46eb9SPierre Jolivet   if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE));
80fdc842d1SBarry Smith #endif
819566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat, m_f, m_c, mx, Mx));
829566063dSJacob Faibussowitsch   PetscCall(ConvertToAIJ(dac->mattype, &mattype));
839566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat, mattype));
849566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat, 2, NULL));
859566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(mat, 2, NULL, 1, NULL));
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
8885afcc9aSBarry Smith   if (!NEWVERSION) {
8947c6ae99SBarry Smith     for (i = i_start; i < i_start + m_f; i++) {
9047c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
91e3fbd1f4SBarry Smith       row = idx_f[i - i_start_ghost];
9247c6ae99SBarry Smith 
9347c6ae99SBarry Smith       i_c = (i / ratio); /* coarse grid node to left of fine grid node */
9408401ef6SPierre Jolivet       PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
959371c9d4SSatish Balay                                           i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT,
969371c9d4SSatish Balay                  i_start, i_c, i_start_ghost_c);
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith       /*
9947c6ae99SBarry Smith        Only include those interpolation points that are truly
10047c6ae99SBarry Smith        nonzero. Note this is very important for final grid lines
10147c6ae99SBarry Smith        in x direction; since they have no right neighbor
10247c6ae99SBarry Smith        */
1036712e2f1SBarry Smith       x        = ((PetscReal)(i - i_c * ratio)) / ((PetscReal)ratio);
10447c6ae99SBarry Smith       nc       = 0;
10547c6ae99SBarry Smith       /* one left and below; or we are right on it */
106e3fbd1f4SBarry Smith       col      = (i_c - i_start_ghost_c);
107e3fbd1f4SBarry Smith       cols[nc] = idx_c[col];
10847c6ae99SBarry Smith       v[nc++]  = -x + 1.0;
10947c6ae99SBarry Smith       /* one right? */
11047c6ae99SBarry Smith       if (i_c * ratio != i) {
111e3fbd1f4SBarry Smith         cols[nc] = idx_c[col + 1];
11247c6ae99SBarry Smith         v[nc++]  = x;
11347c6ae99SBarry Smith       }
1149566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES));
11547c6ae99SBarry Smith     }
116b3bd94feSDave May 
117b3bd94feSDave May   } else {
118b3bd94feSDave May     PetscScalar *xi;
119b3bd94feSDave May     PetscInt     li, nxi, n;
120b3bd94feSDave May     PetscScalar  Ni[2];
121b3bd94feSDave May 
122b3bd94feSDave May     /* compute local coordinate arrays */
123b3bd94feSDave May     nxi = ratio + 1;
1249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nxi, &xi));
1259371c9d4SSatish Balay     for (li = 0; li < nxi; li++) { xi[li] = -1.0 + (PetscScalar)li * (2.0 / (PetscScalar)(nxi - 1)); }
126b3bd94feSDave May 
127b3bd94feSDave May     for (i = i_start; i < i_start + m_f; i++) {
128b3bd94feSDave May       /* convert to local "natural" numbering and then to PETSc global numbering */
129e3fbd1f4SBarry Smith       row = idx_f[(i - i_start_ghost)];
130b3bd94feSDave May 
131b3bd94feSDave May       i_c = (i / ratio); /* coarse grid node to left of fine grid node */
13208401ef6SPierre Jolivet       PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
1339371c9d4SSatish Balay                                           i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT,
1349371c9d4SSatish Balay                  i_start, i_c, i_start_ghost_c);
135b3bd94feSDave May 
136b3bd94feSDave May       /* remainders */
137b3bd94feSDave May       li = i - ratio * (i / ratio);
1388865f1eaSKarl Rupp       if (i == mx - 1) li = nxi - 1;
139b3bd94feSDave May 
140b3bd94feSDave May       /* corners */
141e3fbd1f4SBarry Smith       col     = (i_c - i_start_ghost_c);
142e3fbd1f4SBarry Smith       cols[0] = idx_c[col];
143b3bd94feSDave May       Ni[0]   = 1.0;
144b3bd94feSDave May       if ((li == 0) || (li == nxi - 1)) {
1459566063dSJacob Faibussowitsch         PetscCall(MatSetValue(mat, row, cols[0], Ni[0], INSERT_VALUES));
146b3bd94feSDave May         continue;
147b3bd94feSDave May       }
148b3bd94feSDave May 
149b3bd94feSDave May       /* edges + interior */
150b3bd94feSDave May       /* remainders */
1518865f1eaSKarl Rupp       if (i == mx - 1) i_c--;
152b3bd94feSDave May 
153e3fbd1f4SBarry Smith       col     = (i_c - i_start_ghost_c);
154e3fbd1f4SBarry Smith       cols[0] = idx_c[col]; /* one left and below; or we are right on it */
155e3fbd1f4SBarry Smith       cols[1] = idx_c[col + 1];
156b3bd94feSDave May 
157b3bd94feSDave May       Ni[0] = 0.5 * (1.0 - xi[li]);
158b3bd94feSDave May       Ni[1] = 0.5 * (1.0 + xi[li]);
159b3bd94feSDave May       for (n = 0; n < 2; n++) {
1608865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[n]) < 1.0e-32) cols[n] = -1;
161b3bd94feSDave May       }
1629566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat, 1, &row, 2, cols, Ni, INSERT_VALUES));
163b3bd94feSDave May     }
1649566063dSJacob Faibussowitsch     PetscCall(PetscFree(xi));
165b3bd94feSDave May   }
1669566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
1679566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
1689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
1699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
1709566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(mat, dof, A));
1719566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
17247c6ae99SBarry Smith   PetscFunctionReturn(0);
17347c6ae99SBarry Smith }
17447c6ae99SBarry Smith 
1759371c9d4SSatish Balay PetscErrorCode DMCreateInterpolation_DA_1D_Q0(DM dac, DM daf, Mat *A) {
1768ea3bf28SBarry Smith   PetscInt               i, i_start, m_f, Mx;
1778ea3bf28SBarry Smith   const PetscInt        *idx_f, *idx_c;
178e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f, ltog_c;
1798ea3bf28SBarry Smith   PetscInt               m_ghost, m_ghost_c;
18047c6ae99SBarry Smith   PetscInt               row, col, i_start_ghost, mx, m_c, nc, ratio;
18147c6ae99SBarry Smith   PetscInt               i_c, i_start_c, i_start_ghost_c, cols[2], dof;
18247c6ae99SBarry Smith   PetscScalar            v[2], x;
18347c6ae99SBarry Smith   Mat                    mat;
184bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
185fdc842d1SBarry Smith   MatType                mattype;
18647c6ae99SBarry Smith 
18747c6ae99SBarry Smith   PetscFunctionBegin;
1889566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, NULL, &Mx, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bx, NULL, NULL, NULL));
1899566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, NULL, &mx, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
190bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
19163a3b9bcSJacob Faibussowitsch     PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx);
19247c6ae99SBarry Smith     ratio = mx / Mx;
19363a3b9bcSJacob 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);
19447c6ae99SBarry Smith   } else {
19563a3b9bcSJacob 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);
19647c6ae99SBarry Smith     ratio = (mx - 1) / (Mx - 1);
1971dca8a05SBarry 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);
19847c6ae99SBarry Smith   }
19947c6ae99SBarry Smith 
2009566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf, &i_start, NULL, NULL, &m_f, NULL, NULL));
2019566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, NULL, NULL, &m_ghost, NULL, NULL));
2029566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf, &ltog_f));
2039566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
20447c6ae99SBarry Smith 
2059566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac, &i_start_c, NULL, NULL, &m_c, NULL, NULL));
2069566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, NULL, NULL, &m_ghost_c, NULL, NULL));
2079566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac, &ltog_c));
2089566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
20947c6ae99SBarry Smith 
21047c6ae99SBarry Smith   /* create interpolation matrix */
2119566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dac), &mat));
212fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
213fdc842d1SBarry Smith   /*
214fdc842d1SBarry Smith      Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
215fdc842d1SBarry Smith      we don't want the original unconverted matrix copied to the GPU
216fdc842d1SBarry Smith    */
217*48a46eb9SPierre Jolivet   if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE));
218fdc842d1SBarry Smith #endif
2199566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat, m_f, m_c, mx, Mx));
2209566063dSJacob Faibussowitsch   PetscCall(ConvertToAIJ(dac->mattype, &mattype));
2219566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat, mattype));
2229566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat, 2, NULL));
2239566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(mat, 2, NULL, 0, NULL));
22447c6ae99SBarry Smith 
22547c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
22647c6ae99SBarry Smith   for (i = i_start; i < i_start + m_f; i++) {
22747c6ae99SBarry Smith     /* convert to local "natural" numbering and then to PETSc global numbering */
228e3fbd1f4SBarry Smith     row = idx_f[(i - i_start_ghost)];
22947c6ae99SBarry Smith 
23047c6ae99SBarry Smith     i_c = (i / ratio); /* coarse grid node to left of fine grid node */
23147c6ae99SBarry Smith 
23247c6ae99SBarry Smith     /*
23347c6ae99SBarry Smith          Only include those interpolation points that are truly
23447c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
23547c6ae99SBarry Smith          in x direction; since they have no right neighbor
23647c6ae99SBarry Smith     */
2376712e2f1SBarry Smith     x        = ((PetscReal)(i - i_c * ratio)) / ((PetscReal)ratio);
23847c6ae99SBarry Smith     nc       = 0;
23947c6ae99SBarry Smith     /* one left and below; or we are right on it */
240e3fbd1f4SBarry Smith     col      = (i_c - i_start_ghost_c);
241e3fbd1f4SBarry Smith     cols[nc] = idx_c[col];
24247c6ae99SBarry Smith     v[nc++]  = -x + 1.0;
24347c6ae99SBarry Smith     /* one right? */
24447c6ae99SBarry Smith     if (i_c * ratio != i) {
245e3fbd1f4SBarry Smith       cols[nc] = idx_c[col + 1];
24647c6ae99SBarry Smith       v[nc++]  = x;
24747c6ae99SBarry Smith     }
2489566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES));
24947c6ae99SBarry Smith   }
2509566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
2519566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
2529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
2539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
2549566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(mat, dof, A));
2559566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
2569566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(5.0 * m_f));
25747c6ae99SBarry Smith   PetscFunctionReturn(0);
25847c6ae99SBarry Smith }
25947c6ae99SBarry Smith 
2609371c9d4SSatish Balay PetscErrorCode DMCreateInterpolation_DA_2D_Q1(DM dac, DM daf, Mat *A) {
2618ea3bf28SBarry Smith   PetscInt               i, j, i_start, j_start, m_f, n_f, Mx, My, dof;
2628ea3bf28SBarry Smith   const PetscInt        *idx_c, *idx_f;
263e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f, ltog_c;
2648ea3bf28SBarry Smith   PetscInt               m_ghost, n_ghost, m_ghost_c, n_ghost_c, *dnz, *onz;
26547c6ae99SBarry Smith   PetscInt               row, col, i_start_ghost, j_start_ghost, cols[4], mx, m_c, my, nc, ratioi, ratioj;
26647c6ae99SBarry 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;
26747c6ae99SBarry Smith   PetscMPIInt            size_c, size_f, rank_f;
26847c6ae99SBarry Smith   PetscScalar            v[4], x, y;
26947c6ae99SBarry Smith   Mat                    mat;
270bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx, by;
271fdc842d1SBarry Smith   MatType                mattype;
27247c6ae99SBarry Smith 
27347c6ae99SBarry Smith   PetscFunctionBegin;
2749566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, NULL, NULL, NULL, NULL, NULL, NULL, &bx, &by, NULL, NULL));
2759566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
276bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
27763a3b9bcSJacob Faibussowitsch     PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx);
27847c6ae99SBarry Smith     ratioi = mx / Mx;
27963a3b9bcSJacob 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);
28047c6ae99SBarry Smith   } else {
28163a3b9bcSJacob 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);
28247c6ae99SBarry Smith     ratioi = (mx - 1) / (Mx - 1);
2831dca8a05SBarry 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);
28447c6ae99SBarry Smith   }
285bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) {
28663a3b9bcSJacob Faibussowitsch     PetscCheck(My, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be positive", My);
28747c6ae99SBarry Smith     ratioj = my / My;
28863a3b9bcSJacob 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);
28947c6ae99SBarry Smith   } else {
29063a3b9bcSJacob 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);
29147c6ae99SBarry Smith     ratioj = (my - 1) / (My - 1);
2921dca8a05SBarry 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);
29347c6ae99SBarry Smith   }
29447c6ae99SBarry Smith 
2959566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf, &i_start, &j_start, NULL, &m_f, &n_f, NULL));
2969566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, NULL, &m_ghost, &n_ghost, NULL));
2979566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf, &ltog_f));
2989566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
29947c6ae99SBarry Smith 
3009566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, NULL, &m_c, &n_c, NULL));
3019566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, NULL, &m_ghost_c, &n_ghost_c, NULL));
3029566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac, &ltog_c));
3039566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
30447c6ae99SBarry Smith 
30547c6ae99SBarry Smith   /*
306aa219208SBarry Smith    Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
30747c6ae99SBarry Smith    The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
30847c6ae99SBarry Smith    processors). It's effective length is hence 4 times its normal length, this is
30947c6ae99SBarry Smith    why the col_scale is multiplied by the interpolation matrix column sizes.
31047c6ae99SBarry Smith    sol_shift allows each set of 1/4 processors do its own interpolation using ITS
31147c6ae99SBarry Smith    copy of the coarse vector. A bit of a hack but you do better.
31247c6ae99SBarry Smith 
31347c6ae99SBarry Smith    In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
31447c6ae99SBarry Smith    */
3159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dac), &size_c));
3169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)daf), &size_f));
3179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)daf), &rank_f));
31847c6ae99SBarry Smith   col_scale = size_f / size_c;
31947c6ae99SBarry Smith   col_shift = Mx * My * (rank_f / size_c);
32047c6ae99SBarry Smith 
321d0609cedSBarry Smith   MatPreallocateBegin(PetscObjectComm((PetscObject)daf), m_f * n_f, col_scale * m_c * n_c, dnz, onz);
32247c6ae99SBarry Smith   for (j = j_start; j < j_start + n_f; j++) {
32347c6ae99SBarry Smith     for (i = i_start; i < i_start + m_f; i++) {
32447c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
325e3fbd1f4SBarry Smith       row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
32647c6ae99SBarry Smith 
32747c6ae99SBarry Smith       i_c = (i / ratioi); /* coarse grid node to left of fine grid node */
32847c6ae99SBarry Smith       j_c = (j / ratioj); /* coarse grid node below fine grid node */
32947c6ae99SBarry Smith 
33008401ef6SPierre Jolivet       PetscCheck(j_c >= j_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
3319371c9d4SSatish Balay                                           j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT,
3329371c9d4SSatish Balay                  j_start, j_c, j_start_ghost_c);
33308401ef6SPierre Jolivet       PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
3349371c9d4SSatish Balay                                           i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT,
3359371c9d4SSatish Balay                  i_start, i_c, i_start_ghost_c);
33647c6ae99SBarry Smith 
33747c6ae99SBarry Smith       /*
33847c6ae99SBarry Smith        Only include those interpolation points that are truly
33947c6ae99SBarry Smith        nonzero. Note this is very important for final grid lines
34047c6ae99SBarry Smith        in x and y directions; since they have no right/top neighbors
34147c6ae99SBarry Smith        */
34247c6ae99SBarry Smith       nc         = 0;
34347c6ae99SBarry Smith       /* one left and below; or we are right on it */
344e3fbd1f4SBarry Smith       col        = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c));
345e3fbd1f4SBarry Smith       cols[nc++] = col_shift + idx_c[col];
34647c6ae99SBarry Smith       /* one right and below */
347e3fbd1f4SBarry Smith       if (i_c * ratioi != i) cols[nc++] = col_shift + idx_c[col + 1];
34847c6ae99SBarry Smith       /* one left and above */
349e3fbd1f4SBarry Smith       if (j_c * ratioj != j) cols[nc++] = col_shift + idx_c[col + m_ghost_c];
35047c6ae99SBarry Smith       /* one right and above */
351e3fbd1f4SBarry Smith       if (i_c * ratioi != i && j_c * ratioj != j) cols[nc++] = col_shift + idx_c[col + (m_ghost_c + 1)];
3529566063dSJacob Faibussowitsch       PetscCall(MatPreallocateSet(row, nc, cols, dnz, onz));
35347c6ae99SBarry Smith     }
35447c6ae99SBarry Smith   }
3559566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)daf), &mat));
356fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
357fdc842d1SBarry Smith   /*
358fdc842d1SBarry Smith      Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
359fdc842d1SBarry Smith      we don't want the original unconverted matrix copied to the GPU
360fdc842d1SBarry Smith   */
361*48a46eb9SPierre Jolivet   if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE));
362fdc842d1SBarry Smith #endif
3639566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat, m_f * n_f, col_scale * m_c * n_c, mx * my, col_scale * Mx * My));
3649566063dSJacob Faibussowitsch   PetscCall(ConvertToAIJ(dac->mattype, &mattype));
3659566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat, mattype));
3669566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat, 0, dnz));
3679566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(mat, 0, dnz, 0, onz));
368d0609cedSBarry Smith   MatPreallocateEnd(dnz, onz);
36947c6ae99SBarry Smith 
37047c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
37185afcc9aSBarry Smith   if (!NEWVERSION) {
37247c6ae99SBarry Smith     for (j = j_start; j < j_start + n_f; j++) {
37347c6ae99SBarry Smith       for (i = i_start; i < i_start + m_f; i++) {
37447c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
375e3fbd1f4SBarry Smith         row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
37647c6ae99SBarry Smith 
37747c6ae99SBarry Smith         i_c = (i / ratioi); /* coarse grid node to left of fine grid node */
37847c6ae99SBarry Smith         j_c = (j / ratioj); /* coarse grid node below fine grid node */
37947c6ae99SBarry Smith 
38047c6ae99SBarry Smith         /*
38147c6ae99SBarry Smith          Only include those interpolation points that are truly
38247c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
38347c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
38447c6ae99SBarry Smith          */
3856712e2f1SBarry Smith         x = ((PetscReal)(i - i_c * ratioi)) / ((PetscReal)ratioi);
3866712e2f1SBarry Smith         y = ((PetscReal)(j - j_c * ratioj)) / ((PetscReal)ratioj);
387b3bd94feSDave May 
38847c6ae99SBarry Smith         nc       = 0;
38947c6ae99SBarry Smith         /* one left and below; or we are right on it */
390e3fbd1f4SBarry Smith         col      = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c));
391e3fbd1f4SBarry Smith         cols[nc] = col_shift + idx_c[col];
39247c6ae99SBarry Smith         v[nc++]  = x * y - x - y + 1.0;
39347c6ae99SBarry Smith         /* one right and below */
39447c6ae99SBarry Smith         if (i_c * ratioi != i) {
395e3fbd1f4SBarry Smith           cols[nc] = col_shift + idx_c[col + 1];
39647c6ae99SBarry Smith           v[nc++]  = -x * y + x;
39747c6ae99SBarry Smith         }
39847c6ae99SBarry Smith         /* one left and above */
39947c6ae99SBarry Smith         if (j_c * ratioj != j) {
400e3fbd1f4SBarry Smith           cols[nc] = col_shift + idx_c[col + m_ghost_c];
40147c6ae99SBarry Smith           v[nc++]  = -x * y + y;
40247c6ae99SBarry Smith         }
40347c6ae99SBarry Smith         /* one right and above */
40447c6ae99SBarry Smith         if (j_c * ratioj != j && i_c * ratioi != i) {
405e3fbd1f4SBarry Smith           cols[nc] = col_shift + idx_c[col + (m_ghost_c + 1)];
40647c6ae99SBarry Smith           v[nc++]  = x * y;
40747c6ae99SBarry Smith         }
4089566063dSJacob Faibussowitsch         PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES));
40947c6ae99SBarry Smith       }
41047c6ae99SBarry Smith     }
411b3bd94feSDave May 
412b3bd94feSDave May   } else {
413b3bd94feSDave May     PetscScalar  Ni[4];
414b3bd94feSDave May     PetscScalar *xi, *eta;
415b3bd94feSDave May     PetscInt     li, nxi, lj, neta;
416b3bd94feSDave May 
417b3bd94feSDave May     /* compute local coordinate arrays */
418b3bd94feSDave May     nxi  = ratioi + 1;
419b3bd94feSDave May     neta = ratioj + 1;
4209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nxi, &xi));
4219566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(neta, &eta));
4229371c9d4SSatish Balay     for (li = 0; li < nxi; li++) { xi[li] = -1.0 + (PetscScalar)li * (2.0 / (PetscScalar)(nxi - 1)); }
4239371c9d4SSatish Balay     for (lj = 0; lj < neta; lj++) { eta[lj] = -1.0 + (PetscScalar)lj * (2.0 / (PetscScalar)(neta - 1)); }
424b3bd94feSDave May 
425b3bd94feSDave May     /* loop over local fine grid nodes setting interpolation for those*/
426b3bd94feSDave May     for (j = j_start; j < j_start + n_f; j++) {
427b3bd94feSDave May       for (i = i_start; i < i_start + m_f; i++) {
428b3bd94feSDave May         /* convert to local "natural" numbering and then to PETSc global numbering */
429e3fbd1f4SBarry Smith         row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
430b3bd94feSDave May 
431b3bd94feSDave May         i_c = (i / ratioi); /* coarse grid node to left of fine grid node */
432b3bd94feSDave May         j_c = (j / ratioj); /* coarse grid node below fine grid node */
433b3bd94feSDave May 
434b3bd94feSDave May         /* remainders */
435b3bd94feSDave May         li = i - ratioi * (i / ratioi);
4368865f1eaSKarl Rupp         if (i == mx - 1) li = nxi - 1;
437b3bd94feSDave May         lj = j - ratioj * (j / ratioj);
4388865f1eaSKarl Rupp         if (j == my - 1) lj = neta - 1;
439b3bd94feSDave May 
440b3bd94feSDave May         /* corners */
441e3fbd1f4SBarry Smith         col     = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c));
442e3fbd1f4SBarry Smith         cols[0] = col_shift + idx_c[col]; /* left, below */
443b3bd94feSDave May         Ni[0]   = 1.0;
444b3bd94feSDave May         if ((li == 0) || (li == nxi - 1)) {
445b3bd94feSDave May           if ((lj == 0) || (lj == neta - 1)) {
4469566063dSJacob Faibussowitsch             PetscCall(MatSetValue(mat, row, cols[0], Ni[0], INSERT_VALUES));
447b3bd94feSDave May             continue;
448b3bd94feSDave May           }
449b3bd94feSDave May         }
450b3bd94feSDave May 
451b3bd94feSDave May         /* edges + interior */
452b3bd94feSDave May         /* remainders */
4538865f1eaSKarl Rupp         if (i == mx - 1) i_c--;
4548865f1eaSKarl Rupp         if (j == my - 1) j_c--;
455b3bd94feSDave May 
456e3fbd1f4SBarry Smith         col     = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c));
457e3fbd1f4SBarry Smith         cols[0] = col_shift + idx_c[col];                   /* left, below */
458e3fbd1f4SBarry Smith         cols[1] = col_shift + idx_c[col + 1];               /* right, below */
459e3fbd1f4SBarry Smith         cols[2] = col_shift + idx_c[col + m_ghost_c];       /* left, above */
460e3fbd1f4SBarry Smith         cols[3] = col_shift + idx_c[col + (m_ghost_c + 1)]; /* right, above */
461b3bd94feSDave May 
462b3bd94feSDave May         Ni[0] = 0.25 * (1.0 - xi[li]) * (1.0 - eta[lj]);
463b3bd94feSDave May         Ni[1] = 0.25 * (1.0 + xi[li]) * (1.0 - eta[lj]);
464b3bd94feSDave May         Ni[2] = 0.25 * (1.0 - xi[li]) * (1.0 + eta[lj]);
465b3bd94feSDave May         Ni[3] = 0.25 * (1.0 + xi[li]) * (1.0 + eta[lj]);
466b3bd94feSDave May 
467b3bd94feSDave May         nc = 0;
4688865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[0]) < 1.0e-32) cols[0] = -1;
4698865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[1]) < 1.0e-32) cols[1] = -1;
4708865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[2]) < 1.0e-32) cols[2] = -1;
4718865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[3]) < 1.0e-32) cols[3] = -1;
472b3bd94feSDave May 
4739566063dSJacob Faibussowitsch         PetscCall(MatSetValues(mat, 1, &row, 4, cols, Ni, INSERT_VALUES));
474b3bd94feSDave May       }
475b3bd94feSDave May     }
4769566063dSJacob Faibussowitsch     PetscCall(PetscFree(xi));
4779566063dSJacob Faibussowitsch     PetscCall(PetscFree(eta));
478b3bd94feSDave May   }
4799566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
4809566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
4819566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
4829566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
4839566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(mat, dof, A));
4849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
48547c6ae99SBarry Smith   PetscFunctionReturn(0);
48647c6ae99SBarry Smith }
48747c6ae99SBarry Smith 
48847c6ae99SBarry Smith /*
48947c6ae99SBarry Smith        Contributed by Andrei Draganescu <aidraga@sandia.gov>
49047c6ae99SBarry Smith */
4919371c9d4SSatish Balay PetscErrorCode DMCreateInterpolation_DA_2D_Q0(DM dac, DM daf, Mat *A) {
4928ea3bf28SBarry Smith   PetscInt               i, j, i_start, j_start, m_f, n_f, Mx, My, dof;
4938ea3bf28SBarry Smith   const PetscInt        *idx_c, *idx_f;
494e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f, ltog_c;
4958ea3bf28SBarry Smith   PetscInt               m_ghost, n_ghost, m_ghost_c, n_ghost_c, *dnz, *onz;
49647c6ae99SBarry Smith   PetscInt               row, col, i_start_ghost, j_start_ghost, cols[4], mx, m_c, my, nc, ratioi, ratioj;
49747c6ae99SBarry 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;
49847c6ae99SBarry Smith   PetscMPIInt            size_c, size_f, rank_f;
49947c6ae99SBarry Smith   PetscScalar            v[4];
50047c6ae99SBarry Smith   Mat                    mat;
501bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx, by;
502fdc842d1SBarry Smith   MatType                mattype;
50347c6ae99SBarry Smith 
50447c6ae99SBarry Smith   PetscFunctionBegin;
5059566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, NULL, NULL, NULL, NULL, NULL, NULL, &bx, &by, NULL, NULL));
5069566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
50763a3b9bcSJacob Faibussowitsch   PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx);
50863a3b9bcSJacob Faibussowitsch   PetscCheck(My, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be positive", My);
50947c6ae99SBarry Smith   ratioi = mx / Mx;
51047c6ae99SBarry Smith   ratioj = my / My;
51108401ef6SPierre Jolivet   PetscCheck(ratioi * Mx == mx, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in x");
51208401ef6SPierre Jolivet   PetscCheck(ratioj * My == my, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in y");
51308401ef6SPierre Jolivet   PetscCheck(ratioi == 2, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in x must be 2");
51408401ef6SPierre Jolivet   PetscCheck(ratioj == 2, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in y must be 2");
51547c6ae99SBarry Smith 
5169566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf, &i_start, &j_start, NULL, &m_f, &n_f, NULL));
5179566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, NULL, &m_ghost, &n_ghost, NULL));
5189566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf, &ltog_f));
5199566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
52047c6ae99SBarry Smith 
5219566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, NULL, &m_c, &n_c, NULL));
5229566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, NULL, &m_ghost_c, &n_ghost_c, NULL));
5239566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac, &ltog_c));
5249566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
52547c6ae99SBarry Smith 
52647c6ae99SBarry Smith   /*
527aa219208SBarry Smith      Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
52847c6ae99SBarry Smith      The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
52947c6ae99SBarry Smith      processors). It's effective length is hence 4 times its normal length, this is
53047c6ae99SBarry Smith      why the col_scale is multiplied by the interpolation matrix column sizes.
53147c6ae99SBarry Smith      sol_shift allows each set of 1/4 processors do its own interpolation using ITS
53247c6ae99SBarry Smith      copy of the coarse vector. A bit of a hack but you do better.
53347c6ae99SBarry Smith 
53447c6ae99SBarry Smith      In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
53547c6ae99SBarry Smith   */
5369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dac), &size_c));
5379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)daf), &size_f));
5389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)daf), &rank_f));
53947c6ae99SBarry Smith   col_scale = size_f / size_c;
54047c6ae99SBarry Smith   col_shift = Mx * My * (rank_f / size_c);
54147c6ae99SBarry Smith 
542d0609cedSBarry Smith   MatPreallocateBegin(PetscObjectComm((PetscObject)daf), m_f * n_f, col_scale * m_c * n_c, dnz, onz);
54347c6ae99SBarry Smith   for (j = j_start; j < j_start + n_f; j++) {
54447c6ae99SBarry Smith     for (i = i_start; i < i_start + m_f; i++) {
54547c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
546e3fbd1f4SBarry Smith       row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
54747c6ae99SBarry Smith 
54847c6ae99SBarry Smith       i_c = (i / ratioi); /* coarse grid node to left of fine grid node */
54947c6ae99SBarry Smith       j_c = (j / ratioj); /* coarse grid node below fine grid node */
55047c6ae99SBarry Smith 
55108401ef6SPierre Jolivet       PetscCheck(j_c >= j_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
5529371c9d4SSatish Balay     j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT,
5539371c9d4SSatish Balay                  j_start, j_c, j_start_ghost_c);
55408401ef6SPierre Jolivet       PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
5559371c9d4SSatish Balay     i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT,
5569371c9d4SSatish Balay                  i_start, i_c, i_start_ghost_c);
55747c6ae99SBarry Smith 
55847c6ae99SBarry Smith       /*
55947c6ae99SBarry Smith          Only include those interpolation points that are truly
56047c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
56147c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
56247c6ae99SBarry Smith       */
56347c6ae99SBarry Smith       nc         = 0;
56447c6ae99SBarry Smith       /* one left and below; or we are right on it */
565e3fbd1f4SBarry Smith       col        = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c));
566e3fbd1f4SBarry Smith       cols[nc++] = col_shift + idx_c[col];
5679566063dSJacob Faibussowitsch       PetscCall(MatPreallocateSet(row, nc, cols, dnz, onz));
56847c6ae99SBarry Smith     }
56947c6ae99SBarry Smith   }
5709566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)daf), &mat));
571fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
572fdc842d1SBarry Smith   /*
573fdc842d1SBarry Smith      Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
574fdc842d1SBarry Smith      we don't want the original unconverted matrix copied to the GPU
575fdc842d1SBarry Smith   */
576*48a46eb9SPierre Jolivet   if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE));
577fdc842d1SBarry Smith #endif
5789566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat, m_f * n_f, col_scale * m_c * n_c, mx * my, col_scale * Mx * My));
5799566063dSJacob Faibussowitsch   PetscCall(ConvertToAIJ(dac->mattype, &mattype));
5809566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat, mattype));
5819566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat, 0, dnz));
5829566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(mat, 0, dnz, 0, onz));
583d0609cedSBarry Smith   MatPreallocateEnd(dnz, onz);
58447c6ae99SBarry Smith 
58547c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
58647c6ae99SBarry Smith   for (j = j_start; j < j_start + n_f; j++) {
58747c6ae99SBarry Smith     for (i = i_start; i < i_start + m_f; i++) {
58847c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
589e3fbd1f4SBarry Smith       row = idx_f[(m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
59047c6ae99SBarry Smith 
59147c6ae99SBarry Smith       i_c      = (i / ratioi); /* coarse grid node to left of fine grid node */
59247c6ae99SBarry Smith       j_c      = (j / ratioj); /* coarse grid node below fine grid node */
59347c6ae99SBarry Smith       nc       = 0;
59447c6ae99SBarry Smith       /* one left and below; or we are right on it */
595e3fbd1f4SBarry Smith       col      = (m_ghost_c * (j_c - j_start_ghost_c) + (i_c - i_start_ghost_c));
596e3fbd1f4SBarry Smith       cols[nc] = col_shift + idx_c[col];
59747c6ae99SBarry Smith       v[nc++]  = 1.0;
59847c6ae99SBarry Smith 
5999566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES));
60047c6ae99SBarry Smith     }
60147c6ae99SBarry Smith   }
6029566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
6039566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
6049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
6059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
6069566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(mat, dof, A));
6079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
6089566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(13.0 * m_f * n_f));
60947c6ae99SBarry Smith   PetscFunctionReturn(0);
61047c6ae99SBarry Smith }
61147c6ae99SBarry Smith 
61247c6ae99SBarry Smith /*
61347c6ae99SBarry Smith        Contributed by Jianming Yang <jianming-yang@uiowa.edu>
61447c6ae99SBarry Smith */
6159371c9d4SSatish Balay PetscErrorCode DMCreateInterpolation_DA_3D_Q0(DM dac, DM daf, Mat *A) {
6168ea3bf28SBarry Smith   PetscInt               i, j, l, i_start, j_start, l_start, m_f, n_f, p_f, Mx, My, Mz, dof;
6178ea3bf28SBarry Smith   const PetscInt        *idx_c, *idx_f;
618e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f, ltog_c;
6198ea3bf28SBarry Smith   PetscInt               m_ghost, n_ghost, p_ghost, m_ghost_c, n_ghost_c, p_ghost_c, nc, *dnz, *onz;
62047c6ae99SBarry 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;
62147c6ae99SBarry 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;
62247c6ae99SBarry Smith   PetscMPIInt            size_c, size_f, rank_f;
62347c6ae99SBarry Smith   PetscScalar            v[8];
62447c6ae99SBarry Smith   Mat                    mat;
625bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx, by, bz;
626fdc842d1SBarry Smith   MatType                mattype;
62747c6ae99SBarry Smith 
62847c6ae99SBarry Smith   PetscFunctionBegin;
6299566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, &Mz, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL));
63063a3b9bcSJacob Faibussowitsch   PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx);
63163a3b9bcSJacob Faibussowitsch   PetscCheck(My, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be positive", My);
63263a3b9bcSJacob Faibussowitsch   PetscCheck(Mz, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of z coarse grid points %" PetscInt_FMT " must be positive", Mz);
6339566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, &mz, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
63447c6ae99SBarry Smith   ratioi = mx / Mx;
63547c6ae99SBarry Smith   ratioj = my / My;
63647c6ae99SBarry Smith   ratiol = mz / Mz;
63708401ef6SPierre Jolivet   PetscCheck(ratioi * Mx == mx, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in x");
63808401ef6SPierre Jolivet   PetscCheck(ratioj * My == my, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in y");
63908401ef6SPierre Jolivet   PetscCheck(ratiol * Mz == mz, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Fine grid points must be multiple of coarse grid points in z");
6401dca8a05SBarry Smith   PetscCheck(ratioi == 2 || ratioi == 1, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in x must be 1 or 2");
6411dca8a05SBarry Smith   PetscCheck(ratioj == 2 || ratioj == 1, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in y must be 1 or 2");
6421dca8a05SBarry Smith   PetscCheck(ratiol == 2 || ratiol == 1, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_WRONG, "Coarsening factor in z must be 1 or 2");
64347c6ae99SBarry Smith 
6449566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf, &i_start, &j_start, &l_start, &m_f, &n_f, &p_f));
6459566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, &l_start_ghost, &m_ghost, &n_ghost, &p_ghost));
6469566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf, &ltog_f));
6479566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
64847c6ae99SBarry Smith 
6499566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, &l_start_c, &m_c, &n_c, &p_c));
6509566063dSJacob 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));
6519566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac, &ltog_c));
6529566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
653e3fbd1f4SBarry Smith 
65447c6ae99SBarry Smith   /*
655aa219208SBarry Smith      Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
65647c6ae99SBarry Smith      The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
65747c6ae99SBarry Smith      processors). It's effective length is hence 4 times its normal length, this is
65847c6ae99SBarry Smith      why the col_scale is multiplied by the interpolation matrix column sizes.
65947c6ae99SBarry Smith      sol_shift allows each set of 1/4 processors do its own interpolation using ITS
66047c6ae99SBarry Smith      copy of the coarse vector. A bit of a hack but you do better.
66147c6ae99SBarry Smith 
66247c6ae99SBarry Smith      In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
66347c6ae99SBarry Smith   */
6649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dac), &size_c));
6659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)daf), &size_f));
6669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)daf), &rank_f));
66747c6ae99SBarry Smith   col_scale = size_f / size_c;
66847c6ae99SBarry Smith   col_shift = Mx * My * Mz * (rank_f / size_c);
66947c6ae99SBarry Smith 
670d0609cedSBarry Smith   MatPreallocateBegin(PetscObjectComm((PetscObject)daf), m_f * n_f * p_f, col_scale * m_c * n_c * p_c, dnz, onz);
67147c6ae99SBarry Smith   for (l = l_start; l < l_start + p_f; l++) {
67247c6ae99SBarry Smith     for (j = j_start; j < j_start + n_f; j++) {
67347c6ae99SBarry Smith       for (i = i_start; i < i_start + m_f; i++) {
67447c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
675e3fbd1f4SBarry Smith         row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
67647c6ae99SBarry Smith 
67747c6ae99SBarry Smith         i_c = (i / ratioi); /* coarse grid node to left of fine grid node */
67847c6ae99SBarry Smith         j_c = (j / ratioj); /* coarse grid node below fine grid node */
67947c6ae99SBarry Smith         l_c = (l / ratiol);
68047c6ae99SBarry Smith 
68108401ef6SPierre Jolivet         PetscCheck(l_c >= l_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
6829371c9d4SSatish Balay     l_start %" PetscInt_FMT " l_c %" PetscInt_FMT " l_start_ghost_c %" PetscInt_FMT,
6839371c9d4SSatish Balay                    l_start, l_c, l_start_ghost_c);
68408401ef6SPierre Jolivet         PetscCheck(j_c >= j_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
6859371c9d4SSatish Balay     j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT,
6869371c9d4SSatish Balay                    j_start, j_c, j_start_ghost_c);
68708401ef6SPierre Jolivet         PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
6889371c9d4SSatish Balay     i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT,
6899371c9d4SSatish Balay                    i_start, i_c, i_start_ghost_c);
69047c6ae99SBarry Smith 
69147c6ae99SBarry Smith         /*
69247c6ae99SBarry Smith            Only include those interpolation points that are truly
69347c6ae99SBarry Smith            nonzero. Note this is very important for final grid lines
69447c6ae99SBarry Smith            in x and y directions; since they have no right/top neighbors
69547c6ae99SBarry Smith         */
69647c6ae99SBarry Smith         nc         = 0;
69747c6ae99SBarry Smith         /* one left and below; or we are right on it */
698e3fbd1f4SBarry 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));
699e3fbd1f4SBarry Smith         cols[nc++] = col_shift + idx_c[col];
7009566063dSJacob Faibussowitsch         PetscCall(MatPreallocateSet(row, nc, cols, dnz, onz));
70147c6ae99SBarry Smith       }
70247c6ae99SBarry Smith     }
70347c6ae99SBarry Smith   }
7049566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)daf), &mat));
705fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
706fdc842d1SBarry Smith   /*
707fdc842d1SBarry Smith      Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
708fdc842d1SBarry Smith      we don't want the original unconverted matrix copied to the GPU
709fdc842d1SBarry Smith   */
710*48a46eb9SPierre Jolivet   if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE));
711fdc842d1SBarry Smith #endif
7129566063dSJacob 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));
7139566063dSJacob Faibussowitsch   PetscCall(ConvertToAIJ(dac->mattype, &mattype));
7149566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat, mattype));
7159566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat, 0, dnz));
7169566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(mat, 0, dnz, 0, onz));
717d0609cedSBarry Smith   MatPreallocateEnd(dnz, onz);
71847c6ae99SBarry Smith 
71947c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
72047c6ae99SBarry Smith   for (l = l_start; l < l_start + p_f; l++) {
72147c6ae99SBarry Smith     for (j = j_start; j < j_start + n_f; j++) {
72247c6ae99SBarry Smith       for (i = i_start; i < i_start + m_f; i++) {
72347c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
724e3fbd1f4SBarry Smith         row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
72547c6ae99SBarry Smith 
72647c6ae99SBarry Smith         i_c      = (i / ratioi); /* coarse grid node to left of fine grid node */
72747c6ae99SBarry Smith         j_c      = (j / ratioj); /* coarse grid node below fine grid node */
72847c6ae99SBarry Smith         l_c      = (l / ratiol);
72947c6ae99SBarry Smith         nc       = 0;
73047c6ae99SBarry Smith         /* one left and below; or we are right on it */
731e3fbd1f4SBarry 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));
732e3fbd1f4SBarry Smith         cols[nc] = col_shift + idx_c[col];
73347c6ae99SBarry Smith         v[nc++]  = 1.0;
73447c6ae99SBarry Smith 
7359566063dSJacob Faibussowitsch         PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES));
73647c6ae99SBarry Smith       }
73747c6ae99SBarry Smith     }
73847c6ae99SBarry Smith   }
7399566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
7409566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
7419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
7429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
7439566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(mat, dof, A));
7449566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
7459566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(13.0 * m_f * n_f * p_f));
74647c6ae99SBarry Smith   PetscFunctionReturn(0);
74747c6ae99SBarry Smith }
74847c6ae99SBarry Smith 
7499371c9d4SSatish Balay PetscErrorCode DMCreateInterpolation_DA_3D_Q1(DM dac, DM daf, Mat *A) {
7508ea3bf28SBarry Smith   PetscInt               i, j, i_start, j_start, m_f, n_f, Mx, My, dof, l;
7518ea3bf28SBarry Smith   const PetscInt        *idx_c, *idx_f;
752e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f, ltog_c;
7538ea3bf28SBarry Smith   PetscInt               m_ghost, n_ghost, m_ghost_c, n_ghost_c, Mz, mz;
75447c6ae99SBarry Smith   PetscInt               row, col, i_start_ghost, j_start_ghost, cols[8], mx, m_c, my, nc, ratioi, ratioj, ratiok;
75547c6ae99SBarry Smith   PetscInt               i_c, j_c, i_start_c, j_start_c, n_c, i_start_ghost_c, j_start_ghost_c;
75647c6ae99SBarry Smith   PetscInt               l_start, p_f, l_start_ghost, p_ghost, l_start_c, p_c;
75747c6ae99SBarry Smith   PetscInt               l_start_ghost_c, p_ghost_c, l_c, *dnz, *onz;
75847c6ae99SBarry Smith   PetscScalar            v[8], x, y, z;
75947c6ae99SBarry Smith   Mat                    mat;
760bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx, by, bz;
761fdc842d1SBarry Smith   MatType                mattype;
76247c6ae99SBarry Smith 
76347c6ae99SBarry Smith   PetscFunctionBegin;
7649566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, &Mz, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL));
7659566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, &mz, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
76647c6ae99SBarry Smith   if (mx == Mx) {
76747c6ae99SBarry Smith     ratioi = 1;
768bff4a2f0SMatthew G. Knepley   } else if (bx == DM_BOUNDARY_PERIODIC) {
76963a3b9bcSJacob Faibussowitsch     PetscCheck(Mx, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of x coarse grid points %" PetscInt_FMT " must be positive", Mx);
77047c6ae99SBarry Smith     ratioi = mx / Mx;
77163a3b9bcSJacob 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);
77247c6ae99SBarry Smith   } else {
77363a3b9bcSJacob 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);
77447c6ae99SBarry Smith     ratioi = (mx - 1) / (Mx - 1);
7751dca8a05SBarry 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);
77647c6ae99SBarry Smith   }
77747c6ae99SBarry Smith   if (my == My) {
77847c6ae99SBarry Smith     ratioj = 1;
779bff4a2f0SMatthew G. Knepley   } else if (by == DM_BOUNDARY_PERIODIC) {
78063a3b9bcSJacob Faibussowitsch     PetscCheck(My, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of y coarse grid points %" PetscInt_FMT " must be positive", My);
78147c6ae99SBarry Smith     ratioj = my / My;
78263a3b9bcSJacob 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);
78347c6ae99SBarry Smith   } else {
78463a3b9bcSJacob 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);
78547c6ae99SBarry Smith     ratioj = (my - 1) / (My - 1);
7861dca8a05SBarry 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);
78747c6ae99SBarry Smith   }
78847c6ae99SBarry Smith   if (mz == Mz) {
78947c6ae99SBarry Smith     ratiok = 1;
790bff4a2f0SMatthew G. Knepley   } else if (bz == DM_BOUNDARY_PERIODIC) {
79163a3b9bcSJacob Faibussowitsch     PetscCheck(Mz, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of z coarse grid points %" PetscInt_FMT " must be positive", Mz);
79247c6ae99SBarry Smith     ratiok = mz / Mz;
79363a3b9bcSJacob 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);
79447c6ae99SBarry Smith   } else {
79563a3b9bcSJacob 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);
79647c6ae99SBarry Smith     ratiok = (mz - 1) / (Mz - 1);
7971dca8a05SBarry 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);
79847c6ae99SBarry Smith   }
79947c6ae99SBarry Smith 
8009566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf, &i_start, &j_start, &l_start, &m_f, &n_f, &p_f));
8019566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, &l_start_ghost, &m_ghost, &n_ghost, &p_ghost));
8029566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf, &ltog_f));
8039566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
80447c6ae99SBarry Smith 
8059566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, &l_start_c, &m_c, &n_c, &p_c));
8069566063dSJacob 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));
8079566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac, &ltog_c));
8089566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
80947c6ae99SBarry Smith 
81047c6ae99SBarry Smith   /* create interpolation matrix, determining exact preallocation */
811d0609cedSBarry Smith   MatPreallocateBegin(PetscObjectComm((PetscObject)dac), m_f * n_f * p_f, m_c * n_c * p_c, dnz, onz);
81247c6ae99SBarry Smith   /* loop over local fine grid nodes counting interpolating points */
81347c6ae99SBarry Smith   for (l = l_start; l < l_start + p_f; l++) {
81447c6ae99SBarry Smith     for (j = j_start; j < j_start + n_f; j++) {
81547c6ae99SBarry Smith       for (i = i_start; i < i_start + m_f; i++) {
81647c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
817e3fbd1f4SBarry Smith         row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
81847c6ae99SBarry Smith         i_c = (i / ratioi);
81947c6ae99SBarry Smith         j_c = (j / ratioj);
82047c6ae99SBarry Smith         l_c = (l / ratiok);
82108401ef6SPierre Jolivet         PetscCheck(l_c >= l_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
8229371c9d4SSatish Balay                                             l_start %" PetscInt_FMT " l_c %" PetscInt_FMT " l_start_ghost_c %" PetscInt_FMT,
8239371c9d4SSatish Balay                    l_start, l_c, l_start_ghost_c);
82408401ef6SPierre Jolivet         PetscCheck(j_c >= j_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
8259371c9d4SSatish Balay                                             j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT,
8269371c9d4SSatish Balay                    j_start, j_c, j_start_ghost_c);
82708401ef6SPierre Jolivet         PetscCheck(i_c >= i_start_ghost_c, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
8289371c9d4SSatish Balay                                             i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT,
8299371c9d4SSatish Balay                    i_start, i_c, i_start_ghost_c);
83047c6ae99SBarry Smith 
83147c6ae99SBarry Smith         /*
83247c6ae99SBarry Smith          Only include those interpolation points that are truly
83347c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
83447c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
83547c6ae99SBarry Smith          */
83647c6ae99SBarry Smith         nc         = 0;
837e3fbd1f4SBarry 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));
838e3fbd1f4SBarry Smith         cols[nc++] = idx_c[col];
8399371c9d4SSatish Balay         if (i_c * ratioi != i) { cols[nc++] = idx_c[col + 1]; }
8409371c9d4SSatish Balay         if (j_c * ratioj != j) { cols[nc++] = idx_c[col + m_ghost_c]; }
8419371c9d4SSatish Balay         if (l_c * ratiok != l) { cols[nc++] = idx_c[col + m_ghost_c * n_ghost_c]; }
8429371c9d4SSatish Balay         if (j_c * ratioj != j && i_c * ratioi != i) { cols[nc++] = idx_c[col + (m_ghost_c + 1)]; }
8439371c9d4SSatish Balay         if (j_c * ratioj != j && l_c * ratiok != l) { cols[nc++] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c)]; }
8449371c9d4SSatish Balay         if (i_c * ratioi != i && l_c * ratiok != l) { cols[nc++] = idx_c[col + (m_ghost_c * n_ghost_c + 1)]; }
8459371c9d4SSatish Balay         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)]; }
8469566063dSJacob Faibussowitsch         PetscCall(MatPreallocateSet(row, nc, cols, dnz, onz));
84747c6ae99SBarry Smith       }
84847c6ae99SBarry Smith     }
84947c6ae99SBarry Smith   }
8509566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dac), &mat));
851fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
852fdc842d1SBarry Smith   /*
853fdc842d1SBarry Smith      Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
854fdc842d1SBarry Smith      we don't want the original unconverted matrix copied to the GPU
855fdc842d1SBarry Smith   */
856*48a46eb9SPierre Jolivet   if (dof > 1) PetscCall(MatBindToCPU(mat, PETSC_TRUE));
857fdc842d1SBarry Smith #endif
8589566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat, m_f * n_f * p_f, m_c * n_c * p_c, mx * my * mz, Mx * My * Mz));
8599566063dSJacob Faibussowitsch   PetscCall(ConvertToAIJ(dac->mattype, &mattype));
8609566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat, mattype));
8619566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat, 0, dnz));
8629566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(mat, 0, dnz, 0, onz));
863d0609cedSBarry Smith   MatPreallocateEnd(dnz, onz);
86447c6ae99SBarry Smith 
86547c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
8662adb9dcfSBarry Smith   if (!NEWVERSION) {
86747c6ae99SBarry Smith     for (l = l_start; l < l_start + p_f; l++) {
86847c6ae99SBarry Smith       for (j = j_start; j < j_start + n_f; j++) {
86947c6ae99SBarry Smith         for (i = i_start; i < i_start + m_f; i++) {
87047c6ae99SBarry Smith           /* convert to local "natural" numbering and then to PETSc global numbering */
871e3fbd1f4SBarry Smith           row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
87247c6ae99SBarry Smith 
87347c6ae99SBarry Smith           i_c = (i / ratioi);
87447c6ae99SBarry Smith           j_c = (j / ratioj);
87547c6ae99SBarry Smith           l_c = (l / ratiok);
87647c6ae99SBarry Smith 
87747c6ae99SBarry Smith           /*
87847c6ae99SBarry Smith            Only include those interpolation points that are truly
87947c6ae99SBarry Smith            nonzero. Note this is very important for final grid lines
88047c6ae99SBarry Smith            in x and y directions; since they have no right/top neighbors
88147c6ae99SBarry Smith            */
8826712e2f1SBarry Smith           x = ((PetscReal)(i - i_c * ratioi)) / ((PetscReal)ratioi);
8836712e2f1SBarry Smith           y = ((PetscReal)(j - j_c * ratioj)) / ((PetscReal)ratioj);
8846712e2f1SBarry Smith           z = ((PetscReal)(l - l_c * ratiok)) / ((PetscReal)ratiok);
885b3bd94feSDave May 
88647c6ae99SBarry Smith           nc  = 0;
88747c6ae99SBarry Smith           /* one left and below; or we are right on it */
888e3fbd1f4SBarry 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));
88947c6ae99SBarry Smith 
890e3fbd1f4SBarry Smith           cols[nc] = idx_c[col];
89147c6ae99SBarry Smith           v[nc++]  = .125 * (1. - (2.0 * x - 1.)) * (1. - (2.0 * y - 1.)) * (1. - (2.0 * z - 1.));
89247c6ae99SBarry Smith 
89347c6ae99SBarry Smith           if (i_c * ratioi != i) {
894e3fbd1f4SBarry Smith             cols[nc] = idx_c[col + 1];
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 (j_c * ratioj != j) {
899e3fbd1f4SBarry Smith             cols[nc] = idx_c[col + m_ghost_c];
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 (l_c * ratiok != l) {
904e3fbd1f4SBarry Smith             cols[nc] = idx_c[col + m_ghost_c * n_ghost_c];
90547c6ae99SBarry Smith             v[nc++]  = .125 * (1. - (2.0 * x - 1.)) * (1. - (2.0 * y - 1.)) * (1. + (2.0 * z - 1.));
90647c6ae99SBarry Smith           }
90747c6ae99SBarry Smith 
90847c6ae99SBarry Smith           if (j_c * ratioj != j && i_c * ratioi != i) {
909e3fbd1f4SBarry Smith             cols[nc] = idx_c[col + (m_ghost_c + 1)];
91047c6ae99SBarry Smith             v[nc++]  = .125 * (1. + (2.0 * x - 1.)) * (1. + (2.0 * y - 1.)) * (1. - (2.0 * z - 1.));
91147c6ae99SBarry Smith           }
91247c6ae99SBarry Smith 
91347c6ae99SBarry Smith           if (j_c * ratioj != j && l_c * ratiok != l) {
914e3fbd1f4SBarry Smith             cols[nc] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c)];
91547c6ae99SBarry Smith             v[nc++]  = .125 * (1. - (2.0 * x - 1.)) * (1. + (2.0 * y - 1.)) * (1. + (2.0 * z - 1.));
91647c6ae99SBarry Smith           }
91747c6ae99SBarry Smith 
91847c6ae99SBarry Smith           if (i_c * ratioi != i && l_c * ratiok != l) {
919e3fbd1f4SBarry Smith             cols[nc] = idx_c[col + (m_ghost_c * n_ghost_c + 1)];
92047c6ae99SBarry Smith             v[nc++]  = .125 * (1. + (2.0 * x - 1.)) * (1. - (2.0 * y - 1.)) * (1. + (2.0 * z - 1.));
92147c6ae99SBarry Smith           }
92247c6ae99SBarry Smith 
92347c6ae99SBarry Smith           if (i_c * ratioi != i && l_c * ratiok != l && j_c * ratioj != j) {
924e3fbd1f4SBarry Smith             cols[nc] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c + 1)];
92547c6ae99SBarry Smith             v[nc++]  = .125 * (1. + (2.0 * x - 1.)) * (1. + (2.0 * y - 1.)) * (1. + (2.0 * z - 1.));
92647c6ae99SBarry Smith           }
9279566063dSJacob Faibussowitsch           PetscCall(MatSetValues(mat, 1, &row, nc, cols, v, INSERT_VALUES));
92847c6ae99SBarry Smith         }
92947c6ae99SBarry Smith       }
93047c6ae99SBarry Smith     }
931b3bd94feSDave May 
932b3bd94feSDave May   } else {
933b3bd94feSDave May     PetscScalar *xi, *eta, *zeta;
934b3bd94feSDave May     PetscInt     li, nxi, lj, neta, lk, nzeta, n;
935b3bd94feSDave May     PetscScalar  Ni[8];
936b3bd94feSDave May 
937b3bd94feSDave May     /* compute local coordinate arrays */
938b3bd94feSDave May     nxi   = ratioi + 1;
939b3bd94feSDave May     neta  = ratioj + 1;
940b3bd94feSDave May     nzeta = ratiok + 1;
9419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nxi, &xi));
9429566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(neta, &eta));
9439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nzeta, &zeta));
9448865f1eaSKarl Rupp     for (li = 0; li < nxi; li++) xi[li] = -1.0 + (PetscScalar)li * (2.0 / (PetscScalar)(nxi - 1));
9458865f1eaSKarl Rupp     for (lj = 0; lj < neta; lj++) eta[lj] = -1.0 + (PetscScalar)lj * (2.0 / (PetscScalar)(neta - 1));
9468865f1eaSKarl Rupp     for (lk = 0; lk < nzeta; lk++) zeta[lk] = -1.0 + (PetscScalar)lk * (2.0 / (PetscScalar)(nzeta - 1));
947b3bd94feSDave May 
948b3bd94feSDave May     for (l = l_start; l < l_start + p_f; l++) {
949b3bd94feSDave May       for (j = j_start; j < j_start + n_f; j++) {
950b3bd94feSDave May         for (i = i_start; i < i_start + m_f; i++) {
951b3bd94feSDave May           /* convert to local "natural" numbering and then to PETSc global numbering */
952e3fbd1f4SBarry Smith           row = idx_f[(m_ghost * n_ghost * (l - l_start_ghost) + m_ghost * (j - j_start_ghost) + (i - i_start_ghost))];
953b3bd94feSDave May 
954b3bd94feSDave May           i_c = (i / ratioi);
955b3bd94feSDave May           j_c = (j / ratioj);
956b3bd94feSDave May           l_c = (l / ratiok);
957b3bd94feSDave May 
958b3bd94feSDave May           /* remainders */
959b3bd94feSDave May           li = i - ratioi * (i / ratioi);
9608865f1eaSKarl Rupp           if (i == mx - 1) li = nxi - 1;
961b3bd94feSDave May           lj = j - ratioj * (j / ratioj);
9628865f1eaSKarl Rupp           if (j == my - 1) lj = neta - 1;
963b3bd94feSDave May           lk = l - ratiok * (l / ratiok);
9648865f1eaSKarl Rupp           if (l == mz - 1) lk = nzeta - 1;
965b3bd94feSDave May 
966b3bd94feSDave May           /* corners */
967e3fbd1f4SBarry 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));
968e3fbd1f4SBarry Smith           cols[0] = idx_c[col];
969b3bd94feSDave May           Ni[0]   = 1.0;
970b3bd94feSDave May           if ((li == 0) || (li == nxi - 1)) {
971b3bd94feSDave May             if ((lj == 0) || (lj == neta - 1)) {
972b3bd94feSDave May               if ((lk == 0) || (lk == nzeta - 1)) {
9739566063dSJacob Faibussowitsch                 PetscCall(MatSetValue(mat, row, cols[0], Ni[0], INSERT_VALUES));
974b3bd94feSDave May                 continue;
975b3bd94feSDave May               }
976b3bd94feSDave May             }
977b3bd94feSDave May           }
978b3bd94feSDave May 
979b3bd94feSDave May           /* edges + interior */
980b3bd94feSDave May           /* remainders */
9818865f1eaSKarl Rupp           if (i == mx - 1) i_c--;
9828865f1eaSKarl Rupp           if (j == my - 1) j_c--;
9838865f1eaSKarl Rupp           if (l == mz - 1) l_c--;
984b3bd94feSDave May 
985e3fbd1f4SBarry 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));
986e3fbd1f4SBarry Smith           cols[0] = idx_c[col];                   /* one left and below; or we are right on it */
987e3fbd1f4SBarry Smith           cols[1] = idx_c[col + 1];               /* one right and below */
988e3fbd1f4SBarry Smith           cols[2] = idx_c[col + m_ghost_c];       /* one left and above */
989e3fbd1f4SBarry Smith           cols[3] = idx_c[col + (m_ghost_c + 1)]; /* one right and above */
990b3bd94feSDave May 
991e3fbd1f4SBarry Smith           cols[4] = idx_c[col + m_ghost_c * n_ghost_c];                   /* one left and below and front; or we are right on it */
992e3fbd1f4SBarry Smith           cols[5] = idx_c[col + (m_ghost_c * n_ghost_c + 1)];             /* one right and below, and front */
993e3fbd1f4SBarry Smith           cols[6] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c)];     /* one left and above and front*/
994e3fbd1f4SBarry Smith           cols[7] = idx_c[col + (m_ghost_c * n_ghost_c + m_ghost_c + 1)]; /* one right and above and front */
995b3bd94feSDave May 
996b3bd94feSDave May           Ni[0] = 0.125 * (1.0 - xi[li]) * (1.0 - eta[lj]) * (1.0 - zeta[lk]);
997b3bd94feSDave May           Ni[1] = 0.125 * (1.0 + xi[li]) * (1.0 - eta[lj]) * (1.0 - zeta[lk]);
998b3bd94feSDave May           Ni[2] = 0.125 * (1.0 - xi[li]) * (1.0 + eta[lj]) * (1.0 - zeta[lk]);
999b3bd94feSDave May           Ni[3] = 0.125 * (1.0 + xi[li]) * (1.0 + eta[lj]) * (1.0 - zeta[lk]);
1000b3bd94feSDave May 
1001b3bd94feSDave May           Ni[4] = 0.125 * (1.0 - xi[li]) * (1.0 - eta[lj]) * (1.0 + zeta[lk]);
1002b3bd94feSDave May           Ni[5] = 0.125 * (1.0 + xi[li]) * (1.0 - eta[lj]) * (1.0 + zeta[lk]);
1003b3bd94feSDave May           Ni[6] = 0.125 * (1.0 - xi[li]) * (1.0 + eta[lj]) * (1.0 + zeta[lk]);
1004b3bd94feSDave May           Ni[7] = 0.125 * (1.0 + xi[li]) * (1.0 + eta[lj]) * (1.0 + zeta[lk]);
1005b3bd94feSDave May 
1006b3bd94feSDave May           for (n = 0; n < 8; n++) {
10078865f1eaSKarl Rupp             if (PetscAbsScalar(Ni[n]) < 1.0e-32) cols[n] = -1;
1008b3bd94feSDave May           }
10099566063dSJacob Faibussowitsch           PetscCall(MatSetValues(mat, 1, &row, 8, cols, Ni, INSERT_VALUES));
1010b3bd94feSDave May         }
1011b3bd94feSDave May       }
1012b3bd94feSDave May     }
10139566063dSJacob Faibussowitsch     PetscCall(PetscFree(xi));
10149566063dSJacob Faibussowitsch     PetscCall(PetscFree(eta));
10159566063dSJacob Faibussowitsch     PetscCall(PetscFree(zeta));
1016b3bd94feSDave May   }
10179566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
10189566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
10199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
10209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
102147c6ae99SBarry Smith 
10229566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(mat, dof, A));
10239566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
102447c6ae99SBarry Smith   PetscFunctionReturn(0);
102547c6ae99SBarry Smith }
102647c6ae99SBarry Smith 
10279371c9d4SSatish Balay PetscErrorCode DMCreateInterpolation_DA(DM dac, DM daf, Mat *A, Vec *scale) {
102847c6ae99SBarry Smith   PetscInt        dimc, Mc, Nc, Pc, mc, nc, pc, dofc, sc, dimf, Mf, Nf, Pf, mf, nf, pf, doff, sf;
1029bff4a2f0SMatthew G. Knepley   DMBoundaryType  bxc, byc, bzc, bxf, byf, bzf;
1030aa219208SBarry Smith   DMDAStencilType stc, stf;
103147c6ae99SBarry Smith   DM_DA          *ddc = (DM_DA *)dac->data;
103247c6ae99SBarry Smith 
103347c6ae99SBarry Smith   PetscFunctionBegin;
103447c6ae99SBarry Smith   PetscValidHeaderSpecific(dac, DM_CLASSID, 1);
103547c6ae99SBarry Smith   PetscValidHeaderSpecific(daf, DM_CLASSID, 2);
103647c6ae99SBarry Smith   PetscValidPointer(A, 3);
103747c6ae99SBarry Smith   if (scale) PetscValidPointer(scale, 4);
103847c6ae99SBarry Smith 
10399566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, &dimc, &Mc, &Nc, &Pc, &mc, &nc, &pc, &dofc, &sc, &bxc, &byc, &bzc, &stc));
10409566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, &dimf, &Mf, &Nf, &Pf, &mf, &nf, &pf, &doff, &sf, &bxf, &byf, &bzf, &stf));
104163a3b9bcSJacob Faibussowitsch   PetscCheck(dimc == dimf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Dimensions of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dimc, dimf);
104263a3b9bcSJacob Faibussowitsch   PetscCheck(dofc == doff, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "DOF of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dofc, doff);
104363a3b9bcSJacob Faibussowitsch   PetscCheck(sc == sf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil width of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, sc, sf);
10441dca8a05SBarry Smith   PetscCheck(bxc == bxf && byc == byf && bzc == bzf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Boundary type different in two DMDAs");
104508401ef6SPierre Jolivet   PetscCheck(stc == stf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil type different in two DMDAs");
10461dca8a05SBarry Smith   PetscCheck(Mc >= 2 || Mf <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in x direction");
10471dca8a05SBarry 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");
10481dca8a05SBarry 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");
104947c6ae99SBarry Smith 
1050aa219208SBarry Smith   if (ddc->interptype == DMDA_Q1) {
105147c6ae99SBarry Smith     if (dimc == 1) {
10529566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation_DA_1D_Q1(dac, daf, A));
105347c6ae99SBarry Smith     } else if (dimc == 2) {
10549566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation_DA_2D_Q1(dac, daf, A));
105547c6ae99SBarry Smith     } else if (dimc == 3) {
10569566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation_DA_3D_Q1(dac, daf, A));
105763a3b9bcSJacob 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);
1058aa219208SBarry Smith   } else if (ddc->interptype == DMDA_Q0) {
105947c6ae99SBarry Smith     if (dimc == 1) {
10609566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation_DA_1D_Q0(dac, daf, A));
106147c6ae99SBarry Smith     } else if (dimc == 2) {
10629566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation_DA_2D_Q0(dac, daf, A));
106347c6ae99SBarry Smith     } else if (dimc == 3) {
10649566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation_DA_3D_Q0(dac, daf, A));
106563a3b9bcSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)daf), PETSC_ERR_SUP, "No support for this DMDA dimension %" PetscInt_FMT " for interpolation type %d", dimc, (int)ddc->interptype);
106647c6ae99SBarry Smith   }
10671baa6e33SBarry Smith   if (scale) PetscCall(DMCreateInterpolationScale((DM)dac, (DM)daf, *A, scale));
106847c6ae99SBarry Smith   PetscFunctionReturn(0);
106947c6ae99SBarry Smith }
107047c6ae99SBarry Smith 
10719371c9d4SSatish Balay PetscErrorCode DMCreateInjection_DA_1D(DM dac, DM daf, VecScatter *inject) {
10728ea3bf28SBarry Smith   PetscInt               i, i_start, m_f, Mx, dof;
10738ea3bf28SBarry Smith   const PetscInt        *idx_f;
1074e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f;
10758ea3bf28SBarry Smith   PetscInt               m_ghost, m_ghost_c;
10760bb2b966SJungho Lee   PetscInt               row, i_start_ghost, mx, m_c, nc, ratioi;
10770bb2b966SJungho Lee   PetscInt               i_start_c, i_start_ghost_c;
10780bb2b966SJungho Lee   PetscInt              *cols;
1079bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
10800bb2b966SJungho Lee   Vec                    vecf, vecc;
10810bb2b966SJungho Lee   IS                     isf;
10820bb2b966SJungho Lee 
10830bb2b966SJungho Lee   PetscFunctionBegin;
10849566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, NULL, &Mx, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bx, NULL, NULL, NULL));
10859566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, NULL, &mx, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
1086bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
10870bb2b966SJungho Lee     ratioi = mx / Mx;
108863a3b9bcSJacob 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);
10890bb2b966SJungho Lee   } else {
10900bb2b966SJungho Lee     ratioi = (mx - 1) / (Mx - 1);
10911dca8a05SBarry 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);
10920bb2b966SJungho Lee   }
10930bb2b966SJungho Lee 
10949566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf, &i_start, NULL, NULL, &m_f, NULL, NULL));
10959566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, NULL, NULL, &m_ghost, NULL, NULL));
10969566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf, &ltog_f));
10979566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
10980bb2b966SJungho Lee 
10999566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac, &i_start_c, NULL, NULL, &m_c, NULL, NULL));
11009566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, NULL, NULL, &m_ghost_c, NULL, NULL));
11010bb2b966SJungho Lee 
11020bb2b966SJungho Lee   /* loop over local fine grid nodes setting interpolation for those*/
11030bb2b966SJungho Lee   nc = 0;
11049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m_f, &cols));
11050bb2b966SJungho Lee 
11060bb2b966SJungho Lee   for (i = i_start_c; i < i_start_c + m_c; i++) {
11070bb2b966SJungho Lee     PetscInt i_f = i * ratioi;
11080bb2b966SJungho Lee 
11091dca8a05SBarry Smith     PetscCheck(i_f >= i_start_ghost && i_f < i_start_ghost + m_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\ni_c %" PetscInt_FMT " i_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]", i, i_f, i_start_ghost, i_start_ghost + m_ghost);
11108865f1eaSKarl Rupp 
1111e3fbd1f4SBarry Smith     row        = idx_f[(i_f - i_start_ghost)];
1112e3fbd1f4SBarry Smith     cols[nc++] = row;
11130bb2b966SJungho Lee   }
11140bb2b966SJungho Lee 
11159566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
11169566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)daf), dof, nc, cols, PETSC_OWN_POINTER, &isf));
11179566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dac, &vecc));
11189566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(daf, &vecf));
11199566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(vecf, isf, vecc, NULL, inject));
11209566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dac, &vecc));
11219566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(daf, &vecf));
11229566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isf));
11230bb2b966SJungho Lee   PetscFunctionReturn(0);
11240bb2b966SJungho Lee }
11250bb2b966SJungho Lee 
11269371c9d4SSatish Balay PetscErrorCode DMCreateInjection_DA_2D(DM dac, DM daf, VecScatter *inject) {
11278ea3bf28SBarry Smith   PetscInt               i, j, i_start, j_start, m_f, n_f, Mx, My, dof;
11288ea3bf28SBarry Smith   const PetscInt        *idx_c, *idx_f;
1129e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f, ltog_c;
11308ea3bf28SBarry Smith   PetscInt               m_ghost, n_ghost, m_ghost_c, n_ghost_c;
113147c6ae99SBarry Smith   PetscInt               row, i_start_ghost, j_start_ghost, mx, m_c, my, nc, ratioi, ratioj;
1132076202ddSJed Brown   PetscInt               i_start_c, j_start_c, n_c, i_start_ghost_c, j_start_ghost_c;
113347c6ae99SBarry Smith   PetscInt              *cols;
1134bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx, by;
113547c6ae99SBarry Smith   Vec                    vecf, vecc;
113647c6ae99SBarry Smith   IS                     isf;
113747c6ae99SBarry Smith 
113847c6ae99SBarry Smith   PetscFunctionBegin;
11399566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, NULL, NULL, NULL, NULL, NULL, NULL, &bx, &by, NULL, NULL));
11409566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
1141bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
114247c6ae99SBarry Smith     ratioi = mx / Mx;
114363a3b9bcSJacob 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);
114447c6ae99SBarry Smith   } else {
114547c6ae99SBarry Smith     ratioi = (mx - 1) / (Mx - 1);
11461dca8a05SBarry 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);
114747c6ae99SBarry Smith   }
1148bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) {
114947c6ae99SBarry Smith     ratioj = my / My;
115063a3b9bcSJacob 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);
115147c6ae99SBarry Smith   } else {
115247c6ae99SBarry Smith     ratioj = (my - 1) / (My - 1);
11531dca8a05SBarry 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);
115447c6ae99SBarry Smith   }
115547c6ae99SBarry Smith 
11569566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf, &i_start, &j_start, NULL, &m_f, &n_f, NULL));
11579566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, NULL, &m_ghost, &n_ghost, NULL));
11589566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf, &ltog_f));
11599566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
116047c6ae99SBarry Smith 
11619566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, NULL, &m_c, &n_c, NULL));
11629566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac, &i_start_ghost_c, &j_start_ghost_c, NULL, &m_ghost_c, &n_ghost_c, NULL));
11639566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac, &ltog_c));
11649566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
116547c6ae99SBarry Smith 
116647c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
116747c6ae99SBarry Smith   nc = 0;
11689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_f * m_f, &cols));
1169076202ddSJed Brown   for (j = j_start_c; j < j_start_c + n_c; j++) {
1170076202ddSJed Brown     for (i = i_start_c; i < i_start_c + m_c; i++) {
1171076202ddSJed Brown       PetscInt i_f = i * ratioi, j_f = j * ratioj;
11721dca8a05SBarry Smith       PetscCheck(j_f >= j_start_ghost && j_f < j_start_ghost + n_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
11739371c9d4SSatish Balay     j_c %" PetscInt_FMT " j_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",
11749371c9d4SSatish Balay                  j, j_f, j_start_ghost, j_start_ghost + n_ghost);
11751dca8a05SBarry Smith       PetscCheck(i_f >= i_start_ghost && i_f < i_start_ghost + m_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Processor's coarse DMDA must lie over fine DMDA\n\
11769371c9d4SSatish Balay     i_c %" PetscInt_FMT " i_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",
11779371c9d4SSatish Balay                  i, i_f, i_start_ghost, i_start_ghost + m_ghost);
1178e3fbd1f4SBarry Smith       row        = idx_f[(m_ghost * (j_f - j_start_ghost) + (i_f - i_start_ghost))];
1179e3fbd1f4SBarry Smith       cols[nc++] = row;
118047c6ae99SBarry Smith     }
118147c6ae99SBarry Smith   }
11829566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
11839566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
118447c6ae99SBarry Smith 
11859566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)daf), dof, nc, cols, PETSC_OWN_POINTER, &isf));
11869566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dac, &vecc));
11879566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(daf, &vecf));
11889566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(vecf, isf, vecc, NULL, inject));
11899566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dac, &vecc));
11909566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(daf, &vecf));
11919566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isf));
119247c6ae99SBarry Smith   PetscFunctionReturn(0);
119347c6ae99SBarry Smith }
119447c6ae99SBarry Smith 
11959371c9d4SSatish Balay PetscErrorCode DMCreateInjection_DA_3D(DM dac, DM daf, VecScatter *inject) {
1196b1757049SJed Brown   PetscInt               i, j, k, i_start, j_start, k_start, m_f, n_f, p_f, Mx, My, Mz;
1197b1757049SJed Brown   PetscInt               m_ghost, n_ghost, p_ghost, m_ghost_c, n_ghost_c, p_ghost_c;
1198b1757049SJed Brown   PetscInt               i_start_ghost, j_start_ghost, k_start_ghost;
1199b1757049SJed Brown   PetscInt               mx, my, mz, ratioi, ratioj, ratiok;
1200b1757049SJed Brown   PetscInt               i_start_c, j_start_c, k_start_c;
1201b1757049SJed Brown   PetscInt               m_c, n_c, p_c;
1202b1757049SJed Brown   PetscInt               i_start_ghost_c, j_start_ghost_c, k_start_ghost_c;
1203b1757049SJed Brown   PetscInt               row, nc, dof;
12048ea3bf28SBarry Smith   const PetscInt        *idx_c, *idx_f;
1205e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f, ltog_c;
1206b1757049SJed Brown   PetscInt              *cols;
1207bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx, by, bz;
1208b1757049SJed Brown   Vec                    vecf, vecc;
1209b1757049SJed Brown   IS                     isf;
1210b1757049SJed Brown 
1211b1757049SJed Brown   PetscFunctionBegin;
12129566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, NULL, &Mx, &My, &Mz, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL));
12139566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, NULL, &mx, &my, &mz, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
1214b1757049SJed Brown 
1215bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
1216b1757049SJed Brown     ratioi = mx / Mx;
121763a3b9bcSJacob 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);
1218b1757049SJed Brown   } else {
1219b1757049SJed Brown     ratioi = (mx - 1) / (Mx - 1);
12201dca8a05SBarry 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);
1221b1757049SJed Brown   }
1222bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) {
1223b1757049SJed Brown     ratioj = my / My;
122463a3b9bcSJacob 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);
1225b1757049SJed Brown   } else {
1226b1757049SJed Brown     ratioj = (my - 1) / (My - 1);
12271dca8a05SBarry 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);
1228b1757049SJed Brown   }
1229bff4a2f0SMatthew G. Knepley   if (bz == DM_BOUNDARY_PERIODIC) {
1230b1757049SJed Brown     ratiok = mz / Mz;
123163a3b9bcSJacob 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);
1232b1757049SJed Brown   } else {
1233b1757049SJed Brown     ratiok = (mz - 1) / (Mz - 1);
12341dca8a05SBarry 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);
1235b1757049SJed Brown   }
1236b1757049SJed Brown 
12379566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf, &i_start, &j_start, &k_start, &m_f, &n_f, &p_f));
12389566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, &k_start_ghost, &m_ghost, &n_ghost, &p_ghost));
12399566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf, &ltog_f));
12409566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f, &idx_f));
1241b1757049SJed Brown 
12429566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, &k_start_c, &m_c, &n_c, &p_c));
12439566063dSJacob 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));
12449566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac, &ltog_c));
12459566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c, &idx_c));
1246b1757049SJed Brown 
1247b1757049SJed Brown   /* loop over local fine grid nodes setting interpolation for those*/
1248b1757049SJed Brown   nc = 0;
12499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_f * m_f * p_f, &cols));
1250b1757049SJed Brown   for (k = k_start_c; k < k_start_c + p_c; k++) {
1251b1757049SJed Brown     for (j = j_start_c; j < j_start_c + n_c; j++) {
1252b1757049SJed Brown       for (i = i_start_c; i < i_start_c + m_c; i++) {
1253b1757049SJed Brown         PetscInt i_f = i * ratioi, j_f = j * ratioj, k_f = k * ratiok;
12549371c9d4SSatish Balay         PetscCheck(k_f >= k_start_ghost && k_f < k_start_ghost + p_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
12559371c9d4SSatish Balay                    "Processor's coarse DMDA must lie over fine DMDA  "
12569371c9d4SSatish Balay                    "k_c %" PetscInt_FMT " k_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",
12579371c9d4SSatish Balay                    k, k_f, k_start_ghost, k_start_ghost + p_ghost);
12589371c9d4SSatish Balay         PetscCheck(j_f >= j_start_ghost && j_f < j_start_ghost + n_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
12599371c9d4SSatish Balay                    "Processor's coarse DMDA must lie over fine DMDA  "
12609371c9d4SSatish Balay                    "j_c %" PetscInt_FMT " j_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",
12619371c9d4SSatish Balay                    j, j_f, j_start_ghost, j_start_ghost + n_ghost);
12629371c9d4SSatish Balay         PetscCheck(i_f >= i_start_ghost && i_f < i_start_ghost + m_ghost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
12639371c9d4SSatish Balay                    "Processor's coarse DMDA must lie over fine DMDA  "
12649371c9d4SSatish Balay                    "i_c %" PetscInt_FMT " i_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",
12659371c9d4SSatish Balay                    i, i_f, i_start_ghost, i_start_ghost + m_ghost);
1266e3fbd1f4SBarry 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))];
1267e3fbd1f4SBarry Smith         cols[nc++] = row;
1268b1757049SJed Brown       }
1269b1757049SJed Brown     }
1270b1757049SJed Brown   }
12719566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f, &idx_f));
12729566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c, &idx_c));
1273b1757049SJed Brown 
12749566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)daf), dof, nc, cols, PETSC_OWN_POINTER, &isf));
12759566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dac, &vecc));
12769566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(daf, &vecf));
12779566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(vecf, isf, vecc, NULL, inject));
12789566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dac, &vecc));
12799566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(daf, &vecf));
12809566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isf));
1281b1757049SJed Brown   PetscFunctionReturn(0);
1282b1757049SJed Brown }
1283b1757049SJed Brown 
12849371c9d4SSatish Balay PetscErrorCode DMCreateInjection_DA(DM dac, DM daf, Mat *mat) {
128547c6ae99SBarry Smith   PetscInt        dimc, Mc, Nc, Pc, mc, nc, pc, dofc, sc, dimf, Mf, Nf, Pf, mf, nf, pf, doff, sf;
1286bff4a2f0SMatthew G. Knepley   DMBoundaryType  bxc, byc, bzc, bxf, byf, bzf;
1287aa219208SBarry Smith   DMDAStencilType stc, stf;
128860c16ac1SBarry Smith   VecScatter      inject = NULL;
128947c6ae99SBarry Smith 
129047c6ae99SBarry Smith   PetscFunctionBegin;
129147c6ae99SBarry Smith   PetscValidHeaderSpecific(dac, DM_CLASSID, 1);
129247c6ae99SBarry Smith   PetscValidHeaderSpecific(daf, DM_CLASSID, 2);
12936dbf9973SLawrence Mitchell   PetscValidPointer(mat, 3);
129447c6ae99SBarry Smith 
12959566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, &dimc, &Mc, &Nc, &Pc, &mc, &nc, &pc, &dofc, &sc, &bxc, &byc, &bzc, &stc));
12969566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, &dimf, &Mf, &Nf, &Pf, &mf, &nf, &pf, &doff, &sf, &bxf, &byf, &bzf, &stf));
129763a3b9bcSJacob Faibussowitsch   PetscCheck(dimc == dimf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Dimensions of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dimc, dimf);
129863a3b9bcSJacob Faibussowitsch   PetscCheck(dofc == doff, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "DOF of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dofc, doff);
129963a3b9bcSJacob Faibussowitsch   PetscCheck(sc == sf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil width of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, sc, sf);
13001dca8a05SBarry Smith   PetscCheck(bxc == bxf && byc == byf && bzc == bzf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Boundary type different in two DMDAs");
130108401ef6SPierre Jolivet   PetscCheck(stc == stf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil type different in two DMDAs");
130208401ef6SPierre Jolivet   PetscCheck(Mc >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in x direction");
13031dca8a05SBarry Smith   PetscCheck(dimc <= 1 || Nc >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in y direction");
13041dca8a05SBarry Smith   PetscCheck(dimc <= 2 || Pc >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coarse grid requires at least 2 points in z direction");
130547c6ae99SBarry Smith 
13060bb2b966SJungho Lee   if (dimc == 1) {
13079566063dSJacob Faibussowitsch     PetscCall(DMCreateInjection_DA_1D(dac, daf, &inject));
13080bb2b966SJungho Lee   } else if (dimc == 2) {
13099566063dSJacob Faibussowitsch     PetscCall(DMCreateInjection_DA_2D(dac, daf, &inject));
1310b1757049SJed Brown   } else if (dimc == 3) {
13119566063dSJacob Faibussowitsch     PetscCall(DMCreateInjection_DA_3D(dac, daf, &inject));
131247c6ae99SBarry Smith   }
13139566063dSJacob Faibussowitsch   PetscCall(MatCreateScatter(PetscObjectComm((PetscObject)inject), inject, mat));
13149566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&inject));
131547c6ae99SBarry Smith   PetscFunctionReturn(0);
131647c6ae99SBarry Smith }
131747c6ae99SBarry Smith 
131897779f9aSLisandro Dalcin /*@
131997779f9aSLisandro Dalcin    DMCreateAggregates - Deprecated, see DMDACreateAggregates.
13206c877ef6SSatish Balay 
13216c877ef6SSatish Balay    Level: intermediate
132297779f9aSLisandro Dalcin @*/
13239371c9d4SSatish Balay PetscErrorCode DMCreateAggregates(DM dac, DM daf, Mat *mat) {
1324a5bc1bf3SBarry Smith   return DMDACreateAggregates(dac, daf, mat);
132597779f9aSLisandro Dalcin }
132697779f9aSLisandro Dalcin 
132797779f9aSLisandro Dalcin /*@
132897779f9aSLisandro Dalcin    DMDACreateAggregates - Gets the aggregates that map between
132997779f9aSLisandro Dalcin    grids associated with two DMDAs.
133097779f9aSLisandro Dalcin 
133197779f9aSLisandro Dalcin    Collective on dmc
133297779f9aSLisandro Dalcin 
133397779f9aSLisandro Dalcin    Input Parameters:
133497779f9aSLisandro Dalcin +  dmc - the coarse grid DMDA
133597779f9aSLisandro Dalcin -  dmf - the fine grid DMDA
133697779f9aSLisandro Dalcin 
133797779f9aSLisandro Dalcin    Output Parameters:
133897779f9aSLisandro Dalcin .  rest - the restriction matrix (transpose of the projection matrix)
133997779f9aSLisandro Dalcin 
134097779f9aSLisandro Dalcin    Level: intermediate
134197779f9aSLisandro Dalcin 
134297779f9aSLisandro Dalcin    Note: This routine is not used by PETSc.
134397779f9aSLisandro Dalcin    It is not clear what its use case is and it may be removed in a future release.
134497779f9aSLisandro Dalcin    Users should contact petsc-maint@mcs.anl.gov if they plan to use it.
134597779f9aSLisandro Dalcin 
1346db781477SPatrick Sanan .seealso: `DMRefine()`, `DMCreateInjection()`, `DMCreateInterpolation()`
134797779f9aSLisandro Dalcin @*/
13489371c9d4SSatish Balay PetscErrorCode DMDACreateAggregates(DM dac, DM daf, Mat *rest) {
134947c6ae99SBarry Smith   PetscInt               dimc, Mc, Nc, Pc, mc, nc, pc, dofc, sc;
135047c6ae99SBarry Smith   PetscInt               dimf, Mf, Nf, Pf, mf, nf, pf, doff, sf;
1351bff4a2f0SMatthew G. Knepley   DMBoundaryType         bxc, byc, bzc, bxf, byf, bzf;
1352aa219208SBarry Smith   DMDAStencilType        stc, stf;
135347c6ae99SBarry Smith   PetscInt               i, j, l;
135447c6ae99SBarry Smith   PetscInt               i_start, j_start, l_start, m_f, n_f, p_f;
135547c6ae99SBarry Smith   PetscInt               i_start_ghost, j_start_ghost, l_start_ghost, m_ghost, n_ghost, p_ghost;
13568ea3bf28SBarry Smith   const PetscInt        *idx_f;
135747c6ae99SBarry Smith   PetscInt               i_c, j_c, l_c;
135847c6ae99SBarry Smith   PetscInt               i_start_c, j_start_c, l_start_c, m_c, n_c, p_c;
135947c6ae99SBarry Smith   PetscInt               i_start_ghost_c, j_start_ghost_c, l_start_ghost_c, m_ghost_c, n_ghost_c, p_ghost_c;
13608ea3bf28SBarry Smith   const PetscInt        *idx_c;
136147c6ae99SBarry Smith   PetscInt               d;
136247c6ae99SBarry Smith   PetscInt               a;
136347c6ae99SBarry Smith   PetscInt               max_agg_size;
136447c6ae99SBarry Smith   PetscInt              *fine_nodes;
136547c6ae99SBarry Smith   PetscScalar           *one_vec;
136647c6ae99SBarry Smith   PetscInt               fn_idx;
1367565245c5SBarry Smith   ISLocalToGlobalMapping ltogmf, ltogmc;
136847c6ae99SBarry Smith 
136947c6ae99SBarry Smith   PetscFunctionBegin;
137097779f9aSLisandro Dalcin   PetscValidHeaderSpecificType(dac, DM_CLASSID, 1, DMDA);
137197779f9aSLisandro Dalcin   PetscValidHeaderSpecificType(daf, DM_CLASSID, 2, DMDA);
137247c6ae99SBarry Smith   PetscValidPointer(rest, 3);
137347c6ae99SBarry Smith 
13749566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, &dimc, &Mc, &Nc, &Pc, &mc, &nc, &pc, &dofc, &sc, &bxc, &byc, &bzc, &stc));
13759566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, &dimf, &Mf, &Nf, &Pf, &mf, &nf, &pf, &doff, &sf, &bxf, &byf, &bzf, &stf));
137663a3b9bcSJacob Faibussowitsch   PetscCheck(dimc == dimf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Dimensions of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dimc, dimf);
137763a3b9bcSJacob Faibussowitsch   PetscCheck(dofc == doff, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "DOF of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, dofc, doff);
137863a3b9bcSJacob Faibussowitsch   PetscCheck(sc == sf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil width of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT, sc, sf);
13791dca8a05SBarry Smith   PetscCheck(bxc == bxf && byc == byf && bzc == bzf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Boundary type different in two DMDAs");
138008401ef6SPierre Jolivet   PetscCheck(stc == stf, PetscObjectComm((PetscObject)daf), PETSC_ERR_ARG_INCOMP, "Stencil type different in two DMDAs");
138147c6ae99SBarry Smith 
138263a3b9bcSJacob 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);
138363a3b9bcSJacob 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);
138463a3b9bcSJacob 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);
138547c6ae99SBarry Smith 
138647c6ae99SBarry Smith   if (Pc < 0) Pc = 1;
138747c6ae99SBarry Smith   if (Pf < 0) Pf = 1;
138847c6ae99SBarry Smith   if (Nc < 0) Nc = 1;
138947c6ae99SBarry Smith   if (Nf < 0) Nf = 1;
139047c6ae99SBarry Smith 
13919566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf, &i_start, &j_start, &l_start, &m_f, &n_f, &p_f));
13929566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf, &i_start_ghost, &j_start_ghost, &l_start_ghost, &m_ghost, &n_ghost, &p_ghost));
1393565245c5SBarry Smith 
13949566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf, &ltogmf));
13959566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetIndices(ltogmf, &idx_f));
139647c6ae99SBarry Smith 
13979566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac, &i_start_c, &j_start_c, &l_start_c, &m_c, &n_c, &p_c));
13989566063dSJacob 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));
1399565245c5SBarry Smith 
14009566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac, &ltogmc));
14019566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetIndices(ltogmc, &idx_c));
140247c6ae99SBarry Smith 
140347c6ae99SBarry Smith   /*
140447c6ae99SBarry Smith      Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios
140547c6ae99SBarry Smith      for dimension 1 and 2 respectively.
140647c6ae99SBarry Smith      Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1)
140747c6ae99SBarry Smith      and r_y*j and r_y*(j+1) will be grouped into the same coarse grid agregate.
140847c6ae99SBarry Smith      Each specific dof on the fine grid is mapped to one dof on the coarse grid.
140947c6ae99SBarry Smith   */
141047c6ae99SBarry Smith 
141147c6ae99SBarry Smith   max_agg_size = (Mf / Mc + 1) * (Nf / Nc + 1) * (Pf / Pc + 1);
141247c6ae99SBarry Smith 
141347c6ae99SBarry Smith   /* create the matrix that will contain the restriction operator */
14149371c9d4SSatish 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));
141547c6ae99SBarry Smith 
141647c6ae99SBarry Smith   /* store nodes in the fine grid here */
14179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(max_agg_size, &one_vec, max_agg_size, &fine_nodes));
141847c6ae99SBarry Smith   for (i = 0; i < max_agg_size; i++) one_vec[i] = 1.0;
141947c6ae99SBarry Smith 
142047c6ae99SBarry Smith   /* loop over all coarse nodes */
142147c6ae99SBarry Smith   for (l_c = l_start_c; l_c < l_start_c + p_c; l_c++) {
142247c6ae99SBarry Smith     for (j_c = j_start_c; j_c < j_start_c + n_c; j_c++) {
142347c6ae99SBarry Smith       for (i_c = i_start_c; i_c < i_start_c + m_c; i_c++) {
142447c6ae99SBarry Smith         for (d = 0; d < dofc; d++) {
142547c6ae99SBarry Smith           /* convert to local "natural" numbering and then to PETSc global numbering */
142647c6ae99SBarry 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;
142747c6ae99SBarry Smith 
142847c6ae99SBarry Smith           fn_idx = 0;
142947c6ae99SBarry Smith           /* Corresponding fine points are all points (i_f, j_f, l_f) such that
143047c6ae99SBarry Smith              i_c*Mf/Mc <= i_f < (i_c+1)*Mf/Mc
143147c6ae99SBarry Smith              (same for other dimensions)
143247c6ae99SBarry Smith           */
143347c6ae99SBarry Smith           for (l = l_c * Pf / Pc; l < PetscMin((l_c + 1) * Pf / Pc, Pf); l++) {
143447c6ae99SBarry Smith             for (j = j_c * Nf / Nc; j < PetscMin((j_c + 1) * Nf / Nc, Nf); j++) {
143547c6ae99SBarry Smith               for (i = i_c * Mf / Mc; i < PetscMin((i_c + 1) * Mf / Mc, Mf); i++) {
143647c6ae99SBarry 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;
143747c6ae99SBarry Smith                 fn_idx++;
143847c6ae99SBarry Smith               }
143947c6ae99SBarry Smith             }
144047c6ae99SBarry Smith           }
144147c6ae99SBarry Smith           /* add all these points to one aggregate */
14429566063dSJacob Faibussowitsch           PetscCall(MatSetValues(*rest, 1, &a, fn_idx, fine_nodes, one_vec, INSERT_VALUES));
144347c6ae99SBarry Smith         }
144447c6ae99SBarry Smith       }
144547c6ae99SBarry Smith     }
144647c6ae99SBarry Smith   }
14479566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogmf, &idx_f));
14489566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogmc, &idx_c));
14499566063dSJacob Faibussowitsch   PetscCall(PetscFree2(one_vec, fine_nodes));
14509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*rest, MAT_FINAL_ASSEMBLY));
14519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*rest, MAT_FINAL_ASSEMBLY));
145247c6ae99SBarry Smith   PetscFunctionReturn(0);
145347c6ae99SBarry Smith }
1454