xref: /petsc/src/dm/impls/da/dainterp.c (revision 8ea3bf28dab5b469bf97b96160a139e62b44e343)
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 
144035e84dSBarry Smith #include <petsc-private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith #undef __FUNCT__
17e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolationScale"
1847c6ae99SBarry Smith /*@
19e727c939SJed Brown     DMCreateInterpolationScale - Forms L = R*1/diag(R*1) - L.*v is like a coarse grid average of the
2047c6ae99SBarry Smith       nearby fine grid points.
2147c6ae99SBarry Smith 
2247c6ae99SBarry Smith   Input Parameters:
2347c6ae99SBarry Smith +      dac - DM that defines a coarse mesh
2447c6ae99SBarry Smith .      daf - DM that defines a fine mesh
2547c6ae99SBarry Smith -      mat - the restriction (or interpolation operator) from fine to coarse
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith   Output Parameter:
2847c6ae99SBarry Smith .    scale - the scaled vector
2947c6ae99SBarry Smith 
3047c6ae99SBarry Smith   Level: developer
3147c6ae99SBarry Smith 
32e727c939SJed Brown .seealso: DMCreateInterpolation()
3347c6ae99SBarry Smith 
3447c6ae99SBarry Smith @*/
35e727c939SJed Brown PetscErrorCode  DMCreateInterpolationScale(DM dac,DM daf,Mat mat,Vec *scale)
3647c6ae99SBarry Smith {
3747c6ae99SBarry Smith   PetscErrorCode ierr;
3847c6ae99SBarry Smith   Vec            fine;
3947c6ae99SBarry Smith   PetscScalar    one = 1.0;
4047c6ae99SBarry Smith 
4147c6ae99SBarry Smith   PetscFunctionBegin;
4247c6ae99SBarry Smith   ierr = DMCreateGlobalVector(daf,&fine);CHKERRQ(ierr);
4347c6ae99SBarry Smith   ierr = DMCreateGlobalVector(dac,scale);CHKERRQ(ierr);
4447c6ae99SBarry Smith   ierr = VecSet(fine,one);CHKERRQ(ierr);
4547c6ae99SBarry Smith   ierr = MatRestrict(mat,fine,*scale);CHKERRQ(ierr);
46fcfd50ebSBarry Smith   ierr = VecDestroy(&fine);CHKERRQ(ierr);
4747c6ae99SBarry Smith   ierr = VecReciprocal(*scale);CHKERRQ(ierr);
4847c6ae99SBarry Smith   PetscFunctionReturn(0);
4947c6ae99SBarry Smith }
5047c6ae99SBarry Smith 
5147c6ae99SBarry Smith #undef __FUNCT__
52e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA_1D_Q1"
53e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_1D_Q1(DM dac,DM daf,Mat *A)
5447c6ae99SBarry Smith {
5547c6ae99SBarry Smith   PetscErrorCode   ierr;
56*8ea3bf28SBarry Smith   PetscInt         i,i_start,m_f,Mx;
57*8ea3bf28SBarry Smith   const PetscInt   *idx_f,*idx_c;
58*8ea3bf28SBarry Smith   PetscInt         m_ghost,m_ghost_c;
5947c6ae99SBarry Smith   PetscInt         row,col,i_start_ghost,mx,m_c,nc,ratio;
6047c6ae99SBarry Smith   PetscInt         i_c,i_start_c,i_start_ghost_c,cols[2],dof;
6185afcc9aSBarry Smith   PetscScalar      v[2],x;
6247c6ae99SBarry Smith   Mat              mat;
631321219cSEthan Coon   DMDABoundaryType bx;
6447c6ae99SBarry Smith 
6547c6ae99SBarry Smith   PetscFunctionBegin;
661321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr);
671321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
681321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) {
6947c6ae99SBarry Smith     ratio = mx/Mx;
7047c6ae99SBarry Smith     if (ratio*Mx != mx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: mx/Mx  must be integer: mx %D Mx %D",mx,Mx);
7147c6ae99SBarry Smith   } else {
7247c6ae99SBarry Smith     ratio = (mx-1)/(Mx-1);
7347c6ae99SBarry Smith     if (ratio*(Mx-1) != mx-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %D Mx %D",mx,Mx);
7447c6ae99SBarry Smith   }
7547c6ae99SBarry Smith 
76aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr);
77aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr);
780298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
7947c6ae99SBarry Smith 
80aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr);
81aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr);
820298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
8347c6ae99SBarry Smith 
8447c6ae99SBarry Smith   /* create interpolation matrix */
85ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)dac),&mat);CHKERRQ(ierr);
8647c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr);
8747c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
880298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,2,NULL);CHKERRQ(ierr);
890298fd71SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,2,NULL,1,NULL);CHKERRQ(ierr);
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 */
9647c6ae99SBarry Smith       row = idx_f[dof*(i-i_start_ghost)]/dof;
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith       i_c = (i/ratio);    /* coarse grid node to left of fine grid node */
99aa219208SBarry Smith       if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
10047c6ae99SBarry Smith                                           i_start %D i_c %D i_start_ghost_c %D",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        */
10747c6ae99SBarry Smith       x  = ((double)(i - i_c*ratio))/((double)ratio);
10847c6ae99SBarry Smith       nc = 0;
10947c6ae99SBarry Smith       /* one left and below; or we are right on it */
11047c6ae99SBarry Smith       col      = dof*(i_c-i_start_ghost_c);
11147c6ae99SBarry Smith       cols[nc] = idx_c[col]/dof;
11247c6ae99SBarry Smith       v[nc++]  = -x + 1.0;
11347c6ae99SBarry Smith       /* one right? */
11447c6ae99SBarry Smith       if (i_c*ratio != i) {
11547c6ae99SBarry Smith         cols[nc] = idx_c[col+dof]/dof;
11647c6ae99SBarry Smith         v[nc++]  = x;
11747c6ae99SBarry Smith       }
11847c6ae99SBarry Smith       ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
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;
128b3bd94feSDave May     ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr);
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 */
135b3bd94feSDave May       row = idx_f[dof*(i-i_start_ghost)]/dof;
136b3bd94feSDave May 
137b3bd94feSDave May       i_c = (i/ratio);    /* coarse grid node to left of fine grid node */
138b3bd94feSDave May       if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
139b3bd94feSDave May                                           i_start %D i_c %D i_start_ghost_c %D",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 */
146b3bd94feSDave May       col     = dof*(i_c-i_start_ghost_c);
147b3bd94feSDave May       cols[0] = idx_c[col]/dof;
148b3bd94feSDave May       Ni[0]   = 1.0;
149b3bd94feSDave May       if ((li==0) || (li==nxi-1)) {
150b3bd94feSDave May         ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr);
151b3bd94feSDave May         continue;
152b3bd94feSDave May       }
153b3bd94feSDave May 
154b3bd94feSDave May       /* edges + interior */
155b3bd94feSDave May       /* remainders */
1568865f1eaSKarl Rupp       if (i==mx-1) i_c--;
157b3bd94feSDave May 
158b3bd94feSDave May       col     = dof*(i_c-i_start_ghost_c);
159b3bd94feSDave May       cols[0] = idx_c[col]/dof; /* one left and below; or we are right on it */
160b3bd94feSDave May       cols[1] = idx_c[col+dof]/dof;
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       }
167b3bd94feSDave May       ierr = MatSetValues(mat,1,&row,2,cols,Ni,INSERT_VALUES);CHKERRQ(ierr);
168b3bd94feSDave May     }
169b3bd94feSDave May     ierr = PetscFree(xi);CHKERRQ(ierr);
170b3bd94feSDave May   }
171*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
172*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
17347c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17447c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17547c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
176fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
17747c6ae99SBarry Smith   PetscFunctionReturn(0);
17847c6ae99SBarry Smith }
17947c6ae99SBarry Smith 
18047c6ae99SBarry Smith #undef __FUNCT__
181e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA_1D_Q0"
182e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_1D_Q0(DM dac,DM daf,Mat *A)
18347c6ae99SBarry Smith {
18447c6ae99SBarry Smith   PetscErrorCode   ierr;
185*8ea3bf28SBarry Smith   PetscInt         i,i_start,m_f,Mx;
186*8ea3bf28SBarry Smith   const PetscInt   *idx_f,*idx_c;
187*8ea3bf28SBarry Smith   PetscInt         m_ghost,m_ghost_c;
18847c6ae99SBarry Smith   PetscInt         row,col,i_start_ghost,mx,m_c,nc,ratio;
18947c6ae99SBarry Smith   PetscInt         i_c,i_start_c,i_start_ghost_c,cols[2],dof;
19047c6ae99SBarry Smith   PetscScalar      v[2],x;
19147c6ae99SBarry Smith   Mat              mat;
1921321219cSEthan Coon   DMDABoundaryType bx;
19347c6ae99SBarry Smith 
19447c6ae99SBarry Smith   PetscFunctionBegin;
1951321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr);
1961321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
1971321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) {
19847c6ae99SBarry Smith     ratio = mx/Mx;
19947c6ae99SBarry Smith     if (ratio*Mx != mx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: mx/Mx  must be integer: mx %D Mx %D",mx,Mx);
20047c6ae99SBarry Smith   } else {
20147c6ae99SBarry Smith     ratio = (mx-1)/(Mx-1);
20247c6ae99SBarry Smith     if (ratio*(Mx-1) != mx-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %D Mx %D",mx,Mx);
20347c6ae99SBarry Smith   }
20447c6ae99SBarry Smith 
205aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr);
206aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr);
2070298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
20847c6ae99SBarry Smith 
209aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr);
210aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr);
2110298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
21247c6ae99SBarry Smith 
21347c6ae99SBarry Smith   /* create interpolation matrix */
214ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)dac),&mat);CHKERRQ(ierr);
21547c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr);
21647c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
2170298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,2,NULL);CHKERRQ(ierr);
2180298fd71SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,2,NULL,0,NULL);CHKERRQ(ierr);
21947c6ae99SBarry Smith 
22047c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
22147c6ae99SBarry Smith   for (i=i_start; i<i_start+m_f; i++) {
22247c6ae99SBarry Smith     /* convert to local "natural" numbering and then to PETSc global numbering */
22347c6ae99SBarry Smith     row = idx_f[dof*(i-i_start_ghost)]/dof;
22447c6ae99SBarry Smith 
22547c6ae99SBarry Smith     i_c = (i/ratio);    /* coarse grid node to left of fine grid node */
22647c6ae99SBarry Smith 
22747c6ae99SBarry Smith     /*
22847c6ae99SBarry Smith          Only include those interpolation points that are truly
22947c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
23047c6ae99SBarry Smith          in x direction; since they have no right neighbor
23147c6ae99SBarry Smith     */
23247c6ae99SBarry Smith     x  = ((double)(i - i_c*ratio))/((double)ratio);
23347c6ae99SBarry Smith     nc = 0;
23447c6ae99SBarry Smith     /* one left and below; or we are right on it */
23547c6ae99SBarry Smith     col      = dof*(i_c-i_start_ghost_c);
23647c6ae99SBarry Smith     cols[nc] = idx_c[col]/dof;
23747c6ae99SBarry Smith     v[nc++]  = -x + 1.0;
23847c6ae99SBarry Smith     /* one right? */
23947c6ae99SBarry Smith     if (i_c*ratio != i) {
24047c6ae99SBarry Smith       cols[nc] = idx_c[col+dof]/dof;
24147c6ae99SBarry Smith       v[nc++]  = x;
24247c6ae99SBarry Smith     }
24347c6ae99SBarry Smith     ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
24447c6ae99SBarry Smith   }
245*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
246*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
24747c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
24847c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
24947c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
250fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
25147c6ae99SBarry Smith   ierr = PetscLogFlops(5.0*m_f);CHKERRQ(ierr);
25247c6ae99SBarry Smith   PetscFunctionReturn(0);
25347c6ae99SBarry Smith }
25447c6ae99SBarry Smith 
25547c6ae99SBarry Smith #undef __FUNCT__
256e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA_2D_Q1"
257e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_2D_Q1(DM dac,DM daf,Mat *A)
25847c6ae99SBarry Smith {
25947c6ae99SBarry Smith   PetscErrorCode   ierr;
260*8ea3bf28SBarry Smith   PetscInt         i,j,i_start,j_start,m_f,n_f,Mx,My,dof;
261*8ea3bf28SBarry Smith   const PetscInt   *idx_c,*idx_f;
262*8ea3bf28SBarry Smith   PetscInt         m_ghost,n_ghost,m_ghost_c,n_ghost_c,*dnz,*onz;
26347c6ae99SBarry Smith   PetscInt         row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj;
26447c6ae99SBarry 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;
26547c6ae99SBarry Smith   PetscMPIInt      size_c,size_f,rank_f;
26647c6ae99SBarry Smith   PetscScalar      v[4],x,y;
26747c6ae99SBarry Smith   Mat              mat;
2681321219cSEthan Coon   DMDABoundaryType bx,by;
26947c6ae99SBarry Smith 
27047c6ae99SBarry Smith   PetscFunctionBegin;
2711321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr);
2721321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
2731321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) {
27447c6ae99SBarry Smith     ratioi = mx/Mx;
27547c6ae99SBarry Smith     if (ratioi*Mx != mx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: mx/Mx  must be integer: mx %D Mx %D",mx,Mx);
27647c6ae99SBarry Smith   } else {
27747c6ae99SBarry Smith     ratioi = (mx-1)/(Mx-1);
27847c6ae99SBarry Smith     if (ratioi*(Mx-1) != mx-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %D Mx %D",mx,Mx);
27947c6ae99SBarry Smith   }
2801321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC) {
28147c6ae99SBarry Smith     ratioj = my/My;
28247c6ae99SBarry Smith     if (ratioj*My != my) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: my/My  must be integer: my %D My %D",my,My);
28347c6ae99SBarry Smith   } else {
28447c6ae99SBarry Smith     ratioj = (my-1)/(My-1);
28547c6ae99SBarry Smith     if (ratioj*(My-1) != my-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (my - 1)/(My - 1) must be integer: my %D My %D",my,My);
28647c6ae99SBarry Smith   }
28747c6ae99SBarry Smith 
28847c6ae99SBarry Smith 
289aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr);
290aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr);
2910298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
29247c6ae99SBarry Smith 
293aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr);
294aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr);
2950298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
29647c6ae99SBarry Smith 
29747c6ae99SBarry Smith   /*
298aa219208SBarry Smith    Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
29947c6ae99SBarry Smith    The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
30047c6ae99SBarry Smith    processors). It's effective length is hence 4 times its normal length, this is
30147c6ae99SBarry Smith    why the col_scale is multiplied by the interpolation matrix column sizes.
30247c6ae99SBarry Smith    sol_shift allows each set of 1/4 processors do its own interpolation using ITS
30347c6ae99SBarry Smith    copy of the coarse vector. A bit of a hack but you do better.
30447c6ae99SBarry Smith 
30547c6ae99SBarry Smith    In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
30647c6ae99SBarry Smith    */
307ce94432eSBarry Smith   ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c);CHKERRQ(ierr);
308ce94432eSBarry Smith   ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f);CHKERRQ(ierr);
309ce94432eSBarry Smith   ierr      = MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f);CHKERRQ(ierr);
31047c6ae99SBarry Smith   col_scale = size_f/size_c;
31147c6ae99SBarry Smith   col_shift = Mx*My*(rank_f/size_c);
31247c6ae99SBarry Smith 
313ce94432eSBarry Smith   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)daf),m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr);
31447c6ae99SBarry Smith   for (j=j_start; j<j_start+n_f; j++) {
31547c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
31647c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
31747c6ae99SBarry Smith       row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
31847c6ae99SBarry Smith 
31947c6ae99SBarry Smith       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
32047c6ae99SBarry Smith       j_c = (j/ratioj);    /* coarse grid node below fine grid node */
32147c6ae99SBarry Smith 
322aa219208SBarry Smith       if (j_c < j_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
32347c6ae99SBarry Smith                                           j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c);
324aa219208SBarry Smith       if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
32547c6ae99SBarry Smith                                           i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
32647c6ae99SBarry Smith 
32747c6ae99SBarry Smith       /*
32847c6ae99SBarry Smith        Only include those interpolation points that are truly
32947c6ae99SBarry Smith        nonzero. Note this is very important for final grid lines
33047c6ae99SBarry Smith        in x and y directions; since they have no right/top neighbors
33147c6ae99SBarry Smith        */
33247c6ae99SBarry Smith       nc = 0;
33347c6ae99SBarry Smith       /* one left and below; or we are right on it */
33447c6ae99SBarry Smith       col        = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
33547c6ae99SBarry Smith       cols[nc++] = col_shift + idx_c[col]/dof;
33647c6ae99SBarry Smith       /* one right and below */
3378865f1eaSKarl Rupp       if (i_c*ratioi != i) cols[nc++] = col_shift + idx_c[col+dof]/dof;
33847c6ae99SBarry Smith       /* one left and above */
3398865f1eaSKarl Rupp       if (j_c*ratioj != j) cols[nc++] = col_shift + idx_c[col+m_ghost_c*dof]/dof;
34047c6ae99SBarry Smith       /* one right and above */
3418865f1eaSKarl Rupp       if (i_c*ratioi != i && j_c*ratioj != j) cols[nc++] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof;
34247c6ae99SBarry Smith       ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr);
34347c6ae99SBarry Smith     }
34447c6ae99SBarry Smith   }
345ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)daf),&mat);CHKERRQ(ierr);
34647c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr);
34747c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
34847c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr);
34947c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr);
35047c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
35147c6ae99SBarry Smith 
35247c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
35385afcc9aSBarry Smith   if (!NEWVERSION) {
354b3bd94feSDave May 
35547c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
35647c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
35747c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
35847c6ae99SBarry Smith         row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
35947c6ae99SBarry Smith 
36047c6ae99SBarry Smith         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
36147c6ae99SBarry Smith         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
36247c6ae99SBarry Smith 
36347c6ae99SBarry Smith         /*
36447c6ae99SBarry Smith          Only include those interpolation points that are truly
36547c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
36647c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
36747c6ae99SBarry Smith          */
36847c6ae99SBarry Smith         x = ((double)(i - i_c*ratioi))/((double)ratioi);
36947c6ae99SBarry Smith         y = ((double)(j - j_c*ratioj))/((double)ratioj);
370b3bd94feSDave May 
37147c6ae99SBarry Smith         nc = 0;
37247c6ae99SBarry Smith         /* one left and below; or we are right on it */
37347c6ae99SBarry Smith         col      = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
37447c6ae99SBarry Smith         cols[nc] = col_shift + idx_c[col]/dof;
37547c6ae99SBarry Smith         v[nc++]  = x*y - x - y + 1.0;
37647c6ae99SBarry Smith         /* one right and below */
37747c6ae99SBarry Smith         if (i_c*ratioi != i) {
37847c6ae99SBarry Smith           cols[nc] = col_shift + idx_c[col+dof]/dof;
37947c6ae99SBarry Smith           v[nc++]  = -x*y + x;
38047c6ae99SBarry Smith         }
38147c6ae99SBarry Smith         /* one left and above */
38247c6ae99SBarry Smith         if (j_c*ratioj != j) {
38347c6ae99SBarry Smith           cols[nc] = col_shift + idx_c[col+m_ghost_c*dof]/dof;
38447c6ae99SBarry Smith           v[nc++]  = -x*y + y;
38547c6ae99SBarry Smith         }
38647c6ae99SBarry Smith         /* one right and above */
38747c6ae99SBarry Smith         if (j_c*ratioj != j && i_c*ratioi != i) {
38847c6ae99SBarry Smith           cols[nc] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof;
38947c6ae99SBarry Smith           v[nc++]  = x*y;
39047c6ae99SBarry Smith         }
39147c6ae99SBarry Smith         ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
39247c6ae99SBarry Smith       }
39347c6ae99SBarry Smith     }
394b3bd94feSDave May 
395b3bd94feSDave May   } else {
396b3bd94feSDave May     PetscScalar Ni[4];
397b3bd94feSDave May     PetscScalar *xi,*eta;
398b3bd94feSDave May     PetscInt    li,nxi,lj,neta;
399b3bd94feSDave May 
400b3bd94feSDave May     /* compute local coordinate arrays */
401b3bd94feSDave May     nxi  = ratioi + 1;
402b3bd94feSDave May     neta = ratioj + 1;
403b3bd94feSDave May     ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr);
404b3bd94feSDave May     ierr = PetscMalloc(sizeof(PetscScalar)*neta,&eta);CHKERRQ(ierr);
405b3bd94feSDave May     for (li=0; li<nxi; li++) {
40652218ef2SJed Brown       xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1));
407b3bd94feSDave May     }
408b3bd94feSDave May     for (lj=0; lj<neta; lj++) {
40952218ef2SJed Brown       eta[lj] = -1.0 + (PetscScalar)lj*(2.0/(PetscScalar)(neta-1));
410b3bd94feSDave May     }
411b3bd94feSDave May 
412b3bd94feSDave May     /* loop over local fine grid nodes setting interpolation for those*/
413b3bd94feSDave May     for (j=j_start; j<j_start+n_f; j++) {
414b3bd94feSDave May       for (i=i_start; i<i_start+m_f; i++) {
415b3bd94feSDave May         /* convert to local "natural" numbering and then to PETSc global numbering */
416b3bd94feSDave May         row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
417b3bd94feSDave May 
418b3bd94feSDave May         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
419b3bd94feSDave May         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
420b3bd94feSDave May 
421b3bd94feSDave May         /* remainders */
422b3bd94feSDave May         li = i - ratioi * (i/ratioi);
4238865f1eaSKarl Rupp         if (i==mx-1) li = nxi-1;
424b3bd94feSDave May         lj = j - ratioj * (j/ratioj);
4258865f1eaSKarl Rupp         if (j==my-1) lj = neta-1;
426b3bd94feSDave May 
427b3bd94feSDave May         /* corners */
428b3bd94feSDave May         col     = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
429b3bd94feSDave May         cols[0] = col_shift + idx_c[col]/dof; /* left, below */
430b3bd94feSDave May         Ni[0]   = 1.0;
431b3bd94feSDave May         if ((li==0) || (li==nxi-1)) {
432b3bd94feSDave May           if ((lj==0) || (lj==neta-1)) {
433b3bd94feSDave May             ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr);
434b3bd94feSDave May             continue;
435b3bd94feSDave May           }
436b3bd94feSDave May         }
437b3bd94feSDave May 
438b3bd94feSDave May         /* edges + interior */
439b3bd94feSDave May         /* remainders */
4408865f1eaSKarl Rupp         if (i==mx-1) i_c--;
4418865f1eaSKarl Rupp         if (j==my-1) j_c--;
442b3bd94feSDave May 
443b3bd94feSDave May         col     = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
444b3bd94feSDave May         cols[0] = col_shift + idx_c[col]/dof; /* left, below */
445b3bd94feSDave May         cols[1] = col_shift + idx_c[col+dof]/dof; /* right, below */
446b3bd94feSDave May         cols[2] = col_shift + idx_c[col+m_ghost_c*dof]/dof; /* left, above */
447b3bd94feSDave May         cols[3] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; /* right, above */
448b3bd94feSDave May 
449b3bd94feSDave May         Ni[0] = 0.25*(1.0-xi[li])*(1.0-eta[lj]);
450b3bd94feSDave May         Ni[1] = 0.25*(1.0+xi[li])*(1.0-eta[lj]);
451b3bd94feSDave May         Ni[2] = 0.25*(1.0-xi[li])*(1.0+eta[lj]);
452b3bd94feSDave May         Ni[3] = 0.25*(1.0+xi[li])*(1.0+eta[lj]);
453b3bd94feSDave May 
454b3bd94feSDave May         nc = 0;
4558865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[0])<1.0e-32) cols[0]=-1;
4568865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[1])<1.0e-32) cols[1]=-1;
4578865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[2])<1.0e-32) cols[2]=-1;
4588865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[3])<1.0e-32) cols[3]=-1;
459b3bd94feSDave May 
460b3bd94feSDave May         ierr = MatSetValues(mat,1,&row,4,cols,Ni,INSERT_VALUES);CHKERRQ(ierr);
461b3bd94feSDave May       }
462b3bd94feSDave May     }
463b3bd94feSDave May     ierr = PetscFree(xi);CHKERRQ(ierr);
464b3bd94feSDave May     ierr = PetscFree(eta);CHKERRQ(ierr);
465b3bd94feSDave May   }
466*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
467*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
46847c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
46947c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
47047c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
471fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
47247c6ae99SBarry Smith   PetscFunctionReturn(0);
47347c6ae99SBarry Smith }
47447c6ae99SBarry Smith 
47547c6ae99SBarry Smith /*
47647c6ae99SBarry Smith        Contributed by Andrei Draganescu <aidraga@sandia.gov>
47747c6ae99SBarry Smith */
47847c6ae99SBarry Smith #undef __FUNCT__
479e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA_2D_Q0"
480e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_2D_Q0(DM dac,DM daf,Mat *A)
48147c6ae99SBarry Smith {
48247c6ae99SBarry Smith   PetscErrorCode   ierr;
483*8ea3bf28SBarry Smith   PetscInt         i,j,i_start,j_start,m_f,n_f,Mx,My,dof;
484*8ea3bf28SBarry Smith   const PetscInt   *idx_c,*idx_f;
485*8ea3bf28SBarry Smith   PetscInt         m_ghost,n_ghost,m_ghost_c,n_ghost_c,*dnz,*onz;
48647c6ae99SBarry Smith   PetscInt         row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj;
48747c6ae99SBarry 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;
48847c6ae99SBarry Smith   PetscMPIInt      size_c,size_f,rank_f;
48947c6ae99SBarry Smith   PetscScalar      v[4];
49047c6ae99SBarry Smith   Mat              mat;
4911321219cSEthan Coon   DMDABoundaryType bx,by;
49247c6ae99SBarry Smith 
49347c6ae99SBarry Smith   PetscFunctionBegin;
4941321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr);
4951321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
496ce94432eSBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x");
497ce94432eSBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y");
49847c6ae99SBarry Smith   ratioi = mx/Mx;
49947c6ae99SBarry Smith   ratioj = my/My;
500ce94432eSBarry Smith   if (ratioi*Mx != mx) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in x");
501ce94432eSBarry Smith   if (ratioj*My != my) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in y");
502ce94432eSBarry Smith   if (ratioi != 2) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 2");
503ce94432eSBarry Smith   if (ratioj != 2) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 2");
50447c6ae99SBarry Smith 
505aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr);
506aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr);
5070298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
50847c6ae99SBarry Smith 
509aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr);
510aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr);
5110298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith   /*
514aa219208SBarry Smith      Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
51547c6ae99SBarry Smith      The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
51647c6ae99SBarry Smith      processors). It's effective length is hence 4 times its normal length, this is
51747c6ae99SBarry Smith      why the col_scale is multiplied by the interpolation matrix column sizes.
51847c6ae99SBarry Smith      sol_shift allows each set of 1/4 processors do its own interpolation using ITS
51947c6ae99SBarry Smith      copy of the coarse vector. A bit of a hack but you do better.
52047c6ae99SBarry Smith 
52147c6ae99SBarry Smith      In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
52247c6ae99SBarry Smith   */
523ce94432eSBarry Smith   ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c);CHKERRQ(ierr);
524ce94432eSBarry Smith   ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f);CHKERRQ(ierr);
525ce94432eSBarry Smith   ierr      = MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f);CHKERRQ(ierr);
52647c6ae99SBarry Smith   col_scale = size_f/size_c;
52747c6ae99SBarry Smith   col_shift = Mx*My*(rank_f/size_c);
52847c6ae99SBarry Smith 
529ce94432eSBarry Smith   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)daf),m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr);
53047c6ae99SBarry Smith   for (j=j_start; j<j_start+n_f; j++) {
53147c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
53247c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
53347c6ae99SBarry Smith       row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
53447c6ae99SBarry Smith 
53547c6ae99SBarry Smith       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
53647c6ae99SBarry Smith       j_c = (j/ratioj);    /* coarse grid node below fine grid node */
53747c6ae99SBarry Smith 
538aa219208SBarry Smith       if (j_c < j_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
53947c6ae99SBarry Smith     j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c);
540aa219208SBarry Smith       if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
54147c6ae99SBarry Smith     i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
54247c6ae99SBarry Smith 
54347c6ae99SBarry Smith       /*
54447c6ae99SBarry Smith          Only include those interpolation points that are truly
54547c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
54647c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
54747c6ae99SBarry Smith       */
54847c6ae99SBarry Smith       nc = 0;
54947c6ae99SBarry Smith       /* one left and below; or we are right on it */
55047c6ae99SBarry Smith       col        = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
55147c6ae99SBarry Smith       cols[nc++] = col_shift + idx_c[col]/dof;
55247c6ae99SBarry Smith       ierr       = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr);
55347c6ae99SBarry Smith     }
55447c6ae99SBarry Smith   }
555ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)daf),&mat);CHKERRQ(ierr);
55647c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr);
55747c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
55847c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr);
55947c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr);
56047c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
56147c6ae99SBarry Smith 
56247c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
56347c6ae99SBarry Smith   for (j=j_start; j<j_start+n_f; j++) {
56447c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
56547c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
56647c6ae99SBarry Smith       row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
56747c6ae99SBarry Smith 
56847c6ae99SBarry Smith       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
56947c6ae99SBarry Smith       j_c = (j/ratioj);    /* coarse grid node below fine grid node */
57047c6ae99SBarry Smith       nc  = 0;
57147c6ae99SBarry Smith       /* one left and below; or we are right on it */
57247c6ae99SBarry Smith       col      = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
57347c6ae99SBarry Smith       cols[nc] = col_shift + idx_c[col]/dof;
57447c6ae99SBarry Smith       v[nc++]  = 1.0;
57547c6ae99SBarry Smith 
57647c6ae99SBarry Smith       ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
57747c6ae99SBarry Smith     }
57847c6ae99SBarry Smith   }
579*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
580*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
58147c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
58247c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
58347c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
584fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
58547c6ae99SBarry Smith   ierr = PetscLogFlops(13.0*m_f*n_f);CHKERRQ(ierr);
58647c6ae99SBarry Smith   PetscFunctionReturn(0);
58747c6ae99SBarry Smith }
58847c6ae99SBarry Smith 
58947c6ae99SBarry Smith /*
59047c6ae99SBarry Smith        Contributed by Jianming Yang <jianming-yang@uiowa.edu>
59147c6ae99SBarry Smith */
59247c6ae99SBarry Smith #undef __FUNCT__
593e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA_3D_Q0"
594e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_3D_Q0(DM dac,DM daf,Mat *A)
59547c6ae99SBarry Smith {
59647c6ae99SBarry Smith   PetscErrorCode   ierr;
597*8ea3bf28SBarry Smith   PetscInt         i,j,l,i_start,j_start,l_start,m_f,n_f,p_f,Mx,My,Mz,dof;
598*8ea3bf28SBarry Smith   const PetscInt   *idx_c,*idx_f;
599*8ea3bf28SBarry Smith   PetscInt         m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c,nc,*dnz,*onz;
60047c6ae99SBarry 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;
60147c6ae99SBarry 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;
60247c6ae99SBarry Smith   PetscMPIInt      size_c,size_f,rank_f;
60347c6ae99SBarry Smith   PetscScalar      v[8];
60447c6ae99SBarry Smith   Mat              mat;
6051321219cSEthan Coon   DMDABoundaryType bx,by,bz;
60647c6ae99SBarry Smith 
60747c6ae99SBarry Smith   PetscFunctionBegin;
6081321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr);
6091321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
610ce94432eSBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x");
611ce94432eSBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y");
612ce94432eSBarry Smith   if (bz == DMDA_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in z");
61347c6ae99SBarry Smith   ratioi = mx/Mx;
61447c6ae99SBarry Smith   ratioj = my/My;
61547c6ae99SBarry Smith   ratiol = mz/Mz;
616ce94432eSBarry Smith   if (ratioi*Mx != mx) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in x");
617ce94432eSBarry Smith   if (ratioj*My != my) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in y");
618ce94432eSBarry Smith   if (ratiol*Mz != mz) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in z");
619ce94432eSBarry Smith   if (ratioi != 2 && ratioi != 1) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 1 or 2");
620ce94432eSBarry Smith   if (ratioj != 2 && ratioj != 1) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 1 or 2");
621ce94432eSBarry Smith   if (ratiol != 2 && ratiol != 1) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in z must be 1 or 2");
62247c6ae99SBarry Smith 
623aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr);
624aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr);
6250298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
62647c6ae99SBarry Smith 
627aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr);
628aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,&l_start_ghost_c,&m_ghost_c,&n_ghost_c,&p_ghost_c);CHKERRQ(ierr);
6290298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
63047c6ae99SBarry Smith   /*
631aa219208SBarry Smith      Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
63247c6ae99SBarry Smith      The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
63347c6ae99SBarry Smith      processors). It's effective length is hence 4 times its normal length, this is
63447c6ae99SBarry Smith      why the col_scale is multiplied by the interpolation matrix column sizes.
63547c6ae99SBarry Smith      sol_shift allows each set of 1/4 processors do its own interpolation using ITS
63647c6ae99SBarry Smith      copy of the coarse vector. A bit of a hack but you do better.
63747c6ae99SBarry Smith 
63847c6ae99SBarry Smith      In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
63947c6ae99SBarry Smith   */
640ce94432eSBarry Smith   ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c);CHKERRQ(ierr);
641ce94432eSBarry Smith   ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f);CHKERRQ(ierr);
642ce94432eSBarry Smith   ierr      = MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f);CHKERRQ(ierr);
64347c6ae99SBarry Smith   col_scale = size_f/size_c;
64447c6ae99SBarry Smith   col_shift = Mx*My*Mz*(rank_f/size_c);
64547c6ae99SBarry Smith 
646ce94432eSBarry Smith   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)daf),m_f*n_f*p_f,col_scale*m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr);
64747c6ae99SBarry Smith   for (l=l_start; l<l_start+p_f; l++) {
64847c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
64947c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
65047c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
65147c6ae99SBarry Smith         row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
65247c6ae99SBarry Smith 
65347c6ae99SBarry Smith         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
65447c6ae99SBarry Smith         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
65547c6ae99SBarry Smith         l_c = (l/ratiol);
65647c6ae99SBarry Smith 
657aa219208SBarry Smith         if (l_c < l_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
65847c6ae99SBarry Smith     l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c);
659aa219208SBarry Smith         if (j_c < j_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
66047c6ae99SBarry Smith     j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c);
661aa219208SBarry Smith         if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
66247c6ae99SBarry Smith     i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
66347c6ae99SBarry Smith 
66447c6ae99SBarry Smith         /*
66547c6ae99SBarry Smith            Only include those interpolation points that are truly
66647c6ae99SBarry Smith            nonzero. Note this is very important for final grid lines
66747c6ae99SBarry Smith            in x and y directions; since they have no right/top neighbors
66847c6ae99SBarry Smith         */
66947c6ae99SBarry Smith         nc = 0;
67047c6ae99SBarry Smith         /* one left and below; or we are right on it */
67147c6ae99SBarry Smith         col        = dof*(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));
67247c6ae99SBarry Smith         cols[nc++] = col_shift + idx_c[col]/dof;
67347c6ae99SBarry Smith         ierr       = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr);
67447c6ae99SBarry Smith       }
67547c6ae99SBarry Smith     }
67647c6ae99SBarry Smith   }
677ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)daf),&mat);CHKERRQ(ierr);
67847c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f*n_f*p_f,col_scale*m_c*n_c*p_c,mx*my*mz,col_scale*Mx*My*Mz);CHKERRQ(ierr);
67947c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
68047c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr);
68147c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr);
68247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
68347c6ae99SBarry Smith 
68447c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
68547c6ae99SBarry Smith   for (l=l_start; l<l_start+p_f; l++) {
68647c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
68747c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
68847c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
68947c6ae99SBarry Smith         row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
69047c6ae99SBarry Smith 
69147c6ae99SBarry Smith         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
69247c6ae99SBarry Smith         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
69347c6ae99SBarry Smith         l_c = (l/ratiol);
69447c6ae99SBarry Smith         nc  = 0;
69547c6ae99SBarry Smith         /* one left and below; or we are right on it */
69647c6ae99SBarry Smith         col      = dof*(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));
69747c6ae99SBarry Smith         cols[nc] = col_shift + idx_c[col]/dof;
69847c6ae99SBarry Smith         v[nc++]  = 1.0;
69947c6ae99SBarry Smith 
70047c6ae99SBarry Smith         ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
70147c6ae99SBarry Smith       }
70247c6ae99SBarry Smith     }
70347c6ae99SBarry Smith   }
704*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
705*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
70647c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70747c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70847c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
709fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
71047c6ae99SBarry Smith   ierr = PetscLogFlops(13.0*m_f*n_f*p_f);CHKERRQ(ierr);
71147c6ae99SBarry Smith   PetscFunctionReturn(0);
71247c6ae99SBarry Smith }
71347c6ae99SBarry Smith 
71447c6ae99SBarry Smith #undef __FUNCT__
715e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA_3D_Q1"
716e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_3D_Q1(DM dac,DM daf,Mat *A)
71747c6ae99SBarry Smith {
71847c6ae99SBarry Smith   PetscErrorCode   ierr;
719*8ea3bf28SBarry Smith   PetscInt         i,j,i_start,j_start,m_f,n_f,Mx,My,dof,l;
720*8ea3bf28SBarry Smith   const PetscInt   *idx_c,*idx_f;
721*8ea3bf28SBarry Smith   PetscInt         m_ghost,n_ghost,m_ghost_c,n_ghost_c,Mz,mz;
72247c6ae99SBarry Smith   PetscInt         row,col,i_start_ghost,j_start_ghost,cols[8],mx,m_c,my,nc,ratioi,ratioj,ratiok;
72347c6ae99SBarry Smith   PetscInt         i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c;
72447c6ae99SBarry Smith   PetscInt         l_start,p_f,l_start_ghost,p_ghost,l_start_c,p_c;
72547c6ae99SBarry Smith   PetscInt         l_start_ghost_c,p_ghost_c,l_c,*dnz,*onz;
72647c6ae99SBarry Smith   PetscScalar      v[8],x,y,z;
72747c6ae99SBarry Smith   Mat              mat;
7281321219cSEthan Coon   DMDABoundaryType bx,by,bz;
72947c6ae99SBarry Smith 
73047c6ae99SBarry Smith   PetscFunctionBegin;
7311321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr);
7321321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
73347c6ae99SBarry Smith   if (mx == Mx) {
73447c6ae99SBarry Smith     ratioi = 1;
7351321219cSEthan Coon   } else if (bx == DMDA_BOUNDARY_PERIODIC) {
73647c6ae99SBarry Smith     ratioi = mx/Mx;
73747c6ae99SBarry Smith     if (ratioi*Mx != mx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: mx/Mx  must be integer: mx %D Mx %D",mx,Mx);
73847c6ae99SBarry Smith   } else {
73947c6ae99SBarry Smith     ratioi = (mx-1)/(Mx-1);
74047c6ae99SBarry Smith     if (ratioi*(Mx-1) != mx-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %D Mx %D",mx,Mx);
74147c6ae99SBarry Smith   }
74247c6ae99SBarry Smith   if (my == My) {
74347c6ae99SBarry Smith     ratioj = 1;
7441321219cSEthan Coon   } else if (by == DMDA_BOUNDARY_PERIODIC) {
74547c6ae99SBarry Smith     ratioj = my/My;
74647c6ae99SBarry Smith     if (ratioj*My != my) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: my/My  must be integer: my %D My %D",my,My);
74747c6ae99SBarry Smith   } else {
74847c6ae99SBarry Smith     ratioj = (my-1)/(My-1);
74947c6ae99SBarry Smith     if (ratioj*(My-1) != my-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (my - 1)/(My - 1) must be integer: my %D My %D",my,My);
75047c6ae99SBarry Smith   }
75147c6ae99SBarry Smith   if (mz == Mz) {
75247c6ae99SBarry Smith     ratiok = 1;
7531321219cSEthan Coon   } else if (bz == DMDA_BOUNDARY_PERIODIC) {
75447c6ae99SBarry Smith     ratiok = mz/Mz;
75547c6ae99SBarry Smith     if (ratiok*Mz != mz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: mz/Mz  must be integer: mz %D Mz %D",mz,Mz);
75647c6ae99SBarry Smith   } else {
75747c6ae99SBarry Smith     ratiok = (mz-1)/(Mz-1);
75847c6ae99SBarry Smith     if (ratiok*(Mz-1) != mz-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (mz - 1)/(Mz - 1) must be integer: mz %D Mz %D",mz,Mz);
75947c6ae99SBarry Smith   }
76047c6ae99SBarry Smith 
761aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr);
762aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr);
7630298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
76447c6ae99SBarry Smith 
765aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr);
766aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,&l_start_ghost_c,&m_ghost_c,&n_ghost_c,&p_ghost_c);CHKERRQ(ierr);
7670298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
76847c6ae99SBarry Smith 
76947c6ae99SBarry Smith   /* create interpolation matrix, determining exact preallocation */
770ce94432eSBarry Smith   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)dac),m_f*n_f*p_f,m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr);
77147c6ae99SBarry Smith   /* loop over local fine grid nodes counting interpolating points */
77247c6ae99SBarry Smith   for (l=l_start; l<l_start+p_f; l++) {
77347c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
77447c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
77547c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
77647c6ae99SBarry Smith         row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
77747c6ae99SBarry Smith         i_c = (i/ratioi);
77847c6ae99SBarry Smith         j_c = (j/ratioj);
77947c6ae99SBarry Smith         l_c = (l/ratiok);
780aa219208SBarry Smith         if (l_c < l_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
78147c6ae99SBarry Smith                                             l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c);
782aa219208SBarry Smith         if (j_c < j_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
78347c6ae99SBarry Smith                                             j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c);
784aa219208SBarry Smith         if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
78547c6ae99SBarry Smith                                             i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
78647c6ae99SBarry Smith 
78747c6ae99SBarry Smith         /*
78847c6ae99SBarry Smith          Only include those interpolation points that are truly
78947c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
79047c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
79147c6ae99SBarry Smith          */
79247c6ae99SBarry Smith         nc         = 0;
79347c6ae99SBarry Smith         col        = dof*(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));
79447c6ae99SBarry Smith         cols[nc++] = idx_c[col]/dof;
79547c6ae99SBarry Smith         if (i_c*ratioi != i) {
79647c6ae99SBarry Smith           cols[nc++] = idx_c[col+dof]/dof;
79747c6ae99SBarry Smith         }
79847c6ae99SBarry Smith         if (j_c*ratioj != j) {
79947c6ae99SBarry Smith           cols[nc++] = idx_c[col+m_ghost_c*dof]/dof;
80047c6ae99SBarry Smith         }
80147c6ae99SBarry Smith         if (l_c*ratiok != l) {
80247c6ae99SBarry Smith           cols[nc++] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof;
80347c6ae99SBarry Smith         }
80447c6ae99SBarry Smith         if (j_c*ratioj != j && i_c*ratioi != i) {
80547c6ae99SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c+1)*dof]/dof;
80647c6ae99SBarry Smith         }
80747c6ae99SBarry Smith         if (j_c*ratioj != j && l_c*ratiok != l) {
80847c6ae99SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof;
80947c6ae99SBarry Smith         }
81047c6ae99SBarry Smith         if (i_c*ratioi != i && l_c*ratiok != l) {
81147c6ae99SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof;
81247c6ae99SBarry Smith         }
81347c6ae99SBarry Smith         if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) {
81447c6ae99SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof;
81547c6ae99SBarry Smith         }
81647c6ae99SBarry Smith         ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr);
81747c6ae99SBarry Smith       }
81847c6ae99SBarry Smith     }
81947c6ae99SBarry Smith   }
820ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)dac),&mat);CHKERRQ(ierr);
82147c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f*n_f*p_f,m_c*n_c*p_c,mx*my*mz,Mx*My*Mz);CHKERRQ(ierr);
82247c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
82347c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr);
82447c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr);
82547c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
82647c6ae99SBarry Smith 
82747c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
8282adb9dcfSBarry Smith   if (!NEWVERSION) {
829b3bd94feSDave May 
83047c6ae99SBarry Smith     for (l=l_start; l<l_start+p_f; l++) {
83147c6ae99SBarry Smith       for (j=j_start; j<j_start+n_f; j++) {
83247c6ae99SBarry Smith         for (i=i_start; i<i_start+m_f; i++) {
83347c6ae99SBarry Smith           /* convert to local "natural" numbering and then to PETSc global numbering */
83447c6ae99SBarry Smith           row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
83547c6ae99SBarry Smith 
83647c6ae99SBarry Smith           i_c = (i/ratioi);
83747c6ae99SBarry Smith           j_c = (j/ratioj);
83847c6ae99SBarry Smith           l_c = (l/ratiok);
83947c6ae99SBarry Smith 
84047c6ae99SBarry Smith           /*
84147c6ae99SBarry Smith            Only include those interpolation points that are truly
84247c6ae99SBarry Smith            nonzero. Note this is very important for final grid lines
84347c6ae99SBarry Smith            in x and y directions; since they have no right/top neighbors
84447c6ae99SBarry Smith            */
84547c6ae99SBarry Smith           x = ((double)(i - i_c*ratioi))/((double)ratioi);
84647c6ae99SBarry Smith           y = ((double)(j - j_c*ratioj))/((double)ratioj);
84747c6ae99SBarry Smith           z = ((double)(l - l_c*ratiok))/((double)ratiok);
848b3bd94feSDave May 
84947c6ae99SBarry Smith           nc = 0;
85047c6ae99SBarry Smith           /* one left and below; or we are right on it */
85147c6ae99SBarry Smith           col = dof*(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));
85247c6ae99SBarry Smith 
85347c6ae99SBarry Smith           cols[nc] = idx_c[col]/dof;
85447c6ae99SBarry Smith           v[nc++]  = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.));
85547c6ae99SBarry Smith 
85647c6ae99SBarry Smith           if (i_c*ratioi != i) {
85747c6ae99SBarry Smith             cols[nc] = idx_c[col+dof]/dof;
85847c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.));
85947c6ae99SBarry Smith           }
86047c6ae99SBarry Smith 
86147c6ae99SBarry Smith           if (j_c*ratioj != j) {
86247c6ae99SBarry Smith             cols[nc] = idx_c[col+m_ghost_c*dof]/dof;
86347c6ae99SBarry Smith             v[nc++]  = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.));
86447c6ae99SBarry Smith           }
86547c6ae99SBarry Smith 
86647c6ae99SBarry Smith           if (l_c*ratiok != l) {
86747c6ae99SBarry Smith             cols[nc] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof;
86847c6ae99SBarry Smith             v[nc++]  = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.));
86947c6ae99SBarry Smith           }
87047c6ae99SBarry Smith 
87147c6ae99SBarry Smith           if (j_c*ratioj != j && i_c*ratioi != i) {
87247c6ae99SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c+1)*dof]/dof;
87347c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.));
87447c6ae99SBarry Smith           }
87547c6ae99SBarry Smith 
87647c6ae99SBarry Smith           if (j_c*ratioj != j && l_c*ratiok != l) {
87747c6ae99SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof;
87847c6ae99SBarry Smith             v[nc++]  = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.));
87947c6ae99SBarry Smith           }
88047c6ae99SBarry Smith 
88147c6ae99SBarry Smith           if (i_c*ratioi != i && l_c*ratiok != l) {
88247c6ae99SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof;
88347c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.));
88447c6ae99SBarry Smith           }
88547c6ae99SBarry Smith 
88647c6ae99SBarry Smith           if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) {
88747c6ae99SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof;
88847c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.));
88947c6ae99SBarry Smith           }
89047c6ae99SBarry Smith           ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
89147c6ae99SBarry Smith         }
89247c6ae99SBarry Smith       }
89347c6ae99SBarry Smith     }
894b3bd94feSDave May 
895b3bd94feSDave May   } else {
896b3bd94feSDave May     PetscScalar *xi,*eta,*zeta;
897b3bd94feSDave May     PetscInt    li,nxi,lj,neta,lk,nzeta,n;
898b3bd94feSDave May     PetscScalar Ni[8];
899b3bd94feSDave May 
900b3bd94feSDave May     /* compute local coordinate arrays */
901b3bd94feSDave May     nxi   = ratioi + 1;
902b3bd94feSDave May     neta  = ratioj + 1;
903b3bd94feSDave May     nzeta = ratiok + 1;
904b3bd94feSDave May     ierr  = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr);
905b3bd94feSDave May     ierr  = PetscMalloc(sizeof(PetscScalar)*neta,&eta);CHKERRQ(ierr);
906b3bd94feSDave May     ierr  = PetscMalloc(sizeof(PetscScalar)*nzeta,&zeta);CHKERRQ(ierr);
9078865f1eaSKarl Rupp     for (li=0; li<nxi; li++) xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1));
9088865f1eaSKarl Rupp     for (lj=0; lj<neta; lj++) eta[lj] = -1.0 + (PetscScalar)lj*(2.0/(PetscScalar)(neta-1));
9098865f1eaSKarl Rupp     for (lk=0; lk<nzeta; lk++) zeta[lk] = -1.0 + (PetscScalar)lk*(2.0/(PetscScalar)(nzeta-1));
910b3bd94feSDave May 
911b3bd94feSDave May     for (l=l_start; l<l_start+p_f; l++) {
912b3bd94feSDave May       for (j=j_start; j<j_start+n_f; j++) {
913b3bd94feSDave May         for (i=i_start; i<i_start+m_f; i++) {
914b3bd94feSDave May           /* convert to local "natural" numbering and then to PETSc global numbering */
915b3bd94feSDave May           row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
916b3bd94feSDave May 
917b3bd94feSDave May           i_c = (i/ratioi);
918b3bd94feSDave May           j_c = (j/ratioj);
919b3bd94feSDave May           l_c = (l/ratiok);
920b3bd94feSDave May 
921b3bd94feSDave May           /* remainders */
922b3bd94feSDave May           li = i - ratioi * (i/ratioi);
9238865f1eaSKarl Rupp           if (i==mx-1) li = nxi-1;
924b3bd94feSDave May           lj = j - ratioj * (j/ratioj);
9258865f1eaSKarl Rupp           if (j==my-1) lj = neta-1;
926b3bd94feSDave May           lk = l - ratiok * (l/ratiok);
9278865f1eaSKarl Rupp           if (l==mz-1) lk = nzeta-1;
928b3bd94feSDave May 
929b3bd94feSDave May           /* corners */
930b3bd94feSDave May           col     = dof*(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));
931b3bd94feSDave May           cols[0] = idx_c[col]/dof;
932b3bd94feSDave May           Ni[0]   = 1.0;
933b3bd94feSDave May           if ((li==0) || (li==nxi-1)) {
934b3bd94feSDave May             if ((lj==0) || (lj==neta-1)) {
935b3bd94feSDave May               if ((lk==0) || (lk==nzeta-1)) {
936b3bd94feSDave May                 ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr);
937b3bd94feSDave May                 continue;
938b3bd94feSDave May               }
939b3bd94feSDave May             }
940b3bd94feSDave May           }
941b3bd94feSDave May 
942b3bd94feSDave May           /* edges + interior */
943b3bd94feSDave May           /* remainders */
9448865f1eaSKarl Rupp           if (i==mx-1) i_c--;
9458865f1eaSKarl Rupp           if (j==my-1) j_c--;
9468865f1eaSKarl Rupp           if (l==mz-1) l_c--;
947b3bd94feSDave May 
948b3bd94feSDave May           col     = dof*(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));
949b3bd94feSDave May           cols[0] = idx_c[col]/dof; /* one left and below; or we are right on it */
950b3bd94feSDave May           cols[1] = idx_c[col+dof]/dof; /* one right and below */
951b3bd94feSDave May           cols[2] = idx_c[col+m_ghost_c*dof]/dof;  /* one left and above */
952b3bd94feSDave May           cols[3] = idx_c[col+(m_ghost_c+1)*dof]/dof; /* one right and above */
953b3bd94feSDave May 
954b3bd94feSDave May           cols[4] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; /* one left and below and front; or we are right on it */
955b3bd94feSDave May           cols[5] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; /* one right and below, and front */
956b3bd94feSDave May           cols[6] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof; /* one left and above and front*/
957b3bd94feSDave May           cols[7] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; /* one right and above and front */
958b3bd94feSDave May 
959b3bd94feSDave May           Ni[0] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]);
960b3bd94feSDave May           Ni[1] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]);
961b3bd94feSDave May           Ni[2] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]);
962b3bd94feSDave May           Ni[3] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]);
963b3bd94feSDave May 
964b3bd94feSDave May           Ni[4] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]);
965b3bd94feSDave May           Ni[5] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]);
966b3bd94feSDave May           Ni[6] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]);
967b3bd94feSDave May           Ni[7] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]);
968b3bd94feSDave May 
969b3bd94feSDave May           for (n=0; n<8; n++) {
9708865f1eaSKarl Rupp             if (PetscAbsScalar(Ni[n])<1.0e-32) cols[n]=-1;
971b3bd94feSDave May           }
972b3bd94feSDave May           ierr = MatSetValues(mat,1,&row,8,cols,Ni,INSERT_VALUES);CHKERRQ(ierr);
973b3bd94feSDave May 
974b3bd94feSDave May         }
975b3bd94feSDave May       }
976b3bd94feSDave May     }
977b3bd94feSDave May     ierr = PetscFree(xi);CHKERRQ(ierr);
978b3bd94feSDave May     ierr = PetscFree(eta);CHKERRQ(ierr);
979b3bd94feSDave May     ierr = PetscFree(zeta);CHKERRQ(ierr);
980b3bd94feSDave May   }
981*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
982*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
98347c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
98447c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
98547c6ae99SBarry Smith 
98647c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
987fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
98847c6ae99SBarry Smith   PetscFunctionReturn(0);
98947c6ae99SBarry Smith }
99047c6ae99SBarry Smith 
99147c6ae99SBarry Smith #undef __FUNCT__
992e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA"
993e727c939SJed Brown PetscErrorCode  DMCreateInterpolation_DA(DM dac,DM daf,Mat *A,Vec *scale)
99447c6ae99SBarry Smith {
99547c6ae99SBarry Smith   PetscErrorCode   ierr;
99647c6ae99SBarry Smith   PetscInt         dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
9971321219cSEthan Coon   DMDABoundaryType bxc,byc,bzc,bxf,byf,bzf;
998aa219208SBarry Smith   DMDAStencilType  stc,stf;
99947c6ae99SBarry Smith   DM_DA            *ddc = (DM_DA*)dac->data;
100047c6ae99SBarry Smith 
100147c6ae99SBarry Smith   PetscFunctionBegin;
100247c6ae99SBarry Smith   PetscValidHeaderSpecific(dac,DM_CLASSID,1);
100347c6ae99SBarry Smith   PetscValidHeaderSpecific(daf,DM_CLASSID,2);
100447c6ae99SBarry Smith   PetscValidPointer(A,3);
100547c6ae99SBarry Smith   if (scale) PetscValidPointer(scale,4);
100647c6ae99SBarry Smith 
10071321219cSEthan Coon   ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr);
10081321219cSEthan Coon   ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr);
1009ce94432eSBarry Smith   if (dimc != dimf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr);
1010ce94432eSBarry Smith   if (dofc != doff) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr);
1011ce94432eSBarry Smith   if (sc != sf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr);
1012ce94432eSBarry Smith   if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr);
1013ce94432eSBarry Smith   if (stc != stf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr);
101447c6ae99SBarry Smith   if (Mc < 2 && Mf > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction");
101547c6ae99SBarry Smith   if (dimc > 1 && Nc < 2 && Nf > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in y direction");
101647c6ae99SBarry Smith   if (dimc > 2 && Pc < 2 && Pf > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in z direction");
101747c6ae99SBarry Smith 
1018aa219208SBarry Smith   if (ddc->interptype == DMDA_Q1) {
101947c6ae99SBarry Smith     if (dimc == 1) {
1020e727c939SJed Brown       ierr = DMCreateInterpolation_DA_1D_Q1(dac,daf,A);CHKERRQ(ierr);
102147c6ae99SBarry Smith     } else if (dimc == 2) {
1022e727c939SJed Brown       ierr = DMCreateInterpolation_DA_2D_Q1(dac,daf,A);CHKERRQ(ierr);
102347c6ae99SBarry Smith     } else if (dimc == 3) {
1024e727c939SJed Brown       ierr = DMCreateInterpolation_DA_3D_Q1(dac,daf,A);CHKERRQ(ierr);
1025ce94432eSBarry Smith     } else SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype);
1026aa219208SBarry Smith   } else if (ddc->interptype == DMDA_Q0) {
102747c6ae99SBarry Smith     if (dimc == 1) {
1028e727c939SJed Brown       ierr = DMCreateInterpolation_DA_1D_Q0(dac,daf,A);CHKERRQ(ierr);
102947c6ae99SBarry Smith     } else if (dimc == 2) {
1030e727c939SJed Brown       ierr = DMCreateInterpolation_DA_2D_Q0(dac,daf,A);CHKERRQ(ierr);
103147c6ae99SBarry Smith     } else if (dimc == 3) {
1032e727c939SJed Brown       ierr = DMCreateInterpolation_DA_3D_Q0(dac,daf,A);CHKERRQ(ierr);
1033ce94432eSBarry Smith     } else SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype);
103447c6ae99SBarry Smith   }
103547c6ae99SBarry Smith   if (scale) {
1036e727c939SJed Brown     ierr = DMCreateInterpolationScale((DM)dac,(DM)daf,*A,scale);CHKERRQ(ierr);
103747c6ae99SBarry Smith   }
103847c6ae99SBarry Smith   PetscFunctionReturn(0);
103947c6ae99SBarry Smith }
104047c6ae99SBarry Smith 
104147c6ae99SBarry Smith #undef __FUNCT__
1042e727c939SJed Brown #define __FUNCT__ "DMCreateInjection_DA_1D"
1043e727c939SJed Brown PetscErrorCode DMCreateInjection_DA_1D(DM dac,DM daf,VecScatter *inject)
10440bb2b966SJungho Lee {
10450bb2b966SJungho Lee   PetscErrorCode   ierr;
1046*8ea3bf28SBarry Smith   PetscInt         i,i_start,m_f,Mx,dof;
1047*8ea3bf28SBarry Smith   const PetscInt   *idx_f;
1048*8ea3bf28SBarry Smith   PetscInt         m_ghost,m_ghost_c;
10490bb2b966SJungho Lee   PetscInt         row,i_start_ghost,mx,m_c,nc,ratioi;
10500bb2b966SJungho Lee   PetscInt         i_start_c,i_start_ghost_c;
10510bb2b966SJungho Lee   PetscInt         *cols;
10520bb2b966SJungho Lee   DMDABoundaryType bx;
10530bb2b966SJungho Lee   Vec              vecf,vecc;
10540bb2b966SJungho Lee   IS               isf;
10550bb2b966SJungho Lee 
10560bb2b966SJungho Lee   PetscFunctionBegin;
10570bb2b966SJungho Lee   ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr);
10580bb2b966SJungho Lee   ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
10590bb2b966SJungho Lee   if (bx == DMDA_BOUNDARY_PERIODIC) {
10600bb2b966SJungho Lee     ratioi = mx/Mx;
10610bb2b966SJungho Lee     if (ratioi*Mx != mx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: mx/Mx  must be integer: mx %D Mx %D",mx,Mx);
10620bb2b966SJungho Lee   } else {
10630bb2b966SJungho Lee     ratioi = (mx-1)/(Mx-1);
10640bb2b966SJungho Lee     if (ratioi*(Mx-1) != mx-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %D Mx %D",mx,Mx);
10650bb2b966SJungho Lee   }
10660bb2b966SJungho Lee 
10670bb2b966SJungho Lee   ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr);
10680bb2b966SJungho Lee   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr);
10690298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
10700bb2b966SJungho Lee 
10710bb2b966SJungho Lee   ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr);
10720bb2b966SJungho Lee   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr);
10730bb2b966SJungho Lee 
10740bb2b966SJungho Lee 
10750bb2b966SJungho Lee   /* loop over local fine grid nodes setting interpolation for those*/
10760bb2b966SJungho Lee   nc   = 0;
10770bb2b966SJungho Lee   ierr = PetscMalloc(m_f*sizeof(PetscInt),&cols);CHKERRQ(ierr);
10780bb2b966SJungho Lee 
10790bb2b966SJungho Lee 
10800bb2b966SJungho Lee   for (i=i_start_c; i<i_start_c+m_c; i++) {
10810bb2b966SJungho Lee     PetscInt i_f = i*ratioi;
10820bb2b966SJungho Lee 
10838865f1eaSKarl Rupp     if (i_f < i_start_ghost || i_f >= i_start_ghost+m_ghost) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\ni_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost);
10848865f1eaSKarl Rupp 
10850bb2b966SJungho Lee     row        = idx_f[dof*(i_f-i_start_ghost)];
10860bb2b966SJungho Lee     cols[nc++] = row/dof;
10870bb2b966SJungho Lee   }
10880bb2b966SJungho Lee 
1089*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
10900bb2b966SJungho Lee 
1091ce94432eSBarry Smith   ierr = ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr);
10920bb2b966SJungho Lee   ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr);
10930bb2b966SJungho Lee   ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr);
10940298fd71SBarry Smith   ierr = VecScatterCreate(vecf,isf,vecc,NULL,inject);CHKERRQ(ierr);
10950bb2b966SJungho Lee   ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr);
10960bb2b966SJungho Lee   ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr);
10970bb2b966SJungho Lee   ierr = ISDestroy(&isf);CHKERRQ(ierr);
10980bb2b966SJungho Lee   PetscFunctionReturn(0);
10990bb2b966SJungho Lee }
11000bb2b966SJungho Lee 
11010bb2b966SJungho Lee #undef __FUNCT__
1102e727c939SJed Brown #define __FUNCT__ "DMCreateInjection_DA_2D"
1103e727c939SJed Brown PetscErrorCode DMCreateInjection_DA_2D(DM dac,DM daf,VecScatter *inject)
110447c6ae99SBarry Smith {
110547c6ae99SBarry Smith   PetscErrorCode   ierr;
1106*8ea3bf28SBarry Smith   PetscInt         i,j,i_start,j_start,m_f,n_f,Mx,My,dof;
1107*8ea3bf28SBarry Smith   const PetscInt   *idx_c,*idx_f;
1108*8ea3bf28SBarry Smith   PetscInt         m_ghost,n_ghost,m_ghost_c,n_ghost_c;
110947c6ae99SBarry Smith   PetscInt         row,i_start_ghost,j_start_ghost,mx,m_c,my,nc,ratioi,ratioj;
1110076202ddSJed Brown   PetscInt         i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c;
111147c6ae99SBarry Smith   PetscInt         *cols;
11121321219cSEthan Coon   DMDABoundaryType bx,by;
111347c6ae99SBarry Smith   Vec              vecf,vecc;
111447c6ae99SBarry Smith   IS               isf;
111547c6ae99SBarry Smith 
111647c6ae99SBarry Smith   PetscFunctionBegin;
11171321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr);
11181321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
11191321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) {
112047c6ae99SBarry Smith     ratioi = mx/Mx;
112147c6ae99SBarry Smith     if (ratioi*Mx != mx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: mx/Mx  must be integer: mx %D Mx %D",mx,Mx);
112247c6ae99SBarry Smith   } else {
112347c6ae99SBarry Smith     ratioi = (mx-1)/(Mx-1);
112447c6ae99SBarry Smith     if (ratioi*(Mx-1) != mx-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %D Mx %D",mx,Mx);
112547c6ae99SBarry Smith   }
11261321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC) {
112747c6ae99SBarry Smith     ratioj = my/My;
112847c6ae99SBarry Smith     if (ratioj*My != my) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: my/My  must be integer: my %D My %D",my,My);
112947c6ae99SBarry Smith   } else {
113047c6ae99SBarry Smith     ratioj = (my-1)/(My-1);
113147c6ae99SBarry Smith     if (ratioj*(My-1) != my-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (my - 1)/(My - 1) must be integer: my %D My %D",my,My);
113247c6ae99SBarry Smith   }
113347c6ae99SBarry Smith 
1134aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr);
1135aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr);
11360298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
113747c6ae99SBarry Smith 
1138aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr);
1139aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr);
11400298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
114147c6ae99SBarry Smith 
114247c6ae99SBarry Smith 
114347c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
114447c6ae99SBarry Smith   nc   = 0;
114547c6ae99SBarry Smith   ierr = PetscMalloc(n_f*m_f*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1146076202ddSJed Brown   for (j=j_start_c; j<j_start_c+n_c; j++) {
1147076202ddSJed Brown     for (i=i_start_c; i<i_start_c+m_c; i++) {
1148076202ddSJed Brown       PetscInt i_f = i*ratioi,j_f = j*ratioj;
1149076202ddSJed Brown       if (j_f < j_start_ghost || j_f >= j_start_ghost+n_ghost) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
1150076202ddSJed Brown     j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost);
1151076202ddSJed Brown       if (i_f < i_start_ghost || i_f >= i_start_ghost+m_ghost) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
1152076202ddSJed Brown     i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost);
1153076202ddSJed Brown       row        = idx_f[dof*(m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))];
115447c6ae99SBarry Smith       cols[nc++] = row/dof;
115547c6ae99SBarry Smith     }
115647c6ae99SBarry Smith   }
1157*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
1158*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
115947c6ae99SBarry Smith 
1160ce94432eSBarry Smith   ierr = ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr);
11619a42bb27SBarry Smith   ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr);
11629a42bb27SBarry Smith   ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr);
11630298fd71SBarry Smith   ierr = VecScatterCreate(vecf,isf,vecc,NULL,inject);CHKERRQ(ierr);
11649a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr);
11659a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr);
1166fcfd50ebSBarry Smith   ierr = ISDestroy(&isf);CHKERRQ(ierr);
116747c6ae99SBarry Smith   PetscFunctionReturn(0);
116847c6ae99SBarry Smith }
116947c6ae99SBarry Smith 
117047c6ae99SBarry Smith #undef __FUNCT__
1171e727c939SJed Brown #define __FUNCT__ "DMCreateInjection_DA_3D"
1172e727c939SJed Brown PetscErrorCode DMCreateInjection_DA_3D(DM dac,DM daf,VecScatter *inject)
1173b1757049SJed Brown {
1174b1757049SJed Brown   PetscErrorCode   ierr;
1175b1757049SJed Brown   PetscInt         i,j,k,i_start,j_start,k_start,m_f,n_f,p_f,Mx,My,Mz;
1176b1757049SJed Brown   PetscInt         m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c;
1177b1757049SJed Brown   PetscInt         i_start_ghost,j_start_ghost,k_start_ghost;
1178b1757049SJed Brown   PetscInt         mx,my,mz,ratioi,ratioj,ratiok;
1179b1757049SJed Brown   PetscInt         i_start_c,j_start_c,k_start_c;
1180b1757049SJed Brown   PetscInt         m_c,n_c,p_c;
1181b1757049SJed Brown   PetscInt         i_start_ghost_c,j_start_ghost_c,k_start_ghost_c;
1182b1757049SJed Brown   PetscInt         row,nc,dof;
1183*8ea3bf28SBarry Smith   const PetscInt   *idx_c,*idx_f;
1184b1757049SJed Brown   PetscInt         *cols;
11851321219cSEthan Coon   DMDABoundaryType bx,by,bz;
1186b1757049SJed Brown   Vec              vecf,vecc;
1187b1757049SJed Brown   IS               isf;
1188b1757049SJed Brown 
1189b1757049SJed Brown   PetscFunctionBegin;
11901321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr);
11911321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
1192b1757049SJed Brown 
11931321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) {
1194b1757049SJed Brown     ratioi = mx/Mx;
1195b1757049SJed Brown     if (ratioi*Mx != mx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: mx/Mx  must be integer: mx %D Mx %D",mx,Mx);
1196b1757049SJed Brown   } else {
1197b1757049SJed Brown     ratioi = (mx-1)/(Mx-1);
1198b1757049SJed Brown     if (ratioi*(Mx-1) != mx-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %D Mx %D",mx,Mx);
1199b1757049SJed Brown   }
12001321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC) {
1201b1757049SJed Brown     ratioj = my/My;
1202b1757049SJed Brown     if (ratioj*My != my) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: my/My  must be integer: my %D My %D",my,My);
1203b1757049SJed Brown   } else {
1204b1757049SJed Brown     ratioj = (my-1)/(My-1);
1205b1757049SJed Brown     if (ratioj*(My-1) != my-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (my - 1)/(My - 1) must be integer: my %D My %D",my,My);
1206b1757049SJed Brown   }
12071321219cSEthan Coon   if (bz == DMDA_BOUNDARY_PERIODIC) {
1208b1757049SJed Brown     ratiok = mz/Mz;
1209b1757049SJed Brown     if (ratiok*Mz != mz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: mz/Mz  must be integer: mz %D My %D",mz,Mz);
1210b1757049SJed Brown   } else {
1211b1757049SJed Brown     ratiok = (mz-1)/(Mz-1);
1212b1757049SJed Brown     if (ratiok*(Mz-1) != mz-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (mz - 1)/(Mz - 1) must be integer: mz %D Mz %D",mz,Mz);
1213b1757049SJed Brown   }
1214b1757049SJed Brown 
1215b1757049SJed Brown   ierr = DMDAGetCorners(daf,&i_start,&j_start,&k_start,&m_f,&n_f,&p_f);CHKERRQ(ierr);
1216b1757049SJed Brown   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&k_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr);
12170298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
1218b1757049SJed Brown 
1219b1757049SJed Brown   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&k_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr);
1220b1757049SJed Brown   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,&k_start_ghost_c,&m_ghost_c,&n_ghost_c,&p_ghost_c);CHKERRQ(ierr);
12210298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
1222b1757049SJed Brown 
1223b1757049SJed Brown 
1224b1757049SJed Brown   /* loop over local fine grid nodes setting interpolation for those*/
1225b1757049SJed Brown   nc   = 0;
1226b1757049SJed Brown   ierr = PetscMalloc(n_f*m_f*p_f*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1227b1757049SJed Brown   for (k=k_start_c; k<k_start_c+p_c; k++) {
1228b1757049SJed Brown     for (j=j_start_c; j<j_start_c+n_c; j++) {
1229b1757049SJed Brown       for (i=i_start_c; i<i_start_c+m_c; i++) {
1230b1757049SJed Brown         PetscInt i_f = i*ratioi,j_f = j*ratioj,k_f = k*ratiok;
1231b1757049SJed Brown         if (k_f < k_start_ghost || k_f >= k_start_ghost+p_ghost) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA  "
1232b1757049SJed Brown                                                                           "k_c %D k_f %D fine ghost range [%D,%D]",k,k_f,k_start_ghost,k_start_ghost+p_ghost);
1233b1757049SJed Brown         if (j_f < j_start_ghost || j_f >= j_start_ghost+n_ghost) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA  "
1234b1757049SJed Brown                                                                           "j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost);
1235b1757049SJed Brown         if (i_f < i_start_ghost || i_f >= i_start_ghost+m_ghost) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA  "
1236b1757049SJed Brown                                                                           "i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost);
1237b1757049SJed Brown         row        = idx_f[dof*(m_ghost*n_ghost*(k_f-k_start_ghost) + m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))];
1238b1757049SJed Brown         cols[nc++] = row/dof;
1239b1757049SJed Brown       }
1240b1757049SJed Brown     }
1241b1757049SJed Brown   }
1242*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
1243*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
1244b1757049SJed Brown 
1245ce94432eSBarry Smith   ierr = ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr);
1246b1757049SJed Brown   ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr);
1247b1757049SJed Brown   ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr);
12480298fd71SBarry Smith   ierr = VecScatterCreate(vecf,isf,vecc,NULL,inject);CHKERRQ(ierr);
1249b1757049SJed Brown   ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr);
1250b1757049SJed Brown   ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr);
1251fcfd50ebSBarry Smith   ierr = ISDestroy(&isf);CHKERRQ(ierr);
1252b1757049SJed Brown   PetscFunctionReturn(0);
1253b1757049SJed Brown }
1254b1757049SJed Brown 
1255b1757049SJed Brown #undef __FUNCT__
1256e727c939SJed Brown #define __FUNCT__ "DMCreateInjection_DA"
1257e727c939SJed Brown PetscErrorCode  DMCreateInjection_DA(DM dac,DM daf,VecScatter *inject)
125847c6ae99SBarry Smith {
125947c6ae99SBarry Smith   PetscErrorCode   ierr;
126047c6ae99SBarry Smith   PetscInt         dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
12611321219cSEthan Coon   DMDABoundaryType bxc,byc,bzc,bxf,byf,bzf;
1262aa219208SBarry Smith   DMDAStencilType  stc,stf;
126347c6ae99SBarry Smith 
126447c6ae99SBarry Smith   PetscFunctionBegin;
126547c6ae99SBarry Smith   PetscValidHeaderSpecific(dac,DM_CLASSID,1);
126647c6ae99SBarry Smith   PetscValidHeaderSpecific(daf,DM_CLASSID,2);
126747c6ae99SBarry Smith   PetscValidPointer(inject,3);
126847c6ae99SBarry Smith 
12691321219cSEthan Coon   ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr);
12701321219cSEthan Coon   ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr);
1271ce94432eSBarry Smith   if (dimc != dimf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr);
1272ce94432eSBarry Smith   if (dofc != doff) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr);
1273ce94432eSBarry Smith   if (sc != sf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr);
1274ce94432eSBarry Smith   if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr);
1275ce94432eSBarry Smith   if (stc != stf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr);
127647c6ae99SBarry Smith   if (Mc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction");
127747c6ae99SBarry Smith   if (dimc > 1 && Nc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in y direction");
127847c6ae99SBarry Smith   if (dimc > 2 && Pc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in z direction");
127947c6ae99SBarry Smith 
12800bb2b966SJungho Lee   if (dimc == 1) {
1281e727c939SJed Brown     ierr = DMCreateInjection_DA_1D(dac,daf,inject);CHKERRQ(ierr);
12820bb2b966SJungho Lee   } else if (dimc == 2) {
1283e727c939SJed Brown     ierr = DMCreateInjection_DA_2D(dac,daf,inject);CHKERRQ(ierr);
1284b1757049SJed Brown   } else if (dimc == 3) {
1285e727c939SJed Brown     ierr = DMCreateInjection_DA_3D(dac,daf,inject);CHKERRQ(ierr);
128647c6ae99SBarry Smith   }
128747c6ae99SBarry Smith   PetscFunctionReturn(0);
128847c6ae99SBarry Smith }
128947c6ae99SBarry Smith 
129047c6ae99SBarry Smith #undef __FUNCT__
1291e727c939SJed Brown #define __FUNCT__ "DMCreateAggregates_DA"
1292e727c939SJed Brown PetscErrorCode  DMCreateAggregates_DA(DM dac,DM daf,Mat *rest)
129347c6ae99SBarry Smith {
129447c6ae99SBarry Smith   PetscErrorCode   ierr;
129547c6ae99SBarry Smith   PetscInt         dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc;
129647c6ae99SBarry Smith   PetscInt         dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
12971321219cSEthan Coon   DMDABoundaryType bxc,byc,bzc,bxf,byf,bzf;
1298aa219208SBarry Smith   DMDAStencilType  stc,stf;
129947c6ae99SBarry Smith   PetscInt         i,j,l;
130047c6ae99SBarry Smith   PetscInt         i_start,j_start,l_start, m_f,n_f,p_f;
130147c6ae99SBarry Smith   PetscInt         i_start_ghost,j_start_ghost,l_start_ghost,m_ghost,n_ghost,p_ghost;
1302*8ea3bf28SBarry Smith   const PetscInt   *idx_f;
130347c6ae99SBarry Smith   PetscInt         i_c,j_c,l_c;
130447c6ae99SBarry Smith   PetscInt         i_start_c,j_start_c,l_start_c, m_c,n_c,p_c;
130547c6ae99SBarry Smith   PetscInt         i_start_ghost_c,j_start_ghost_c,l_start_ghost_c,m_ghost_c,n_ghost_c,p_ghost_c;
1306*8ea3bf28SBarry Smith   const PetscInt   *idx_c;
130747c6ae99SBarry Smith   PetscInt         d;
130847c6ae99SBarry Smith   PetscInt         a;
130947c6ae99SBarry Smith   PetscInt         max_agg_size;
131047c6ae99SBarry Smith   PetscInt         *fine_nodes;
131147c6ae99SBarry Smith   PetscScalar      *one_vec;
131247c6ae99SBarry Smith   PetscInt         fn_idx;
131347c6ae99SBarry Smith 
131447c6ae99SBarry Smith   PetscFunctionBegin;
131547c6ae99SBarry Smith   PetscValidHeaderSpecific(dac,DM_CLASSID,1);
131647c6ae99SBarry Smith   PetscValidHeaderSpecific(daf,DM_CLASSID,2);
131747c6ae99SBarry Smith   PetscValidPointer(rest,3);
131847c6ae99SBarry Smith 
13191321219cSEthan Coon   ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr);
13201321219cSEthan Coon   ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr);
1321ce94432eSBarry Smith   if (dimc != dimf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr);
1322ce94432eSBarry Smith   if (dofc != doff) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr);
1323ce94432eSBarry Smith   if (sc != sf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr);
1324ce94432eSBarry Smith   if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr);
1325ce94432eSBarry Smith   if (stc != stf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr);
132647c6ae99SBarry Smith 
132747c6ae99SBarry Smith   if (Mf < Mc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Coarse grid has more points than fine grid, Mc %D, Mf %D", Mc, Mf);
132847c6ae99SBarry Smith   if (Nf < Nc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Coarse grid has more points than fine grid, Nc %D, Nf %D", Nc, Nf);
132947c6ae99SBarry Smith   if (Pf < Pc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Coarse grid has more points than fine grid, Pc %D, Pf %D", Pc, Pf);
133047c6ae99SBarry Smith 
133147c6ae99SBarry Smith   if (Pc < 0) Pc = 1;
133247c6ae99SBarry Smith   if (Pf < 0) Pf = 1;
133347c6ae99SBarry Smith   if (Nc < 0) Nc = 1;
133447c6ae99SBarry Smith   if (Nf < 0) Nf = 1;
133547c6ae99SBarry Smith 
1336aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr);
1337aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr);
13380298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
133947c6ae99SBarry Smith 
1340aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr);
1341aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,&l_start_ghost_c,&m_ghost_c,&n_ghost_c,&p_ghost_c);CHKERRQ(ierr);
13420298fd71SBarry Smith   ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
134347c6ae99SBarry Smith 
134447c6ae99SBarry Smith   /*
134547c6ae99SBarry Smith      Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios
134647c6ae99SBarry Smith      for dimension 1 and 2 respectively.
134747c6ae99SBarry Smith      Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1)
134847c6ae99SBarry Smith      and r_y*j and r_y*(j+1) will be grouped into the same coarse grid agregate.
134947c6ae99SBarry Smith      Each specific dof on the fine grid is mapped to one dof on the coarse grid.
135047c6ae99SBarry Smith   */
135147c6ae99SBarry Smith 
135247c6ae99SBarry Smith   max_agg_size = (Mf/Mc+1)*(Nf/Nc+1)*(Pf/Pc+1);
135347c6ae99SBarry Smith 
135447c6ae99SBarry Smith   /* create the matrix that will contain the restriction operator */
1355ce94432eSBarry Smith   ierr = 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,
13560298fd71SBarry Smith                       max_agg_size, NULL, max_agg_size, NULL, rest);CHKERRQ(ierr);
135747c6ae99SBarry Smith 
135847c6ae99SBarry Smith   /* store nodes in the fine grid here */
135947c6ae99SBarry Smith   ierr = PetscMalloc2(max_agg_size,PetscScalar, &one_vec,max_agg_size,PetscInt, &fine_nodes);CHKERRQ(ierr);
136047c6ae99SBarry Smith   for (i=0; i<max_agg_size; i++) one_vec[i] = 1.0;
136147c6ae99SBarry Smith 
136247c6ae99SBarry Smith   /* loop over all coarse nodes */
136347c6ae99SBarry Smith   for (l_c=l_start_c; l_c<l_start_c+p_c; l_c++) {
136447c6ae99SBarry Smith     for (j_c=j_start_c; j_c<j_start_c+n_c; j_c++) {
136547c6ae99SBarry Smith       for (i_c=i_start_c; i_c<i_start_c+m_c; i_c++) {
136647c6ae99SBarry Smith         for (d=0; d<dofc; d++) {
136747c6ae99SBarry Smith           /* convert to local "natural" numbering and then to PETSc global numbering */
136847c6ae99SBarry 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;
136947c6ae99SBarry Smith 
137047c6ae99SBarry Smith           fn_idx = 0;
137147c6ae99SBarry Smith           /* Corresponding fine points are all points (i_f, j_f, l_f) such that
137247c6ae99SBarry Smith              i_c*Mf/Mc <= i_f < (i_c+1)*Mf/Mc
137347c6ae99SBarry Smith              (same for other dimensions)
137447c6ae99SBarry Smith           */
137547c6ae99SBarry Smith           for (l=l_c*Pf/Pc; l<PetscMin((l_c+1)*Pf/Pc,Pf); l++) {
137647c6ae99SBarry Smith             for (j=j_c*Nf/Nc; j<PetscMin((j_c+1)*Nf/Nc,Nf); j++) {
137747c6ae99SBarry Smith               for (i=i_c*Mf/Mc; i<PetscMin((i_c+1)*Mf/Mc,Mf); i++) {
137847c6ae99SBarry 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;
137947c6ae99SBarry Smith                 fn_idx++;
138047c6ae99SBarry Smith               }
138147c6ae99SBarry Smith             }
138247c6ae99SBarry Smith           }
138347c6ae99SBarry Smith           /* add all these points to one aggregate */
138447c6ae99SBarry Smith           ierr = MatSetValues(*rest, 1, &a, fn_idx, fine_nodes, one_vec, INSERT_VALUES);CHKERRQ(ierr);
138547c6ae99SBarry Smith         }
138647c6ae99SBarry Smith       }
138747c6ae99SBarry Smith     }
138847c6ae99SBarry Smith   }
1389*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr);
1390*8ea3bf28SBarry Smith   ierr = DMDARestoreGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr);
139147c6ae99SBarry Smith   ierr = PetscFree2(one_vec,fine_nodes);CHKERRQ(ierr);
139247c6ae99SBarry Smith   ierr = MatAssemblyBegin(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
139347c6ae99SBarry Smith   ierr = MatAssemblyEnd(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
139447c6ae99SBarry Smith   PetscFunctionReturn(0);
139547c6ae99SBarry Smith }
1396