xref: /petsc/src/dm/impls/da/dainterp.c (revision 0bb2b9666eb3df51ad3fb835297fd983308a762b)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
3aa219208SBarry Smith   Code for interpolating between grids represented by DMDAs
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
6c6db04a5SJed Brown #include <private/daimpl.h>    /*I   "petscdmda.h"   I*/
7c6db04a5SJed Brown #include <petscpcmg.h>
847c6ae99SBarry Smith 
947c6ae99SBarry Smith #undef __FUNCT__
1047c6ae99SBarry Smith #define __FUNCT__ "DMGetInterpolationScale"
1147c6ae99SBarry Smith /*@
1247c6ae99SBarry Smith     DMGetInterpolationScale - Forms L = R*1/diag(R*1) - L.*v is like a coarse grid average of the
1347c6ae99SBarry Smith       nearby fine grid points.
1447c6ae99SBarry Smith 
1547c6ae99SBarry Smith   Input Parameters:
1647c6ae99SBarry Smith +      dac - DM that defines a coarse mesh
1747c6ae99SBarry Smith .      daf - DM that defines a fine mesh
1847c6ae99SBarry Smith -      mat - the restriction (or interpolation operator) from fine to coarse
1947c6ae99SBarry Smith 
2047c6ae99SBarry Smith   Output Parameter:
2147c6ae99SBarry Smith .    scale - the scaled vector
2247c6ae99SBarry Smith 
2347c6ae99SBarry Smith   Level: developer
2447c6ae99SBarry Smith 
2547c6ae99SBarry Smith .seealso: DMGetInterpolation()
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith @*/
287087cfbeSBarry Smith PetscErrorCode  DMGetInterpolationScale(DM dac,DM daf,Mat mat,Vec *scale)
2947c6ae99SBarry Smith {
3047c6ae99SBarry Smith   PetscErrorCode ierr;
3147c6ae99SBarry Smith   Vec            fine;
3247c6ae99SBarry Smith   PetscScalar    one = 1.0;
3347c6ae99SBarry Smith 
3447c6ae99SBarry Smith   PetscFunctionBegin;
3547c6ae99SBarry Smith   ierr = DMCreateGlobalVector(daf,&fine);CHKERRQ(ierr);
3647c6ae99SBarry Smith   ierr = DMCreateGlobalVector(dac,scale);CHKERRQ(ierr);
3747c6ae99SBarry Smith   ierr = VecSet(fine,one);CHKERRQ(ierr);
3847c6ae99SBarry Smith   ierr = MatRestrict(mat,fine,*scale);CHKERRQ(ierr);
39fcfd50ebSBarry Smith   ierr = VecDestroy(&fine);CHKERRQ(ierr);
4047c6ae99SBarry Smith   ierr = VecReciprocal(*scale);CHKERRQ(ierr);
4147c6ae99SBarry Smith   PetscFunctionReturn(0);
4247c6ae99SBarry Smith }
4347c6ae99SBarry Smith 
4447c6ae99SBarry Smith #undef __FUNCT__
4594013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_1D_Q1"
4694013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_1D_Q1(DM dac,DM daf,Mat *A)
4747c6ae99SBarry Smith {
4847c6ae99SBarry Smith   PetscErrorCode   ierr;
4947c6ae99SBarry Smith   PetscInt         i,i_start,m_f,Mx,*idx_f;
5047c6ae99SBarry Smith   PetscInt         m_ghost,*idx_c,m_ghost_c;
5147c6ae99SBarry Smith   PetscInt         row,col,i_start_ghost,mx,m_c,nc,ratio;
5247c6ae99SBarry Smith   PetscInt         i_c,i_start_c,i_start_ghost_c,cols[2],dof;
5347c6ae99SBarry Smith   PetscScalar      v[2],x,*coors = 0,*ccoors;
5447c6ae99SBarry Smith   Mat              mat;
551321219cSEthan Coon   DMDABoundaryType bx;
5647c6ae99SBarry Smith   Vec              vcoors,cvcoors;
5747c6ae99SBarry Smith   DM_DA            *ddc = (DM_DA*)dac->data, *ddf = (DM_DA*)daf->data;
5847c6ae99SBarry Smith 
5947c6ae99SBarry Smith   PetscFunctionBegin;
601321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr);
611321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
621321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC){
6347c6ae99SBarry Smith     ratio = mx/Mx;
6447c6ae99SBarry 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);
6547c6ae99SBarry Smith   } else {
6647c6ae99SBarry Smith     ratio = (mx-1)/(Mx-1);
6747c6ae99SBarry 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);
6847c6ae99SBarry Smith   }
6947c6ae99SBarry Smith 
70aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr);
71aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr);
72aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr);
7347c6ae99SBarry Smith 
74aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr);
75aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr);
76aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr);
7747c6ae99SBarry Smith 
7847c6ae99SBarry Smith   /* create interpolation matrix */
7947c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr);
8047c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr);
8147c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
8247c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,2,PETSC_NULL);CHKERRQ(ierr);
8347c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,2,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
8447c6ae99SBarry Smith 
85aa219208SBarry Smith   ierr = DMDAGetCoordinates(daf,&vcoors);CHKERRQ(ierr);
8647c6ae99SBarry Smith   if (vcoors) {
87aa219208SBarry Smith     ierr = DMDAGetGhostedCoordinates(dac,&cvcoors);CHKERRQ(ierr);
88aa219208SBarry Smith     ierr = DMDAVecGetArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr);
89aa219208SBarry Smith     ierr = DMDAVecGetArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr);
9047c6ae99SBarry Smith   }
9147c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
92b3bd94feSDave May   if (!vcoors) {
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);
143b3bd94feSDave May       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 */
156b3bd94feSDave May       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++) {
165b3bd94feSDave May         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   }
17147c6ae99SBarry Smith   if (vcoors) {
172aa219208SBarry Smith     ierr = DMDAVecRestoreArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr);
173aa219208SBarry Smith     ierr = DMDAVecRestoreArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr);
17447c6ae99SBarry Smith   }
17547c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17647c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17747c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
178fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
17947c6ae99SBarry Smith   PetscFunctionReturn(0);
18047c6ae99SBarry Smith }
18147c6ae99SBarry Smith 
18247c6ae99SBarry Smith #undef __FUNCT__
18394013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_1D_Q0"
18494013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_1D_Q0(DM dac,DM daf,Mat *A)
18547c6ae99SBarry Smith {
18647c6ae99SBarry Smith   PetscErrorCode   ierr;
18747c6ae99SBarry Smith   PetscInt         i,i_start,m_f,Mx,*idx_f;
18847c6ae99SBarry Smith   PetscInt         m_ghost,*idx_c,m_ghost_c;
18947c6ae99SBarry Smith   PetscInt         row,col,i_start_ghost,mx,m_c,nc,ratio;
19047c6ae99SBarry Smith   PetscInt         i_c,i_start_c,i_start_ghost_c,cols[2],dof;
19147c6ae99SBarry Smith   PetscScalar      v[2],x;
19247c6ae99SBarry Smith   Mat              mat;
1931321219cSEthan Coon   DMDABoundaryType bx;
19447c6ae99SBarry Smith 
19547c6ae99SBarry Smith   PetscFunctionBegin;
1961321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr);
1971321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
1981321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC){
19947c6ae99SBarry Smith     ratio = mx/Mx;
20047c6ae99SBarry 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);
20147c6ae99SBarry Smith   } else {
20247c6ae99SBarry Smith     ratio = (mx-1)/(Mx-1);
20347c6ae99SBarry 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);
20447c6ae99SBarry Smith   }
20547c6ae99SBarry Smith 
206aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr);
207aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr);
208aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr);
20947c6ae99SBarry Smith 
210aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr);
211aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr);
212aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr);
21347c6ae99SBarry Smith 
21447c6ae99SBarry Smith   /* create interpolation matrix */
21547c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr);
21647c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr);
21747c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
21847c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,2,PETSC_NULL);CHKERRQ(ierr);
21947c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,2,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
22047c6ae99SBarry Smith 
22147c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
22247c6ae99SBarry Smith   for (i=i_start; i<i_start+m_f; i++) {
22347c6ae99SBarry Smith     /* convert to local "natural" numbering and then to PETSc global numbering */
22447c6ae99SBarry Smith     row    = idx_f[dof*(i-i_start_ghost)]/dof;
22547c6ae99SBarry Smith 
22647c6ae99SBarry Smith     i_c = (i/ratio);    /* coarse grid node to left of fine grid node */
22747c6ae99SBarry Smith 
22847c6ae99SBarry Smith     /*
22947c6ae99SBarry Smith          Only include those interpolation points that are truly
23047c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
23147c6ae99SBarry Smith          in x direction; since they have no right neighbor
23247c6ae99SBarry Smith     */
23347c6ae99SBarry Smith     x  = ((double)(i - i_c*ratio))/((double)ratio);
23447c6ae99SBarry Smith     nc = 0;
23547c6ae99SBarry Smith       /* one left and below; or we are right on it */
23647c6ae99SBarry Smith     col      = dof*(i_c-i_start_ghost_c);
23747c6ae99SBarry Smith     cols[nc] = idx_c[col]/dof;
23847c6ae99SBarry Smith     v[nc++]  = - x + 1.0;
23947c6ae99SBarry Smith     /* one right? */
24047c6ae99SBarry Smith     if (i_c*ratio != i) {
24147c6ae99SBarry Smith       cols[nc] = idx_c[col+dof]/dof;
24247c6ae99SBarry Smith       v[nc++]  = x;
24347c6ae99SBarry Smith     }
24447c6ae99SBarry Smith     ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
24547c6ae99SBarry Smith   }
24647c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
24747c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
24847c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
249fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
25047c6ae99SBarry Smith   ierr = PetscLogFlops(5.0*m_f);CHKERRQ(ierr);
25147c6ae99SBarry Smith   PetscFunctionReturn(0);
25247c6ae99SBarry Smith }
25347c6ae99SBarry Smith 
25447c6ae99SBarry Smith #undef __FUNCT__
25594013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_2D_Q1"
25694013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_2D_Q1(DM dac,DM daf,Mat *A)
25747c6ae99SBarry Smith {
25847c6ae99SBarry Smith   PetscErrorCode   ierr;
25947c6ae99SBarry Smith   PetscInt         i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof;
26047c6ae99SBarry Smith   PetscInt         m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz;
26147c6ae99SBarry Smith   PetscInt         row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj;
26247c6ae99SBarry 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;
26347c6ae99SBarry Smith   PetscMPIInt      size_c,size_f,rank_f;
26447c6ae99SBarry Smith   PetscScalar      v[4],x,y;
26547c6ae99SBarry Smith   Mat              mat;
2661321219cSEthan Coon   DMDABoundaryType bx,by;
267aa219208SBarry Smith   DMDACoor2d       **coors = 0,**ccoors;
26847c6ae99SBarry Smith   Vec              vcoors,cvcoors;
26947c6ae99SBarry Smith   DM_DA            *ddc = (DM_DA*)dac->data, *ddf = (DM_DA*)daf->data;
27047c6ae99SBarry Smith 
27147c6ae99SBarry Smith   PetscFunctionBegin;
2721321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr);
2731321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
2741321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC){
27547c6ae99SBarry Smith     ratioi = mx/Mx;
27647c6ae99SBarry 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);
27747c6ae99SBarry Smith   } else {
27847c6ae99SBarry Smith     ratioi = (mx-1)/(Mx-1);
27947c6ae99SBarry 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);
28047c6ae99SBarry Smith   }
2811321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC){
28247c6ae99SBarry Smith     ratioj = my/My;
28347c6ae99SBarry 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);
28447c6ae99SBarry Smith   } else {
28547c6ae99SBarry Smith     ratioj = (my-1)/(My-1);
28647c6ae99SBarry 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);
28747c6ae99SBarry Smith   }
28847c6ae99SBarry Smith 
28947c6ae99SBarry Smith 
290aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr);
291aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr);
292aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr);
29347c6ae99SBarry Smith 
294aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr);
295aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr);
296aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr);
29747c6ae99SBarry Smith 
29847c6ae99SBarry Smith   /*
299aa219208SBarry Smith    Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
30047c6ae99SBarry Smith    The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
30147c6ae99SBarry Smith    processors). It's effective length is hence 4 times its normal length, this is
30247c6ae99SBarry Smith    why the col_scale is multiplied by the interpolation matrix column sizes.
30347c6ae99SBarry Smith    sol_shift allows each set of 1/4 processors do its own interpolation using ITS
30447c6ae99SBarry Smith    copy of the coarse vector. A bit of a hack but you do better.
30547c6ae99SBarry Smith 
30647c6ae99SBarry Smith    In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
30747c6ae99SBarry Smith    */
30847c6ae99SBarry Smith   ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr);
30947c6ae99SBarry Smith   ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr);
31047c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr);
31147c6ae99SBarry Smith   col_scale = size_f/size_c;
31247c6ae99SBarry Smith   col_shift = Mx*My*(rank_f/size_c);
31347c6ae99SBarry Smith 
31447c6ae99SBarry Smith   ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr);
31547c6ae99SBarry Smith   for (j=j_start; j<j_start+n_f; j++) {
31647c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
31747c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
31847c6ae99SBarry Smith       row    = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
31947c6ae99SBarry Smith 
32047c6ae99SBarry Smith       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
32147c6ae99SBarry Smith       j_c = (j/ratioj);    /* coarse grid node below fine grid node */
32247c6ae99SBarry Smith 
323aa219208SBarry 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\
32447c6ae99SBarry Smith                                           j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c);
325aa219208SBarry 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\
32647c6ae99SBarry Smith                                           i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
32747c6ae99SBarry Smith 
32847c6ae99SBarry Smith       /*
32947c6ae99SBarry Smith        Only include those interpolation points that are truly
33047c6ae99SBarry Smith        nonzero. Note this is very important for final grid lines
33147c6ae99SBarry Smith        in x and y directions; since they have no right/top neighbors
33247c6ae99SBarry Smith        */
33347c6ae99SBarry Smith       nc = 0;
33447c6ae99SBarry Smith       /* one left and below; or we are right on it */
33547c6ae99SBarry Smith       col        = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
33647c6ae99SBarry Smith       cols[nc++] = col_shift + idx_c[col]/dof;
33747c6ae99SBarry Smith       /* one right and below */
33847c6ae99SBarry Smith       if (i_c*ratioi != i) {
33947c6ae99SBarry Smith         cols[nc++] = col_shift + idx_c[col+dof]/dof;
34047c6ae99SBarry Smith       }
34147c6ae99SBarry Smith       /* one left and above */
34247c6ae99SBarry Smith       if (j_c*ratioj != j) {
34347c6ae99SBarry Smith         cols[nc++] = col_shift + idx_c[col+m_ghost_c*dof]/dof;
34447c6ae99SBarry Smith       }
34547c6ae99SBarry Smith       /* one right and above */
3460c82216cSJed Brown       if (i_c*ratioi != i && j_c*ratioj != j) {
34747c6ae99SBarry Smith         cols[nc++] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof;
34847c6ae99SBarry Smith       }
34947c6ae99SBarry Smith       ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr);
35047c6ae99SBarry Smith     }
35147c6ae99SBarry Smith   }
35247c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr);
35347c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr);
35447c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
35547c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr);
35647c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr);
35747c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
35847c6ae99SBarry Smith 
359aa219208SBarry Smith   ierr = DMDAGetCoordinates(daf,&vcoors);CHKERRQ(ierr);
36047c6ae99SBarry Smith   if (vcoors) {
361aa219208SBarry Smith     ierr = DMDAGetGhostedCoordinates(dac,&cvcoors);CHKERRQ(ierr);
362aa219208SBarry Smith     ierr = DMDAVecGetArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr);
363aa219208SBarry Smith     ierr = DMDAVecGetArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr);
36447c6ae99SBarry Smith   }
36547c6ae99SBarry Smith 
36647c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
367b3bd94feSDave May   if (!vcoors) {
368b3bd94feSDave May 
36947c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
37047c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
37147c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
37247c6ae99SBarry Smith         row    = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
37347c6ae99SBarry Smith 
37447c6ae99SBarry Smith         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
37547c6ae99SBarry Smith         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
37647c6ae99SBarry Smith 
37747c6ae99SBarry Smith         /*
37847c6ae99SBarry Smith          Only include those interpolation points that are truly
37947c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
38047c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
38147c6ae99SBarry Smith          */
38247c6ae99SBarry Smith         x  = ((double)(i - i_c*ratioi))/((double)ratioi);
38347c6ae99SBarry Smith         y  = ((double)(j - j_c*ratioj))/((double)ratioj);
384b3bd94feSDave May 
38547c6ae99SBarry Smith         nc = 0;
38647c6ae99SBarry Smith         /* one left and below; or we are right on it */
38747c6ae99SBarry Smith         col      = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
38847c6ae99SBarry Smith         cols[nc] = col_shift + idx_c[col]/dof;
38947c6ae99SBarry Smith         v[nc++]  = x*y - x - y + 1.0;
39047c6ae99SBarry Smith         /* one right and below */
39147c6ae99SBarry Smith         if (i_c*ratioi != i) {
39247c6ae99SBarry Smith           cols[nc] = col_shift + idx_c[col+dof]/dof;
39347c6ae99SBarry Smith           v[nc++]  = -x*y + x;
39447c6ae99SBarry Smith         }
39547c6ae99SBarry Smith         /* one left and above */
39647c6ae99SBarry Smith         if (j_c*ratioj != j) {
39747c6ae99SBarry Smith           cols[nc] = col_shift + idx_c[col+m_ghost_c*dof]/dof;
39847c6ae99SBarry Smith           v[nc++]  = -x*y + y;
39947c6ae99SBarry Smith         }
40047c6ae99SBarry Smith         /* one right and above */
40147c6ae99SBarry Smith         if (j_c*ratioj != j && i_c*ratioi != i) {
40247c6ae99SBarry Smith           cols[nc] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof;
40347c6ae99SBarry Smith           v[nc++]  = x*y;
40447c6ae99SBarry Smith         }
40547c6ae99SBarry Smith         ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
40647c6ae99SBarry Smith       }
40747c6ae99SBarry Smith     }
408b3bd94feSDave May 
409b3bd94feSDave May   } else {
410b3bd94feSDave May     PetscScalar    Ni[4];
411b3bd94feSDave May     PetscScalar    *xi,*eta;
412b3bd94feSDave May     PetscInt       li,nxi,lj,neta;
413b3bd94feSDave May 
414b3bd94feSDave May     /* compute local coordinate arrays */
415b3bd94feSDave May     nxi  = ratioi + 1;
416b3bd94feSDave May     neta = ratioj + 1;
417b3bd94feSDave May     ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr);
418b3bd94feSDave May     ierr = PetscMalloc(sizeof(PetscScalar)*neta,&eta);CHKERRQ(ierr);
419b3bd94feSDave May     for (li=0; li<nxi; li++) {
42052218ef2SJed Brown       xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1));
421b3bd94feSDave May     }
422b3bd94feSDave May     for (lj=0; lj<neta; lj++) {
42352218ef2SJed Brown       eta[lj] = -1.0 + (PetscScalar)lj*(2.0/(PetscScalar)(neta-1));
424b3bd94feSDave May     }
425b3bd94feSDave May 
426b3bd94feSDave May     /* loop over local fine grid nodes setting interpolation for those*/
427b3bd94feSDave May     for (j=j_start; j<j_start+n_f; j++) {
428b3bd94feSDave May       for (i=i_start; i<i_start+m_f; i++) {
429b3bd94feSDave May         /* convert to local "natural" numbering and then to PETSc global numbering */
430b3bd94feSDave May         row    = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
431b3bd94feSDave May 
432b3bd94feSDave May         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
433b3bd94feSDave May         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
434b3bd94feSDave May 
435b3bd94feSDave May         /* remainders */
436b3bd94feSDave May         li = i - ratioi * (i/ratioi);
437b3bd94feSDave May         if (i==mx-1){ li = nxi-1; }
438b3bd94feSDave May         lj = j - ratioj * (j/ratioj);
439b3bd94feSDave May         if (j==my-1){ lj = neta-1; }
440b3bd94feSDave May 
441b3bd94feSDave May         /* corners */
442b3bd94feSDave May         col     = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
443b3bd94feSDave May         cols[0] = col_shift + idx_c[col]/dof; /* left, below */
444b3bd94feSDave May         Ni[0]   = 1.0;
445b3bd94feSDave May         if ( (li==0) || (li==nxi-1) ) {
446b3bd94feSDave May           if ( (lj==0) || (lj==neta-1) ) {
447b3bd94feSDave May             ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr);
448b3bd94feSDave May             continue;
449b3bd94feSDave May           }
450b3bd94feSDave May         }
451b3bd94feSDave May 
452b3bd94feSDave May         /* edges + interior */
453b3bd94feSDave May         /* remainders */
454b3bd94feSDave May         if (i==mx-1){ i_c--; }
455b3bd94feSDave May         if (j==my-1){ j_c--; }
456b3bd94feSDave May 
457b3bd94feSDave May         col     = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
458b3bd94feSDave May         cols[0] = col_shift + idx_c[col]/dof; /* left, below */
459b3bd94feSDave May         cols[1] = col_shift + idx_c[col+dof]/dof; /* right, below */
460b3bd94feSDave May         cols[2] = col_shift + idx_c[col+m_ghost_c*dof]/dof; /* left, above */
461b3bd94feSDave May         cols[3] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; /* right, above */
462b3bd94feSDave May 
463b3bd94feSDave May         Ni[0] = 0.25*(1.0-xi[li])*(1.0-eta[lj]);
464b3bd94feSDave May         Ni[1] = 0.25*(1.0+xi[li])*(1.0-eta[lj]);
465b3bd94feSDave May         Ni[2] = 0.25*(1.0-xi[li])*(1.0+eta[lj]);
466b3bd94feSDave May         Ni[3] = 0.25*(1.0+xi[li])*(1.0+eta[lj]);
467b3bd94feSDave May 
468b3bd94feSDave May         nc = 0;
469b3bd94feSDave May         if( PetscAbsScalar(Ni[0])<1.0e-32) { cols[0]=-1; }
470b3bd94feSDave May         if( PetscAbsScalar(Ni[1])<1.0e-32) { cols[1]=-1; }
471b3bd94feSDave May         if( PetscAbsScalar(Ni[2])<1.0e-32) { cols[2]=-1; }
472b3bd94feSDave May         if( PetscAbsScalar(Ni[3])<1.0e-32) { cols[3]=-1; }
473b3bd94feSDave May 
474b3bd94feSDave May         ierr = MatSetValues(mat,1,&row,4,cols,Ni,INSERT_VALUES);CHKERRQ(ierr);
475b3bd94feSDave May       }
476b3bd94feSDave May     }
477b3bd94feSDave May     ierr = PetscFree(xi);CHKERRQ(ierr);
478b3bd94feSDave May     ierr = PetscFree(eta);CHKERRQ(ierr);
479b3bd94feSDave May   }
48047c6ae99SBarry Smith   if (vcoors) {
481aa219208SBarry Smith     ierr = DMDAVecRestoreArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr);
482aa219208SBarry Smith     ierr = DMDAVecRestoreArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr);
48347c6ae99SBarry Smith   }
48447c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48547c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48647c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
487fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
48847c6ae99SBarry Smith   PetscFunctionReturn(0);
48947c6ae99SBarry Smith }
49047c6ae99SBarry Smith 
49147c6ae99SBarry Smith /*
49247c6ae99SBarry Smith        Contributed by Andrei Draganescu <aidraga@sandia.gov>
49347c6ae99SBarry Smith */
49447c6ae99SBarry Smith #undef __FUNCT__
49594013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_2D_Q0"
49694013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_2D_Q0(DM dac,DM daf,Mat *A)
49747c6ae99SBarry Smith {
49847c6ae99SBarry Smith   PetscErrorCode   ierr;
49947c6ae99SBarry Smith   PetscInt         i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof;
50047c6ae99SBarry Smith   PetscInt         m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz;
50147c6ae99SBarry Smith   PetscInt         row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj;
50247c6ae99SBarry Smith   PetscInt         i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c,col_shift,col_scale;
50347c6ae99SBarry Smith   PetscMPIInt      size_c,size_f,rank_f;
50447c6ae99SBarry Smith   PetscScalar      v[4];
50547c6ae99SBarry Smith   Mat              mat;
5061321219cSEthan Coon   DMDABoundaryType bx,by;
50747c6ae99SBarry Smith 
50847c6ae99SBarry Smith   PetscFunctionBegin;
5091321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr);
5101321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
5111321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x");
5121321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y");
51347c6ae99SBarry Smith   ratioi = mx/Mx;
51447c6ae99SBarry Smith   ratioj = my/My;
51547c6ae99SBarry Smith   if (ratioi*Mx != mx) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in x");
51647c6ae99SBarry Smith   if (ratioj*My != my) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in y");
51747c6ae99SBarry Smith   if (ratioi != 2) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 2");
51847c6ae99SBarry Smith   if (ratioj != 2) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 2");
51947c6ae99SBarry Smith 
520aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr);
521aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr);
522aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr);
52347c6ae99SBarry Smith 
524aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr);
525aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr);
526aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr);
52747c6ae99SBarry Smith 
52847c6ae99SBarry Smith   /*
529aa219208SBarry Smith      Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
53047c6ae99SBarry Smith      The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
53147c6ae99SBarry Smith      processors). It's effective length is hence 4 times its normal length, this is
53247c6ae99SBarry Smith      why the col_scale is multiplied by the interpolation matrix column sizes.
53347c6ae99SBarry Smith      sol_shift allows each set of 1/4 processors do its own interpolation using ITS
53447c6ae99SBarry Smith      copy of the coarse vector. A bit of a hack but you do better.
53547c6ae99SBarry Smith 
53647c6ae99SBarry Smith      In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
53747c6ae99SBarry Smith   */
53847c6ae99SBarry Smith   ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr);
53947c6ae99SBarry Smith   ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr);
54047c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr);
54147c6ae99SBarry Smith   col_scale = size_f/size_c;
54247c6ae99SBarry Smith   col_shift = Mx*My*(rank_f/size_c);
54347c6ae99SBarry Smith 
54447c6ae99SBarry Smith   ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr);
54547c6ae99SBarry Smith   for (j=j_start; j<j_start+n_f; j++) {
54647c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
54747c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
54847c6ae99SBarry Smith       row    = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
54947c6ae99SBarry Smith 
55047c6ae99SBarry Smith       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
55147c6ae99SBarry Smith       j_c = (j/ratioj);    /* coarse grid node below fine grid node */
55247c6ae99SBarry Smith 
553aa219208SBarry 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\
55447c6ae99SBarry Smith     j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c);
555aa219208SBarry 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\
55647c6ae99SBarry Smith     i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
55747c6ae99SBarry Smith 
55847c6ae99SBarry Smith       /*
55947c6ae99SBarry Smith          Only include those interpolation points that are truly
56047c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
56147c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
56247c6ae99SBarry Smith       */
56347c6ae99SBarry Smith       nc = 0;
56447c6ae99SBarry Smith       /* one left and below; or we are right on it */
56547c6ae99SBarry Smith       col        = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
56647c6ae99SBarry Smith       cols[nc++] = col_shift + idx_c[col]/dof;
56747c6ae99SBarry Smith       ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr);
56847c6ae99SBarry Smith     }
56947c6ae99SBarry Smith   }
57047c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr);
57147c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr);
57247c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
57347c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr);
57447c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr);
57547c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
57647c6ae99SBarry Smith 
57747c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
57847c6ae99SBarry Smith   for (j=j_start; j<j_start+n_f; j++) {
57947c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
58047c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
58147c6ae99SBarry Smith       row    = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
58247c6ae99SBarry Smith 
58347c6ae99SBarry Smith       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
58447c6ae99SBarry Smith       j_c = (j/ratioj);    /* coarse grid node below fine grid node */
58547c6ae99SBarry Smith       nc = 0;
58647c6ae99SBarry Smith       /* one left and below; or we are right on it */
58747c6ae99SBarry Smith       col      = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
58847c6ae99SBarry Smith       cols[nc] = col_shift + idx_c[col]/dof;
58947c6ae99SBarry Smith       v[nc++]  = 1.0;
59047c6ae99SBarry Smith 
59147c6ae99SBarry Smith       ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
59247c6ae99SBarry Smith     }
59347c6ae99SBarry Smith   }
59447c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59547c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59647c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
597fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
59847c6ae99SBarry Smith   ierr = PetscLogFlops(13.0*m_f*n_f);CHKERRQ(ierr);
59947c6ae99SBarry Smith   PetscFunctionReturn(0);
60047c6ae99SBarry Smith }
60147c6ae99SBarry Smith 
60247c6ae99SBarry Smith /*
60347c6ae99SBarry Smith        Contributed by Jianming Yang <jianming-yang@uiowa.edu>
60447c6ae99SBarry Smith */
60547c6ae99SBarry Smith #undef __FUNCT__
60694013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_3D_Q0"
60794013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_3D_Q0(DM dac,DM daf,Mat *A)
60847c6ae99SBarry Smith {
60947c6ae99SBarry Smith   PetscErrorCode   ierr;
61047c6ae99SBarry Smith   PetscInt         i,j,l,i_start,j_start,l_start,m_f,n_f,p_f,Mx,My,Mz,*idx_f,dof;
61147c6ae99SBarry Smith   PetscInt         m_ghost,n_ghost,p_ghost,*idx_c,m_ghost_c,n_ghost_c,p_ghost_c,nc,*dnz,*onz;
61247c6ae99SBarry 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;
61347c6ae99SBarry 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;
61447c6ae99SBarry Smith   PetscMPIInt      size_c,size_f,rank_f;
61547c6ae99SBarry Smith   PetscScalar      v[8];
61647c6ae99SBarry Smith   Mat              mat;
6171321219cSEthan Coon   DMDABoundaryType bx,by,bz;
61847c6ae99SBarry Smith 
61947c6ae99SBarry Smith   PetscFunctionBegin;
6201321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr);
6211321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
6221321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x");
6231321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y");
6241321219cSEthan Coon   if (bz == DMDA_BOUNDARY_PERIODIC) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in z");
62547c6ae99SBarry Smith   ratioi = mx/Mx;
62647c6ae99SBarry Smith   ratioj = my/My;
62747c6ae99SBarry Smith   ratiol = mz/Mz;
62847c6ae99SBarry Smith   if (ratioi*Mx != mx) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in x");
62947c6ae99SBarry Smith   if (ratioj*My != my) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in y");
63047c6ae99SBarry Smith   if (ratiol*Mz != mz) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in z");
63147c6ae99SBarry Smith   if (ratioi != 2 && ratioi != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 1 or 2");
63247c6ae99SBarry Smith   if (ratioj != 2 && ratioj != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 1 or 2");
63347c6ae99SBarry Smith   if (ratiol != 2 && ratiol != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in z must be 1 or 2");
63447c6ae99SBarry Smith 
635aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr);
636aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr);
637aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr);
63847c6ae99SBarry Smith 
639aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr);
640aa219208SBarry 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);
641aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr);
64247c6ae99SBarry Smith   /*
643aa219208SBarry Smith      Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
64447c6ae99SBarry Smith      The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
64547c6ae99SBarry Smith      processors). It's effective length is hence 4 times its normal length, this is
64647c6ae99SBarry Smith      why the col_scale is multiplied by the interpolation matrix column sizes.
64747c6ae99SBarry Smith      sol_shift allows each set of 1/4 processors do its own interpolation using ITS
64847c6ae99SBarry Smith      copy of the coarse vector. A bit of a hack but you do better.
64947c6ae99SBarry Smith 
65047c6ae99SBarry Smith      In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
65147c6ae99SBarry Smith   */
65247c6ae99SBarry Smith   ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr);
65347c6ae99SBarry Smith   ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr);
65447c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr);
65547c6ae99SBarry Smith   col_scale = size_f/size_c;
65647c6ae99SBarry Smith   col_shift = Mx*My*Mz*(rank_f/size_c);
65747c6ae99SBarry Smith 
65847c6ae99SBarry Smith   ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f*p_f,col_scale*m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr);
65947c6ae99SBarry Smith   for (l=l_start; l<l_start+p_f; l++) {
66047c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
66147c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
66247c6ae99SBarry Smith 	/* convert to local "natural" numbering and then to PETSc global numbering */
66347c6ae99SBarry Smith 	row    = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
66447c6ae99SBarry Smith 
66547c6ae99SBarry Smith 	i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
66647c6ae99SBarry Smith 	j_c = (j/ratioj);    /* coarse grid node below fine grid node */
66747c6ae99SBarry Smith 	l_c = (l/ratiol);
66847c6ae99SBarry Smith 
669aa219208SBarry 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\
67047c6ae99SBarry Smith     l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c);
671aa219208SBarry 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\
67247c6ae99SBarry Smith     j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c);
673aa219208SBarry 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\
67447c6ae99SBarry Smith     i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
67547c6ae99SBarry Smith 
67647c6ae99SBarry Smith 	/*
67747c6ae99SBarry Smith 	   Only include those interpolation points that are truly
67847c6ae99SBarry Smith 	   nonzero. Note this is very important for final grid lines
67947c6ae99SBarry Smith 	   in x and y directions; since they have no right/top neighbors
68047c6ae99SBarry Smith 	*/
68147c6ae99SBarry Smith 	nc = 0;
68247c6ae99SBarry Smith 	/* one left and below; or we are right on it */
68347c6ae99SBarry 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));
68447c6ae99SBarry Smith 	cols[nc++] = col_shift + idx_c[col]/dof;
68547c6ae99SBarry Smith 	ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr);
68647c6ae99SBarry Smith       }
68747c6ae99SBarry Smith     }
68847c6ae99SBarry Smith   }
68947c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr);
69047c6ae99SBarry 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);
69147c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
69247c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr);
69347c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr);
69447c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
69547c6ae99SBarry Smith 
69647c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
69747c6ae99SBarry Smith   for (l=l_start; l<l_start+p_f; l++) {
69847c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
69947c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
70047c6ae99SBarry Smith 	/* convert to local "natural" numbering and then to PETSc global numbering */
70147c6ae99SBarry Smith 	row    = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
70247c6ae99SBarry Smith 
70347c6ae99SBarry Smith 	i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
70447c6ae99SBarry Smith 	j_c = (j/ratioj);    /* coarse grid node below fine grid node */
70547c6ae99SBarry Smith 	l_c = (l/ratiol);
70647c6ae99SBarry Smith 	nc = 0;
70747c6ae99SBarry Smith 	/* one left and below; or we are right on it */
70847c6ae99SBarry 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));
70947c6ae99SBarry Smith 	cols[nc] = col_shift + idx_c[col]/dof;
71047c6ae99SBarry Smith 	v[nc++]  = 1.0;
71147c6ae99SBarry Smith 
71247c6ae99SBarry Smith 	ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
71347c6ae99SBarry Smith       }
71447c6ae99SBarry Smith     }
71547c6ae99SBarry Smith   }
71647c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71747c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71847c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
719fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
72047c6ae99SBarry Smith   ierr = PetscLogFlops(13.0*m_f*n_f*p_f);CHKERRQ(ierr);
72147c6ae99SBarry Smith   PetscFunctionReturn(0);
72247c6ae99SBarry Smith }
72347c6ae99SBarry Smith 
72447c6ae99SBarry Smith #undef __FUNCT__
72594013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_3D_Q1"
72694013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_3D_Q1(DM dac,DM daf,Mat *A)
72747c6ae99SBarry Smith {
72847c6ae99SBarry Smith   PetscErrorCode   ierr;
72947c6ae99SBarry Smith   PetscInt         i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof,l;
73047c6ae99SBarry Smith   PetscInt         m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,Mz,mz;
73147c6ae99SBarry Smith   PetscInt         row,col,i_start_ghost,j_start_ghost,cols[8],mx,m_c,my,nc,ratioi,ratioj,ratiok;
73247c6ae99SBarry Smith   PetscInt         i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c;
73347c6ae99SBarry Smith   PetscInt         l_start,p_f,l_start_ghost,p_ghost,l_start_c,p_c;
73447c6ae99SBarry Smith   PetscInt         l_start_ghost_c,p_ghost_c,l_c,*dnz,*onz;
73547c6ae99SBarry Smith   PetscScalar      v[8],x,y,z;
73647c6ae99SBarry Smith   Mat              mat;
7371321219cSEthan Coon   DMDABoundaryType bx,by,bz;
738aa219208SBarry Smith   DMDACoor3d       ***coors = 0,***ccoors;
73947c6ae99SBarry Smith   Vec              vcoors,cvcoors;
74047c6ae99SBarry Smith   DM_DA            *ddc = (DM_DA*)dac->data, *ddf = (DM_DA*)daf->data;
74147c6ae99SBarry Smith 
74247c6ae99SBarry Smith   PetscFunctionBegin;
7431321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr);
7441321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
74547c6ae99SBarry Smith   if (mx == Mx) {
74647c6ae99SBarry Smith     ratioi = 1;
7471321219cSEthan Coon   } else if (bx == DMDA_BOUNDARY_PERIODIC) {
74847c6ae99SBarry Smith     ratioi = mx/Mx;
74947c6ae99SBarry 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);
75047c6ae99SBarry Smith   } else {
75147c6ae99SBarry Smith     ratioi = (mx-1)/(Mx-1);
75247c6ae99SBarry 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);
75347c6ae99SBarry Smith   }
75447c6ae99SBarry Smith   if (my == My) {
75547c6ae99SBarry Smith     ratioj = 1;
7561321219cSEthan Coon   } else if (by == DMDA_BOUNDARY_PERIODIC) {
75747c6ae99SBarry Smith     ratioj = my/My;
75847c6ae99SBarry 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);
75947c6ae99SBarry Smith   } else {
76047c6ae99SBarry Smith     ratioj = (my-1)/(My-1);
76147c6ae99SBarry 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);
76247c6ae99SBarry Smith   }
76347c6ae99SBarry Smith   if (mz == Mz) {
76447c6ae99SBarry Smith     ratiok = 1;
7651321219cSEthan Coon   } else if (bz == DMDA_BOUNDARY_PERIODIC) {
76647c6ae99SBarry Smith     ratiok = mz/Mz;
76747c6ae99SBarry 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);
76847c6ae99SBarry Smith   } else {
76947c6ae99SBarry Smith     ratiok = (mz-1)/(Mz-1);
77047c6ae99SBarry 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);
77147c6ae99SBarry Smith   }
77247c6ae99SBarry Smith 
773aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr);
774aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr);
775aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr);
77647c6ae99SBarry Smith 
777aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr);
778aa219208SBarry 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);
779aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr);
78047c6ae99SBarry Smith 
78147c6ae99SBarry Smith   /* create interpolation matrix, determining exact preallocation */
78247c6ae99SBarry Smith   ierr = MatPreallocateInitialize(((PetscObject)dac)->comm,m_f*n_f*p_f,m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr);
78347c6ae99SBarry Smith   /* loop over local fine grid nodes counting interpolating points */
78447c6ae99SBarry Smith   for (l=l_start; l<l_start+p_f; l++) {
78547c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
78647c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
78747c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
78847c6ae99SBarry Smith         row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
78947c6ae99SBarry Smith         i_c = (i/ratioi);
79047c6ae99SBarry Smith         j_c = (j/ratioj);
79147c6ae99SBarry Smith         l_c = (l/ratiok);
792aa219208SBarry 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\
79347c6ae99SBarry Smith                                             l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c);
794aa219208SBarry 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\
79547c6ae99SBarry Smith                                             j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c);
796aa219208SBarry 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\
79747c6ae99SBarry Smith                                             i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
79847c6ae99SBarry Smith 
79947c6ae99SBarry Smith         /*
80047c6ae99SBarry Smith          Only include those interpolation points that are truly
80147c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
80247c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
80347c6ae99SBarry Smith          */
80447c6ae99SBarry Smith         nc       = 0;
80547c6ae99SBarry 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));
80647c6ae99SBarry Smith         cols[nc++] = idx_c[col]/dof;
80747c6ae99SBarry Smith         if (i_c*ratioi != i) {
80847c6ae99SBarry Smith           cols[nc++] = idx_c[col+dof]/dof;
80947c6ae99SBarry Smith         }
81047c6ae99SBarry Smith         if (j_c*ratioj != j) {
81147c6ae99SBarry Smith           cols[nc++] = idx_c[col+m_ghost_c*dof]/dof;
81247c6ae99SBarry Smith         }
81347c6ae99SBarry Smith         if (l_c*ratiok != l) {
81447c6ae99SBarry Smith           cols[nc++] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof;
81547c6ae99SBarry Smith         }
81647c6ae99SBarry Smith         if (j_c*ratioj != j && i_c*ratioi != i) {
81747c6ae99SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c+1)*dof]/dof;
81847c6ae99SBarry Smith         }
81947c6ae99SBarry Smith         if (j_c*ratioj != j && l_c*ratiok != l) {
82047c6ae99SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof;
82147c6ae99SBarry Smith         }
82247c6ae99SBarry Smith         if (i_c*ratioi != i && l_c*ratiok != l) {
82347c6ae99SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof;
82447c6ae99SBarry Smith         }
82547c6ae99SBarry Smith         if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) {
82647c6ae99SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof;
82747c6ae99SBarry Smith         }
82847c6ae99SBarry Smith         ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr);
82947c6ae99SBarry Smith       }
83047c6ae99SBarry Smith     }
83147c6ae99SBarry Smith   }
83247c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr);
83347c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f*n_f*p_f,m_c*n_c*p_c,mx*my*mz,Mx*My*Mz);CHKERRQ(ierr);
83447c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
83547c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr);
83647c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr);
83747c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
83847c6ae99SBarry Smith 
839aa219208SBarry Smith   ierr = DMDAGetCoordinates(daf,&vcoors);CHKERRQ(ierr);
84047c6ae99SBarry Smith   if (vcoors) {
841aa219208SBarry Smith     ierr = DMDAGetGhostedCoordinates(dac,&cvcoors);CHKERRQ(ierr);
842aa219208SBarry Smith     ierr = DMDAVecGetArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr);
843aa219208SBarry Smith     ierr = DMDAVecGetArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr);
84447c6ae99SBarry Smith   }
84547c6ae99SBarry Smith 
84647c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
847b3bd94feSDave May   if (!vcoors) {
848b3bd94feSDave May 
84947c6ae99SBarry Smith     for (l=l_start; l<l_start+p_f; l++) {
85047c6ae99SBarry Smith       for (j=j_start; j<j_start+n_f; j++) {
85147c6ae99SBarry Smith         for (i=i_start; i<i_start+m_f; i++) {
85247c6ae99SBarry Smith           /* convert to local "natural" numbering and then to PETSc global numbering */
85347c6ae99SBarry Smith           row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
85447c6ae99SBarry Smith 
85547c6ae99SBarry Smith           i_c = (i/ratioi);
85647c6ae99SBarry Smith           j_c = (j/ratioj);
85747c6ae99SBarry Smith           l_c = (l/ratiok);
85847c6ae99SBarry Smith 
85947c6ae99SBarry Smith         /*
86047c6ae99SBarry Smith          Only include those interpolation points that are truly
86147c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
86247c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
86347c6ae99SBarry Smith          */
86447c6ae99SBarry Smith           x  = ((double)(i - i_c*ratioi))/((double)ratioi);
86547c6ae99SBarry Smith           y  = ((double)(j - j_c*ratioj))/((double)ratioj);
86647c6ae99SBarry Smith           z  = ((double)(l - l_c*ratiok))/((double)ratiok);
867b3bd94feSDave May 
86847c6ae99SBarry Smith           nc = 0;
86947c6ae99SBarry Smith           /* one left and below; or we are right on it */
87047c6ae99SBarry 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));
87147c6ae99SBarry Smith 
87247c6ae99SBarry Smith           cols[nc] = idx_c[col]/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           if (i_c*ratioi != i) {
87647c6ae99SBarry Smith             cols[nc] = idx_c[col+dof]/dof;
87747c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.));
87847c6ae99SBarry Smith           }
87947c6ae99SBarry Smith 
88047c6ae99SBarry Smith           if (j_c*ratioj != j) {
88147c6ae99SBarry Smith             cols[nc] = idx_c[col+m_ghost_c*dof]/dof;
88247c6ae99SBarry Smith             v[nc++]  = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.));
88347c6ae99SBarry Smith           }
88447c6ae99SBarry Smith 
88547c6ae99SBarry Smith           if (l_c*ratiok != l) {
88647c6ae99SBarry Smith             cols[nc] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof;
88747c6ae99SBarry Smith             v[nc++]  = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.));
88847c6ae99SBarry Smith           }
88947c6ae99SBarry Smith 
89047c6ae99SBarry Smith           if (j_c*ratioj != j && i_c*ratioi != i) {
89147c6ae99SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c+1)*dof]/dof;
89247c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.));
89347c6ae99SBarry Smith           }
89447c6ae99SBarry Smith 
89547c6ae99SBarry Smith           if (j_c*ratioj != j && l_c*ratiok != l) {
89647c6ae99SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof;
89747c6ae99SBarry Smith             v[nc++]  = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.));
89847c6ae99SBarry Smith           }
89947c6ae99SBarry Smith 
90047c6ae99SBarry Smith           if (i_c*ratioi != i && l_c*ratiok != l) {
90147c6ae99SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof;
90247c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.));
90347c6ae99SBarry Smith           }
90447c6ae99SBarry Smith 
90547c6ae99SBarry Smith           if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) {
90647c6ae99SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof;
90747c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.));
90847c6ae99SBarry Smith           }
90947c6ae99SBarry Smith           ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
91047c6ae99SBarry Smith         }
91147c6ae99SBarry Smith       }
91247c6ae99SBarry Smith     }
913b3bd94feSDave May 
914b3bd94feSDave May   } else {
915b3bd94feSDave May     PetscScalar    *xi,*eta,*zeta;
916b3bd94feSDave May     PetscInt       li,nxi,lj,neta,lk,nzeta,n;
917b3bd94feSDave May     PetscScalar    Ni[8];
918b3bd94feSDave May 
919b3bd94feSDave May     /* compute local coordinate arrays */
920b3bd94feSDave May     nxi   = ratioi + 1;
921b3bd94feSDave May     neta  = ratioj + 1;
922b3bd94feSDave May     nzeta = ratiok + 1;
923b3bd94feSDave May     ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr);
924b3bd94feSDave May     ierr = PetscMalloc(sizeof(PetscScalar)*neta,&eta);CHKERRQ(ierr);
925b3bd94feSDave May     ierr = PetscMalloc(sizeof(PetscScalar)*nzeta,&zeta);CHKERRQ(ierr);
926b3bd94feSDave May     for (li=0; li<nxi; li++) {
92752218ef2SJed Brown       xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1));
928b3bd94feSDave May     }
929b3bd94feSDave May     for (lj=0; lj<neta; lj++) {
93052218ef2SJed Brown       eta[lj] = -1.0 + (PetscScalar)lj*(2.0/(PetscScalar)(neta-1));
931b3bd94feSDave May     }
932b3bd94feSDave May     for (lk=0; lk<nzeta; lk++) {
93352218ef2SJed Brown       zeta[lk] = -1.0 + (PetscScalar)lk*(2.0/(PetscScalar)(nzeta-1));
934b3bd94feSDave May     }
935b3bd94feSDave May 
936b3bd94feSDave May     for (l=l_start; l<l_start+p_f; l++) {
937b3bd94feSDave May       for (j=j_start; j<j_start+n_f; j++) {
938b3bd94feSDave May         for (i=i_start; i<i_start+m_f; i++) {
939b3bd94feSDave May           /* convert to local "natural" numbering and then to PETSc global numbering */
940b3bd94feSDave May           row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
941b3bd94feSDave May 
942b3bd94feSDave May           i_c = (i/ratioi);
943b3bd94feSDave May           j_c = (j/ratioj);
944b3bd94feSDave May           l_c = (l/ratiok);
945b3bd94feSDave May 
946b3bd94feSDave May           /* remainders */
947b3bd94feSDave May           li = i - ratioi * (i/ratioi);
948b3bd94feSDave May           if (i==mx-1){ li = nxi-1; }
949b3bd94feSDave May           lj = j - ratioj * (j/ratioj);
950b3bd94feSDave May           if (j==my-1){ lj = neta-1; }
951b3bd94feSDave May           lk = l - ratiok * (l/ratiok);
952b3bd94feSDave May           if (l==mz-1){ lk = nzeta-1; }
953b3bd94feSDave May 
954b3bd94feSDave May           /* corners */
955b3bd94feSDave 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));
956b3bd94feSDave May           cols[0] = idx_c[col]/dof;
957b3bd94feSDave May           Ni[0]   = 1.0;
958b3bd94feSDave May           if ( (li==0) || (li==nxi-1) ) {
959b3bd94feSDave May             if ( (lj==0) || (lj==neta-1) ) {
960b3bd94feSDave May               if ( (lk==0) || (lk==nzeta-1) ) {
961b3bd94feSDave May                 ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr);
962b3bd94feSDave May                 continue;
963b3bd94feSDave May               }
964b3bd94feSDave May             }
965b3bd94feSDave May           }
966b3bd94feSDave May 
967b3bd94feSDave May           /* edges + interior */
968b3bd94feSDave May           /* remainders */
969b3bd94feSDave May           if (i==mx-1){ i_c--; }
970b3bd94feSDave May           if (j==my-1){ j_c--; }
971b3bd94feSDave May           if (l==mz-1){ l_c--; }
972b3bd94feSDave May 
973b3bd94feSDave 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));
974b3bd94feSDave May           cols[0] = idx_c[col]/dof; /* one left and below; or we are right on it */
975b3bd94feSDave May           cols[1] = idx_c[col+dof]/dof; /* one right and below */
976b3bd94feSDave May           cols[2] = idx_c[col+m_ghost_c*dof]/dof;  /* one left and above */
977b3bd94feSDave May           cols[3] = idx_c[col+(m_ghost_c+1)*dof]/dof; /* one right and above */
978b3bd94feSDave May 
979b3bd94feSDave 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 */
980b3bd94feSDave May           cols[5] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; /* one right and below, and front */
981b3bd94feSDave May           cols[6] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof;/* one left and above and front*/
982b3bd94feSDave May           cols[7] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; /* one right and above and front */
983b3bd94feSDave May 
984b3bd94feSDave May           Ni[0] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]);
985b3bd94feSDave May           Ni[1] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]);
986b3bd94feSDave May           Ni[2] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]);
987b3bd94feSDave May           Ni[3] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]);
988b3bd94feSDave May 
989b3bd94feSDave May           Ni[4] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]);
990b3bd94feSDave May           Ni[5] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]);
991b3bd94feSDave May           Ni[6] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]);
992b3bd94feSDave May           Ni[7] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]);
993b3bd94feSDave May 
994b3bd94feSDave May           for (n=0; n<8; n++) {
995b3bd94feSDave May             if( PetscAbsScalar(Ni[n])<1.0e-32) { cols[n]=-1; }
996b3bd94feSDave May           }
997b3bd94feSDave May           ierr = MatSetValues(mat,1,&row,8,cols,Ni,INSERT_VALUES);CHKERRQ(ierr);
998b3bd94feSDave May 
999b3bd94feSDave May         }
1000b3bd94feSDave May       }
1001b3bd94feSDave May     }
1002b3bd94feSDave May     ierr = PetscFree(xi);CHKERRQ(ierr);
1003b3bd94feSDave May     ierr = PetscFree(eta);CHKERRQ(ierr);
1004b3bd94feSDave May     ierr = PetscFree(zeta);CHKERRQ(ierr);
1005b3bd94feSDave May   }
1006b3bd94feSDave May 
100747c6ae99SBarry Smith   if (vcoors) {
1008aa219208SBarry Smith     ierr = DMDAVecRestoreArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr);
1009aa219208SBarry Smith     ierr = DMDAVecRestoreArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr);
101047c6ae99SBarry Smith   }
101147c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101247c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101347c6ae99SBarry Smith 
101447c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
1015fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
101647c6ae99SBarry Smith   PetscFunctionReturn(0);
101747c6ae99SBarry Smith }
101847c6ae99SBarry Smith 
101947c6ae99SBarry Smith #undef __FUNCT__
102094013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA"
10217087cfbeSBarry Smith PetscErrorCode  DMGetInterpolation_DA(DM dac,DM daf,Mat *A,Vec *scale)
102247c6ae99SBarry Smith {
102347c6ae99SBarry Smith   PetscErrorCode   ierr;
102447c6ae99SBarry Smith   PetscInt         dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
10251321219cSEthan Coon   DMDABoundaryType bxc,byc,bzc,bxf,byf,bzf;
1026aa219208SBarry Smith   DMDAStencilType  stc,stf;
102747c6ae99SBarry Smith   DM_DA            *ddc = (DM_DA*)dac->data;
102847c6ae99SBarry Smith 
102947c6ae99SBarry Smith   PetscFunctionBegin;
103047c6ae99SBarry Smith   PetscValidHeaderSpecific(dac,DM_CLASSID,1);
103147c6ae99SBarry Smith   PetscValidHeaderSpecific(daf,DM_CLASSID,2);
103247c6ae99SBarry Smith   PetscValidPointer(A,3);
103347c6ae99SBarry Smith   if (scale) PetscValidPointer(scale,4);
103447c6ae99SBarry Smith 
10351321219cSEthan Coon   ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr);
10361321219cSEthan Coon   ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr);
1037aa219208SBarry Smith   if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr);
1038aa219208SBarry Smith   if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr);
1039aa219208SBarry Smith   if (sc != sf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr);
10401321219cSEthan Coon   if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr);
1041aa219208SBarry Smith   if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr);
104247c6ae99SBarry Smith   if (Mc < 2 && Mf > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction");
104347c6ae99SBarry 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");
104447c6ae99SBarry 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");
104547c6ae99SBarry Smith 
1046aa219208SBarry Smith   if (ddc->interptype == DMDA_Q1){
104747c6ae99SBarry Smith     if (dimc == 1){
104894013140SBarry Smith       ierr = DMGetInterpolation_DA_1D_Q1(dac,daf,A);CHKERRQ(ierr);
104947c6ae99SBarry Smith     } else if (dimc == 2){
105094013140SBarry Smith       ierr = DMGetInterpolation_DA_2D_Q1(dac,daf,A);CHKERRQ(ierr);
105147c6ae99SBarry Smith     } else if (dimc == 3){
105294013140SBarry Smith       ierr = DMGetInterpolation_DA_3D_Q1(dac,daf,A);CHKERRQ(ierr);
105347c6ae99SBarry Smith     } else {
1054aa219208SBarry Smith       SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype);
105547c6ae99SBarry Smith     }
1056aa219208SBarry Smith   } else if (ddc->interptype == DMDA_Q0){
105747c6ae99SBarry Smith     if (dimc == 1){
105894013140SBarry Smith       ierr = DMGetInterpolation_DA_1D_Q0(dac,daf,A);CHKERRQ(ierr);
105947c6ae99SBarry Smith     } else if (dimc == 2){
106094013140SBarry Smith        ierr = DMGetInterpolation_DA_2D_Q0(dac,daf,A);CHKERRQ(ierr);
106147c6ae99SBarry Smith     } else if (dimc == 3){
106294013140SBarry Smith        ierr = DMGetInterpolation_DA_3D_Q0(dac,daf,A);CHKERRQ(ierr);
106347c6ae99SBarry Smith     } else {
1064aa219208SBarry Smith       SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype);
106547c6ae99SBarry Smith     }
106647c6ae99SBarry Smith   }
106747c6ae99SBarry Smith   if (scale) {
106847c6ae99SBarry Smith     ierr = DMGetInterpolationScale((DM)dac,(DM)daf,*A,scale);CHKERRQ(ierr);
106947c6ae99SBarry Smith   }
107047c6ae99SBarry Smith   PetscFunctionReturn(0);
107147c6ae99SBarry Smith }
107247c6ae99SBarry Smith 
107347c6ae99SBarry Smith #undef __FUNCT__
1074*0bb2b966SJungho Lee #define __FUNCT__ "DMGetInjection_DA_1D"
1075*0bb2b966SJungho Lee PetscErrorCode DMGetInjection_DA_1D(DM dac,DM daf,VecScatter *inject)
1076*0bb2b966SJungho Lee {
1077*0bb2b966SJungho Lee     PetscErrorCode   ierr;
1078*0bb2b966SJungho Lee     PetscInt         i,i_start,m_f,Mx,*idx_f,dof;
1079*0bb2b966SJungho Lee     PetscInt         m_ghost,*idx_c,m_ghost_c;
1080*0bb2b966SJungho Lee     PetscInt         row,i_start_ghost,mx,m_c,nc,ratioi;
1081*0bb2b966SJungho Lee     PetscInt         i_start_c,i_start_ghost_c;
1082*0bb2b966SJungho Lee     PetscInt         *cols;
1083*0bb2b966SJungho Lee     DMDABoundaryType bx;
1084*0bb2b966SJungho Lee     Vec              vecf,vecc;
1085*0bb2b966SJungho Lee     IS               isf;
1086*0bb2b966SJungho Lee 
1087*0bb2b966SJungho Lee     PetscFunctionBegin;
1088*0bb2b966SJungho Lee     ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr);
1089*0bb2b966SJungho Lee     ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
1090*0bb2b966SJungho Lee     if (bx == DMDA_BOUNDARY_PERIODIC) {
1091*0bb2b966SJungho Lee         ratioi = mx/Mx;
1092*0bb2b966SJungho 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);
1093*0bb2b966SJungho Lee     } else {
1094*0bb2b966SJungho Lee         ratioi = (mx-1)/(Mx-1);
1095*0bb2b966SJungho 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);
1096*0bb2b966SJungho Lee     }
1097*0bb2b966SJungho Lee 
1098*0bb2b966SJungho Lee     ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr);
1099*0bb2b966SJungho Lee     ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr);
1100*0bb2b966SJungho Lee     ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr);
1101*0bb2b966SJungho Lee 
1102*0bb2b966SJungho Lee     ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr);
1103*0bb2b966SJungho Lee     ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr);
1104*0bb2b966SJungho Lee     ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr);
1105*0bb2b966SJungho Lee 
1106*0bb2b966SJungho Lee 
1107*0bb2b966SJungho Lee     /* loop over local fine grid nodes setting interpolation for those*/
1108*0bb2b966SJungho Lee     nc = 0;
1109*0bb2b966SJungho Lee     ierr = PetscMalloc(m_f*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1110*0bb2b966SJungho Lee 
1111*0bb2b966SJungho Lee 
1112*0bb2b966SJungho Lee     for (i=i_start_c; i<i_start_c+m_c; i++) {
1113*0bb2b966SJungho Lee         PetscInt i_f = i*ratioi;
1114*0bb2b966SJungho Lee 
1115*0bb2b966SJungho Lee            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\
1116*0bb2b966SJungho Lee  i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost);
1117*0bb2b966SJungho Lee             row = idx_f[dof*(i_f-i_start_ghost)];
1118*0bb2b966SJungho Lee             cols[nc++] = row/dof;
1119*0bb2b966SJungho Lee     }
1120*0bb2b966SJungho Lee 
1121*0bb2b966SJungho Lee 
1122*0bb2b966SJungho Lee     ierr = ISCreateBlock(((PetscObject)daf)->comm,dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr);
1123*0bb2b966SJungho Lee     ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr);
1124*0bb2b966SJungho Lee     ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr);
1125*0bb2b966SJungho Lee     ierr = VecScatterCreate(vecf,isf,vecc,PETSC_NULL,inject);CHKERRQ(ierr);
1126*0bb2b966SJungho Lee     ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr);
1127*0bb2b966SJungho Lee     ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr);
1128*0bb2b966SJungho Lee     ierr = ISDestroy(&isf);CHKERRQ(ierr);
1129*0bb2b966SJungho Lee     PetscFunctionReturn(0);
1130*0bb2b966SJungho Lee }
1131*0bb2b966SJungho Lee 
1132*0bb2b966SJungho Lee #undef __FUNCT__
113394013140SBarry Smith #define __FUNCT__ "DMGetInjection_DA_2D"
113494013140SBarry Smith PetscErrorCode DMGetInjection_DA_2D(DM dac,DM daf,VecScatter *inject)
113547c6ae99SBarry Smith {
113647c6ae99SBarry Smith   PetscErrorCode   ierr;
113747c6ae99SBarry Smith   PetscInt         i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof;
113847c6ae99SBarry Smith   PetscInt         m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c;
113947c6ae99SBarry Smith   PetscInt         row,i_start_ghost,j_start_ghost,mx,m_c,my,nc,ratioi,ratioj;
1140076202ddSJed Brown   PetscInt         i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c;
114147c6ae99SBarry Smith   PetscInt         *cols;
11421321219cSEthan Coon   DMDABoundaryType bx,by;
114347c6ae99SBarry Smith   Vec              vecf,vecc;
114447c6ae99SBarry Smith   IS               isf;
114547c6ae99SBarry Smith 
114647c6ae99SBarry Smith   PetscFunctionBegin;
11471321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr);
11481321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
11491321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) {
115047c6ae99SBarry Smith     ratioi = mx/Mx;
115147c6ae99SBarry 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);
115247c6ae99SBarry Smith   } else {
115347c6ae99SBarry Smith     ratioi = (mx-1)/(Mx-1);
115447c6ae99SBarry 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);
115547c6ae99SBarry Smith   }
11561321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC) {
115747c6ae99SBarry Smith     ratioj = my/My;
115847c6ae99SBarry 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);
115947c6ae99SBarry Smith   } else {
116047c6ae99SBarry Smith     ratioj = (my-1)/(My-1);
116147c6ae99SBarry 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);
116247c6ae99SBarry Smith   }
116347c6ae99SBarry Smith 
1164aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr);
1165aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr);
1166aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr);
116747c6ae99SBarry Smith 
1168aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr);
1169aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr);
1170aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr);
117147c6ae99SBarry Smith 
117247c6ae99SBarry Smith 
117347c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
117447c6ae99SBarry Smith   nc = 0;
117547c6ae99SBarry Smith   ierr = PetscMalloc(n_f*m_f*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1176076202ddSJed Brown   for (j=j_start_c; j<j_start_c+n_c; j++) {
1177076202ddSJed Brown     for (i=i_start_c; i<i_start_c+m_c; i++) {
1178076202ddSJed Brown       PetscInt i_f = i*ratioi,j_f = j*ratioj;
1179076202ddSJed 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\
1180076202ddSJed Brown     j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost);
1181076202ddSJed 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\
1182076202ddSJed Brown     i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost);
1183076202ddSJed Brown       row = idx_f[dof*(m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))];
118447c6ae99SBarry Smith       cols[nc++] = row/dof;
118547c6ae99SBarry Smith     }
118647c6ae99SBarry Smith   }
118747c6ae99SBarry Smith 
118847c6ae99SBarry Smith   ierr = ISCreateBlock(((PetscObject)daf)->comm,dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr);
11899a42bb27SBarry Smith   ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr);
11909a42bb27SBarry Smith   ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr);
119147c6ae99SBarry Smith   ierr = VecScatterCreate(vecf,isf,vecc,PETSC_NULL,inject);CHKERRQ(ierr);
11929a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr);
11939a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr);
1194fcfd50ebSBarry Smith   ierr = ISDestroy(&isf);CHKERRQ(ierr);
119547c6ae99SBarry Smith   PetscFunctionReturn(0);
119647c6ae99SBarry Smith }
119747c6ae99SBarry Smith 
119847c6ae99SBarry Smith #undef __FUNCT__
1199b1757049SJed Brown #define __FUNCT__ "DMGetInjection_DA_3D"
1200b1757049SJed Brown PetscErrorCode DMGetInjection_DA_3D(DM dac,DM daf,VecScatter *inject)
1201b1757049SJed Brown {
1202b1757049SJed Brown   PetscErrorCode   ierr;
1203b1757049SJed Brown   PetscInt         i,j,k,i_start,j_start,k_start,m_f,n_f,p_f,Mx,My,Mz;
1204b1757049SJed Brown   PetscInt         m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c;
1205b1757049SJed Brown   PetscInt         i_start_ghost,j_start_ghost,k_start_ghost;
1206b1757049SJed Brown   PetscInt         mx,my,mz,ratioi,ratioj,ratiok;
1207b1757049SJed Brown   PetscInt         i_start_c,j_start_c,k_start_c;
1208b1757049SJed Brown   PetscInt         m_c,n_c,p_c;
1209b1757049SJed Brown   PetscInt         i_start_ghost_c,j_start_ghost_c,k_start_ghost_c;
1210b1757049SJed Brown   PetscInt         row,nc,dof;
1211b1757049SJed Brown   PetscInt         *idx_c,*idx_f;
1212b1757049SJed Brown   PetscInt         *cols;
12131321219cSEthan Coon   DMDABoundaryType bx,by,bz;
1214b1757049SJed Brown   Vec              vecf,vecc;
1215b1757049SJed Brown   IS               isf;
1216b1757049SJed Brown 
1217b1757049SJed Brown   PetscFunctionBegin;
12181321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr);
12191321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
1220b1757049SJed Brown 
12211321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC){
1222b1757049SJed Brown     ratioi = mx/Mx;
1223b1757049SJed 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);
1224b1757049SJed Brown   } else {
1225b1757049SJed Brown     ratioi = (mx-1)/(Mx-1);
1226b1757049SJed 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);
1227b1757049SJed Brown   }
12281321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC){
1229b1757049SJed Brown     ratioj = my/My;
1230b1757049SJed 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);
1231b1757049SJed Brown   } else {
1232b1757049SJed Brown     ratioj = (my-1)/(My-1);
1233b1757049SJed 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);
1234b1757049SJed Brown   }
12351321219cSEthan Coon   if (bz == DMDA_BOUNDARY_PERIODIC){
1236b1757049SJed Brown     ratiok = mz/Mz;
1237b1757049SJed 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);
1238b1757049SJed Brown   } else {
1239b1757049SJed Brown     ratiok = (mz-1)/(Mz-1);
1240b1757049SJed 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);
1241b1757049SJed Brown   }
1242b1757049SJed Brown 
1243b1757049SJed Brown   ierr = DMDAGetCorners(daf,&i_start,&j_start,&k_start,&m_f,&n_f,&p_f);CHKERRQ(ierr);
1244b1757049SJed Brown   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&k_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr);
1245b1757049SJed Brown   ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr);
1246b1757049SJed Brown 
1247b1757049SJed Brown   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&k_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr);
1248b1757049SJed 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);
1249b1757049SJed Brown   ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr);
1250b1757049SJed Brown 
1251b1757049SJed Brown 
1252b1757049SJed Brown   /* loop over local fine grid nodes setting interpolation for those*/
1253b1757049SJed Brown   nc = 0;
1254b1757049SJed Brown   ierr = PetscMalloc(n_f*m_f*p_f*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1255b1757049SJed Brown   for (k=k_start_c; k<k_start_c+p_c; k++) {
1256b1757049SJed Brown     for (j=j_start_c; j<j_start_c+n_c; j++) {
1257b1757049SJed Brown       for (i=i_start_c; i<i_start_c+m_c; i++) {
1258b1757049SJed Brown         PetscInt i_f = i*ratioi,j_f = j*ratioj,k_f = k*ratiok;
1259b1757049SJed 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  "
1260b1757049SJed Brown                                                                           "k_c %D k_f %D fine ghost range [%D,%D]",k,k_f,k_start_ghost,k_start_ghost+p_ghost);
1261b1757049SJed 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  "
1262b1757049SJed Brown                                                                           "j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost);
1263b1757049SJed 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  "
1264b1757049SJed Brown                                                                           "i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost);
1265b1757049SJed 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))];
1266b1757049SJed Brown         cols[nc++] = row/dof;
1267b1757049SJed Brown       }
1268b1757049SJed Brown     }
1269b1757049SJed Brown   }
1270b1757049SJed Brown 
1271b1757049SJed Brown   ierr = ISCreateBlock(((PetscObject)daf)->comm,dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr);
1272b1757049SJed Brown   ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr);
1273b1757049SJed Brown   ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr);
1274b1757049SJed Brown   ierr = VecScatterCreate(vecf,isf,vecc,PETSC_NULL,inject);CHKERRQ(ierr);
1275b1757049SJed Brown   ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr);
1276b1757049SJed Brown   ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr);
1277fcfd50ebSBarry Smith   ierr = ISDestroy(&isf);CHKERRQ(ierr);
1278b1757049SJed Brown   PetscFunctionReturn(0);
1279b1757049SJed Brown }
1280b1757049SJed Brown 
1281b1757049SJed Brown #undef __FUNCT__
128294013140SBarry Smith #define __FUNCT__ "DMGetInjection_DA"
12837087cfbeSBarry Smith PetscErrorCode  DMGetInjection_DA(DM dac,DM daf,VecScatter *inject)
128447c6ae99SBarry Smith {
128547c6ae99SBarry Smith   PetscErrorCode   ierr;
128647c6ae99SBarry Smith   PetscInt         dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
12871321219cSEthan Coon   DMDABoundaryType bxc,byc,bzc,bxf,byf,bzf;
1288aa219208SBarry Smith   DMDAStencilType  stc,stf;
128947c6ae99SBarry Smith 
129047c6ae99SBarry Smith   PetscFunctionBegin;
129147c6ae99SBarry Smith   PetscValidHeaderSpecific(dac,DM_CLASSID,1);
129247c6ae99SBarry Smith   PetscValidHeaderSpecific(daf,DM_CLASSID,2);
129347c6ae99SBarry Smith   PetscValidPointer(inject,3);
129447c6ae99SBarry Smith 
12951321219cSEthan Coon   ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr);
12961321219cSEthan Coon   ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr);
1297aa219208SBarry Smith   if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr);
1298aa219208SBarry Smith   if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr);
1299aa219208SBarry Smith   if (sc != sf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr);
13001321219cSEthan Coon   if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr);
1301aa219208SBarry Smith   if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr);
130247c6ae99SBarry Smith   if (Mc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction");
130347c6ae99SBarry Smith   if (dimc > 1 && Nc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in y direction");
130447c6ae99SBarry Smith   if (dimc > 2 && Pc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in z direction");
130547c6ae99SBarry Smith 
1306*0bb2b966SJungho Lee   if (dimc == 1){
1307*0bb2b966SJungho Lee     ierr = DMGetInjection_DA_1D(dac,daf,inject);CHKERRQ(ierr);
1308*0bb2b966SJungho Lee   } else if (dimc == 2) {
130994013140SBarry Smith     ierr = DMGetInjection_DA_2D(dac,daf,inject);CHKERRQ(ierr);
1310b1757049SJed Brown   } else if (dimc == 3) {
1311b1757049SJed Brown     ierr = DMGetInjection_DA_3D(dac,daf,inject);CHKERRQ(ierr);
131247c6ae99SBarry Smith   }
131347c6ae99SBarry Smith   PetscFunctionReturn(0);
131447c6ae99SBarry Smith }
131547c6ae99SBarry Smith 
131647c6ae99SBarry Smith #undef __FUNCT__
131794013140SBarry Smith #define __FUNCT__ "DMGetAggregates_DA"
13187087cfbeSBarry Smith PetscErrorCode  DMGetAggregates_DA(DM dac,DM daf,Mat *rest)
131947c6ae99SBarry Smith {
132047c6ae99SBarry Smith   PetscErrorCode   ierr;
132147c6ae99SBarry Smith   PetscInt         dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc;
132247c6ae99SBarry Smith   PetscInt         dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
13231321219cSEthan Coon   DMDABoundaryType bxc,byc,bzc,bxf,byf,bzf;
1324aa219208SBarry Smith   DMDAStencilType  stc,stf;
132547c6ae99SBarry Smith   PetscInt         i,j,l;
132647c6ae99SBarry Smith   PetscInt         i_start,j_start,l_start, m_f,n_f,p_f;
132747c6ae99SBarry Smith   PetscInt         i_start_ghost,j_start_ghost,l_start_ghost,m_ghost,n_ghost,p_ghost;
132847c6ae99SBarry Smith   PetscInt         *idx_f;
132947c6ae99SBarry Smith   PetscInt         i_c,j_c,l_c;
133047c6ae99SBarry Smith   PetscInt         i_start_c,j_start_c,l_start_c, m_c,n_c,p_c;
133147c6ae99SBarry Smith   PetscInt         i_start_ghost_c,j_start_ghost_c,l_start_ghost_c,m_ghost_c,n_ghost_c,p_ghost_c;
133247c6ae99SBarry Smith   PetscInt         *idx_c;
133347c6ae99SBarry Smith   PetscInt         d;
133447c6ae99SBarry Smith   PetscInt         a;
133547c6ae99SBarry Smith   PetscInt         max_agg_size;
133647c6ae99SBarry Smith   PetscInt         *fine_nodes;
133747c6ae99SBarry Smith   PetscScalar      *one_vec;
133847c6ae99SBarry Smith   PetscInt         fn_idx;
133947c6ae99SBarry Smith 
134047c6ae99SBarry Smith   PetscFunctionBegin;
134147c6ae99SBarry Smith   PetscValidHeaderSpecific(dac,DM_CLASSID,1);
134247c6ae99SBarry Smith   PetscValidHeaderSpecific(daf,DM_CLASSID,2);
134347c6ae99SBarry Smith   PetscValidPointer(rest,3);
134447c6ae99SBarry Smith 
13451321219cSEthan Coon   ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr);
13461321219cSEthan Coon   ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr);
1347aa219208SBarry Smith   if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr);
1348aa219208SBarry Smith   if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr);
1349aa219208SBarry Smith   if (sc != sf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr);
13501321219cSEthan Coon   if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr);
1351aa219208SBarry Smith   if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr);
135247c6ae99SBarry Smith 
135347c6ae99SBarry 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);
135447c6ae99SBarry 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);
135547c6ae99SBarry 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);
135647c6ae99SBarry Smith 
135747c6ae99SBarry Smith   if (Pc < 0) Pc = 1;
135847c6ae99SBarry Smith   if (Pf < 0) Pf = 1;
135947c6ae99SBarry Smith   if (Nc < 0) Nc = 1;
136047c6ae99SBarry Smith   if (Nf < 0) Nf = 1;
136147c6ae99SBarry Smith 
1362aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr);
1363aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr);
1364aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr);
136547c6ae99SBarry Smith 
1366aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr);
1367aa219208SBarry 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);
1368aa219208SBarry Smith   ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr);
136947c6ae99SBarry Smith 
137047c6ae99SBarry Smith   /*
137147c6ae99SBarry Smith      Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios
137247c6ae99SBarry Smith      for dimension 1 and 2 respectively.
137347c6ae99SBarry Smith      Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1)
137447c6ae99SBarry Smith      and r_y*j and r_y*(j+1) will be grouped into the same coarse grid agregate.
137547c6ae99SBarry Smith      Each specific dof on the fine grid is mapped to one dof on the coarse grid.
137647c6ae99SBarry Smith   */
137747c6ae99SBarry Smith 
137847c6ae99SBarry Smith   max_agg_size = (Mf/Mc+1)*(Nf/Nc+1)*(Pf/Pc+1);
137947c6ae99SBarry Smith 
138047c6ae99SBarry Smith   /* create the matrix that will contain the restriction operator */
138147c6ae99SBarry Smith   ierr = MatCreateMPIAIJ( ((PetscObject)daf)->comm, m_c*n_c*p_c*dofc, m_f*n_f*p_f*doff, Mc*Nc*Pc*dofc, Mf*Nf*Pf*doff,
138247c6ae99SBarry Smith 			  max_agg_size, PETSC_NULL, max_agg_size, PETSC_NULL, rest);CHKERRQ(ierr);
138347c6ae99SBarry Smith 
138447c6ae99SBarry Smith   /* store nodes in the fine grid here */
138547c6ae99SBarry Smith   ierr = PetscMalloc2(max_agg_size,PetscScalar, &one_vec,max_agg_size,PetscInt, &fine_nodes);CHKERRQ(ierr);
138647c6ae99SBarry Smith   for(i=0; i<max_agg_size; i++) one_vec[i] = 1.0;
138747c6ae99SBarry Smith 
138847c6ae99SBarry Smith   /* loop over all coarse nodes */
138947c6ae99SBarry Smith   for (l_c=l_start_c; l_c<l_start_c+p_c; l_c++) {
139047c6ae99SBarry Smith     for (j_c=j_start_c; j_c<j_start_c+n_c; j_c++) {
139147c6ae99SBarry Smith       for (i_c=i_start_c; i_c<i_start_c+m_c; i_c++) {
139247c6ae99SBarry Smith 	for(d=0; d<dofc; d++) {
139347c6ae99SBarry Smith 	  /* convert to local "natural" numbering and then to PETSc global numbering */
139447c6ae99SBarry 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;
139547c6ae99SBarry Smith 
139647c6ae99SBarry Smith 	  fn_idx = 0;
139747c6ae99SBarry Smith 	  /* Corresponding fine points are all points (i_f, j_f, l_f) such that
139847c6ae99SBarry Smith 	     i_c*Mf/Mc <= i_f < (i_c+1)*Mf/Mc
139947c6ae99SBarry Smith 	     (same for other dimensions)
140047c6ae99SBarry Smith 	  */
140147c6ae99SBarry Smith 	  for (l=l_c*Pf/Pc; l<PetscMin((l_c+1)*Pf/Pc,Pf); l++) {
140247c6ae99SBarry Smith 	    for (j=j_c*Nf/Nc; j<PetscMin((j_c+1)*Nf/Nc,Nf); j++) {
140347c6ae99SBarry Smith 	      for (i=i_c*Mf/Mc; i<PetscMin((i_c+1)*Mf/Mc,Mf); i++) {
140447c6ae99SBarry 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;
140547c6ae99SBarry Smith 		fn_idx++;
140647c6ae99SBarry Smith 	      }
140747c6ae99SBarry Smith 	    }
140847c6ae99SBarry Smith 	  }
140947c6ae99SBarry Smith 	  /* add all these points to one aggregate */
141047c6ae99SBarry Smith 	  ierr = MatSetValues(*rest, 1, &a, fn_idx, fine_nodes, one_vec, INSERT_VALUES);CHKERRQ(ierr);
141147c6ae99SBarry Smith 	}
141247c6ae99SBarry Smith       }
141347c6ae99SBarry Smith     }
141447c6ae99SBarry Smith   }
141547c6ae99SBarry Smith   ierr = PetscFree2(one_vec,fine_nodes);CHKERRQ(ierr);
141647c6ae99SBarry Smith   ierr = MatAssemblyBegin(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
141747c6ae99SBarry Smith   ierr = MatAssemblyEnd(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
141847c6ae99SBarry Smith   PetscFunctionReturn(0);
141947c6ae99SBarry Smith }
1420