xref: /petsc/src/dm/impls/da/dainterp.c (revision 1dca8a0504492127e77eac64bc165d7372dd6d63)
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 */
22fdc842d1SBarry Smith static PetscErrorCode ConvertToAIJ(MatType intype,MatType *outtype)
23fdc842d1SBarry Smith {
24fdc842d1SBarry Smith   PetscInt       i;
25fdc842d1SBarry Smith   char           const *types[3] = {MATAIJ,MATSEQAIJ,MATMPIAIJ};
26fdc842d1SBarry Smith   PetscBool      flg;
27fdc842d1SBarry Smith 
28fdc842d1SBarry Smith   PetscFunctionBegin;
2933a4d587SStefano Zampini   *outtype = MATAIJ;
30fdc842d1SBarry Smith   for (i=0; i<3; i++) {
319566063dSJacob Faibussowitsch     PetscCall(PetscStrbeginswith(intype,types[i],&flg));
32fdc842d1SBarry Smith     if (flg) {
33fdc842d1SBarry Smith       *outtype = intype;
34fdc842d1SBarry Smith       PetscFunctionReturn(0);
35fdc842d1SBarry Smith     }
36fdc842d1SBarry Smith   }
37fdc842d1SBarry Smith   PetscFunctionReturn(0);
38fdc842d1SBarry Smith }
39fdc842d1SBarry Smith 
40e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_1D_Q1(DM dac,DM daf,Mat *A)
4147c6ae99SBarry Smith {
428ea3bf28SBarry Smith   PetscInt               i,i_start,m_f,Mx;
438ea3bf28SBarry Smith   const PetscInt         *idx_f,*idx_c;
448ea3bf28SBarry Smith   PetscInt               m_ghost,m_ghost_c;
4547c6ae99SBarry Smith   PetscInt               row,col,i_start_ghost,mx,m_c,nc,ratio;
4647c6ae99SBarry Smith   PetscInt               i_c,i_start_c,i_start_ghost_c,cols[2],dof;
4785afcc9aSBarry Smith   PetscScalar            v[2],x;
4847c6ae99SBarry Smith   Mat                    mat;
49bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
50e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
51fdc842d1SBarry Smith   MatType                mattype;
5247c6ae99SBarry Smith 
5347c6ae99SBarry Smith   PetscFunctionBegin;
549566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,NULL,&Mx,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bx,NULL,NULL,NULL));
559566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf,NULL,&mx,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
56bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
5747c6ae99SBarry Smith     ratio = mx/Mx;
5863a3b9bcSJacob Faibussowitsch     PetscCheck(ratio*Mx == mx,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: mx/Mx  must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT,mx,Mx);
5947c6ae99SBarry Smith   } else {
6047c6ae99SBarry Smith     ratio = (mx-1)/(Mx-1);
61*1dca8a05SBarry Smith     PetscCheck(ratio*(Mx-1) == mx-1,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %" PetscInt_FMT " Mx %" PetscInt_FMT,mx,Mx);
6247c6ae99SBarry Smith   }
6347c6ae99SBarry Smith 
649566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf,&i_start,NULL,NULL,&m_f,NULL,NULL));
659566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf,&i_start_ghost,NULL,NULL,&m_ghost,NULL,NULL));
669566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf,&ltog_f));
679566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f));
6847c6ae99SBarry Smith 
699566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac,&i_start_c,NULL,NULL,&m_c,NULL,NULL));
709566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac,&i_start_ghost_c,NULL,NULL,&m_ghost_c,NULL,NULL));
719566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac,&ltog_c));
729566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c));
7347c6ae99SBarry Smith 
7447c6ae99SBarry Smith   /* create interpolation matrix */
759566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dac),&mat));
76fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
77fdc842d1SBarry Smith   /*
78fdc842d1SBarry Smith      Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
79fdc842d1SBarry Smith      we don't want the original unconverted matrix copied to the GPU
80fdc842d1SBarry Smith    */
81fdc842d1SBarry Smith   if (dof > 1) {
829566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(mat,PETSC_TRUE));
83fdc842d1SBarry Smith   }
84fdc842d1SBarry Smith   #endif
859566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat,m_f,m_c,mx,Mx));
869566063dSJacob Faibussowitsch   PetscCall(ConvertToAIJ(dac->mattype,&mattype));
879566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat,mattype));
889566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat,2,NULL));
899566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(mat,2,NULL,1,NULL));
9047c6ae99SBarry Smith 
9147c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
9285afcc9aSBarry Smith   if (!NEWVERSION) {
93b3bd94feSDave May 
9447c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
9547c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
96e3fbd1f4SBarry Smith       row = idx_f[i-i_start_ghost];
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith       i_c = (i/ratio);    /* coarse grid node to left of fine grid node */
9908401ef6SPierre 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\
10063a3b9bcSJacob Faibussowitsch                                           i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT,i_start,i_c,i_start_ghost_c);
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith       /*
10347c6ae99SBarry Smith        Only include those interpolation points that are truly
10447c6ae99SBarry Smith        nonzero. Note this is very important for final grid lines
10547c6ae99SBarry Smith        in x direction; since they have no right neighbor
10647c6ae99SBarry Smith        */
1076712e2f1SBarry Smith       x  = ((PetscReal)(i - i_c*ratio))/((PetscReal)ratio);
10847c6ae99SBarry Smith       nc = 0;
10947c6ae99SBarry Smith       /* one left and below; or we are right on it */
110e3fbd1f4SBarry Smith       col      = (i_c-i_start_ghost_c);
111e3fbd1f4SBarry Smith       cols[nc] = idx_c[col];
11247c6ae99SBarry Smith       v[nc++]  = -x + 1.0;
11347c6ae99SBarry Smith       /* one right? */
11447c6ae99SBarry Smith       if (i_c*ratio != i) {
115e3fbd1f4SBarry Smith         cols[nc] = idx_c[col+1];
11647c6ae99SBarry Smith         v[nc++]  = x;
11747c6ae99SBarry Smith       }
1189566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES));
11947c6ae99SBarry Smith     }
120b3bd94feSDave May 
121b3bd94feSDave May   } else {
122b3bd94feSDave May     PetscScalar *xi;
123b3bd94feSDave May     PetscInt    li,nxi,n;
124b3bd94feSDave May     PetscScalar Ni[2];
125b3bd94feSDave May 
126b3bd94feSDave May     /* compute local coordinate arrays */
127b3bd94feSDave May     nxi  = ratio + 1;
1289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nxi,&xi));
129b3bd94feSDave May     for (li=0; li<nxi; li++) {
13052218ef2SJed Brown       xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1));
131b3bd94feSDave May     }
132b3bd94feSDave May 
133b3bd94feSDave May     for (i=i_start; i<i_start+m_f; i++) {
134b3bd94feSDave May       /* convert to local "natural" numbering and then to PETSc global numbering */
135e3fbd1f4SBarry Smith       row = idx_f[(i-i_start_ghost)];
136b3bd94feSDave May 
137b3bd94feSDave May       i_c = (i/ratio);    /* coarse grid node to left of fine grid node */
13808401ef6SPierre 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\
13963a3b9bcSJacob Faibussowitsch                                           i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT,i_start,i_c,i_start_ghost_c);
140b3bd94feSDave May 
141b3bd94feSDave May       /* remainders */
142b3bd94feSDave May       li = i - ratio * (i/ratio);
1438865f1eaSKarl Rupp       if (i==mx-1) li = nxi-1;
144b3bd94feSDave May 
145b3bd94feSDave May       /* corners */
146e3fbd1f4SBarry Smith       col     = (i_c-i_start_ghost_c);
147e3fbd1f4SBarry Smith       cols[0] = idx_c[col];
148b3bd94feSDave May       Ni[0]   = 1.0;
149b3bd94feSDave May       if ((li==0) || (li==nxi-1)) {
1509566063dSJacob Faibussowitsch         PetscCall(MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES));
151b3bd94feSDave May         continue;
152b3bd94feSDave May       }
153b3bd94feSDave May 
154b3bd94feSDave May       /* edges + interior */
155b3bd94feSDave May       /* remainders */
1568865f1eaSKarl Rupp       if (i==mx-1) i_c--;
157b3bd94feSDave May 
158e3fbd1f4SBarry Smith       col     = (i_c-i_start_ghost_c);
159e3fbd1f4SBarry Smith       cols[0] = idx_c[col]; /* one left and below; or we are right on it */
160e3fbd1f4SBarry Smith       cols[1] = idx_c[col+1];
161b3bd94feSDave May 
162b3bd94feSDave May       Ni[0] = 0.5*(1.0-xi[li]);
163b3bd94feSDave May       Ni[1] = 0.5*(1.0+xi[li]);
164b3bd94feSDave May       for (n=0; n<2; n++) {
1658865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[n])<1.0e-32) cols[n]=-1;
166b3bd94feSDave May       }
1679566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat,1,&row,2,cols,Ni,INSERT_VALUES));
168b3bd94feSDave May     }
1699566063dSJacob Faibussowitsch     PetscCall(PetscFree(xi));
170b3bd94feSDave May   }
1719566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f));
1729566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c));
1739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
1749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
1759566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(mat,dof,A));
1769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
17747c6ae99SBarry Smith   PetscFunctionReturn(0);
17847c6ae99SBarry Smith }
17947c6ae99SBarry Smith 
180e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_1D_Q0(DM dac,DM daf,Mat *A)
18147c6ae99SBarry Smith {
1828ea3bf28SBarry Smith   PetscInt               i,i_start,m_f,Mx;
1838ea3bf28SBarry Smith   const PetscInt         *idx_f,*idx_c;
184e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
1858ea3bf28SBarry Smith   PetscInt               m_ghost,m_ghost_c;
18647c6ae99SBarry Smith   PetscInt               row,col,i_start_ghost,mx,m_c,nc,ratio;
18747c6ae99SBarry Smith   PetscInt               i_c,i_start_c,i_start_ghost_c,cols[2],dof;
18847c6ae99SBarry Smith   PetscScalar            v[2],x;
18947c6ae99SBarry Smith   Mat                    mat;
190bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
191fdc842d1SBarry Smith   MatType                mattype;
19247c6ae99SBarry Smith 
19347c6ae99SBarry Smith   PetscFunctionBegin;
1949566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,NULL,&Mx,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bx,NULL,NULL,NULL));
1959566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf,NULL,&mx,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
196bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
19763a3b9bcSJacob Faibussowitsch     PetscCheck(Mx,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %" PetscInt_FMT " must be positive",Mx);
19847c6ae99SBarry Smith     ratio = mx/Mx;
19963a3b9bcSJacob 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);
20047c6ae99SBarry Smith   } else {
20163a3b9bcSJacob 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);
20247c6ae99SBarry Smith     ratio = (mx-1)/(Mx-1);
203*1dca8a05SBarry 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);
20447c6ae99SBarry Smith   }
20547c6ae99SBarry Smith 
2069566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf,&i_start,NULL,NULL,&m_f,NULL,NULL));
2079566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf,&i_start_ghost,NULL,NULL,&m_ghost,NULL,NULL));
2089566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf,&ltog_f));
2099566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f));
21047c6ae99SBarry Smith 
2119566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac,&i_start_c,NULL,NULL,&m_c,NULL,NULL));
2129566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac,&i_start_ghost_c,NULL,NULL,&m_ghost_c,NULL,NULL));
2139566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac,&ltog_c));
2149566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c));
21547c6ae99SBarry Smith 
21647c6ae99SBarry Smith   /* create interpolation matrix */
2179566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dac),&mat));
218fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
219fdc842d1SBarry Smith   /*
220fdc842d1SBarry Smith      Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
221fdc842d1SBarry Smith      we don't want the original unconverted matrix copied to the GPU
222fdc842d1SBarry Smith    */
223fdc842d1SBarry Smith   if (dof > 1) {
2249566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(mat,PETSC_TRUE));
225fdc842d1SBarry Smith   }
226fdc842d1SBarry Smith   #endif
2279566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat,m_f,m_c,mx,Mx));
2289566063dSJacob Faibussowitsch   PetscCall(ConvertToAIJ(dac->mattype,&mattype));
2299566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat,mattype));
2309566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat,2,NULL));
2319566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(mat,2,NULL,0,NULL));
23247c6ae99SBarry Smith 
23347c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
23447c6ae99SBarry Smith   for (i=i_start; i<i_start+m_f; i++) {
23547c6ae99SBarry Smith     /* convert to local "natural" numbering and then to PETSc global numbering */
236e3fbd1f4SBarry Smith     row = idx_f[(i-i_start_ghost)];
23747c6ae99SBarry Smith 
23847c6ae99SBarry Smith     i_c = (i/ratio);    /* coarse grid node to left of fine grid node */
23947c6ae99SBarry Smith 
24047c6ae99SBarry Smith     /*
24147c6ae99SBarry Smith          Only include those interpolation points that are truly
24247c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
24347c6ae99SBarry Smith          in x direction; since they have no right neighbor
24447c6ae99SBarry Smith     */
2456712e2f1SBarry Smith     x  = ((PetscReal)(i - i_c*ratio))/((PetscReal)ratio);
24647c6ae99SBarry Smith     nc = 0;
24747c6ae99SBarry Smith     /* one left and below; or we are right on it */
248e3fbd1f4SBarry Smith     col      = (i_c-i_start_ghost_c);
249e3fbd1f4SBarry Smith     cols[nc] = idx_c[col];
25047c6ae99SBarry Smith     v[nc++]  = -x + 1.0;
25147c6ae99SBarry Smith     /* one right? */
25247c6ae99SBarry Smith     if (i_c*ratio != i) {
253e3fbd1f4SBarry Smith       cols[nc] = idx_c[col+1];
25447c6ae99SBarry Smith       v[nc++]  = x;
25547c6ae99SBarry Smith     }
2569566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES));
25747c6ae99SBarry Smith   }
2589566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f));
2599566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c));
2609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
2619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
2629566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(mat,dof,A));
2639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
2649566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(5.0*m_f));
26547c6ae99SBarry Smith   PetscFunctionReturn(0);
26647c6ae99SBarry Smith }
26747c6ae99SBarry Smith 
268e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_2D_Q1(DM dac,DM daf,Mat *A)
26947c6ae99SBarry Smith {
2708ea3bf28SBarry Smith   PetscInt               i,j,i_start,j_start,m_f,n_f,Mx,My,dof;
2718ea3bf28SBarry Smith   const PetscInt         *idx_c,*idx_f;
272e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
2738ea3bf28SBarry Smith   PetscInt               m_ghost,n_ghost,m_ghost_c,n_ghost_c,*dnz,*onz;
27447c6ae99SBarry Smith   PetscInt               row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj;
27547c6ae99SBarry 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;
27647c6ae99SBarry Smith   PetscMPIInt            size_c,size_f,rank_f;
27747c6ae99SBarry Smith   PetscScalar            v[4],x,y;
27847c6ae99SBarry Smith   Mat                    mat;
279bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
280fdc842d1SBarry Smith   MatType                mattype;
28147c6ae99SBarry Smith 
28247c6ae99SBarry Smith   PetscFunctionBegin;
2839566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,NULL,&Mx,&My,NULL,NULL,NULL,NULL,NULL,NULL,&bx,&by,NULL,NULL));
2849566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf,NULL,&mx,&my,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
285bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
28663a3b9bcSJacob Faibussowitsch     PetscCheck(Mx,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %" PetscInt_FMT " must be positive",Mx);
28747c6ae99SBarry Smith     ratioi = mx/Mx;
28863a3b9bcSJacob 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);
28947c6ae99SBarry Smith   } else {
29063a3b9bcSJacob 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);
29147c6ae99SBarry Smith     ratioi = (mx-1)/(Mx-1);
292*1dca8a05SBarry 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);
29347c6ae99SBarry Smith   }
294bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) {
29563a3b9bcSJacob Faibussowitsch     PetscCheck(My,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of y coarse grid points %" PetscInt_FMT " must be positive",My);
29647c6ae99SBarry Smith     ratioj = my/My;
29763a3b9bcSJacob 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);
29847c6ae99SBarry Smith   } else {
29963a3b9bcSJacob 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);
30047c6ae99SBarry Smith     ratioj = (my-1)/(My-1);
301*1dca8a05SBarry 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);
30247c6ae99SBarry Smith   }
30347c6ae99SBarry Smith 
3049566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf,&i_start,&j_start,NULL,&m_f,&n_f,NULL));
3059566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,NULL,&m_ghost,&n_ghost,NULL));
3069566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf,&ltog_f));
3079566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f));
30847c6ae99SBarry Smith 
3099566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac,&i_start_c,&j_start_c,NULL,&m_c,&n_c,NULL));
3109566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,NULL,&m_ghost_c,&n_ghost_c,NULL));
3119566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac,&ltog_c));
3129566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c));
31347c6ae99SBarry Smith 
31447c6ae99SBarry Smith   /*
315aa219208SBarry Smith    Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
31647c6ae99SBarry Smith    The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
31747c6ae99SBarry Smith    processors). It's effective length is hence 4 times its normal length, this is
31847c6ae99SBarry Smith    why the col_scale is multiplied by the interpolation matrix column sizes.
31947c6ae99SBarry Smith    sol_shift allows each set of 1/4 processors do its own interpolation using ITS
32047c6ae99SBarry Smith    copy of the coarse vector. A bit of a hack but you do better.
32147c6ae99SBarry Smith 
32247c6ae99SBarry Smith    In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
32347c6ae99SBarry Smith    */
3249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c));
3259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f));
3269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f));
32747c6ae99SBarry Smith   col_scale = size_f/size_c;
32847c6ae99SBarry Smith   col_shift = Mx*My*(rank_f/size_c);
32947c6ae99SBarry Smith 
330d0609cedSBarry Smith   MatPreallocateBegin(PetscObjectComm((PetscObject)daf),m_f*n_f,col_scale*m_c*n_c,dnz,onz);
33147c6ae99SBarry Smith   for (j=j_start; j<j_start+n_f; j++) {
33247c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
33347c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
334e3fbd1f4SBarry Smith       row = idx_f[(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
33547c6ae99SBarry Smith 
33647c6ae99SBarry Smith       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
33747c6ae99SBarry Smith       j_c = (j/ratioj);    /* coarse grid node below fine grid node */
33847c6ae99SBarry Smith 
33908401ef6SPierre 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\
34063a3b9bcSJacob Faibussowitsch                                           j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT,j_start,j_c,j_start_ghost_c);
34108401ef6SPierre 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\
34263a3b9bcSJacob Faibussowitsch                                           i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT,i_start,i_c,i_start_ghost_c);
34347c6ae99SBarry Smith 
34447c6ae99SBarry Smith       /*
34547c6ae99SBarry Smith        Only include those interpolation points that are truly
34647c6ae99SBarry Smith        nonzero. Note this is very important for final grid lines
34747c6ae99SBarry Smith        in x and y directions; since they have no right/top neighbors
34847c6ae99SBarry Smith        */
34947c6ae99SBarry Smith       nc = 0;
35047c6ae99SBarry Smith       /* one left and below; or we are right on it */
351e3fbd1f4SBarry Smith       col        = (m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
352e3fbd1f4SBarry Smith       cols[nc++] = col_shift + idx_c[col];
35347c6ae99SBarry Smith       /* one right and below */
354e3fbd1f4SBarry Smith       if (i_c*ratioi != i) cols[nc++] = col_shift + idx_c[col+1];
35547c6ae99SBarry Smith       /* one left and above */
356e3fbd1f4SBarry Smith       if (j_c*ratioj != j) cols[nc++] = col_shift + idx_c[col+m_ghost_c];
35747c6ae99SBarry Smith       /* one right and above */
358e3fbd1f4SBarry Smith       if (i_c*ratioi != i && j_c*ratioj != j) cols[nc++] = col_shift + idx_c[col+(m_ghost_c+1)];
3599566063dSJacob Faibussowitsch       PetscCall(MatPreallocateSet(row,nc,cols,dnz,onz));
36047c6ae99SBarry Smith     }
36147c6ae99SBarry Smith   }
3629566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)daf),&mat));
363fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
364fdc842d1SBarry Smith   /*
365fdc842d1SBarry Smith      Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
366fdc842d1SBarry Smith      we don't want the original unconverted matrix copied to the GPU
367fdc842d1SBarry Smith   */
368fdc842d1SBarry Smith   if (dof > 1) {
3699566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(mat,PETSC_TRUE));
370fdc842d1SBarry Smith   }
371fdc842d1SBarry Smith #endif
3729566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My));
3739566063dSJacob Faibussowitsch   PetscCall(ConvertToAIJ(dac->mattype,&mattype));
3749566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat,mattype));
3759566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat,0,dnz));
3769566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(mat,0,dnz,0,onz));
377d0609cedSBarry Smith   MatPreallocateEnd(dnz,onz);
37847c6ae99SBarry Smith 
37947c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
38085afcc9aSBarry Smith   if (!NEWVERSION) {
381b3bd94feSDave May 
38247c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
38347c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
38447c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
385e3fbd1f4SBarry Smith         row = idx_f[(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
38647c6ae99SBarry Smith 
38747c6ae99SBarry Smith         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
38847c6ae99SBarry Smith         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
38947c6ae99SBarry Smith 
39047c6ae99SBarry Smith         /*
39147c6ae99SBarry Smith          Only include those interpolation points that are truly
39247c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
39347c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
39447c6ae99SBarry Smith          */
3956712e2f1SBarry Smith         x = ((PetscReal)(i - i_c*ratioi))/((PetscReal)ratioi);
3966712e2f1SBarry Smith         y = ((PetscReal)(j - j_c*ratioj))/((PetscReal)ratioj);
397b3bd94feSDave May 
39847c6ae99SBarry Smith         nc = 0;
39947c6ae99SBarry Smith         /* one left and below; or we are right on it */
400e3fbd1f4SBarry Smith         col      = (m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
401e3fbd1f4SBarry Smith         cols[nc] = col_shift + idx_c[col];
40247c6ae99SBarry Smith         v[nc++]  = x*y - x - y + 1.0;
40347c6ae99SBarry Smith         /* one right and below */
40447c6ae99SBarry Smith         if (i_c*ratioi != i) {
405e3fbd1f4SBarry Smith           cols[nc] = col_shift + idx_c[col+1];
40647c6ae99SBarry Smith           v[nc++]  = -x*y + x;
40747c6ae99SBarry Smith         }
40847c6ae99SBarry Smith         /* one left and above */
40947c6ae99SBarry Smith         if (j_c*ratioj != j) {
410e3fbd1f4SBarry Smith           cols[nc] = col_shift + idx_c[col+m_ghost_c];
41147c6ae99SBarry Smith           v[nc++]  = -x*y + y;
41247c6ae99SBarry Smith         }
41347c6ae99SBarry Smith         /* one right and above */
41447c6ae99SBarry Smith         if (j_c*ratioj != j && i_c*ratioi != i) {
415e3fbd1f4SBarry Smith           cols[nc] = col_shift + idx_c[col+(m_ghost_c+1)];
41647c6ae99SBarry Smith           v[nc++]  = x*y;
41747c6ae99SBarry Smith         }
4189566063dSJacob Faibussowitsch         PetscCall(MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES));
41947c6ae99SBarry Smith       }
42047c6ae99SBarry Smith     }
421b3bd94feSDave May 
422b3bd94feSDave May   } else {
423b3bd94feSDave May     PetscScalar Ni[4];
424b3bd94feSDave May     PetscScalar *xi,*eta;
425b3bd94feSDave May     PetscInt    li,nxi,lj,neta;
426b3bd94feSDave May 
427b3bd94feSDave May     /* compute local coordinate arrays */
428b3bd94feSDave May     nxi  = ratioi + 1;
429b3bd94feSDave May     neta = ratioj + 1;
4309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nxi,&xi));
4319566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(neta,&eta));
432b3bd94feSDave May     for (li=0; li<nxi; li++) {
43352218ef2SJed Brown       xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1));
434b3bd94feSDave May     }
435b3bd94feSDave May     for (lj=0; lj<neta; lj++) {
43652218ef2SJed Brown       eta[lj] = -1.0 + (PetscScalar)lj*(2.0/(PetscScalar)(neta-1));
437b3bd94feSDave May     }
438b3bd94feSDave May 
439b3bd94feSDave May     /* loop over local fine grid nodes setting interpolation for those*/
440b3bd94feSDave May     for (j=j_start; j<j_start+n_f; j++) {
441b3bd94feSDave May       for (i=i_start; i<i_start+m_f; i++) {
442b3bd94feSDave May         /* convert to local "natural" numbering and then to PETSc global numbering */
443e3fbd1f4SBarry Smith         row = idx_f[(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
444b3bd94feSDave May 
445b3bd94feSDave May         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
446b3bd94feSDave May         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
447b3bd94feSDave May 
448b3bd94feSDave May         /* remainders */
449b3bd94feSDave May         li = i - ratioi * (i/ratioi);
4508865f1eaSKarl Rupp         if (i==mx-1) li = nxi-1;
451b3bd94feSDave May         lj = j - ratioj * (j/ratioj);
4528865f1eaSKarl Rupp         if (j==my-1) lj = neta-1;
453b3bd94feSDave May 
454b3bd94feSDave May         /* corners */
455e3fbd1f4SBarry Smith         col     = (m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
456e3fbd1f4SBarry Smith         cols[0] = col_shift + idx_c[col]; /* left, below */
457b3bd94feSDave May         Ni[0]   = 1.0;
458b3bd94feSDave May         if ((li==0) || (li==nxi-1)) {
459b3bd94feSDave May           if ((lj==0) || (lj==neta-1)) {
4609566063dSJacob Faibussowitsch             PetscCall(MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES));
461b3bd94feSDave May             continue;
462b3bd94feSDave May           }
463b3bd94feSDave May         }
464b3bd94feSDave May 
465b3bd94feSDave May         /* edges + interior */
466b3bd94feSDave May         /* remainders */
4678865f1eaSKarl Rupp         if (i==mx-1) i_c--;
4688865f1eaSKarl Rupp         if (j==my-1) j_c--;
469b3bd94feSDave May 
470e3fbd1f4SBarry Smith         col     = (m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
471e3fbd1f4SBarry Smith         cols[0] = col_shift + idx_c[col]; /* left, below */
472e3fbd1f4SBarry Smith         cols[1] = col_shift + idx_c[col+1]; /* right, below */
473e3fbd1f4SBarry Smith         cols[2] = col_shift + idx_c[col+m_ghost_c]; /* left, above */
474e3fbd1f4SBarry Smith         cols[3] = col_shift + idx_c[col+(m_ghost_c+1)]; /* right, above */
475b3bd94feSDave May 
476b3bd94feSDave May         Ni[0] = 0.25*(1.0-xi[li])*(1.0-eta[lj]);
477b3bd94feSDave May         Ni[1] = 0.25*(1.0+xi[li])*(1.0-eta[lj]);
478b3bd94feSDave May         Ni[2] = 0.25*(1.0-xi[li])*(1.0+eta[lj]);
479b3bd94feSDave May         Ni[3] = 0.25*(1.0+xi[li])*(1.0+eta[lj]);
480b3bd94feSDave May 
481b3bd94feSDave May         nc = 0;
4828865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[0])<1.0e-32) cols[0]=-1;
4838865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[1])<1.0e-32) cols[1]=-1;
4848865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[2])<1.0e-32) cols[2]=-1;
4858865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[3])<1.0e-32) cols[3]=-1;
486b3bd94feSDave May 
4879566063dSJacob Faibussowitsch         PetscCall(MatSetValues(mat,1,&row,4,cols,Ni,INSERT_VALUES));
488b3bd94feSDave May       }
489b3bd94feSDave May     }
4909566063dSJacob Faibussowitsch     PetscCall(PetscFree(xi));
4919566063dSJacob Faibussowitsch     PetscCall(PetscFree(eta));
492b3bd94feSDave May   }
4939566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f));
4949566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c));
4959566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
4969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
4979566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(mat,dof,A));
4989566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
49947c6ae99SBarry Smith   PetscFunctionReturn(0);
50047c6ae99SBarry Smith }
50147c6ae99SBarry Smith 
50247c6ae99SBarry Smith /*
50347c6ae99SBarry Smith        Contributed by Andrei Draganescu <aidraga@sandia.gov>
50447c6ae99SBarry Smith */
505e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_2D_Q0(DM dac,DM daf,Mat *A)
50647c6ae99SBarry Smith {
5078ea3bf28SBarry Smith   PetscInt               i,j,i_start,j_start,m_f,n_f,Mx,My,dof;
5088ea3bf28SBarry Smith   const PetscInt         *idx_c,*idx_f;
509e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
5108ea3bf28SBarry Smith   PetscInt               m_ghost,n_ghost,m_ghost_c,n_ghost_c,*dnz,*onz;
51147c6ae99SBarry Smith   PetscInt               row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj;
51247c6ae99SBarry 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;
51347c6ae99SBarry Smith   PetscMPIInt            size_c,size_f,rank_f;
51447c6ae99SBarry Smith   PetscScalar            v[4];
51547c6ae99SBarry Smith   Mat                    mat;
516bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
517fdc842d1SBarry Smith   MatType                mattype;
51847c6ae99SBarry Smith 
51947c6ae99SBarry Smith   PetscFunctionBegin;
5209566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,NULL,&Mx,&My,NULL,NULL,NULL,NULL,NULL,NULL,&bx,&by,NULL,NULL));
5219566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf,NULL,&mx,&my,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
52263a3b9bcSJacob Faibussowitsch   PetscCheck(Mx,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %" PetscInt_FMT " must be positive",Mx);
52363a3b9bcSJacob Faibussowitsch   PetscCheck(My,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of y coarse grid points %" PetscInt_FMT " must be positive",My);
52447c6ae99SBarry Smith   ratioi = mx/Mx;
52547c6ae99SBarry Smith   ratioj = my/My;
52608401ef6SPierre Jolivet   PetscCheck(ratioi*Mx == mx,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in x");
52708401ef6SPierre Jolivet   PetscCheck(ratioj*My == my,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in y");
52808401ef6SPierre Jolivet   PetscCheck(ratioi == 2,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 2");
52908401ef6SPierre Jolivet   PetscCheck(ratioj == 2,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 2");
53047c6ae99SBarry Smith 
5319566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf,&i_start,&j_start,NULL,&m_f,&n_f,NULL));
5329566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,NULL,&m_ghost,&n_ghost,NULL));
5339566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf,&ltog_f));
5349566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f));
53547c6ae99SBarry Smith 
5369566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac,&i_start_c,&j_start_c,NULL,&m_c,&n_c,NULL));
5379566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,NULL,&m_ghost_c,&n_ghost_c,NULL));
5389566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac,&ltog_c));
5399566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c));
54047c6ae99SBarry Smith 
54147c6ae99SBarry Smith   /*
542aa219208SBarry Smith      Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
54347c6ae99SBarry Smith      The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
54447c6ae99SBarry Smith      processors). It's effective length is hence 4 times its normal length, this is
54547c6ae99SBarry Smith      why the col_scale is multiplied by the interpolation matrix column sizes.
54647c6ae99SBarry Smith      sol_shift allows each set of 1/4 processors do its own interpolation using ITS
54747c6ae99SBarry Smith      copy of the coarse vector. A bit of a hack but you do better.
54847c6ae99SBarry Smith 
54947c6ae99SBarry Smith      In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
55047c6ae99SBarry Smith   */
5519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c));
5529566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f));
5539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f));
55447c6ae99SBarry Smith   col_scale = size_f/size_c;
55547c6ae99SBarry Smith   col_shift = Mx*My*(rank_f/size_c);
55647c6ae99SBarry Smith 
557d0609cedSBarry Smith   MatPreallocateBegin(PetscObjectComm((PetscObject)daf),m_f*n_f,col_scale*m_c*n_c,dnz,onz);
55847c6ae99SBarry Smith   for (j=j_start; j<j_start+n_f; j++) {
55947c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
56047c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
561e3fbd1f4SBarry Smith       row = idx_f[(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
56247c6ae99SBarry Smith 
56347c6ae99SBarry Smith       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
56447c6ae99SBarry Smith       j_c = (j/ratioj);    /* coarse grid node below fine grid node */
56547c6ae99SBarry Smith 
56608401ef6SPierre 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\
56763a3b9bcSJacob Faibussowitsch     j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT,j_start,j_c,j_start_ghost_c);
56808401ef6SPierre 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\
56963a3b9bcSJacob Faibussowitsch     i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT,i_start,i_c,i_start_ghost_c);
57047c6ae99SBarry Smith 
57147c6ae99SBarry Smith       /*
57247c6ae99SBarry Smith          Only include those interpolation points that are truly
57347c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
57447c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
57547c6ae99SBarry Smith       */
57647c6ae99SBarry Smith       nc = 0;
57747c6ae99SBarry Smith       /* one left and below; or we are right on it */
578e3fbd1f4SBarry Smith       col        = (m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
579e3fbd1f4SBarry Smith       cols[nc++] = col_shift + idx_c[col];
5809566063dSJacob Faibussowitsch       PetscCall(MatPreallocateSet(row,nc,cols,dnz,onz));
58147c6ae99SBarry Smith     }
58247c6ae99SBarry Smith   }
5839566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)daf),&mat));
584fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
585fdc842d1SBarry Smith   /*
586fdc842d1SBarry Smith      Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
587fdc842d1SBarry Smith      we don't want the original unconverted matrix copied to the GPU
588fdc842d1SBarry Smith   */
589fdc842d1SBarry Smith   if (dof > 1) {
5909566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(mat,PETSC_TRUE));
591fdc842d1SBarry Smith   }
592fdc842d1SBarry Smith   #endif
5939566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My));
5949566063dSJacob Faibussowitsch   PetscCall(ConvertToAIJ(dac->mattype,&mattype));
5959566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat,mattype));
5969566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat,0,dnz));
5979566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(mat,0,dnz,0,onz));
598d0609cedSBarry Smith   MatPreallocateEnd(dnz,onz);
59947c6ae99SBarry Smith 
60047c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
60147c6ae99SBarry Smith   for (j=j_start; j<j_start+n_f; j++) {
60247c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
60347c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
604e3fbd1f4SBarry Smith       row = idx_f[(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
60547c6ae99SBarry Smith 
60647c6ae99SBarry Smith       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
60747c6ae99SBarry Smith       j_c = (j/ratioj);    /* coarse grid node below fine grid node */
60847c6ae99SBarry Smith       nc  = 0;
60947c6ae99SBarry Smith       /* one left and below; or we are right on it */
610e3fbd1f4SBarry Smith       col      = (m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
611e3fbd1f4SBarry Smith       cols[nc] = col_shift + idx_c[col];
61247c6ae99SBarry Smith       v[nc++]  = 1.0;
61347c6ae99SBarry Smith 
6149566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES));
61547c6ae99SBarry Smith     }
61647c6ae99SBarry Smith   }
6179566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f));
6189566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c));
6199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
6209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
6219566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(mat,dof,A));
6229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
6239566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(13.0*m_f*n_f));
62447c6ae99SBarry Smith   PetscFunctionReturn(0);
62547c6ae99SBarry Smith }
62647c6ae99SBarry Smith 
62747c6ae99SBarry Smith /*
62847c6ae99SBarry Smith        Contributed by Jianming Yang <jianming-yang@uiowa.edu>
62947c6ae99SBarry Smith */
630e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_3D_Q0(DM dac,DM daf,Mat *A)
63147c6ae99SBarry Smith {
6328ea3bf28SBarry Smith   PetscInt               i,j,l,i_start,j_start,l_start,m_f,n_f,p_f,Mx,My,Mz,dof;
6338ea3bf28SBarry Smith   const PetscInt         *idx_c,*idx_f;
634e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
6358ea3bf28SBarry Smith   PetscInt               m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c,nc,*dnz,*onz;
63647c6ae99SBarry 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;
63747c6ae99SBarry 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;
63847c6ae99SBarry Smith   PetscMPIInt            size_c,size_f,rank_f;
63947c6ae99SBarry Smith   PetscScalar            v[8];
64047c6ae99SBarry Smith   Mat                    mat;
641bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
642fdc842d1SBarry Smith   MatType                mattype;
64347c6ae99SBarry Smith 
64447c6ae99SBarry Smith   PetscFunctionBegin;
6459566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,NULL,&Mx,&My,&Mz,NULL,NULL,NULL,NULL,NULL,&bx,&by,&bz,NULL));
64663a3b9bcSJacob Faibussowitsch   PetscCheck(Mx,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %" PetscInt_FMT " must be positive",Mx);
64763a3b9bcSJacob Faibussowitsch   PetscCheck(My,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of y coarse grid points %" PetscInt_FMT " must be positive",My);
64863a3b9bcSJacob Faibussowitsch   PetscCheck(Mz,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of z coarse grid points %" PetscInt_FMT " must be positive",Mz);
6499566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf,NULL,&mx,&my,&mz,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
65047c6ae99SBarry Smith   ratioi = mx/Mx;
65147c6ae99SBarry Smith   ratioj = my/My;
65247c6ae99SBarry Smith   ratiol = mz/Mz;
65308401ef6SPierre Jolivet   PetscCheck(ratioi*Mx == mx,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in x");
65408401ef6SPierre Jolivet   PetscCheck(ratioj*My == my,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in y");
65508401ef6SPierre Jolivet   PetscCheck(ratiol*Mz == mz,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in z");
656*1dca8a05SBarry Smith   PetscCheck(ratioi == 2 || ratioi == 1,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 1 or 2");
657*1dca8a05SBarry Smith   PetscCheck(ratioj == 2 || ratioj == 1,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 1 or 2");
658*1dca8a05SBarry Smith   PetscCheck(ratiol == 2 || ratiol == 1,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in z must be 1 or 2");
65947c6ae99SBarry Smith 
6609566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f));
6619566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost));
6629566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf,&ltog_f));
6639566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f));
66447c6ae99SBarry Smith 
6659566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c));
6669566063dSJacob 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));
6679566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac,&ltog_c));
6689566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c));
669e3fbd1f4SBarry Smith 
67047c6ae99SBarry Smith   /*
671aa219208SBarry Smith      Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
67247c6ae99SBarry Smith      The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
67347c6ae99SBarry Smith      processors). It's effective length is hence 4 times its normal length, this is
67447c6ae99SBarry Smith      why the col_scale is multiplied by the interpolation matrix column sizes.
67547c6ae99SBarry Smith      sol_shift allows each set of 1/4 processors do its own interpolation using ITS
67647c6ae99SBarry Smith      copy of the coarse vector. A bit of a hack but you do better.
67747c6ae99SBarry Smith 
67847c6ae99SBarry Smith      In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
67947c6ae99SBarry Smith   */
6809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c));
6819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f));
6829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f));
68347c6ae99SBarry Smith   col_scale = size_f/size_c;
68447c6ae99SBarry Smith   col_shift = Mx*My*Mz*(rank_f/size_c);
68547c6ae99SBarry Smith 
686d0609cedSBarry Smith   MatPreallocateBegin(PetscObjectComm((PetscObject)daf),m_f*n_f*p_f,col_scale*m_c*n_c*p_c,dnz,onz);
68747c6ae99SBarry Smith   for (l=l_start; l<l_start+p_f; l++) {
68847c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
68947c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
69047c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
691e3fbd1f4SBarry Smith         row = idx_f[(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
69247c6ae99SBarry Smith 
69347c6ae99SBarry Smith         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
69447c6ae99SBarry Smith         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
69547c6ae99SBarry Smith         l_c = (l/ratiol);
69647c6ae99SBarry Smith 
69708401ef6SPierre 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\
69863a3b9bcSJacob Faibussowitsch     l_start %" PetscInt_FMT " l_c %" PetscInt_FMT " l_start_ghost_c %" PetscInt_FMT,l_start,l_c,l_start_ghost_c);
69908401ef6SPierre 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\
70063a3b9bcSJacob Faibussowitsch     j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT,j_start,j_c,j_start_ghost_c);
70108401ef6SPierre 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\
70263a3b9bcSJacob Faibussowitsch     i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT,i_start,i_c,i_start_ghost_c);
70347c6ae99SBarry Smith 
70447c6ae99SBarry Smith         /*
70547c6ae99SBarry Smith            Only include those interpolation points that are truly
70647c6ae99SBarry Smith            nonzero. Note this is very important for final grid lines
70747c6ae99SBarry Smith            in x and y directions; since they have no right/top neighbors
70847c6ae99SBarry Smith         */
70947c6ae99SBarry Smith         nc = 0;
71047c6ae99SBarry Smith         /* one left and below; or we are right on it */
711e3fbd1f4SBarry 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));
712e3fbd1f4SBarry Smith         cols[nc++] = col_shift + idx_c[col];
7139566063dSJacob Faibussowitsch         PetscCall(MatPreallocateSet(row,nc,cols,dnz,onz));
71447c6ae99SBarry Smith       }
71547c6ae99SBarry Smith     }
71647c6ae99SBarry Smith   }
7179566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)daf),&mat));
718fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
719fdc842d1SBarry Smith   /*
720fdc842d1SBarry Smith      Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
721fdc842d1SBarry Smith      we don't want the original unconverted matrix copied to the GPU
722fdc842d1SBarry Smith   */
723fdc842d1SBarry Smith   if (dof > 1) {
7249566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(mat,PETSC_TRUE));
725fdc842d1SBarry Smith   }
726fdc842d1SBarry Smith   #endif
7279566063dSJacob 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));
7289566063dSJacob Faibussowitsch   PetscCall(ConvertToAIJ(dac->mattype,&mattype));
7299566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat,mattype));
7309566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat,0,dnz));
7319566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(mat,0,dnz,0,onz));
732d0609cedSBarry Smith   MatPreallocateEnd(dnz,onz);
73347c6ae99SBarry Smith 
73447c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
73547c6ae99SBarry Smith   for (l=l_start; l<l_start+p_f; l++) {
73647c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
73747c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
73847c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
739e3fbd1f4SBarry Smith         row = idx_f[(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
74047c6ae99SBarry Smith 
74147c6ae99SBarry Smith         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
74247c6ae99SBarry Smith         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
74347c6ae99SBarry Smith         l_c = (l/ratiol);
74447c6ae99SBarry Smith         nc  = 0;
74547c6ae99SBarry Smith         /* one left and below; or we are right on it */
746e3fbd1f4SBarry 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));
747e3fbd1f4SBarry Smith         cols[nc] = col_shift + idx_c[col];
74847c6ae99SBarry Smith         v[nc++]  = 1.0;
74947c6ae99SBarry Smith 
7509566063dSJacob Faibussowitsch         PetscCall(MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES));
75147c6ae99SBarry Smith       }
75247c6ae99SBarry Smith     }
75347c6ae99SBarry Smith   }
7549566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f));
7559566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c));
7569566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
7579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
7589566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(mat,dof,A));
7599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
7609566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(13.0*m_f*n_f*p_f));
76147c6ae99SBarry Smith   PetscFunctionReturn(0);
76247c6ae99SBarry Smith }
76347c6ae99SBarry Smith 
764e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_3D_Q1(DM dac,DM daf,Mat *A)
76547c6ae99SBarry Smith {
7668ea3bf28SBarry Smith   PetscInt               i,j,i_start,j_start,m_f,n_f,Mx,My,dof,l;
7678ea3bf28SBarry Smith   const PetscInt         *idx_c,*idx_f;
768e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
7698ea3bf28SBarry Smith   PetscInt               m_ghost,n_ghost,m_ghost_c,n_ghost_c,Mz,mz;
77047c6ae99SBarry Smith   PetscInt               row,col,i_start_ghost,j_start_ghost,cols[8],mx,m_c,my,nc,ratioi,ratioj,ratiok;
77147c6ae99SBarry Smith   PetscInt               i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c;
77247c6ae99SBarry Smith   PetscInt               l_start,p_f,l_start_ghost,p_ghost,l_start_c,p_c;
77347c6ae99SBarry Smith   PetscInt               l_start_ghost_c,p_ghost_c,l_c,*dnz,*onz;
77447c6ae99SBarry Smith   PetscScalar            v[8],x,y,z;
77547c6ae99SBarry Smith   Mat                    mat;
776bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
777fdc842d1SBarry Smith   MatType                mattype;
77847c6ae99SBarry Smith 
77947c6ae99SBarry Smith   PetscFunctionBegin;
7809566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,NULL,&Mx,&My,&Mz,NULL,NULL,NULL,NULL,NULL,&bx,&by,&bz,NULL));
7819566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf,NULL,&mx,&my,&mz,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
78247c6ae99SBarry Smith   if (mx == Mx) {
78347c6ae99SBarry Smith     ratioi = 1;
784bff4a2f0SMatthew G. Knepley   } else if (bx == DM_BOUNDARY_PERIODIC) {
78563a3b9bcSJacob Faibussowitsch     PetscCheck(Mx,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %" PetscInt_FMT " must be positive",Mx);
78647c6ae99SBarry Smith     ratioi = mx/Mx;
78763a3b9bcSJacob 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);
78847c6ae99SBarry Smith   } else {
78963a3b9bcSJacob 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);
79047c6ae99SBarry Smith     ratioi = (mx-1)/(Mx-1);
791*1dca8a05SBarry 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);
79247c6ae99SBarry Smith   }
79347c6ae99SBarry Smith   if (my == My) {
79447c6ae99SBarry Smith     ratioj = 1;
795bff4a2f0SMatthew G. Knepley   } else if (by == DM_BOUNDARY_PERIODIC) {
79663a3b9bcSJacob Faibussowitsch     PetscCheck(My,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of y coarse grid points %" PetscInt_FMT " must be positive",My);
79747c6ae99SBarry Smith     ratioj = my/My;
79863a3b9bcSJacob 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);
79947c6ae99SBarry Smith   } else {
80063a3b9bcSJacob 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);
80147c6ae99SBarry Smith     ratioj = (my-1)/(My-1);
802*1dca8a05SBarry 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);
80347c6ae99SBarry Smith   }
80447c6ae99SBarry Smith   if (mz == Mz) {
80547c6ae99SBarry Smith     ratiok = 1;
806bff4a2f0SMatthew G. Knepley   } else if (bz == DM_BOUNDARY_PERIODIC) {
80763a3b9bcSJacob Faibussowitsch     PetscCheck(Mz,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of z coarse grid points %" PetscInt_FMT " must be positive",Mz);
80847c6ae99SBarry Smith     ratiok = mz/Mz;
80963a3b9bcSJacob 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);
81047c6ae99SBarry Smith   } else {
81163a3b9bcSJacob 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);
81247c6ae99SBarry Smith     ratiok = (mz-1)/(Mz-1);
813*1dca8a05SBarry 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);
81447c6ae99SBarry Smith   }
81547c6ae99SBarry Smith 
8169566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f));
8179566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost));
8189566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf,&ltog_f));
8199566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f));
82047c6ae99SBarry Smith 
8219566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c));
8229566063dSJacob 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));
8239566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac,&ltog_c));
8249566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c));
82547c6ae99SBarry Smith 
82647c6ae99SBarry Smith   /* create interpolation matrix, determining exact preallocation */
827d0609cedSBarry Smith   MatPreallocateBegin(PetscObjectComm((PetscObject)dac),m_f*n_f*p_f,m_c*n_c*p_c,dnz,onz);
82847c6ae99SBarry Smith   /* loop over local fine grid nodes counting interpolating points */
82947c6ae99SBarry Smith   for (l=l_start; l<l_start+p_f; l++) {
83047c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
83147c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
83247c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
833e3fbd1f4SBarry Smith         row = idx_f[(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
83447c6ae99SBarry Smith         i_c = (i/ratioi);
83547c6ae99SBarry Smith         j_c = (j/ratioj);
83647c6ae99SBarry Smith         l_c = (l/ratiok);
83708401ef6SPierre 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\
83863a3b9bcSJacob Faibussowitsch                                             l_start %" PetscInt_FMT " l_c %" PetscInt_FMT " l_start_ghost_c %" PetscInt_FMT,l_start,l_c,l_start_ghost_c);
83908401ef6SPierre 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\
84063a3b9bcSJacob Faibussowitsch                                             j_start %" PetscInt_FMT " j_c %" PetscInt_FMT " j_start_ghost_c %" PetscInt_FMT,j_start,j_c,j_start_ghost_c);
84108401ef6SPierre 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\
84263a3b9bcSJacob Faibussowitsch                                             i_start %" PetscInt_FMT " i_c %" PetscInt_FMT " i_start_ghost_c %" PetscInt_FMT,i_start,i_c,i_start_ghost_c);
84347c6ae99SBarry Smith 
84447c6ae99SBarry Smith         /*
84547c6ae99SBarry Smith          Only include those interpolation points that are truly
84647c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
84747c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
84847c6ae99SBarry Smith          */
84947c6ae99SBarry Smith         nc         = 0;
850e3fbd1f4SBarry 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));
851e3fbd1f4SBarry Smith         cols[nc++] = idx_c[col];
85247c6ae99SBarry Smith         if (i_c*ratioi != i) {
853e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+1];
85447c6ae99SBarry Smith         }
85547c6ae99SBarry Smith         if (j_c*ratioj != j) {
856e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+m_ghost_c];
85747c6ae99SBarry Smith         }
85847c6ae99SBarry Smith         if (l_c*ratiok != l) {
859e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+m_ghost_c*n_ghost_c];
86047c6ae99SBarry Smith         }
86147c6ae99SBarry Smith         if (j_c*ratioj != j && i_c*ratioi != i) {
862e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c+1)];
86347c6ae99SBarry Smith         }
86447c6ae99SBarry Smith         if (j_c*ratioj != j && l_c*ratiok != l) {
865e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)];
86647c6ae99SBarry Smith         }
86747c6ae99SBarry Smith         if (i_c*ratioi != i && l_c*ratiok != l) {
868e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+1)];
86947c6ae99SBarry Smith         }
87047c6ae99SBarry Smith         if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) {
871e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)];
87247c6ae99SBarry Smith         }
8739566063dSJacob Faibussowitsch         PetscCall(MatPreallocateSet(row,nc,cols,dnz,onz));
87447c6ae99SBarry Smith       }
87547c6ae99SBarry Smith     }
87647c6ae99SBarry Smith   }
8779566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dac),&mat));
878fdc842d1SBarry Smith #if defined(PETSC_HAVE_CUDA)
879fdc842d1SBarry Smith   /*
880fdc842d1SBarry Smith      Temporary hack: Since the MAIJ matrix must be converted to AIJ before being used by the GPU
881fdc842d1SBarry Smith      we don't want the original unconverted matrix copied to the GPU
882fdc842d1SBarry Smith   */
883fdc842d1SBarry Smith   if (dof > 1) {
8849566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(mat,PETSC_TRUE));
885fdc842d1SBarry Smith   }
886fdc842d1SBarry Smith   #endif
8879566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat,m_f*n_f*p_f,m_c*n_c*p_c,mx*my*mz,Mx*My*Mz));
8889566063dSJacob Faibussowitsch   PetscCall(ConvertToAIJ(dac->mattype,&mattype));
8899566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat,mattype));
8909566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat,0,dnz));
8919566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(mat,0,dnz,0,onz));
892d0609cedSBarry Smith   MatPreallocateEnd(dnz,onz);
89347c6ae99SBarry Smith 
89447c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
8952adb9dcfSBarry Smith   if (!NEWVERSION) {
896b3bd94feSDave May 
89747c6ae99SBarry Smith     for (l=l_start; l<l_start+p_f; l++) {
89847c6ae99SBarry Smith       for (j=j_start; j<j_start+n_f; j++) {
89947c6ae99SBarry Smith         for (i=i_start; i<i_start+m_f; i++) {
90047c6ae99SBarry Smith           /* convert to local "natural" numbering and then to PETSc global numbering */
901e3fbd1f4SBarry Smith           row = idx_f[(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
90247c6ae99SBarry Smith 
90347c6ae99SBarry Smith           i_c = (i/ratioi);
90447c6ae99SBarry Smith           j_c = (j/ratioj);
90547c6ae99SBarry Smith           l_c = (l/ratiok);
90647c6ae99SBarry Smith 
90747c6ae99SBarry Smith           /*
90847c6ae99SBarry Smith            Only include those interpolation points that are truly
90947c6ae99SBarry Smith            nonzero. Note this is very important for final grid lines
91047c6ae99SBarry Smith            in x and y directions; since they have no right/top neighbors
91147c6ae99SBarry Smith            */
9126712e2f1SBarry Smith           x = ((PetscReal)(i - i_c*ratioi))/((PetscReal)ratioi);
9136712e2f1SBarry Smith           y = ((PetscReal)(j - j_c*ratioj))/((PetscReal)ratioj);
9146712e2f1SBarry Smith           z = ((PetscReal)(l - l_c*ratiok))/((PetscReal)ratiok);
915b3bd94feSDave May 
91647c6ae99SBarry Smith           nc = 0;
91747c6ae99SBarry Smith           /* one left and below; or we are right on it */
918e3fbd1f4SBarry 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));
91947c6ae99SBarry Smith 
920e3fbd1f4SBarry Smith           cols[nc] = idx_c[col];
92147c6ae99SBarry Smith           v[nc++]  = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.));
92247c6ae99SBarry Smith 
92347c6ae99SBarry Smith           if (i_c*ratioi != i) {
924e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+1];
92547c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.));
92647c6ae99SBarry Smith           }
92747c6ae99SBarry Smith 
92847c6ae99SBarry Smith           if (j_c*ratioj != j) {
929e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+m_ghost_c];
93047c6ae99SBarry Smith             v[nc++]  = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.));
93147c6ae99SBarry Smith           }
93247c6ae99SBarry Smith 
93347c6ae99SBarry Smith           if (l_c*ratiok != l) {
934e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+m_ghost_c*n_ghost_c];
93547c6ae99SBarry Smith             v[nc++]  = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.));
93647c6ae99SBarry Smith           }
93747c6ae99SBarry Smith 
93847c6ae99SBarry Smith           if (j_c*ratioj != j && i_c*ratioi != i) {
939e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c+1)];
94047c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.));
94147c6ae99SBarry Smith           }
94247c6ae99SBarry Smith 
94347c6ae99SBarry Smith           if (j_c*ratioj != j && l_c*ratiok != l) {
944e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)];
94547c6ae99SBarry Smith             v[nc++]  = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.));
94647c6ae99SBarry Smith           }
94747c6ae99SBarry Smith 
94847c6ae99SBarry Smith           if (i_c*ratioi != i && l_c*ratiok != l) {
949e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+1)];
95047c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.));
95147c6ae99SBarry Smith           }
95247c6ae99SBarry Smith 
95347c6ae99SBarry Smith           if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) {
954e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)];
95547c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.));
95647c6ae99SBarry Smith           }
9579566063dSJacob Faibussowitsch           PetscCall(MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES));
95847c6ae99SBarry Smith         }
95947c6ae99SBarry Smith       }
96047c6ae99SBarry Smith     }
961b3bd94feSDave May 
962b3bd94feSDave May   } else {
963b3bd94feSDave May     PetscScalar *xi,*eta,*zeta;
964b3bd94feSDave May     PetscInt    li,nxi,lj,neta,lk,nzeta,n;
965b3bd94feSDave May     PetscScalar Ni[8];
966b3bd94feSDave May 
967b3bd94feSDave May     /* compute local coordinate arrays */
968b3bd94feSDave May     nxi   = ratioi + 1;
969b3bd94feSDave May     neta  = ratioj + 1;
970b3bd94feSDave May     nzeta = ratiok + 1;
9719566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nxi,&xi));
9729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(neta,&eta));
9739566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nzeta,&zeta));
9748865f1eaSKarl Rupp     for (li=0; li<nxi; li++) xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1));
9758865f1eaSKarl Rupp     for (lj=0; lj<neta; lj++) eta[lj] = -1.0 + (PetscScalar)lj*(2.0/(PetscScalar)(neta-1));
9768865f1eaSKarl Rupp     for (lk=0; lk<nzeta; lk++) zeta[lk] = -1.0 + (PetscScalar)lk*(2.0/(PetscScalar)(nzeta-1));
977b3bd94feSDave May 
978b3bd94feSDave May     for (l=l_start; l<l_start+p_f; l++) {
979b3bd94feSDave May       for (j=j_start; j<j_start+n_f; j++) {
980b3bd94feSDave May         for (i=i_start; i<i_start+m_f; i++) {
981b3bd94feSDave May           /* convert to local "natural" numbering and then to PETSc global numbering */
982e3fbd1f4SBarry Smith           row = idx_f[(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
983b3bd94feSDave May 
984b3bd94feSDave May           i_c = (i/ratioi);
985b3bd94feSDave May           j_c = (j/ratioj);
986b3bd94feSDave May           l_c = (l/ratiok);
987b3bd94feSDave May 
988b3bd94feSDave May           /* remainders */
989b3bd94feSDave May           li = i - ratioi * (i/ratioi);
9908865f1eaSKarl Rupp           if (i==mx-1) li = nxi-1;
991b3bd94feSDave May           lj = j - ratioj * (j/ratioj);
9928865f1eaSKarl Rupp           if (j==my-1) lj = neta-1;
993b3bd94feSDave May           lk = l - ratiok * (l/ratiok);
9948865f1eaSKarl Rupp           if (l==mz-1) lk = nzeta-1;
995b3bd94feSDave May 
996b3bd94feSDave May           /* corners */
997e3fbd1f4SBarry 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));
998e3fbd1f4SBarry Smith           cols[0] = idx_c[col];
999b3bd94feSDave May           Ni[0]   = 1.0;
1000b3bd94feSDave May           if ((li==0) || (li==nxi-1)) {
1001b3bd94feSDave May             if ((lj==0) || (lj==neta-1)) {
1002b3bd94feSDave May               if ((lk==0) || (lk==nzeta-1)) {
10039566063dSJacob Faibussowitsch                 PetscCall(MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES));
1004b3bd94feSDave May                 continue;
1005b3bd94feSDave May               }
1006b3bd94feSDave May             }
1007b3bd94feSDave May           }
1008b3bd94feSDave May 
1009b3bd94feSDave May           /* edges + interior */
1010b3bd94feSDave May           /* remainders */
10118865f1eaSKarl Rupp           if (i==mx-1) i_c--;
10128865f1eaSKarl Rupp           if (j==my-1) j_c--;
10138865f1eaSKarl Rupp           if (l==mz-1) l_c--;
1014b3bd94feSDave May 
1015e3fbd1f4SBarry 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));
1016e3fbd1f4SBarry Smith           cols[0] = idx_c[col]; /* one left and below; or we are right on it */
1017e3fbd1f4SBarry Smith           cols[1] = idx_c[col+1]; /* one right and below */
1018e3fbd1f4SBarry Smith           cols[2] = idx_c[col+m_ghost_c];  /* one left and above */
1019e3fbd1f4SBarry Smith           cols[3] = idx_c[col+(m_ghost_c+1)]; /* one right and above */
1020b3bd94feSDave May 
1021e3fbd1f4SBarry Smith           cols[4] = idx_c[col+m_ghost_c*n_ghost_c]; /* one left and below and front; or we are right on it */
1022e3fbd1f4SBarry Smith           cols[5] = idx_c[col+(m_ghost_c*n_ghost_c+1)]; /* one right and below, and front */
1023e3fbd1f4SBarry Smith           cols[6] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)]; /* one left and above and front*/
1024e3fbd1f4SBarry Smith           cols[7] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)]; /* one right and above and front */
1025b3bd94feSDave May 
1026b3bd94feSDave May           Ni[0] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]);
1027b3bd94feSDave May           Ni[1] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]);
1028b3bd94feSDave May           Ni[2] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]);
1029b3bd94feSDave May           Ni[3] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]);
1030b3bd94feSDave May 
1031b3bd94feSDave May           Ni[4] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]);
1032b3bd94feSDave May           Ni[5] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]);
1033b3bd94feSDave May           Ni[6] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]);
1034b3bd94feSDave May           Ni[7] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]);
1035b3bd94feSDave May 
1036b3bd94feSDave May           for (n=0; n<8; n++) {
10378865f1eaSKarl Rupp             if (PetscAbsScalar(Ni[n])<1.0e-32) cols[n]=-1;
1038b3bd94feSDave May           }
10399566063dSJacob Faibussowitsch           PetscCall(MatSetValues(mat,1,&row,8,cols,Ni,INSERT_VALUES));
1040b3bd94feSDave May 
1041b3bd94feSDave May         }
1042b3bd94feSDave May       }
1043b3bd94feSDave May     }
10449566063dSJacob Faibussowitsch     PetscCall(PetscFree(xi));
10459566063dSJacob Faibussowitsch     PetscCall(PetscFree(eta));
10469566063dSJacob Faibussowitsch     PetscCall(PetscFree(zeta));
1047b3bd94feSDave May   }
10489566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f));
10499566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c));
10509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
10519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
105247c6ae99SBarry Smith 
10539566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(mat,dof,A));
10549566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
105547c6ae99SBarry Smith   PetscFunctionReturn(0);
105647c6ae99SBarry Smith }
105747c6ae99SBarry Smith 
1058e727c939SJed Brown PetscErrorCode  DMCreateInterpolation_DA(DM dac,DM daf,Mat *A,Vec *scale)
105947c6ae99SBarry Smith {
106047c6ae99SBarry Smith   PetscInt         dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
1061bff4a2f0SMatthew G. Knepley   DMBoundaryType   bxc,byc,bzc,bxf,byf,bzf;
1062aa219208SBarry Smith   DMDAStencilType  stc,stf;
106347c6ae99SBarry Smith   DM_DA            *ddc = (DM_DA*)dac->data;
106447c6ae99SBarry Smith 
106547c6ae99SBarry Smith   PetscFunctionBegin;
106647c6ae99SBarry Smith   PetscValidHeaderSpecific(dac,DM_CLASSID,1);
106747c6ae99SBarry Smith   PetscValidHeaderSpecific(daf,DM_CLASSID,2);
106847c6ae99SBarry Smith   PetscValidPointer(A,3);
106947c6ae99SBarry Smith   if (scale) PetscValidPointer(scale,4);
107047c6ae99SBarry Smith 
10719566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc));
10729566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf));
107363a3b9bcSJacob Faibussowitsch   PetscCheck(dimc == dimf,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT,dimc,dimf);
107463a3b9bcSJacob Faibussowitsch   PetscCheck(dofc == doff,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT,dofc,doff);
107563a3b9bcSJacob Faibussowitsch   PetscCheck(sc == sf,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT,sc,sf);
1076*1dca8a05SBarry Smith   PetscCheck(bxc == bxf && byc == byf && bzc == bzf,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");
107708401ef6SPierre Jolivet   PetscCheck(stc == stf,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");
1078*1dca8a05SBarry Smith   PetscCheck(Mc >= 2 || Mf <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction");
1079*1dca8a05SBarry 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");
1080*1dca8a05SBarry 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");
108147c6ae99SBarry Smith 
1082aa219208SBarry Smith   if (ddc->interptype == DMDA_Q1) {
108347c6ae99SBarry Smith     if (dimc == 1) {
10849566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation_DA_1D_Q1(dac,daf,A));
108547c6ae99SBarry Smith     } else if (dimc == 2) {
10869566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation_DA_2D_Q1(dac,daf,A));
108747c6ae99SBarry Smith     } else if (dimc == 3) {
10889566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation_DA_3D_Q1(dac,daf,A));
108963a3b9bcSJacob 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);
1090aa219208SBarry Smith   } else if (ddc->interptype == DMDA_Q0) {
109147c6ae99SBarry Smith     if (dimc == 1) {
10929566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation_DA_1D_Q0(dac,daf,A));
109347c6ae99SBarry Smith     } else if (dimc == 2) {
10949566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation_DA_2D_Q0(dac,daf,A));
109547c6ae99SBarry Smith     } else if (dimc == 3) {
10969566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation_DA_3D_Q0(dac,daf,A));
109763a3b9bcSJacob 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);
109847c6ae99SBarry Smith   }
109947c6ae99SBarry Smith   if (scale) {
11009566063dSJacob Faibussowitsch     PetscCall(DMCreateInterpolationScale((DM)dac,(DM)daf,*A,scale));
110147c6ae99SBarry Smith   }
110247c6ae99SBarry Smith   PetscFunctionReturn(0);
110347c6ae99SBarry Smith }
110447c6ae99SBarry Smith 
1105e727c939SJed Brown PetscErrorCode DMCreateInjection_DA_1D(DM dac,DM daf,VecScatter *inject)
11060bb2b966SJungho Lee {
11078ea3bf28SBarry Smith   PetscInt               i,i_start,m_f,Mx,dof;
11088ea3bf28SBarry Smith   const PetscInt         *idx_f;
1109e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f;
11108ea3bf28SBarry Smith   PetscInt               m_ghost,m_ghost_c;
11110bb2b966SJungho Lee   PetscInt               row,i_start_ghost,mx,m_c,nc,ratioi;
11120bb2b966SJungho Lee   PetscInt               i_start_c,i_start_ghost_c;
11130bb2b966SJungho Lee   PetscInt               *cols;
1114bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
11150bb2b966SJungho Lee   Vec                    vecf,vecc;
11160bb2b966SJungho Lee   IS                     isf;
11170bb2b966SJungho Lee 
11180bb2b966SJungho Lee   PetscFunctionBegin;
11199566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,NULL,&Mx,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bx,NULL,NULL,NULL));
11209566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf,NULL,&mx,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
1121bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
11220bb2b966SJungho Lee     ratioi = mx/Mx;
112363a3b9bcSJacob 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);
11240bb2b966SJungho Lee   } else {
11250bb2b966SJungho Lee     ratioi = (mx-1)/(Mx-1);
1126*1dca8a05SBarry 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);
11270bb2b966SJungho Lee   }
11280bb2b966SJungho Lee 
11299566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf,&i_start,NULL,NULL,&m_f,NULL,NULL));
11309566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf,&i_start_ghost,NULL,NULL,&m_ghost,NULL,NULL));
11319566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf,&ltog_f));
11329566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f));
11330bb2b966SJungho Lee 
11349566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac,&i_start_c,NULL,NULL,&m_c,NULL,NULL));
11359566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac,&i_start_ghost_c,NULL,NULL,&m_ghost_c,NULL,NULL));
11360bb2b966SJungho Lee 
11370bb2b966SJungho Lee   /* loop over local fine grid nodes setting interpolation for those*/
11380bb2b966SJungho Lee   nc   = 0;
11399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m_f,&cols));
11400bb2b966SJungho Lee 
11410bb2b966SJungho Lee   for (i=i_start_c; i<i_start_c+m_c; i++) {
11420bb2b966SJungho Lee     PetscInt i_f = i*ratioi;
11430bb2b966SJungho Lee 
1144*1dca8a05SBarry 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);
11458865f1eaSKarl Rupp 
1146e3fbd1f4SBarry Smith     row        = idx_f[(i_f-i_start_ghost)];
1147e3fbd1f4SBarry Smith     cols[nc++] = row;
11480bb2b966SJungho Lee   }
11490bb2b966SJungho Lee 
11509566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f));
11519566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf));
11529566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dac,&vecc));
11539566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(daf,&vecf));
11549566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(vecf,isf,vecc,NULL,inject));
11559566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dac,&vecc));
11569566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(daf,&vecf));
11579566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isf));
11580bb2b966SJungho Lee   PetscFunctionReturn(0);
11590bb2b966SJungho Lee }
11600bb2b966SJungho Lee 
1161e727c939SJed Brown PetscErrorCode DMCreateInjection_DA_2D(DM dac,DM daf,VecScatter *inject)
116247c6ae99SBarry Smith {
11638ea3bf28SBarry Smith   PetscInt               i,j,i_start,j_start,m_f,n_f,Mx,My,dof;
11648ea3bf28SBarry Smith   const PetscInt         *idx_c,*idx_f;
1165e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
11668ea3bf28SBarry Smith   PetscInt               m_ghost,n_ghost,m_ghost_c,n_ghost_c;
116747c6ae99SBarry Smith   PetscInt               row,i_start_ghost,j_start_ghost,mx,m_c,my,nc,ratioi,ratioj;
1168076202ddSJed Brown   PetscInt               i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c;
116947c6ae99SBarry Smith   PetscInt               *cols;
1170bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
117147c6ae99SBarry Smith   Vec                    vecf,vecc;
117247c6ae99SBarry Smith   IS                     isf;
117347c6ae99SBarry Smith 
117447c6ae99SBarry Smith   PetscFunctionBegin;
11759566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,NULL,&Mx,&My,NULL,NULL,NULL,NULL,NULL,NULL,&bx,&by,NULL,NULL));
11769566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf,NULL,&mx,&my,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
1177bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
117847c6ae99SBarry Smith     ratioi = mx/Mx;
117963a3b9bcSJacob 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);
118047c6ae99SBarry Smith   } else {
118147c6ae99SBarry Smith     ratioi = (mx-1)/(Mx-1);
1182*1dca8a05SBarry 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);
118347c6ae99SBarry Smith   }
1184bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) {
118547c6ae99SBarry Smith     ratioj = my/My;
118663a3b9bcSJacob 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);
118747c6ae99SBarry Smith   } else {
118847c6ae99SBarry Smith     ratioj = (my-1)/(My-1);
1189*1dca8a05SBarry 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);
119047c6ae99SBarry Smith   }
119147c6ae99SBarry Smith 
11929566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf,&i_start,&j_start,NULL,&m_f,&n_f,NULL));
11939566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,NULL,&m_ghost,&n_ghost,NULL));
11949566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf,&ltog_f));
11959566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f));
119647c6ae99SBarry Smith 
11979566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac,&i_start_c,&j_start_c,NULL,&m_c,&n_c,NULL));
11989566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,NULL,&m_ghost_c,&n_ghost_c,NULL));
11999566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac,&ltog_c));
12009566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c));
120147c6ae99SBarry Smith 
120247c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
120347c6ae99SBarry Smith   nc   = 0;
12049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_f*m_f,&cols));
1205076202ddSJed Brown   for (j=j_start_c; j<j_start_c+n_c; j++) {
1206076202ddSJed Brown     for (i=i_start_c; i<i_start_c+m_c; i++) {
1207076202ddSJed Brown       PetscInt i_f = i*ratioi,j_f = j*ratioj;
1208*1dca8a05SBarry 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\
120963a3b9bcSJacob Faibussowitsch     j_c %" PetscInt_FMT " j_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",j,j_f,j_start_ghost,j_start_ghost+n_ghost);
1210*1dca8a05SBarry 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\
121163a3b9bcSJacob Faibussowitsch     i_c %" PetscInt_FMT " i_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",i,i_f,i_start_ghost,i_start_ghost+m_ghost);
1212e3fbd1f4SBarry Smith       row        = idx_f[(m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))];
1213e3fbd1f4SBarry Smith       cols[nc++] = row;
121447c6ae99SBarry Smith     }
121547c6ae99SBarry Smith   }
12169566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f));
12179566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c));
121847c6ae99SBarry Smith 
12199566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf));
12209566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dac,&vecc));
12219566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(daf,&vecf));
12229566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(vecf,isf,vecc,NULL,inject));
12239566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dac,&vecc));
12249566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(daf,&vecf));
12259566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isf));
122647c6ae99SBarry Smith   PetscFunctionReturn(0);
122747c6ae99SBarry Smith }
122847c6ae99SBarry Smith 
1229e727c939SJed Brown PetscErrorCode DMCreateInjection_DA_3D(DM dac,DM daf,VecScatter *inject)
1230b1757049SJed Brown {
1231b1757049SJed Brown   PetscInt               i,j,k,i_start,j_start,k_start,m_f,n_f,p_f,Mx,My,Mz;
1232b1757049SJed Brown   PetscInt               m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c;
1233b1757049SJed Brown   PetscInt               i_start_ghost,j_start_ghost,k_start_ghost;
1234b1757049SJed Brown   PetscInt               mx,my,mz,ratioi,ratioj,ratiok;
1235b1757049SJed Brown   PetscInt               i_start_c,j_start_c,k_start_c;
1236b1757049SJed Brown   PetscInt               m_c,n_c,p_c;
1237b1757049SJed Brown   PetscInt               i_start_ghost_c,j_start_ghost_c,k_start_ghost_c;
1238b1757049SJed Brown   PetscInt               row,nc,dof;
12398ea3bf28SBarry Smith   const PetscInt         *idx_c,*idx_f;
1240e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
1241b1757049SJed Brown   PetscInt               *cols;
1242bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
1243b1757049SJed Brown   Vec                    vecf,vecc;
1244b1757049SJed Brown   IS                     isf;
1245b1757049SJed Brown 
1246b1757049SJed Brown   PetscFunctionBegin;
12479566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,NULL,&Mx,&My,&Mz,NULL,NULL,NULL,NULL,NULL,&bx,&by,&bz,NULL));
12489566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf,NULL,&mx,&my,&mz,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
1249b1757049SJed Brown 
1250bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
1251b1757049SJed Brown     ratioi = mx/Mx;
125263a3b9bcSJacob 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);
1253b1757049SJed Brown   } else {
1254b1757049SJed Brown     ratioi = (mx-1)/(Mx-1);
1255*1dca8a05SBarry 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);
1256b1757049SJed Brown   }
1257bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) {
1258b1757049SJed Brown     ratioj = my/My;
125963a3b9bcSJacob 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);
1260b1757049SJed Brown   } else {
1261b1757049SJed Brown     ratioj = (my-1)/(My-1);
1262*1dca8a05SBarry 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);
1263b1757049SJed Brown   }
1264bff4a2f0SMatthew G. Knepley   if (bz == DM_BOUNDARY_PERIODIC) {
1265b1757049SJed Brown     ratiok = mz/Mz;
126663a3b9bcSJacob 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);
1267b1757049SJed Brown   } else {
1268b1757049SJed Brown     ratiok = (mz-1)/(Mz-1);
1269*1dca8a05SBarry 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);
1270b1757049SJed Brown   }
1271b1757049SJed Brown 
12729566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf,&i_start,&j_start,&k_start,&m_f,&n_f,&p_f));
12739566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&k_start_ghost,&m_ghost,&n_ghost,&p_ghost));
12749566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf,&ltog_f));
12759566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f));
1276b1757049SJed Brown 
12779566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac,&i_start_c,&j_start_c,&k_start_c,&m_c,&n_c,&p_c));
12789566063dSJacob 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));
12799566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac,&ltog_c));
12809566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c));
1281b1757049SJed Brown 
1282b1757049SJed Brown   /* loop over local fine grid nodes setting interpolation for those*/
1283b1757049SJed Brown   nc   = 0;
12849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_f*m_f*p_f,&cols));
1285b1757049SJed Brown   for (k=k_start_c; k<k_start_c+p_c; k++) {
1286b1757049SJed Brown     for (j=j_start_c; j<j_start_c+n_c; j++) {
1287b1757049SJed Brown       for (i=i_start_c; i<i_start_c+m_c; i++) {
1288b1757049SJed Brown         PetscInt i_f = i*ratioi,j_f = j*ratioj,k_f = k*ratiok;
1289*1dca8a05SBarry Smith         PetscCheck(k_f >= k_start_ghost && k_f < k_start_ghost+p_ghost,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA  "
129063a3b9bcSJacob Faibussowitsch                                                                           "k_c %" PetscInt_FMT " k_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",k,k_f,k_start_ghost,k_start_ghost+p_ghost);
1291*1dca8a05SBarry 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  "
129263a3b9bcSJacob Faibussowitsch                                                                           "j_c %" PetscInt_FMT " j_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",j,j_f,j_start_ghost,j_start_ghost+n_ghost);
1293*1dca8a05SBarry 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  "
129463a3b9bcSJacob Faibussowitsch                                                                           "i_c %" PetscInt_FMT " i_f %" PetscInt_FMT " fine ghost range [%" PetscInt_FMT ",%" PetscInt_FMT "]",i,i_f,i_start_ghost,i_start_ghost+m_ghost);
1295e3fbd1f4SBarry 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))];
1296e3fbd1f4SBarry Smith         cols[nc++] = row;
1297b1757049SJed Brown       }
1298b1757049SJed Brown     }
1299b1757049SJed Brown   }
13009566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f));
13019566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c));
1302b1757049SJed Brown 
13039566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf));
13049566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dac,&vecc));
13059566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(daf,&vecf));
13069566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(vecf,isf,vecc,NULL,inject));
13079566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dac,&vecc));
13089566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(daf,&vecf));
13099566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isf));
1310b1757049SJed Brown   PetscFunctionReturn(0);
1311b1757049SJed Brown }
1312b1757049SJed Brown 
13136dbf9973SLawrence Mitchell PetscErrorCode  DMCreateInjection_DA(DM dac,DM daf,Mat *mat)
131447c6ae99SBarry Smith {
131547c6ae99SBarry Smith   PetscInt        dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
1316bff4a2f0SMatthew G. Knepley   DMBoundaryType  bxc,byc,bzc,bxf,byf,bzf;
1317aa219208SBarry Smith   DMDAStencilType stc,stf;
131860c16ac1SBarry Smith   VecScatter      inject = NULL;
131947c6ae99SBarry Smith 
132047c6ae99SBarry Smith   PetscFunctionBegin;
132147c6ae99SBarry Smith   PetscValidHeaderSpecific(dac,DM_CLASSID,1);
132247c6ae99SBarry Smith   PetscValidHeaderSpecific(daf,DM_CLASSID,2);
13236dbf9973SLawrence Mitchell   PetscValidPointer(mat,3);
132447c6ae99SBarry Smith 
13259566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc));
13269566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf));
132763a3b9bcSJacob Faibussowitsch   PetscCheck(dimc == dimf,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT,dimc,dimf);
132863a3b9bcSJacob Faibussowitsch   PetscCheck(dofc == doff,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT,dofc,doff);
132963a3b9bcSJacob Faibussowitsch   PetscCheck(sc == sf,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT,sc,sf);
1330*1dca8a05SBarry Smith   PetscCheck(bxc == bxf && byc == byf && bzc == bzf,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");
133108401ef6SPierre Jolivet   PetscCheck(stc == stf,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");
133208401ef6SPierre Jolivet   PetscCheck(Mc >= 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction");
1333*1dca8a05SBarry Smith   PetscCheck(dimc <= 1 || Nc >= 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in y direction");
1334*1dca8a05SBarry Smith   PetscCheck(dimc <= 2 || Pc >= 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in z direction");
133547c6ae99SBarry Smith 
13360bb2b966SJungho Lee   if (dimc == 1) {
13379566063dSJacob Faibussowitsch     PetscCall(DMCreateInjection_DA_1D(dac,daf,&inject));
13380bb2b966SJungho Lee   } else if (dimc == 2) {
13399566063dSJacob Faibussowitsch     PetscCall(DMCreateInjection_DA_2D(dac,daf,&inject));
1340b1757049SJed Brown   } else if (dimc == 3) {
13419566063dSJacob Faibussowitsch     PetscCall(DMCreateInjection_DA_3D(dac,daf,&inject));
134247c6ae99SBarry Smith   }
13439566063dSJacob Faibussowitsch   PetscCall(MatCreateScatter(PetscObjectComm((PetscObject)inject), inject, mat));
13449566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&inject));
134547c6ae99SBarry Smith   PetscFunctionReturn(0);
134647c6ae99SBarry Smith }
134747c6ae99SBarry Smith 
134897779f9aSLisandro Dalcin /*@
134997779f9aSLisandro Dalcin    DMCreateAggregates - Deprecated, see DMDACreateAggregates.
13506c877ef6SSatish Balay 
13516c877ef6SSatish Balay    Level: intermediate
135297779f9aSLisandro Dalcin @*/
1353a5bc1bf3SBarry Smith PetscErrorCode DMCreateAggregates(DM dac,DM daf,Mat *mat)
135497779f9aSLisandro Dalcin {
1355a5bc1bf3SBarry Smith   return DMDACreateAggregates(dac,daf,mat);
135697779f9aSLisandro Dalcin }
135797779f9aSLisandro Dalcin 
135897779f9aSLisandro Dalcin /*@
135997779f9aSLisandro Dalcin    DMDACreateAggregates - Gets the aggregates that map between
136097779f9aSLisandro Dalcin    grids associated with two DMDAs.
136197779f9aSLisandro Dalcin 
136297779f9aSLisandro Dalcin    Collective on dmc
136397779f9aSLisandro Dalcin 
136497779f9aSLisandro Dalcin    Input Parameters:
136597779f9aSLisandro Dalcin +  dmc - the coarse grid DMDA
136697779f9aSLisandro Dalcin -  dmf - the fine grid DMDA
136797779f9aSLisandro Dalcin 
136897779f9aSLisandro Dalcin    Output Parameters:
136997779f9aSLisandro Dalcin .  rest - the restriction matrix (transpose of the projection matrix)
137097779f9aSLisandro Dalcin 
137197779f9aSLisandro Dalcin    Level: intermediate
137297779f9aSLisandro Dalcin 
137397779f9aSLisandro Dalcin    Note: This routine is not used by PETSc.
137497779f9aSLisandro Dalcin    It is not clear what its use case is and it may be removed in a future release.
137597779f9aSLisandro Dalcin    Users should contact petsc-maint@mcs.anl.gov if they plan to use it.
137697779f9aSLisandro Dalcin 
137797779f9aSLisandro Dalcin .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
137897779f9aSLisandro Dalcin @*/
137997779f9aSLisandro Dalcin PetscErrorCode DMDACreateAggregates(DM dac,DM daf,Mat *rest)
138047c6ae99SBarry Smith {
138147c6ae99SBarry Smith   PetscInt               dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc;
138247c6ae99SBarry Smith   PetscInt               dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
1383bff4a2f0SMatthew G. Knepley   DMBoundaryType         bxc,byc,bzc,bxf,byf,bzf;
1384aa219208SBarry Smith   DMDAStencilType        stc,stf;
138547c6ae99SBarry Smith   PetscInt               i,j,l;
138647c6ae99SBarry Smith   PetscInt               i_start,j_start,l_start, m_f,n_f,p_f;
138747c6ae99SBarry Smith   PetscInt               i_start_ghost,j_start_ghost,l_start_ghost,m_ghost,n_ghost,p_ghost;
13888ea3bf28SBarry Smith   const PetscInt         *idx_f;
138947c6ae99SBarry Smith   PetscInt               i_c,j_c,l_c;
139047c6ae99SBarry Smith   PetscInt               i_start_c,j_start_c,l_start_c, m_c,n_c,p_c;
139147c6ae99SBarry Smith   PetscInt               i_start_ghost_c,j_start_ghost_c,l_start_ghost_c,m_ghost_c,n_ghost_c,p_ghost_c;
13928ea3bf28SBarry Smith   const PetscInt         *idx_c;
139347c6ae99SBarry Smith   PetscInt               d;
139447c6ae99SBarry Smith   PetscInt               a;
139547c6ae99SBarry Smith   PetscInt               max_agg_size;
139647c6ae99SBarry Smith   PetscInt               *fine_nodes;
139747c6ae99SBarry Smith   PetscScalar            *one_vec;
139847c6ae99SBarry Smith   PetscInt               fn_idx;
1399565245c5SBarry Smith   ISLocalToGlobalMapping ltogmf,ltogmc;
140047c6ae99SBarry Smith 
140147c6ae99SBarry Smith   PetscFunctionBegin;
140297779f9aSLisandro Dalcin   PetscValidHeaderSpecificType(dac,DM_CLASSID,1,DMDA);
140397779f9aSLisandro Dalcin   PetscValidHeaderSpecificType(daf,DM_CLASSID,2,DMDA);
140447c6ae99SBarry Smith   PetscValidPointer(rest,3);
140547c6ae99SBarry Smith 
14069566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc));
14079566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf));
140863a3b9bcSJacob Faibussowitsch   PetscCheck(dimc == dimf,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT,dimc,dimf);
140963a3b9bcSJacob Faibussowitsch   PetscCheck(dofc == doff,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT,dofc,doff);
141063a3b9bcSJacob Faibussowitsch   PetscCheck(sc == sf,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %" PetscInt_FMT " %" PetscInt_FMT,sc,sf);
1411*1dca8a05SBarry Smith   PetscCheck(bxc == bxf && byc == byf && bzc == bzf,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");
141208401ef6SPierre Jolivet   PetscCheck(stc == stf,PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");
141347c6ae99SBarry Smith 
141463a3b9bcSJacob 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);
141563a3b9bcSJacob 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);
141663a3b9bcSJacob 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);
141747c6ae99SBarry Smith 
141847c6ae99SBarry Smith   if (Pc < 0) Pc = 1;
141947c6ae99SBarry Smith   if (Pf < 0) Pf = 1;
142047c6ae99SBarry Smith   if (Nc < 0) Nc = 1;
142147c6ae99SBarry Smith   if (Nf < 0) Nf = 1;
142247c6ae99SBarry Smith 
14239566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f));
14249566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost));
1425565245c5SBarry Smith 
14269566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(daf,&ltogmf));
14279566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetIndices(ltogmf,&idx_f));
142847c6ae99SBarry Smith 
14299566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c));
14309566063dSJacob 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));
1431565245c5SBarry Smith 
14329566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dac,&ltogmc));
14339566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetIndices(ltogmc,&idx_c));
143447c6ae99SBarry Smith 
143547c6ae99SBarry Smith   /*
143647c6ae99SBarry Smith      Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios
143747c6ae99SBarry Smith      for dimension 1 and 2 respectively.
143847c6ae99SBarry Smith      Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1)
143947c6ae99SBarry Smith      and r_y*j and r_y*(j+1) will be grouped into the same coarse grid agregate.
144047c6ae99SBarry Smith      Each specific dof on the fine grid is mapped to one dof on the coarse grid.
144147c6ae99SBarry Smith   */
144247c6ae99SBarry Smith 
144347c6ae99SBarry Smith   max_agg_size = (Mf/Mc+1)*(Nf/Nc+1)*(Pf/Pc+1);
144447c6ae99SBarry Smith 
144547c6ae99SBarry Smith   /* create the matrix that will contain the restriction operator */
1446d0609cedSBarry Smith   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,
1447d0609cedSBarry Smith                          max_agg_size, NULL, max_agg_size, NULL, rest));
144847c6ae99SBarry Smith 
144947c6ae99SBarry Smith   /* store nodes in the fine grid here */
14509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(max_agg_size, &one_vec,max_agg_size, &fine_nodes));
145147c6ae99SBarry Smith   for (i=0; i<max_agg_size; i++) one_vec[i] = 1.0;
145247c6ae99SBarry Smith 
145347c6ae99SBarry Smith   /* loop over all coarse nodes */
145447c6ae99SBarry Smith   for (l_c=l_start_c; l_c<l_start_c+p_c; l_c++) {
145547c6ae99SBarry Smith     for (j_c=j_start_c; j_c<j_start_c+n_c; j_c++) {
145647c6ae99SBarry Smith       for (i_c=i_start_c; i_c<i_start_c+m_c; i_c++) {
145747c6ae99SBarry Smith         for (d=0; d<dofc; d++) {
145847c6ae99SBarry Smith           /* convert to local "natural" numbering and then to PETSc global numbering */
145947c6ae99SBarry 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;
146047c6ae99SBarry Smith 
146147c6ae99SBarry Smith           fn_idx = 0;
146247c6ae99SBarry Smith           /* Corresponding fine points are all points (i_f, j_f, l_f) such that
146347c6ae99SBarry Smith              i_c*Mf/Mc <= i_f < (i_c+1)*Mf/Mc
146447c6ae99SBarry Smith              (same for other dimensions)
146547c6ae99SBarry Smith           */
146647c6ae99SBarry Smith           for (l=l_c*Pf/Pc; l<PetscMin((l_c+1)*Pf/Pc,Pf); l++) {
146747c6ae99SBarry Smith             for (j=j_c*Nf/Nc; j<PetscMin((j_c+1)*Nf/Nc,Nf); j++) {
146847c6ae99SBarry Smith               for (i=i_c*Mf/Mc; i<PetscMin((i_c+1)*Mf/Mc,Mf); i++) {
146947c6ae99SBarry 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;
147047c6ae99SBarry Smith                 fn_idx++;
147147c6ae99SBarry Smith               }
147247c6ae99SBarry Smith             }
147347c6ae99SBarry Smith           }
147447c6ae99SBarry Smith           /* add all these points to one aggregate */
14759566063dSJacob Faibussowitsch           PetscCall(MatSetValues(*rest, 1, &a, fn_idx, fine_nodes, one_vec, INSERT_VALUES));
147647c6ae99SBarry Smith         }
147747c6ae99SBarry Smith       }
147847c6ae99SBarry Smith     }
147947c6ae99SBarry Smith   }
14809566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogmf,&idx_f));
14819566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogmc,&idx_c));
14829566063dSJacob Faibussowitsch   PetscCall(PetscFree2(one_vec,fine_nodes));
14839566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*rest, MAT_FINAL_ASSEMBLY));
14849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*rest, MAT_FINAL_ASSEMBLY));
148547c6ae99SBarry Smith   PetscFunctionReturn(0);
148647c6ae99SBarry Smith }
1487