xref: /petsc/src/dm/impls/da/dainterp.c (revision 3a5864876ba7ac693351187792d7d4a4cbdb72a9)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
3aa219208SBarry Smith   Code for interpolating between grids represented by DMDAs
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
6d52bd9f3SBarry Smith /*
7d52bd9f3SBarry Smith       For linear elements there are two branches of code to compute the interpolation. They should compute the same results but may not. The "new version" does
8d52bd9f3SBarry Smith    not work for periodic domains, the old does. Change NEWVERSION to 1 to compile in the new version. Eventually when we are sure the two produce identical results
92adb9dcfSBarry Smith    we will remove/merge the new version. Based on current tests, these both produce the same results. We are leaving NEWVERSION for now in the code since some
102adb9dcfSBarry Smith    consider it cleaner, but old version is turned on since it handles periodic case.
11d52bd9f3SBarry Smith */
129314d695SBarry Smith #define NEWVERSION 0
1385afcc9aSBarry Smith 
14af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith /*@
17e727c939SJed Brown     DMCreateInterpolationScale - Forms L = R*1/diag(R*1) - L.*v is like a coarse grid average of the
1847c6ae99SBarry Smith       nearby fine grid points.
1947c6ae99SBarry Smith 
2047c6ae99SBarry Smith   Input Parameters:
2147c6ae99SBarry Smith +      dac - DM that defines a coarse mesh
2247c6ae99SBarry Smith .      daf - DM that defines a fine mesh
2347c6ae99SBarry Smith -      mat - the restriction (or interpolation operator) from fine to coarse
2447c6ae99SBarry Smith 
2547c6ae99SBarry Smith   Output Parameter:
2647c6ae99SBarry Smith .    scale - the scaled vector
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith   Level: developer
2947c6ae99SBarry Smith 
30e727c939SJed Brown .seealso: DMCreateInterpolation()
3147c6ae99SBarry Smith 
3247c6ae99SBarry Smith @*/
33e727c939SJed Brown PetscErrorCode  DMCreateInterpolationScale(DM dac,DM daf,Mat mat,Vec *scale)
3447c6ae99SBarry Smith {
3547c6ae99SBarry Smith   PetscErrorCode ierr;
3647c6ae99SBarry Smith   Vec            fine;
3747c6ae99SBarry Smith   PetscScalar    one = 1.0;
3847c6ae99SBarry Smith 
3947c6ae99SBarry Smith   PetscFunctionBegin;
4047c6ae99SBarry Smith   ierr = DMCreateGlobalVector(daf,&fine);CHKERRQ(ierr);
4147c6ae99SBarry Smith   ierr = DMCreateGlobalVector(dac,scale);CHKERRQ(ierr);
4247c6ae99SBarry Smith   ierr = VecSet(fine,one);CHKERRQ(ierr);
4347c6ae99SBarry Smith   ierr = MatRestrict(mat,fine,*scale);CHKERRQ(ierr);
44fcfd50ebSBarry Smith   ierr = VecDestroy(&fine);CHKERRQ(ierr);
4547c6ae99SBarry Smith   ierr = VecReciprocal(*scale);CHKERRQ(ierr);
4647c6ae99SBarry Smith   PetscFunctionReturn(0);
4747c6ae99SBarry Smith }
4847c6ae99SBarry Smith 
49e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_1D_Q1(DM dac,DM daf,Mat *A)
5047c6ae99SBarry Smith {
5147c6ae99SBarry Smith   PetscErrorCode         ierr;
528ea3bf28SBarry Smith   PetscInt               i,i_start,m_f,Mx;
538ea3bf28SBarry Smith   const PetscInt         *idx_f,*idx_c;
548ea3bf28SBarry Smith   PetscInt               m_ghost,m_ghost_c;
5547c6ae99SBarry Smith   PetscInt               row,col,i_start_ghost,mx,m_c,nc,ratio;
5647c6ae99SBarry Smith   PetscInt               i_c,i_start_c,i_start_ghost_c,cols[2],dof;
5785afcc9aSBarry Smith   PetscScalar            v[2],x;
5847c6ae99SBarry Smith   Mat                    mat;
59bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
60e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
61e3fbd1f4SBarry Smith 
6247c6ae99SBarry Smith 
6347c6ae99SBarry Smith   PetscFunctionBegin;
641321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr);
651321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
66bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
6747c6ae99SBarry Smith     ratio = mx/Mx;
6847c6ae99SBarry 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);
6947c6ae99SBarry Smith   } else {
7047c6ae99SBarry Smith     ratio = (mx-1)/(Mx-1);
7147c6ae99SBarry 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);
7247c6ae99SBarry Smith   }
7347c6ae99SBarry Smith 
74aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr);
75aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr);
7645b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(daf,&ltog_f);CHKERRQ(ierr);
7745b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
7847c6ae99SBarry Smith 
79aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr);
80aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr);
8145b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(dac,&ltog_c);CHKERRQ(ierr);
8245b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
8347c6ae99SBarry Smith 
8447c6ae99SBarry Smith   /* create interpolation matrix */
85ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)dac),&mat);CHKERRQ(ierr);
8647c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr);
8747c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
880298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,2,NULL);CHKERRQ(ierr);
890298fd71SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,2,NULL,1,NULL);CHKERRQ(ierr);
9047c6ae99SBarry Smith 
9147c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
9285afcc9aSBarry Smith   if (!NEWVERSION) {
93b3bd94feSDave May 
9447c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
9547c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
96e3fbd1f4SBarry Smith       row = idx_f[i-i_start_ghost];
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith       i_c = (i/ratio);    /* coarse grid node to left of fine grid node */
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        */
1076712e2f1SBarry Smith       x  = ((PetscReal)(i - i_c*ratio))/((PetscReal)ratio);
10847c6ae99SBarry Smith       nc = 0;
10947c6ae99SBarry Smith       /* one left and below; or we are right on it */
110e3fbd1f4SBarry Smith       col      = (i_c-i_start_ghost_c);
111e3fbd1f4SBarry Smith       cols[nc] = idx_c[col];
11247c6ae99SBarry Smith       v[nc++]  = -x + 1.0;
11347c6ae99SBarry Smith       /* one right? */
11447c6ae99SBarry Smith       if (i_c*ratio != i) {
115e3fbd1f4SBarry Smith         cols[nc] = idx_c[col+1];
11647c6ae99SBarry Smith         v[nc++]  = x;
11747c6ae99SBarry Smith       }
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;
128854ce69bSBarry Smith     ierr = PetscMalloc1(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 */
135e3fbd1f4SBarry Smith       row = idx_f[(i-i_start_ghost)];
136b3bd94feSDave May 
137b3bd94feSDave May       i_c = (i/ratio);    /* coarse grid node to left of fine grid node */
138b3bd94feSDave May       if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\
139b3bd94feSDave May                                           i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
140b3bd94feSDave May 
141b3bd94feSDave May       /* remainders */
142b3bd94feSDave May       li = i - ratio * (i/ratio);
1438865f1eaSKarl Rupp       if (i==mx-1) li = nxi-1;
144b3bd94feSDave May 
145b3bd94feSDave May       /* corners */
146e3fbd1f4SBarry Smith       col     = (i_c-i_start_ghost_c);
147e3fbd1f4SBarry Smith       cols[0] = idx_c[col];
148b3bd94feSDave May       Ni[0]   = 1.0;
149b3bd94feSDave May       if ((li==0) || (li==nxi-1)) {
150b3bd94feSDave May         ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr);
151b3bd94feSDave May         continue;
152b3bd94feSDave May       }
153b3bd94feSDave May 
154b3bd94feSDave May       /* edges + interior */
155b3bd94feSDave May       /* remainders */
1568865f1eaSKarl Rupp       if (i==mx-1) i_c--;
157b3bd94feSDave May 
158e3fbd1f4SBarry Smith       col     = (i_c-i_start_ghost_c);
159e3fbd1f4SBarry Smith       cols[0] = idx_c[col]; /* one left and below; or we are right on it */
160e3fbd1f4SBarry Smith       cols[1] = idx_c[col+1];
161b3bd94feSDave May 
162b3bd94feSDave May       Ni[0] = 0.5*(1.0-xi[li]);
163b3bd94feSDave May       Ni[1] = 0.5*(1.0+xi[li]);
164b3bd94feSDave May       for (n=0; n<2; n++) {
1658865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[n])<1.0e-32) cols[n]=-1;
166b3bd94feSDave May       }
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   }
17145b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
17245b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
17347c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17447c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17547c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
176fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
17747c6ae99SBarry Smith   PetscFunctionReturn(0);
17847c6ae99SBarry Smith }
17947c6ae99SBarry Smith 
180e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_1D_Q0(DM dac,DM daf,Mat *A)
18147c6ae99SBarry Smith {
18247c6ae99SBarry Smith   PetscErrorCode         ierr;
1838ea3bf28SBarry Smith   PetscInt               i,i_start,m_f,Mx;
1848ea3bf28SBarry Smith   const PetscInt         *idx_f,*idx_c;
185e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
1868ea3bf28SBarry Smith   PetscInt               m_ghost,m_ghost_c;
18747c6ae99SBarry Smith   PetscInt               row,col,i_start_ghost,mx,m_c,nc,ratio;
18847c6ae99SBarry Smith   PetscInt               i_c,i_start_c,i_start_ghost_c,cols[2],dof;
18947c6ae99SBarry Smith   PetscScalar            v[2],x;
19047c6ae99SBarry Smith   Mat                    mat;
191bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
19247c6ae99SBarry Smith 
19347c6ae99SBarry Smith   PetscFunctionBegin;
1941321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr);
1951321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
196bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
197*3a586487SStefano Zampini     if (!Mx) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %D must be positive",Mx);
19847c6ae99SBarry Smith     ratio = mx/Mx;
19947c6ae99SBarry Smith     if (ratio*Mx != mx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: mx/Mx  must be integer: mx %D Mx %D",mx,Mx);
20047c6ae99SBarry Smith   } else {
201*3a586487SStefano Zampini     if (Mx < 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %D must be greater than 1",Mx);
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);
20845b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(daf,&ltog_f);CHKERRQ(ierr);
20945b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
21047c6ae99SBarry Smith 
211aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr);
212aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr);
21345b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(dac,&ltog_c);CHKERRQ(ierr);
21445b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
21547c6ae99SBarry Smith 
21647c6ae99SBarry Smith   /* create interpolation matrix */
217ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)dac),&mat);CHKERRQ(ierr);
21847c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr);
21947c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
2200298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,2,NULL);CHKERRQ(ierr);
2210298fd71SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,2,NULL,0,NULL);CHKERRQ(ierr);
22247c6ae99SBarry Smith 
22347c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
22447c6ae99SBarry Smith   for (i=i_start; i<i_start+m_f; i++) {
22547c6ae99SBarry Smith     /* convert to local "natural" numbering and then to PETSc global numbering */
226e3fbd1f4SBarry Smith     row = idx_f[(i-i_start_ghost)];
22747c6ae99SBarry Smith 
22847c6ae99SBarry Smith     i_c = (i/ratio);    /* coarse grid node to left of fine grid node */
22947c6ae99SBarry Smith 
23047c6ae99SBarry Smith     /*
23147c6ae99SBarry Smith          Only include those interpolation points that are truly
23247c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
23347c6ae99SBarry Smith          in x direction; since they have no right neighbor
23447c6ae99SBarry Smith     */
2356712e2f1SBarry Smith     x  = ((PetscReal)(i - i_c*ratio))/((PetscReal)ratio);
23647c6ae99SBarry Smith     nc = 0;
23747c6ae99SBarry Smith     /* one left and below; or we are right on it */
238e3fbd1f4SBarry Smith     col      = (i_c-i_start_ghost_c);
239e3fbd1f4SBarry Smith     cols[nc] = idx_c[col];
24047c6ae99SBarry Smith     v[nc++]  = -x + 1.0;
24147c6ae99SBarry Smith     /* one right? */
24247c6ae99SBarry Smith     if (i_c*ratio != i) {
243e3fbd1f4SBarry Smith       cols[nc] = idx_c[col+1];
24447c6ae99SBarry Smith       v[nc++]  = x;
24547c6ae99SBarry Smith     }
24647c6ae99SBarry Smith     ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
24747c6ae99SBarry Smith   }
24845b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
24945b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
25047c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
25147c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
25247c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
253fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
25447c6ae99SBarry Smith   ierr = PetscLogFlops(5.0*m_f);CHKERRQ(ierr);
25547c6ae99SBarry Smith   PetscFunctionReturn(0);
25647c6ae99SBarry Smith }
25747c6ae99SBarry Smith 
258e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_2D_Q1(DM dac,DM daf,Mat *A)
25947c6ae99SBarry Smith {
26047c6ae99SBarry Smith   PetscErrorCode         ierr;
2618ea3bf28SBarry Smith   PetscInt               i,j,i_start,j_start,m_f,n_f,Mx,My,dof;
2628ea3bf28SBarry Smith   const PetscInt         *idx_c,*idx_f;
263e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
2648ea3bf28SBarry Smith   PetscInt               m_ghost,n_ghost,m_ghost_c,n_ghost_c,*dnz,*onz;
26547c6ae99SBarry Smith   PetscInt               row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj;
26647c6ae99SBarry Smith   PetscInt               i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c,col_shift,col_scale;
26747c6ae99SBarry Smith   PetscMPIInt            size_c,size_f,rank_f;
26847c6ae99SBarry Smith   PetscScalar            v[4],x,y;
26947c6ae99SBarry Smith   Mat                    mat;
270bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
27147c6ae99SBarry Smith 
27247c6ae99SBarry Smith   PetscFunctionBegin;
2731321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr);
2741321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
275bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
276*3a586487SStefano Zampini     if (!Mx) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %D must be positive",Mx);
27747c6ae99SBarry Smith     ratioi = mx/Mx;
27847c6ae99SBarry 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);
27947c6ae99SBarry Smith   } else {
280*3a586487SStefano Zampini     if (Mx < 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %D must be greater than 1",Mx);
28147c6ae99SBarry Smith     ratioi = (mx-1)/(Mx-1);
28247c6ae99SBarry 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);
28347c6ae99SBarry Smith   }
284bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) {
285*3a586487SStefano Zampini     if (!My) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of y coarse grid points %D must be positive",My);
28647c6ae99SBarry Smith     ratioj = my/My;
28747c6ae99SBarry 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);
28847c6ae99SBarry Smith   } else {
289*3a586487SStefano Zampini     if (My < 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of y coarse grid points %D must be greater than 1",My);
29047c6ae99SBarry Smith     ratioj = (my-1)/(My-1);
29147c6ae99SBarry 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);
29247c6ae99SBarry Smith   }
29347c6ae99SBarry Smith 
294aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr);
295aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr);
29645b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(daf,&ltog_f);CHKERRQ(ierr);
29745b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
29847c6ae99SBarry Smith 
299aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr);
300aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr);
30145b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(dac,&ltog_c);CHKERRQ(ierr);
30245b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
30347c6ae99SBarry Smith 
30447c6ae99SBarry Smith   /*
305aa219208SBarry Smith    Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
30647c6ae99SBarry Smith    The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
30747c6ae99SBarry Smith    processors). It's effective length is hence 4 times its normal length, this is
30847c6ae99SBarry Smith    why the col_scale is multiplied by the interpolation matrix column sizes.
30947c6ae99SBarry Smith    sol_shift allows each set of 1/4 processors do its own interpolation using ITS
31047c6ae99SBarry Smith    copy of the coarse vector. A bit of a hack but you do better.
31147c6ae99SBarry Smith 
31247c6ae99SBarry Smith    In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
31347c6ae99SBarry Smith    */
314ce94432eSBarry Smith   ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c);CHKERRQ(ierr);
315ce94432eSBarry Smith   ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f);CHKERRQ(ierr);
316ce94432eSBarry Smith   ierr      = MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f);CHKERRQ(ierr);
31747c6ae99SBarry Smith   col_scale = size_f/size_c;
31847c6ae99SBarry Smith   col_shift = Mx*My*(rank_f/size_c);
31947c6ae99SBarry Smith 
320ce94432eSBarry Smith   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)daf),m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr);
32147c6ae99SBarry Smith   for (j=j_start; j<j_start+n_f; j++) {
32247c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
32347c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
324e3fbd1f4SBarry Smith       row = idx_f[(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
32547c6ae99SBarry Smith 
32647c6ae99SBarry Smith       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
32747c6ae99SBarry Smith       j_c = (j/ratioj);    /* coarse grid node below fine grid node */
32847c6ae99SBarry Smith 
329aa219208SBarry 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\
33047c6ae99SBarry Smith                                           j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c);
331aa219208SBarry 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\
33247c6ae99SBarry Smith                                           i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
33347c6ae99SBarry Smith 
33447c6ae99SBarry Smith       /*
33547c6ae99SBarry Smith        Only include those interpolation points that are truly
33647c6ae99SBarry Smith        nonzero. Note this is very important for final grid lines
33747c6ae99SBarry Smith        in x and y directions; since they have no right/top neighbors
33847c6ae99SBarry Smith        */
33947c6ae99SBarry Smith       nc = 0;
34047c6ae99SBarry Smith       /* one left and below; or we are right on it */
341e3fbd1f4SBarry Smith       col        = (m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
342e3fbd1f4SBarry Smith       cols[nc++] = col_shift + idx_c[col];
34347c6ae99SBarry Smith       /* one right and below */
344e3fbd1f4SBarry Smith       if (i_c*ratioi != i) cols[nc++] = col_shift + idx_c[col+1];
34547c6ae99SBarry Smith       /* one left and above */
346e3fbd1f4SBarry Smith       if (j_c*ratioj != j) cols[nc++] = col_shift + idx_c[col+m_ghost_c];
34747c6ae99SBarry Smith       /* one right and above */
348e3fbd1f4SBarry Smith       if (i_c*ratioi != i && j_c*ratioj != j) cols[nc++] = col_shift + idx_c[col+(m_ghost_c+1)];
34947c6ae99SBarry Smith       ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr);
35047c6ae99SBarry Smith     }
35147c6ae99SBarry Smith   }
352ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)daf),&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 
35947c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
36085afcc9aSBarry Smith   if (!NEWVERSION) {
361b3bd94feSDave May 
36247c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
36347c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
36447c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
365e3fbd1f4SBarry Smith         row = idx_f[(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
36847c6ae99SBarry Smith         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
36947c6ae99SBarry Smith 
37047c6ae99SBarry Smith         /*
37147c6ae99SBarry Smith          Only include those interpolation points that are truly
37247c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
37347c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
37447c6ae99SBarry Smith          */
3756712e2f1SBarry Smith         x = ((PetscReal)(i - i_c*ratioi))/((PetscReal)ratioi);
3766712e2f1SBarry Smith         y = ((PetscReal)(j - j_c*ratioj))/((PetscReal)ratioj);
377b3bd94feSDave May 
37847c6ae99SBarry Smith         nc = 0;
37947c6ae99SBarry Smith         /* one left and below; or we are right on it */
380e3fbd1f4SBarry Smith         col      = (m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
381e3fbd1f4SBarry Smith         cols[nc] = col_shift + idx_c[col];
38247c6ae99SBarry Smith         v[nc++]  = x*y - x - y + 1.0;
38347c6ae99SBarry Smith         /* one right and below */
38447c6ae99SBarry Smith         if (i_c*ratioi != i) {
385e3fbd1f4SBarry Smith           cols[nc] = col_shift + idx_c[col+1];
38647c6ae99SBarry Smith           v[nc++]  = -x*y + x;
38747c6ae99SBarry Smith         }
38847c6ae99SBarry Smith         /* one left and above */
38947c6ae99SBarry Smith         if (j_c*ratioj != j) {
390e3fbd1f4SBarry Smith           cols[nc] = col_shift + idx_c[col+m_ghost_c];
39147c6ae99SBarry Smith           v[nc++]  = -x*y + y;
39247c6ae99SBarry Smith         }
39347c6ae99SBarry Smith         /* one right and above */
39447c6ae99SBarry Smith         if (j_c*ratioj != j && i_c*ratioi != i) {
395e3fbd1f4SBarry Smith           cols[nc] = col_shift + idx_c[col+(m_ghost_c+1)];
39647c6ae99SBarry Smith           v[nc++]  = x*y;
39747c6ae99SBarry Smith         }
39847c6ae99SBarry Smith         ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
39947c6ae99SBarry Smith       }
40047c6ae99SBarry Smith     }
401b3bd94feSDave May 
402b3bd94feSDave May   } else {
403b3bd94feSDave May     PetscScalar Ni[4];
404b3bd94feSDave May     PetscScalar *xi,*eta;
405b3bd94feSDave May     PetscInt    li,nxi,lj,neta;
406b3bd94feSDave May 
407b3bd94feSDave May     /* compute local coordinate arrays */
408b3bd94feSDave May     nxi  = ratioi + 1;
409b3bd94feSDave May     neta = ratioj + 1;
410854ce69bSBarry Smith     ierr = PetscMalloc1(nxi,&xi);CHKERRQ(ierr);
411854ce69bSBarry Smith     ierr = PetscMalloc1(neta,&eta);CHKERRQ(ierr);
412b3bd94feSDave May     for (li=0; li<nxi; li++) {
41352218ef2SJed Brown       xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1));
414b3bd94feSDave May     }
415b3bd94feSDave May     for (lj=0; lj<neta; lj++) {
41652218ef2SJed Brown       eta[lj] = -1.0 + (PetscScalar)lj*(2.0/(PetscScalar)(neta-1));
417b3bd94feSDave May     }
418b3bd94feSDave May 
419b3bd94feSDave May     /* loop over local fine grid nodes setting interpolation for those*/
420b3bd94feSDave May     for (j=j_start; j<j_start+n_f; j++) {
421b3bd94feSDave May       for (i=i_start; i<i_start+m_f; i++) {
422b3bd94feSDave May         /* convert to local "natural" numbering and then to PETSc global numbering */
423e3fbd1f4SBarry Smith         row = idx_f[(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
424b3bd94feSDave May 
425b3bd94feSDave May         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
426b3bd94feSDave May         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
427b3bd94feSDave May 
428b3bd94feSDave May         /* remainders */
429b3bd94feSDave May         li = i - ratioi * (i/ratioi);
4308865f1eaSKarl Rupp         if (i==mx-1) li = nxi-1;
431b3bd94feSDave May         lj = j - ratioj * (j/ratioj);
4328865f1eaSKarl Rupp         if (j==my-1) lj = neta-1;
433b3bd94feSDave May 
434b3bd94feSDave May         /* corners */
435e3fbd1f4SBarry Smith         col     = (m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
436e3fbd1f4SBarry Smith         cols[0] = col_shift + idx_c[col]; /* left, below */
437b3bd94feSDave May         Ni[0]   = 1.0;
438b3bd94feSDave May         if ((li==0) || (li==nxi-1)) {
439b3bd94feSDave May           if ((lj==0) || (lj==neta-1)) {
440b3bd94feSDave May             ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr);
441b3bd94feSDave May             continue;
442b3bd94feSDave May           }
443b3bd94feSDave May         }
444b3bd94feSDave May 
445b3bd94feSDave May         /* edges + interior */
446b3bd94feSDave May         /* remainders */
4478865f1eaSKarl Rupp         if (i==mx-1) i_c--;
4488865f1eaSKarl Rupp         if (j==my-1) j_c--;
449b3bd94feSDave May 
450e3fbd1f4SBarry Smith         col     = (m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
451e3fbd1f4SBarry Smith         cols[0] = col_shift + idx_c[col]; /* left, below */
452e3fbd1f4SBarry Smith         cols[1] = col_shift + idx_c[col+1]; /* right, below */
453e3fbd1f4SBarry Smith         cols[2] = col_shift + idx_c[col+m_ghost_c]; /* left, above */
454e3fbd1f4SBarry Smith         cols[3] = col_shift + idx_c[col+(m_ghost_c+1)]; /* right, above */
455b3bd94feSDave May 
456b3bd94feSDave May         Ni[0] = 0.25*(1.0-xi[li])*(1.0-eta[lj]);
457b3bd94feSDave May         Ni[1] = 0.25*(1.0+xi[li])*(1.0-eta[lj]);
458b3bd94feSDave May         Ni[2] = 0.25*(1.0-xi[li])*(1.0+eta[lj]);
459b3bd94feSDave May         Ni[3] = 0.25*(1.0+xi[li])*(1.0+eta[lj]);
460b3bd94feSDave May 
461b3bd94feSDave May         nc = 0;
4628865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[0])<1.0e-32) cols[0]=-1;
4638865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[1])<1.0e-32) cols[1]=-1;
4648865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[2])<1.0e-32) cols[2]=-1;
4658865f1eaSKarl Rupp         if (PetscAbsScalar(Ni[3])<1.0e-32) cols[3]=-1;
466b3bd94feSDave May 
467b3bd94feSDave May         ierr = MatSetValues(mat,1,&row,4,cols,Ni,INSERT_VALUES);CHKERRQ(ierr);
468b3bd94feSDave May       }
469b3bd94feSDave May     }
470b3bd94feSDave May     ierr = PetscFree(xi);CHKERRQ(ierr);
471b3bd94feSDave May     ierr = PetscFree(eta);CHKERRQ(ierr);
472b3bd94feSDave May   }
47345b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
47445b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
47547c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
47647c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
47747c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
478fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
47947c6ae99SBarry Smith   PetscFunctionReturn(0);
48047c6ae99SBarry Smith }
48147c6ae99SBarry Smith 
48247c6ae99SBarry Smith /*
48347c6ae99SBarry Smith        Contributed by Andrei Draganescu <aidraga@sandia.gov>
48447c6ae99SBarry Smith */
485e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_2D_Q0(DM dac,DM daf,Mat *A)
48647c6ae99SBarry Smith {
48747c6ae99SBarry Smith   PetscErrorCode         ierr;
4888ea3bf28SBarry Smith   PetscInt               i,j,i_start,j_start,m_f,n_f,Mx,My,dof;
4898ea3bf28SBarry Smith   const PetscInt         *idx_c,*idx_f;
490e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
4918ea3bf28SBarry Smith   PetscInt               m_ghost,n_ghost,m_ghost_c,n_ghost_c,*dnz,*onz;
49247c6ae99SBarry Smith   PetscInt               row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj;
49347c6ae99SBarry 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;
49447c6ae99SBarry Smith   PetscMPIInt            size_c,size_f,rank_f;
49547c6ae99SBarry Smith   PetscScalar            v[4];
49647c6ae99SBarry Smith   Mat                    mat;
497bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith   PetscFunctionBegin;
5001321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr);
5011321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
502*3a586487SStefano Zampini   if (!Mx) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %D must be positive",Mx);
503*3a586487SStefano Zampini   if (!My) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of y coarse grid points %D must be positive",My);
50447c6ae99SBarry Smith   ratioi = mx/Mx;
50547c6ae99SBarry Smith   ratioj = my/My;
506ce94432eSBarry Smith   if (ratioi*Mx != mx) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in x");
507ce94432eSBarry Smith   if (ratioj*My != my) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in y");
508ce94432eSBarry Smith   if (ratioi != 2) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 2");
509ce94432eSBarry Smith   if (ratioj != 2) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 2");
51047c6ae99SBarry Smith 
511aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr);
512aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr);
51345b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(daf,&ltog_f);CHKERRQ(ierr);
51445b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
51547c6ae99SBarry Smith 
516aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr);
517aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr);
51845b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(dac,&ltog_c);CHKERRQ(ierr);
51945b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
52047c6ae99SBarry Smith 
52147c6ae99SBarry Smith   /*
522aa219208SBarry Smith      Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
52347c6ae99SBarry Smith      The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
52447c6ae99SBarry Smith      processors). It's effective length is hence 4 times its normal length, this is
52547c6ae99SBarry Smith      why the col_scale is multiplied by the interpolation matrix column sizes.
52647c6ae99SBarry Smith      sol_shift allows each set of 1/4 processors do its own interpolation using ITS
52747c6ae99SBarry Smith      copy of the coarse vector. A bit of a hack but you do better.
52847c6ae99SBarry Smith 
52947c6ae99SBarry Smith      In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
53047c6ae99SBarry Smith   */
531ce94432eSBarry Smith   ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c);CHKERRQ(ierr);
532ce94432eSBarry Smith   ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f);CHKERRQ(ierr);
533ce94432eSBarry Smith   ierr      = MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f);CHKERRQ(ierr);
53447c6ae99SBarry Smith   col_scale = size_f/size_c;
53547c6ae99SBarry Smith   col_shift = Mx*My*(rank_f/size_c);
53647c6ae99SBarry Smith 
537ce94432eSBarry Smith   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)daf),m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr);
53847c6ae99SBarry Smith   for (j=j_start; j<j_start+n_f; j++) {
53947c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
54047c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
541e3fbd1f4SBarry Smith       row = idx_f[(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
54247c6ae99SBarry Smith 
54347c6ae99SBarry Smith       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
54447c6ae99SBarry Smith       j_c = (j/ratioj);    /* coarse grid node below fine grid node */
54547c6ae99SBarry Smith 
546aa219208SBarry 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\
54747c6ae99SBarry Smith     j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c);
548aa219208SBarry 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\
54947c6ae99SBarry Smith     i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
55047c6ae99SBarry Smith 
55147c6ae99SBarry Smith       /*
55247c6ae99SBarry Smith          Only include those interpolation points that are truly
55347c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
55447c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
55547c6ae99SBarry Smith       */
55647c6ae99SBarry Smith       nc = 0;
55747c6ae99SBarry Smith       /* one left and below; or we are right on it */
558e3fbd1f4SBarry Smith       col        = (m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
559e3fbd1f4SBarry Smith       cols[nc++] = col_shift + idx_c[col];
56047c6ae99SBarry Smith       ierr       = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr);
56147c6ae99SBarry Smith     }
56247c6ae99SBarry Smith   }
563ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)daf),&mat);CHKERRQ(ierr);
56447c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr);
56547c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
56647c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr);
56747c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr);
56847c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
56947c6ae99SBarry Smith 
57047c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
57147c6ae99SBarry Smith   for (j=j_start; j<j_start+n_f; j++) {
57247c6ae99SBarry Smith     for (i=i_start; i<i_start+m_f; i++) {
57347c6ae99SBarry Smith       /* convert to local "natural" numbering and then to PETSc global numbering */
574e3fbd1f4SBarry Smith       row = idx_f[(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
57547c6ae99SBarry Smith 
57647c6ae99SBarry Smith       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
57747c6ae99SBarry Smith       j_c = (j/ratioj);    /* coarse grid node below fine grid node */
57847c6ae99SBarry Smith       nc  = 0;
57947c6ae99SBarry Smith       /* one left and below; or we are right on it */
580e3fbd1f4SBarry Smith       col      = (m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
581e3fbd1f4SBarry Smith       cols[nc] = col_shift + idx_c[col];
58247c6ae99SBarry Smith       v[nc++]  = 1.0;
58347c6ae99SBarry Smith 
58447c6ae99SBarry Smith       ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
58547c6ae99SBarry Smith     }
58647c6ae99SBarry Smith   }
58745b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
58845b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
58947c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59047c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59147c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
592fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
59347c6ae99SBarry Smith   ierr = PetscLogFlops(13.0*m_f*n_f);CHKERRQ(ierr);
59447c6ae99SBarry Smith   PetscFunctionReturn(0);
59547c6ae99SBarry Smith }
59647c6ae99SBarry Smith 
59747c6ae99SBarry Smith /*
59847c6ae99SBarry Smith        Contributed by Jianming Yang <jianming-yang@uiowa.edu>
59947c6ae99SBarry Smith */
600e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_3D_Q0(DM dac,DM daf,Mat *A)
60147c6ae99SBarry Smith {
60247c6ae99SBarry Smith   PetscErrorCode         ierr;
6038ea3bf28SBarry Smith   PetscInt               i,j,l,i_start,j_start,l_start,m_f,n_f,p_f,Mx,My,Mz,dof;
6048ea3bf28SBarry Smith   const PetscInt         *idx_c,*idx_f;
605e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
6068ea3bf28SBarry Smith   PetscInt               m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c,nc,*dnz,*onz;
60747c6ae99SBarry 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;
60847c6ae99SBarry 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;
60947c6ae99SBarry Smith   PetscMPIInt            size_c,size_f,rank_f;
61047c6ae99SBarry Smith   PetscScalar            v[8];
61147c6ae99SBarry Smith   Mat                    mat;
612bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
61347c6ae99SBarry Smith 
61447c6ae99SBarry Smith   PetscFunctionBegin;
6151321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr);
616*3a586487SStefano Zampini   if (!Mx) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %D must be positive",Mx);
617*3a586487SStefano Zampini   if (!My) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of y coarse grid points %D must be positive",My);
618*3a586487SStefano Zampini   if (!Mz) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of z coarse grid points %D must be positive",Mz);
6191321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
62047c6ae99SBarry Smith   ratioi = mx/Mx;
62147c6ae99SBarry Smith   ratioj = my/My;
62247c6ae99SBarry Smith   ratiol = mz/Mz;
623ce94432eSBarry Smith   if (ratioi*Mx != mx) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in x");
624ce94432eSBarry Smith   if (ratioj*My != my) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in y");
625ce94432eSBarry Smith   if (ratiol*Mz != mz) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in z");
626ce94432eSBarry Smith   if (ratioi != 2 && ratioi != 1) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 1 or 2");
627ce94432eSBarry Smith   if (ratioj != 2 && ratioj != 1) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 1 or 2");
628ce94432eSBarry Smith   if (ratiol != 2 && ratiol != 1) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in z must be 1 or 2");
62947c6ae99SBarry Smith 
630aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr);
631aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr);
63245b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(daf,&ltog_f);CHKERRQ(ierr);
63345b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
63447c6ae99SBarry Smith 
635aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr);
636aa219208SBarry 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);
63745b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(dac,&ltog_c);CHKERRQ(ierr);
63845b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
639e3fbd1f4SBarry Smith 
64047c6ae99SBarry Smith   /*
641aa219208SBarry Smith      Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA.
64247c6ae99SBarry Smith      The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the
64347c6ae99SBarry Smith      processors). It's effective length is hence 4 times its normal length, this is
64447c6ae99SBarry Smith      why the col_scale is multiplied by the interpolation matrix column sizes.
64547c6ae99SBarry Smith      sol_shift allows each set of 1/4 processors do its own interpolation using ITS
64647c6ae99SBarry Smith      copy of the coarse vector. A bit of a hack but you do better.
64747c6ae99SBarry Smith 
64847c6ae99SBarry Smith      In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
64947c6ae99SBarry Smith   */
650ce94432eSBarry Smith   ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c);CHKERRQ(ierr);
651ce94432eSBarry Smith   ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f);CHKERRQ(ierr);
652ce94432eSBarry Smith   ierr      = MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f);CHKERRQ(ierr);
65347c6ae99SBarry Smith   col_scale = size_f/size_c;
65447c6ae99SBarry Smith   col_shift = Mx*My*Mz*(rank_f/size_c);
65547c6ae99SBarry Smith 
656ce94432eSBarry Smith   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)daf),m_f*n_f*p_f,col_scale*m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr);
65747c6ae99SBarry Smith   for (l=l_start; l<l_start+p_f; l++) {
65847c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
65947c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
66047c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
661e3fbd1f4SBarry Smith         row = idx_f[(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
66247c6ae99SBarry Smith 
66347c6ae99SBarry Smith         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
66447c6ae99SBarry Smith         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
66547c6ae99SBarry Smith         l_c = (l/ratiol);
66647c6ae99SBarry Smith 
667aa219208SBarry 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\
66847c6ae99SBarry Smith     l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c);
669aa219208SBarry 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\
67047c6ae99SBarry Smith     j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c);
671aa219208SBarry 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\
67247c6ae99SBarry Smith     i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
67347c6ae99SBarry Smith 
67447c6ae99SBarry Smith         /*
67547c6ae99SBarry Smith            Only include those interpolation points that are truly
67647c6ae99SBarry Smith            nonzero. Note this is very important for final grid lines
67747c6ae99SBarry Smith            in x and y directions; since they have no right/top neighbors
67847c6ae99SBarry Smith         */
67947c6ae99SBarry Smith         nc = 0;
68047c6ae99SBarry Smith         /* one left and below; or we are right on it */
681e3fbd1f4SBarry Smith         col        = (m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c) + m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
682e3fbd1f4SBarry Smith         cols[nc++] = col_shift + idx_c[col];
68347c6ae99SBarry Smith         ierr       = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr);
68447c6ae99SBarry Smith       }
68547c6ae99SBarry Smith     }
68647c6ae99SBarry Smith   }
687ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)daf),&mat);CHKERRQ(ierr);
68847c6ae99SBarry 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);
68947c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
69047c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr);
69147c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr);
69247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
69347c6ae99SBarry Smith 
69447c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
69547c6ae99SBarry Smith   for (l=l_start; l<l_start+p_f; l++) {
69647c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
69747c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
69847c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
699e3fbd1f4SBarry Smith         row = idx_f[(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
70047c6ae99SBarry Smith 
70147c6ae99SBarry Smith         i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
70247c6ae99SBarry Smith         j_c = (j/ratioj);    /* coarse grid node below fine grid node */
70347c6ae99SBarry Smith         l_c = (l/ratiol);
70447c6ae99SBarry Smith         nc  = 0;
70547c6ae99SBarry Smith         /* one left and below; or we are right on it */
706e3fbd1f4SBarry Smith         col      = (m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c) + m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
707e3fbd1f4SBarry Smith         cols[nc] = col_shift + idx_c[col];
70847c6ae99SBarry Smith         v[nc++]  = 1.0;
70947c6ae99SBarry Smith 
71047c6ae99SBarry Smith         ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
71147c6ae99SBarry Smith       }
71247c6ae99SBarry Smith     }
71347c6ae99SBarry Smith   }
71445b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
71545b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
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 
724e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_3D_Q1(DM dac,DM daf,Mat *A)
72547c6ae99SBarry Smith {
72647c6ae99SBarry Smith   PetscErrorCode         ierr;
7278ea3bf28SBarry Smith   PetscInt               i,j,i_start,j_start,m_f,n_f,Mx,My,dof,l;
7288ea3bf28SBarry Smith   const PetscInt         *idx_c,*idx_f;
729e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
7308ea3bf28SBarry Smith   PetscInt               m_ghost,n_ghost,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;
737bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
73847c6ae99SBarry Smith 
73947c6ae99SBarry Smith   PetscFunctionBegin;
7401321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr);
7411321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
74247c6ae99SBarry Smith   if (mx == Mx) {
74347c6ae99SBarry Smith     ratioi = 1;
744bff4a2f0SMatthew G. Knepley   } else if (bx == DM_BOUNDARY_PERIODIC) {
745*3a586487SStefano Zampini     if (!Mx) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %D must be positive",Mx);
74647c6ae99SBarry Smith     ratioi = mx/Mx;
74747c6ae99SBarry 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);
74847c6ae99SBarry Smith   } else {
749*3a586487SStefano Zampini     if (Mx < 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of x coarse grid points %D must be greater than 1",Mx);
75047c6ae99SBarry Smith     ratioi = (mx-1)/(Mx-1);
75147c6ae99SBarry 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);
75247c6ae99SBarry Smith   }
75347c6ae99SBarry Smith   if (my == My) {
75447c6ae99SBarry Smith     ratioj = 1;
755bff4a2f0SMatthew G. Knepley   } else if (by == DM_BOUNDARY_PERIODIC) {
756*3a586487SStefano Zampini     if (!My) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of y coarse grid points %D must be positive",My);
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 {
760*3a586487SStefano Zampini     if (My < 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of y coarse grid points %D must be greater than 1",My);
76147c6ae99SBarry Smith     ratioj = (my-1)/(My-1);
76247c6ae99SBarry 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);
76347c6ae99SBarry Smith   }
76447c6ae99SBarry Smith   if (mz == Mz) {
76547c6ae99SBarry Smith     ratiok = 1;
766bff4a2f0SMatthew G. Knepley   } else if (bz == DM_BOUNDARY_PERIODIC) {
767*3a586487SStefano Zampini     if (!Mz) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of z coarse grid points %D must be positive",Mz);
76847c6ae99SBarry Smith     ratiok = mz/Mz;
76947c6ae99SBarry 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);
77047c6ae99SBarry Smith   } else {
771*3a586487SStefano Zampini     if (Mz < 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of z coarse grid points %D must be greater than 1",Mz);
77247c6ae99SBarry Smith     ratiok = (mz-1)/(Mz-1);
77347c6ae99SBarry 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);
77447c6ae99SBarry Smith   }
77547c6ae99SBarry Smith 
776aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr);
777aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr);
77845b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(daf,&ltog_f);CHKERRQ(ierr);
77945b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
78047c6ae99SBarry Smith 
781aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr);
782aa219208SBarry 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);
78345b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(dac,&ltog_c);CHKERRQ(ierr);
78445b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
78547c6ae99SBarry Smith 
78647c6ae99SBarry Smith   /* create interpolation matrix, determining exact preallocation */
787ce94432eSBarry Smith   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)dac),m_f*n_f*p_f,m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr);
78847c6ae99SBarry Smith   /* loop over local fine grid nodes counting interpolating points */
78947c6ae99SBarry Smith   for (l=l_start; l<l_start+p_f; l++) {
79047c6ae99SBarry Smith     for (j=j_start; j<j_start+n_f; j++) {
79147c6ae99SBarry Smith       for (i=i_start; i<i_start+m_f; i++) {
79247c6ae99SBarry Smith         /* convert to local "natural" numbering and then to PETSc global numbering */
793e3fbd1f4SBarry Smith         row = idx_f[(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
79447c6ae99SBarry Smith         i_c = (i/ratioi);
79547c6ae99SBarry Smith         j_c = (j/ratioj);
79647c6ae99SBarry Smith         l_c = (l/ratiok);
797aa219208SBarry 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\
79847c6ae99SBarry Smith                                             l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c);
799aa219208SBarry 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\
80047c6ae99SBarry Smith                                             j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c);
801aa219208SBarry 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\
80247c6ae99SBarry Smith                                             i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c);
80347c6ae99SBarry Smith 
80447c6ae99SBarry Smith         /*
80547c6ae99SBarry Smith          Only include those interpolation points that are truly
80647c6ae99SBarry Smith          nonzero. Note this is very important for final grid lines
80747c6ae99SBarry Smith          in x and y directions; since they have no right/top neighbors
80847c6ae99SBarry Smith          */
80947c6ae99SBarry Smith         nc         = 0;
810e3fbd1f4SBarry Smith         col        = (m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c) + m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
811e3fbd1f4SBarry Smith         cols[nc++] = idx_c[col];
81247c6ae99SBarry Smith         if (i_c*ratioi != i) {
813e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+1];
81447c6ae99SBarry Smith         }
81547c6ae99SBarry Smith         if (j_c*ratioj != j) {
816e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+m_ghost_c];
81747c6ae99SBarry Smith         }
81847c6ae99SBarry Smith         if (l_c*ratiok != l) {
819e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+m_ghost_c*n_ghost_c];
82047c6ae99SBarry Smith         }
82147c6ae99SBarry Smith         if (j_c*ratioj != j && i_c*ratioi != i) {
822e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c+1)];
82347c6ae99SBarry Smith         }
82447c6ae99SBarry Smith         if (j_c*ratioj != j && l_c*ratiok != l) {
825e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)];
82647c6ae99SBarry Smith         }
82747c6ae99SBarry Smith         if (i_c*ratioi != i && l_c*ratiok != l) {
828e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+1)];
82947c6ae99SBarry Smith         }
83047c6ae99SBarry Smith         if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) {
831e3fbd1f4SBarry Smith           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)];
83247c6ae99SBarry Smith         }
83347c6ae99SBarry Smith         ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr);
83447c6ae99SBarry Smith       }
83547c6ae99SBarry Smith     }
83647c6ae99SBarry Smith   }
837ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)dac),&mat);CHKERRQ(ierr);
83847c6ae99SBarry Smith   ierr = MatSetSizes(mat,m_f*n_f*p_f,m_c*n_c*p_c,mx*my*mz,Mx*My*Mz);CHKERRQ(ierr);
83947c6ae99SBarry Smith   ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr);
84047c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr);
84147c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr);
84247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
84347c6ae99SBarry Smith 
84447c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
8452adb9dcfSBarry Smith   if (!NEWVERSION) {
846b3bd94feSDave May 
84747c6ae99SBarry Smith     for (l=l_start; l<l_start+p_f; l++) {
84847c6ae99SBarry Smith       for (j=j_start; j<j_start+n_f; j++) {
84947c6ae99SBarry Smith         for (i=i_start; i<i_start+m_f; i++) {
85047c6ae99SBarry Smith           /* convert to local "natural" numbering and then to PETSc global numbering */
851e3fbd1f4SBarry Smith           row = idx_f[(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
85247c6ae99SBarry Smith 
85347c6ae99SBarry Smith           i_c = (i/ratioi);
85447c6ae99SBarry Smith           j_c = (j/ratioj);
85547c6ae99SBarry Smith           l_c = (l/ratiok);
85647c6ae99SBarry Smith 
85747c6ae99SBarry Smith           /*
85847c6ae99SBarry Smith            Only include those interpolation points that are truly
85947c6ae99SBarry Smith            nonzero. Note this is very important for final grid lines
86047c6ae99SBarry Smith            in x and y directions; since they have no right/top neighbors
86147c6ae99SBarry Smith            */
8626712e2f1SBarry Smith           x = ((PetscReal)(i - i_c*ratioi))/((PetscReal)ratioi);
8636712e2f1SBarry Smith           y = ((PetscReal)(j - j_c*ratioj))/((PetscReal)ratioj);
8646712e2f1SBarry Smith           z = ((PetscReal)(l - l_c*ratiok))/((PetscReal)ratiok);
865b3bd94feSDave May 
86647c6ae99SBarry Smith           nc = 0;
86747c6ae99SBarry Smith           /* one left and below; or we are right on it */
868e3fbd1f4SBarry Smith           col = (m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c)+m_ghost_c*(j_c-j_start_ghost_c)+(i_c-i_start_ghost_c));
86947c6ae99SBarry Smith 
870e3fbd1f4SBarry Smith           cols[nc] = idx_c[col];
87147c6ae99SBarry Smith           v[nc++]  = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.));
87247c6ae99SBarry Smith 
87347c6ae99SBarry Smith           if (i_c*ratioi != i) {
874e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+1];
87547c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.));
87647c6ae99SBarry Smith           }
87747c6ae99SBarry Smith 
87847c6ae99SBarry Smith           if (j_c*ratioj != j) {
879e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+m_ghost_c];
88047c6ae99SBarry Smith             v[nc++]  = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.));
88147c6ae99SBarry Smith           }
88247c6ae99SBarry Smith 
88347c6ae99SBarry Smith           if (l_c*ratiok != l) {
884e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+m_ghost_c*n_ghost_c];
88547c6ae99SBarry Smith             v[nc++]  = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.));
88647c6ae99SBarry Smith           }
88747c6ae99SBarry Smith 
88847c6ae99SBarry Smith           if (j_c*ratioj != j && i_c*ratioi != i) {
889e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c+1)];
89047c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.));
89147c6ae99SBarry Smith           }
89247c6ae99SBarry Smith 
89347c6ae99SBarry Smith           if (j_c*ratioj != j && l_c*ratiok != l) {
894e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)];
89547c6ae99SBarry Smith             v[nc++]  = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.));
89647c6ae99SBarry Smith           }
89747c6ae99SBarry Smith 
89847c6ae99SBarry Smith           if (i_c*ratioi != i && l_c*ratiok != l) {
899e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+1)];
90047c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.));
90147c6ae99SBarry Smith           }
90247c6ae99SBarry Smith 
90347c6ae99SBarry Smith           if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) {
904e3fbd1f4SBarry Smith             cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)];
90547c6ae99SBarry Smith             v[nc++]  = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.));
90647c6ae99SBarry Smith           }
90747c6ae99SBarry Smith           ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr);
90847c6ae99SBarry Smith         }
90947c6ae99SBarry Smith       }
91047c6ae99SBarry Smith     }
911b3bd94feSDave May 
912b3bd94feSDave May   } else {
913b3bd94feSDave May     PetscScalar *xi,*eta,*zeta;
914b3bd94feSDave May     PetscInt    li,nxi,lj,neta,lk,nzeta,n;
915b3bd94feSDave May     PetscScalar Ni[8];
916b3bd94feSDave May 
917b3bd94feSDave May     /* compute local coordinate arrays */
918b3bd94feSDave May     nxi   = ratioi + 1;
919b3bd94feSDave May     neta  = ratioj + 1;
920b3bd94feSDave May     nzeta = ratiok + 1;
921854ce69bSBarry Smith     ierr  = PetscMalloc1(nxi,&xi);CHKERRQ(ierr);
922854ce69bSBarry Smith     ierr  = PetscMalloc1(neta,&eta);CHKERRQ(ierr);
923854ce69bSBarry Smith     ierr  = PetscMalloc1(nzeta,&zeta);CHKERRQ(ierr);
9248865f1eaSKarl Rupp     for (li=0; li<nxi; li++) xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1));
9258865f1eaSKarl Rupp     for (lj=0; lj<neta; lj++) eta[lj] = -1.0 + (PetscScalar)lj*(2.0/(PetscScalar)(neta-1));
9268865f1eaSKarl Rupp     for (lk=0; lk<nzeta; lk++) zeta[lk] = -1.0 + (PetscScalar)lk*(2.0/(PetscScalar)(nzeta-1));
927b3bd94feSDave May 
928b3bd94feSDave May     for (l=l_start; l<l_start+p_f; l++) {
929b3bd94feSDave May       for (j=j_start; j<j_start+n_f; j++) {
930b3bd94feSDave May         for (i=i_start; i<i_start+m_f; i++) {
931b3bd94feSDave May           /* convert to local "natural" numbering and then to PETSc global numbering */
932e3fbd1f4SBarry Smith           row = idx_f[(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))];
933b3bd94feSDave May 
934b3bd94feSDave May           i_c = (i/ratioi);
935b3bd94feSDave May           j_c = (j/ratioj);
936b3bd94feSDave May           l_c = (l/ratiok);
937b3bd94feSDave May 
938b3bd94feSDave May           /* remainders */
939b3bd94feSDave May           li = i - ratioi * (i/ratioi);
9408865f1eaSKarl Rupp           if (i==mx-1) li = nxi-1;
941b3bd94feSDave May           lj = j - ratioj * (j/ratioj);
9428865f1eaSKarl Rupp           if (j==my-1) lj = neta-1;
943b3bd94feSDave May           lk = l - ratiok * (l/ratiok);
9448865f1eaSKarl Rupp           if (l==mz-1) lk = nzeta-1;
945b3bd94feSDave May 
946b3bd94feSDave May           /* corners */
947e3fbd1f4SBarry Smith           col     = (m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c)+m_ghost_c*(j_c-j_start_ghost_c)+(i_c-i_start_ghost_c));
948e3fbd1f4SBarry Smith           cols[0] = idx_c[col];
949b3bd94feSDave May           Ni[0]   = 1.0;
950b3bd94feSDave May           if ((li==0) || (li==nxi-1)) {
951b3bd94feSDave May             if ((lj==0) || (lj==neta-1)) {
952b3bd94feSDave May               if ((lk==0) || (lk==nzeta-1)) {
953b3bd94feSDave May                 ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr);
954b3bd94feSDave May                 continue;
955b3bd94feSDave May               }
956b3bd94feSDave May             }
957b3bd94feSDave May           }
958b3bd94feSDave May 
959b3bd94feSDave May           /* edges + interior */
960b3bd94feSDave May           /* remainders */
9618865f1eaSKarl Rupp           if (i==mx-1) i_c--;
9628865f1eaSKarl Rupp           if (j==my-1) j_c--;
9638865f1eaSKarl Rupp           if (l==mz-1) l_c--;
964b3bd94feSDave May 
965e3fbd1f4SBarry Smith           col     = (m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c) + m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
966e3fbd1f4SBarry Smith           cols[0] = idx_c[col]; /* one left and below; or we are right on it */
967e3fbd1f4SBarry Smith           cols[1] = idx_c[col+1]; /* one right and below */
968e3fbd1f4SBarry Smith           cols[2] = idx_c[col+m_ghost_c];  /* one left and above */
969e3fbd1f4SBarry Smith           cols[3] = idx_c[col+(m_ghost_c+1)]; /* one right and above */
970b3bd94feSDave May 
971e3fbd1f4SBarry Smith           cols[4] = idx_c[col+m_ghost_c*n_ghost_c]; /* one left and below and front; or we are right on it */
972e3fbd1f4SBarry Smith           cols[5] = idx_c[col+(m_ghost_c*n_ghost_c+1)]; /* one right and below, and front */
973e3fbd1f4SBarry Smith           cols[6] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)]; /* one left and above and front*/
974e3fbd1f4SBarry Smith           cols[7] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)]; /* one right and above and front */
975b3bd94feSDave May 
976b3bd94feSDave May           Ni[0] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]);
977b3bd94feSDave May           Ni[1] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]);
978b3bd94feSDave May           Ni[2] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]);
979b3bd94feSDave May           Ni[3] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]);
980b3bd94feSDave May 
981b3bd94feSDave May           Ni[4] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]);
982b3bd94feSDave May           Ni[5] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]);
983b3bd94feSDave May           Ni[6] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]);
984b3bd94feSDave May           Ni[7] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]);
985b3bd94feSDave May 
986b3bd94feSDave May           for (n=0; n<8; n++) {
9878865f1eaSKarl Rupp             if (PetscAbsScalar(Ni[n])<1.0e-32) cols[n]=-1;
988b3bd94feSDave May           }
989b3bd94feSDave May           ierr = MatSetValues(mat,1,&row,8,cols,Ni,INSERT_VALUES);CHKERRQ(ierr);
990b3bd94feSDave May 
991b3bd94feSDave May         }
992b3bd94feSDave May       }
993b3bd94feSDave May     }
994b3bd94feSDave May     ierr = PetscFree(xi);CHKERRQ(ierr);
995b3bd94feSDave May     ierr = PetscFree(eta);CHKERRQ(ierr);
996b3bd94feSDave May     ierr = PetscFree(zeta);CHKERRQ(ierr);
997b3bd94feSDave May   }
99845b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
99945b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
100047c6ae99SBarry Smith   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
100147c6ae99SBarry Smith   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
100247c6ae99SBarry Smith 
100347c6ae99SBarry Smith   ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr);
1004fcfd50ebSBarry Smith   ierr = MatDestroy(&mat);CHKERRQ(ierr);
100547c6ae99SBarry Smith   PetscFunctionReturn(0);
100647c6ae99SBarry Smith }
100747c6ae99SBarry Smith 
1008e727c939SJed Brown PetscErrorCode  DMCreateInterpolation_DA(DM dac,DM daf,Mat *A,Vec *scale)
100947c6ae99SBarry Smith {
101047c6ae99SBarry Smith   PetscErrorCode   ierr;
101147c6ae99SBarry Smith   PetscInt         dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
1012bff4a2f0SMatthew G. Knepley   DMBoundaryType   bxc,byc,bzc,bxf,byf,bzf;
1013aa219208SBarry Smith   DMDAStencilType  stc,stf;
101447c6ae99SBarry Smith   DM_DA            *ddc = (DM_DA*)dac->data;
101547c6ae99SBarry Smith 
101647c6ae99SBarry Smith   PetscFunctionBegin;
101747c6ae99SBarry Smith   PetscValidHeaderSpecific(dac,DM_CLASSID,1);
101847c6ae99SBarry Smith   PetscValidHeaderSpecific(daf,DM_CLASSID,2);
101947c6ae99SBarry Smith   PetscValidPointer(A,3);
102047c6ae99SBarry Smith   if (scale) PetscValidPointer(scale,4);
102147c6ae99SBarry Smith 
10221321219cSEthan Coon   ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr);
10231321219cSEthan Coon   ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr);
102413903a91SSatish Balay   if (dimc != dimf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);
102513903a91SSatish Balay   if (dofc != doff) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);
102613903a91SSatish Balay   if (sc != sf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);
102713903a91SSatish Balay   if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");
102813903a91SSatish Balay   if (stc != stf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");
102947c6ae99SBarry Smith   if (Mc < 2 && Mf > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction");
103047c6ae99SBarry 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");
103147c6ae99SBarry 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");
103247c6ae99SBarry Smith 
1033aa219208SBarry Smith   if (ddc->interptype == DMDA_Q1) {
103447c6ae99SBarry Smith     if (dimc == 1) {
1035e727c939SJed Brown       ierr = DMCreateInterpolation_DA_1D_Q1(dac,daf,A);CHKERRQ(ierr);
103647c6ae99SBarry Smith     } else if (dimc == 2) {
1037e727c939SJed Brown       ierr = DMCreateInterpolation_DA_2D_Q1(dac,daf,A);CHKERRQ(ierr);
103847c6ae99SBarry Smith     } else if (dimc == 3) {
1039e727c939SJed Brown       ierr = DMCreateInterpolation_DA_3D_Q1(dac,daf,A);CHKERRQ(ierr);
1040ce94432eSBarry Smith     } else SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype);
1041aa219208SBarry Smith   } else if (ddc->interptype == DMDA_Q0) {
104247c6ae99SBarry Smith     if (dimc == 1) {
1043e727c939SJed Brown       ierr = DMCreateInterpolation_DA_1D_Q0(dac,daf,A);CHKERRQ(ierr);
104447c6ae99SBarry Smith     } else if (dimc == 2) {
1045e727c939SJed Brown       ierr = DMCreateInterpolation_DA_2D_Q0(dac,daf,A);CHKERRQ(ierr);
104647c6ae99SBarry Smith     } else if (dimc == 3) {
1047e727c939SJed Brown       ierr = DMCreateInterpolation_DA_3D_Q0(dac,daf,A);CHKERRQ(ierr);
1048ce94432eSBarry Smith     } else SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype);
104947c6ae99SBarry Smith   }
105047c6ae99SBarry Smith   if (scale) {
1051e727c939SJed Brown     ierr = DMCreateInterpolationScale((DM)dac,(DM)daf,*A,scale);CHKERRQ(ierr);
105247c6ae99SBarry Smith   }
105347c6ae99SBarry Smith   PetscFunctionReturn(0);
105447c6ae99SBarry Smith }
105547c6ae99SBarry Smith 
1056e727c939SJed Brown PetscErrorCode DMCreateInjection_DA_1D(DM dac,DM daf,VecScatter *inject)
10570bb2b966SJungho Lee {
10580bb2b966SJungho Lee   PetscErrorCode         ierr;
10598ea3bf28SBarry Smith   PetscInt               i,i_start,m_f,Mx,dof;
10608ea3bf28SBarry Smith   const PetscInt         *idx_f;
1061e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f;
10628ea3bf28SBarry Smith   PetscInt               m_ghost,m_ghost_c;
10630bb2b966SJungho Lee   PetscInt               row,i_start_ghost,mx,m_c,nc,ratioi;
10640bb2b966SJungho Lee   PetscInt               i_start_c,i_start_ghost_c;
10650bb2b966SJungho Lee   PetscInt               *cols;
1066bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
10670bb2b966SJungho Lee   Vec                    vecf,vecc;
10680bb2b966SJungho Lee   IS                     isf;
10690bb2b966SJungho Lee 
10700bb2b966SJungho Lee   PetscFunctionBegin;
10710bb2b966SJungho Lee   ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr);
10720bb2b966SJungho Lee   ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
1073bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
10740bb2b966SJungho Lee     ratioi = mx/Mx;
10750bb2b966SJungho 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);
10760bb2b966SJungho Lee   } else {
10770bb2b966SJungho Lee     ratioi = (mx-1)/(Mx-1);
10780bb2b966SJungho 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);
10790bb2b966SJungho Lee   }
10800bb2b966SJungho Lee 
10810bb2b966SJungho Lee   ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr);
10820bb2b966SJungho Lee   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr);
108345b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(daf,&ltog_f);CHKERRQ(ierr);
108445b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
10850bb2b966SJungho Lee 
10860bb2b966SJungho Lee   ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr);
10870bb2b966SJungho Lee   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr);
10880bb2b966SJungho Lee 
10890bb2b966SJungho Lee 
10900bb2b966SJungho Lee   /* loop over local fine grid nodes setting interpolation for those*/
10910bb2b966SJungho Lee   nc   = 0;
1092785e854fSJed Brown   ierr = PetscMalloc1(m_f,&cols);CHKERRQ(ierr);
10930bb2b966SJungho Lee 
10940bb2b966SJungho Lee 
10950bb2b966SJungho Lee   for (i=i_start_c; i<i_start_c+m_c; i++) {
10960bb2b966SJungho Lee     PetscInt i_f = i*ratioi;
10970bb2b966SJungho Lee 
10988865f1eaSKarl Rupp     if (i_f < i_start_ghost || i_f >= i_start_ghost+m_ghost) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\ni_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost);
10998865f1eaSKarl Rupp 
1100e3fbd1f4SBarry Smith     row        = idx_f[(i_f-i_start_ghost)];
1101e3fbd1f4SBarry Smith     cols[nc++] = row;
11020bb2b966SJungho Lee   }
11030bb2b966SJungho Lee 
110445b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
1105ce94432eSBarry Smith   ierr = ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr);
11060bb2b966SJungho Lee   ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr);
11070bb2b966SJungho Lee   ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr);
11080298fd71SBarry Smith   ierr = VecScatterCreate(vecf,isf,vecc,NULL,inject);CHKERRQ(ierr);
11090bb2b966SJungho Lee   ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr);
11100bb2b966SJungho Lee   ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr);
11110bb2b966SJungho Lee   ierr = ISDestroy(&isf);CHKERRQ(ierr);
11120bb2b966SJungho Lee   PetscFunctionReturn(0);
11130bb2b966SJungho Lee }
11140bb2b966SJungho Lee 
1115e727c939SJed Brown PetscErrorCode DMCreateInjection_DA_2D(DM dac,DM daf,VecScatter *inject)
111647c6ae99SBarry Smith {
111747c6ae99SBarry Smith   PetscErrorCode         ierr;
11188ea3bf28SBarry Smith   PetscInt               i,j,i_start,j_start,m_f,n_f,Mx,My,dof;
11198ea3bf28SBarry Smith   const PetscInt         *idx_c,*idx_f;
1120e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
11218ea3bf28SBarry Smith   PetscInt               m_ghost,n_ghost,m_ghost_c,n_ghost_c;
112247c6ae99SBarry Smith   PetscInt               row,i_start_ghost,j_start_ghost,mx,m_c,my,nc,ratioi,ratioj;
1123076202ddSJed Brown   PetscInt               i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c;
112447c6ae99SBarry Smith   PetscInt               *cols;
1125bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
112647c6ae99SBarry Smith   Vec                    vecf,vecc;
112747c6ae99SBarry Smith   IS                     isf;
112847c6ae99SBarry Smith 
112947c6ae99SBarry Smith   PetscFunctionBegin;
11301321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr);
11311321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
1132bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
113347c6ae99SBarry Smith     ratioi = mx/Mx;
113447c6ae99SBarry 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);
113547c6ae99SBarry Smith   } else {
113647c6ae99SBarry Smith     ratioi = (mx-1)/(Mx-1);
113747c6ae99SBarry 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);
113847c6ae99SBarry Smith   }
1139bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) {
114047c6ae99SBarry Smith     ratioj = my/My;
114147c6ae99SBarry 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);
114247c6ae99SBarry Smith   } else {
114347c6ae99SBarry Smith     ratioj = (my-1)/(My-1);
114447c6ae99SBarry 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);
114547c6ae99SBarry Smith   }
114647c6ae99SBarry Smith 
1147aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr);
1148aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr);
114945b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(daf,&ltog_f);CHKERRQ(ierr);
115045b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
115147c6ae99SBarry Smith 
1152aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr);
1153aa219208SBarry Smith   ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr);
115445b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(dac,&ltog_c);CHKERRQ(ierr);
115545b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
115647c6ae99SBarry Smith 
115747c6ae99SBarry Smith   /* loop over local fine grid nodes setting interpolation for those*/
115847c6ae99SBarry Smith   nc   = 0;
1159785e854fSJed Brown   ierr = PetscMalloc1(n_f*m_f,&cols);CHKERRQ(ierr);
1160076202ddSJed Brown   for (j=j_start_c; j<j_start_c+n_c; j++) {
1161076202ddSJed Brown     for (i=i_start_c; i<i_start_c+m_c; i++) {
1162076202ddSJed Brown       PetscInt i_f = i*ratioi,j_f = j*ratioj;
1163076202ddSJed 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\
1164076202ddSJed Brown     j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost);
1165076202ddSJed 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\
1166076202ddSJed Brown     i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost);
1167e3fbd1f4SBarry Smith       row        = idx_f[(m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))];
1168e3fbd1f4SBarry Smith       cols[nc++] = row;
116947c6ae99SBarry Smith     }
117047c6ae99SBarry Smith   }
117145b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
117245b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
117347c6ae99SBarry Smith 
1174ce94432eSBarry Smith   ierr = ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr);
11759a42bb27SBarry Smith   ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr);
11769a42bb27SBarry Smith   ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr);
11770298fd71SBarry Smith   ierr = VecScatterCreate(vecf,isf,vecc,NULL,inject);CHKERRQ(ierr);
11789a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr);
11799a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr);
1180fcfd50ebSBarry Smith   ierr = ISDestroy(&isf);CHKERRQ(ierr);
118147c6ae99SBarry Smith   PetscFunctionReturn(0);
118247c6ae99SBarry Smith }
118347c6ae99SBarry Smith 
1184e727c939SJed Brown PetscErrorCode DMCreateInjection_DA_3D(DM dac,DM daf,VecScatter *inject)
1185b1757049SJed Brown {
1186b1757049SJed Brown   PetscErrorCode         ierr;
1187b1757049SJed Brown   PetscInt               i,j,k,i_start,j_start,k_start,m_f,n_f,p_f,Mx,My,Mz;
1188b1757049SJed Brown   PetscInt               m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c;
1189b1757049SJed Brown   PetscInt               i_start_ghost,j_start_ghost,k_start_ghost;
1190b1757049SJed Brown   PetscInt               mx,my,mz,ratioi,ratioj,ratiok;
1191b1757049SJed Brown   PetscInt               i_start_c,j_start_c,k_start_c;
1192b1757049SJed Brown   PetscInt               m_c,n_c,p_c;
1193b1757049SJed Brown   PetscInt               i_start_ghost_c,j_start_ghost_c,k_start_ghost_c;
1194b1757049SJed Brown   PetscInt               row,nc,dof;
11958ea3bf28SBarry Smith   const PetscInt         *idx_c,*idx_f;
1196e3fbd1f4SBarry Smith   ISLocalToGlobalMapping ltog_f,ltog_c;
1197b1757049SJed Brown   PetscInt               *cols;
1198bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
1199b1757049SJed Brown   Vec                    vecf,vecc;
1200b1757049SJed Brown   IS                     isf;
1201b1757049SJed Brown 
1202b1757049SJed Brown   PetscFunctionBegin;
12031321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr);
12041321219cSEthan Coon   ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
1205b1757049SJed Brown 
1206bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) {
1207b1757049SJed Brown     ratioi = mx/Mx;
1208b1757049SJed 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);
1209b1757049SJed Brown   } else {
1210b1757049SJed Brown     ratioi = (mx-1)/(Mx-1);
1211b1757049SJed 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);
1212b1757049SJed Brown   }
1213bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) {
1214b1757049SJed Brown     ratioj = my/My;
1215b1757049SJed 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);
1216b1757049SJed Brown   } else {
1217b1757049SJed Brown     ratioj = (my-1)/(My-1);
1218b1757049SJed 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);
1219b1757049SJed Brown   }
1220bff4a2f0SMatthew G. Knepley   if (bz == DM_BOUNDARY_PERIODIC) {
1221b1757049SJed Brown     ratiok = mz/Mz;
1222b1757049SJed 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);
1223b1757049SJed Brown   } else {
1224b1757049SJed Brown     ratiok = (mz-1)/(Mz-1);
1225b1757049SJed 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);
1226b1757049SJed Brown   }
1227b1757049SJed Brown 
1228b1757049SJed Brown   ierr = DMDAGetCorners(daf,&i_start,&j_start,&k_start,&m_f,&n_f,&p_f);CHKERRQ(ierr);
1229b1757049SJed Brown   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&k_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr);
123045b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(daf,&ltog_f);CHKERRQ(ierr);
123145b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
1232b1757049SJed Brown 
1233b1757049SJed Brown   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&k_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr);
1234b1757049SJed 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);
123545b6f7e9SBarry Smith   ierr = DMGetLocalToGlobalMapping(dac,&ltog_c);CHKERRQ(ierr);
123645b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
1237b1757049SJed Brown 
1238b1757049SJed Brown 
1239b1757049SJed Brown   /* loop over local fine grid nodes setting interpolation for those*/
1240b1757049SJed Brown   nc   = 0;
1241785e854fSJed Brown   ierr = PetscMalloc1(n_f*m_f*p_f,&cols);CHKERRQ(ierr);
1242b1757049SJed Brown   for (k=k_start_c; k<k_start_c+p_c; k++) {
1243b1757049SJed Brown     for (j=j_start_c; j<j_start_c+n_c; j++) {
1244b1757049SJed Brown       for (i=i_start_c; i<i_start_c+m_c; i++) {
1245b1757049SJed Brown         PetscInt i_f = i*ratioi,j_f = j*ratioj,k_f = k*ratiok;
1246b1757049SJed 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  "
1247b1757049SJed Brown                                                                           "k_c %D k_f %D fine ghost range [%D,%D]",k,k_f,k_start_ghost,k_start_ghost+p_ghost);
1248b1757049SJed 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  "
1249b1757049SJed Brown                                                                           "j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost);
1250b1757049SJed 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  "
1251b1757049SJed Brown                                                                           "i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost);
1252e3fbd1f4SBarry Smith         row        = idx_f[(m_ghost*n_ghost*(k_f-k_start_ghost) + m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))];
1253e3fbd1f4SBarry Smith         cols[nc++] = row;
1254b1757049SJed Brown       }
1255b1757049SJed Brown     }
1256b1757049SJed Brown   }
125745b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr);
125845b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr);
1259b1757049SJed Brown 
1260ce94432eSBarry Smith   ierr = ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr);
1261b1757049SJed Brown   ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr);
1262b1757049SJed Brown   ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr);
12630298fd71SBarry Smith   ierr = VecScatterCreate(vecf,isf,vecc,NULL,inject);CHKERRQ(ierr);
1264b1757049SJed Brown   ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr);
1265b1757049SJed Brown   ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr);
1266fcfd50ebSBarry Smith   ierr = ISDestroy(&isf);CHKERRQ(ierr);
1267b1757049SJed Brown   PetscFunctionReturn(0);
1268b1757049SJed Brown }
1269b1757049SJed Brown 
12706dbf9973SLawrence Mitchell PetscErrorCode  DMCreateInjection_DA(DM dac,DM daf,Mat *mat)
127147c6ae99SBarry Smith {
127247c6ae99SBarry Smith   PetscErrorCode  ierr;
127347c6ae99SBarry Smith   PetscInt        dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
1274bff4a2f0SMatthew G. Knepley   DMBoundaryType  bxc,byc,bzc,bxf,byf,bzf;
1275aa219208SBarry Smith   DMDAStencilType stc,stf;
127660c16ac1SBarry Smith   VecScatter      inject = NULL;
127747c6ae99SBarry Smith 
127847c6ae99SBarry Smith   PetscFunctionBegin;
127947c6ae99SBarry Smith   PetscValidHeaderSpecific(dac,DM_CLASSID,1);
128047c6ae99SBarry Smith   PetscValidHeaderSpecific(daf,DM_CLASSID,2);
12816dbf9973SLawrence Mitchell   PetscValidPointer(mat,3);
128247c6ae99SBarry Smith 
12831321219cSEthan Coon   ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr);
12841321219cSEthan Coon   ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr);
128513903a91SSatish Balay   if (dimc != dimf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);
128613903a91SSatish Balay   if (dofc != doff) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);
128713903a91SSatish Balay   if (sc != sf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);
128813903a91SSatish Balay   if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");
128913903a91SSatish Balay   if (stc != stf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");
129047c6ae99SBarry Smith   if (Mc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction");
129147c6ae99SBarry Smith   if (dimc > 1 && Nc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in y direction");
129247c6ae99SBarry Smith   if (dimc > 2 && Pc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in z direction");
129347c6ae99SBarry Smith 
12940bb2b966SJungho Lee   if (dimc == 1) {
12956dbf9973SLawrence Mitchell     ierr = DMCreateInjection_DA_1D(dac,daf,&inject);CHKERRQ(ierr);
12960bb2b966SJungho Lee   } else if (dimc == 2) {
12976dbf9973SLawrence Mitchell     ierr = DMCreateInjection_DA_2D(dac,daf,&inject);CHKERRQ(ierr);
1298b1757049SJed Brown   } else if (dimc == 3) {
12996dbf9973SLawrence Mitchell     ierr = DMCreateInjection_DA_3D(dac,daf,&inject);CHKERRQ(ierr);
130047c6ae99SBarry Smith   }
13016dbf9973SLawrence Mitchell   ierr = MatCreateScatter(PetscObjectComm((PetscObject)inject), inject, mat);CHKERRQ(ierr);
13026dbf9973SLawrence Mitchell   ierr = VecScatterDestroy(&inject);CHKERRQ(ierr);
130347c6ae99SBarry Smith   PetscFunctionReturn(0);
130447c6ae99SBarry Smith }
130547c6ae99SBarry Smith 
1306e727c939SJed Brown PetscErrorCode  DMCreateAggregates_DA(DM dac,DM daf,Mat *rest)
130747c6ae99SBarry Smith {
130847c6ae99SBarry Smith   PetscErrorCode         ierr;
130947c6ae99SBarry Smith   PetscInt               dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc;
131047c6ae99SBarry Smith   PetscInt               dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
1311bff4a2f0SMatthew G. Knepley   DMBoundaryType         bxc,byc,bzc,bxf,byf,bzf;
1312aa219208SBarry Smith   DMDAStencilType        stc,stf;
131347c6ae99SBarry Smith   PetscInt               i,j,l;
131447c6ae99SBarry Smith   PetscInt               i_start,j_start,l_start, m_f,n_f,p_f;
131547c6ae99SBarry Smith   PetscInt               i_start_ghost,j_start_ghost,l_start_ghost,m_ghost,n_ghost,p_ghost;
13168ea3bf28SBarry Smith   const PetscInt         *idx_f;
131747c6ae99SBarry Smith   PetscInt               i_c,j_c,l_c;
131847c6ae99SBarry Smith   PetscInt               i_start_c,j_start_c,l_start_c, m_c,n_c,p_c;
131947c6ae99SBarry Smith   PetscInt               i_start_ghost_c,j_start_ghost_c,l_start_ghost_c,m_ghost_c,n_ghost_c,p_ghost_c;
13208ea3bf28SBarry Smith   const PetscInt         *idx_c;
132147c6ae99SBarry Smith   PetscInt               d;
132247c6ae99SBarry Smith   PetscInt               a;
132347c6ae99SBarry Smith   PetscInt               max_agg_size;
132447c6ae99SBarry Smith   PetscInt               *fine_nodes;
132547c6ae99SBarry Smith   PetscScalar            *one_vec;
132647c6ae99SBarry Smith   PetscInt               fn_idx;
1327565245c5SBarry Smith   ISLocalToGlobalMapping ltogmf,ltogmc;
132847c6ae99SBarry Smith 
132947c6ae99SBarry Smith   PetscFunctionBegin;
133047c6ae99SBarry Smith   PetscValidHeaderSpecific(dac,DM_CLASSID,1);
133147c6ae99SBarry Smith   PetscValidHeaderSpecific(daf,DM_CLASSID,2);
133247c6ae99SBarry Smith   PetscValidPointer(rest,3);
133347c6ae99SBarry Smith 
13341321219cSEthan Coon   ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr);
13351321219cSEthan Coon   ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr);
133613903a91SSatish Balay   if (dimc != dimf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);
133713903a91SSatish Balay   if (dofc != doff) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);
133813903a91SSatish Balay   if (sc != sf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);
133913903a91SSatish Balay   if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");
134013903a91SSatish Balay   if (stc != stf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");
134147c6ae99SBarry Smith 
134247c6ae99SBarry 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);
134347c6ae99SBarry 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);
134447c6ae99SBarry 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);
134547c6ae99SBarry Smith 
134647c6ae99SBarry Smith   if (Pc < 0) Pc = 1;
134747c6ae99SBarry Smith   if (Pf < 0) Pf = 1;
134847c6ae99SBarry Smith   if (Nc < 0) Nc = 1;
134947c6ae99SBarry Smith   if (Nf < 0) Nf = 1;
135047c6ae99SBarry Smith 
1351aa219208SBarry Smith   ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr);
1352aa219208SBarry Smith   ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr);
1353565245c5SBarry Smith 
1354565245c5SBarry Smith   ierr = DMGetLocalToGlobalMapping(daf,&ltogmf);CHKERRQ(ierr);
1355565245c5SBarry Smith   ierr = ISLocalToGlobalMappingGetIndices(ltogmf,&idx_f);CHKERRQ(ierr);
135647c6ae99SBarry Smith 
1357aa219208SBarry Smith   ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr);
1358aa219208SBarry 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);
1359565245c5SBarry Smith 
1360565245c5SBarry Smith   ierr = DMGetLocalToGlobalMapping(dac,&ltogmc);CHKERRQ(ierr);
1361565245c5SBarry Smith   ierr = ISLocalToGlobalMappingGetIndices(ltogmc,&idx_c);CHKERRQ(ierr);
136247c6ae99SBarry Smith 
136347c6ae99SBarry Smith   /*
136447c6ae99SBarry Smith      Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios
136547c6ae99SBarry Smith      for dimension 1 and 2 respectively.
136647c6ae99SBarry Smith      Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1)
136747c6ae99SBarry Smith      and r_y*j and r_y*(j+1) will be grouped into the same coarse grid agregate.
136847c6ae99SBarry Smith      Each specific dof on the fine grid is mapped to one dof on the coarse grid.
136947c6ae99SBarry Smith   */
137047c6ae99SBarry Smith 
137147c6ae99SBarry Smith   max_agg_size = (Mf/Mc+1)*(Nf/Nc+1)*(Pf/Pc+1);
137247c6ae99SBarry Smith 
137347c6ae99SBarry Smith   /* create the matrix that will contain the restriction operator */
1374ce94432eSBarry Smith   ierr = MatCreateAIJ(PetscObjectComm((PetscObject)daf), m_c*n_c*p_c*dofc, m_f*n_f*p_f*doff, Mc*Nc*Pc*dofc, Mf*Nf*Pf*doff,
13750298fd71SBarry Smith                       max_agg_size, NULL, max_agg_size, NULL, rest);CHKERRQ(ierr);
137647c6ae99SBarry Smith 
137747c6ae99SBarry Smith   /* store nodes in the fine grid here */
1378dcca6d9dSJed Brown   ierr = PetscMalloc2(max_agg_size, &one_vec,max_agg_size, &fine_nodes);CHKERRQ(ierr);
137947c6ae99SBarry Smith   for (i=0; i<max_agg_size; i++) one_vec[i] = 1.0;
138047c6ae99SBarry Smith 
138147c6ae99SBarry Smith   /* loop over all coarse nodes */
138247c6ae99SBarry Smith   for (l_c=l_start_c; l_c<l_start_c+p_c; l_c++) {
138347c6ae99SBarry Smith     for (j_c=j_start_c; j_c<j_start_c+n_c; j_c++) {
138447c6ae99SBarry Smith       for (i_c=i_start_c; i_c<i_start_c+m_c; i_c++) {
138547c6ae99SBarry Smith         for (d=0; d<dofc; d++) {
138647c6ae99SBarry Smith           /* convert to local "natural" numbering and then to PETSc global numbering */
138747c6ae99SBarry 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;
138847c6ae99SBarry Smith 
138947c6ae99SBarry Smith           fn_idx = 0;
139047c6ae99SBarry Smith           /* Corresponding fine points are all points (i_f, j_f, l_f) such that
139147c6ae99SBarry Smith              i_c*Mf/Mc <= i_f < (i_c+1)*Mf/Mc
139247c6ae99SBarry Smith              (same for other dimensions)
139347c6ae99SBarry Smith           */
139447c6ae99SBarry Smith           for (l=l_c*Pf/Pc; l<PetscMin((l_c+1)*Pf/Pc,Pf); l++) {
139547c6ae99SBarry Smith             for (j=j_c*Nf/Nc; j<PetscMin((j_c+1)*Nf/Nc,Nf); j++) {
139647c6ae99SBarry Smith               for (i=i_c*Mf/Mc; i<PetscMin((i_c+1)*Mf/Mc,Mf); i++) {
139747c6ae99SBarry 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;
139847c6ae99SBarry Smith                 fn_idx++;
139947c6ae99SBarry Smith               }
140047c6ae99SBarry Smith             }
140147c6ae99SBarry Smith           }
140247c6ae99SBarry Smith           /* add all these points to one aggregate */
140347c6ae99SBarry Smith           ierr = MatSetValues(*rest, 1, &a, fn_idx, fine_nodes, one_vec, INSERT_VALUES);CHKERRQ(ierr);
140447c6ae99SBarry Smith         }
140547c6ae99SBarry Smith       }
140647c6ae99SBarry Smith     }
140747c6ae99SBarry Smith   }
1408565245c5SBarry Smith   ierr = ISLocalToGlobalMappingRestoreIndices(ltogmf,&idx_f);CHKERRQ(ierr);
1409565245c5SBarry Smith   ierr = ISLocalToGlobalMappingRestoreIndices(ltogmc,&idx_c);CHKERRQ(ierr);
141047c6ae99SBarry Smith   ierr = PetscFree2(one_vec,fine_nodes);CHKERRQ(ierr);
141147c6ae99SBarry Smith   ierr = MatAssemblyBegin(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
141247c6ae99SBarry Smith   ierr = MatAssemblyEnd(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
141347c6ae99SBarry Smith   PetscFunctionReturn(0);
141447c6ae99SBarry Smith }
1415