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 14b45d2f2cSJed Brown #include <petsc-private/daimpl.h> /*I "petscdmda.h" I*/ 1547c6ae99SBarry Smith 1647c6ae99SBarry Smith #undef __FUNCT__ 17e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolationScale" 1847c6ae99SBarry Smith /*@ 19e727c939SJed Brown DMCreateInterpolationScale - Forms L = R*1/diag(R*1) - L.*v is like a coarse grid average of the 2047c6ae99SBarry Smith nearby fine grid points. 2147c6ae99SBarry Smith 2247c6ae99SBarry Smith Input Parameters: 2347c6ae99SBarry Smith + dac - DM that defines a coarse mesh 2447c6ae99SBarry Smith . daf - DM that defines a fine mesh 2547c6ae99SBarry Smith - mat - the restriction (or interpolation operator) from fine to coarse 2647c6ae99SBarry Smith 2747c6ae99SBarry Smith Output Parameter: 2847c6ae99SBarry Smith . scale - the scaled vector 2947c6ae99SBarry Smith 3047c6ae99SBarry Smith Level: developer 3147c6ae99SBarry Smith 32e727c939SJed Brown .seealso: DMCreateInterpolation() 3347c6ae99SBarry Smith 3447c6ae99SBarry Smith @*/ 35e727c939SJed Brown PetscErrorCode DMCreateInterpolationScale(DM dac,DM daf,Mat mat,Vec *scale) 3647c6ae99SBarry Smith { 3747c6ae99SBarry Smith PetscErrorCode ierr; 3847c6ae99SBarry Smith Vec fine; 3947c6ae99SBarry Smith PetscScalar one = 1.0; 4047c6ae99SBarry Smith 4147c6ae99SBarry Smith PetscFunctionBegin; 4247c6ae99SBarry Smith ierr = DMCreateGlobalVector(daf,&fine);CHKERRQ(ierr); 4347c6ae99SBarry Smith ierr = DMCreateGlobalVector(dac,scale);CHKERRQ(ierr); 4447c6ae99SBarry Smith ierr = VecSet(fine,one);CHKERRQ(ierr); 4547c6ae99SBarry Smith ierr = MatRestrict(mat,fine,*scale);CHKERRQ(ierr); 46fcfd50ebSBarry Smith ierr = VecDestroy(&fine);CHKERRQ(ierr); 4747c6ae99SBarry Smith ierr = VecReciprocal(*scale);CHKERRQ(ierr); 4847c6ae99SBarry Smith PetscFunctionReturn(0); 4947c6ae99SBarry Smith } 5047c6ae99SBarry Smith 5147c6ae99SBarry Smith #undef __FUNCT__ 52e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA_1D_Q1" 53e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_1D_Q1(DM dac,DM daf,Mat *A) 5447c6ae99SBarry Smith { 5547c6ae99SBarry Smith PetscErrorCode ierr; 5647c6ae99SBarry Smith PetscInt i,i_start,m_f,Mx,*idx_f; 5747c6ae99SBarry Smith PetscInt m_ghost,*idx_c,m_ghost_c; 5847c6ae99SBarry Smith PetscInt row,col,i_start_ghost,mx,m_c,nc,ratio; 5947c6ae99SBarry Smith PetscInt i_c,i_start_c,i_start_ghost_c,cols[2],dof; 6085afcc9aSBarry Smith PetscScalar v[2],x; 6147c6ae99SBarry Smith Mat mat; 621321219cSEthan Coon DMDABoundaryType bx; 6347c6ae99SBarry Smith 6447c6ae99SBarry Smith PetscFunctionBegin; 651321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr); 661321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 671321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC) { 6847c6ae99SBarry Smith ratio = mx/Mx; 6947c6ae99SBarry 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); 7047c6ae99SBarry Smith } else { 7147c6ae99SBarry Smith ratio = (mx-1)/(Mx-1); 7247c6ae99SBarry 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); 7347c6ae99SBarry Smith } 7447c6ae99SBarry Smith 75aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); 76aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); 770298fd71SBarry Smith ierr = DMDAGetGlobalIndices(daf,NULL,&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); 810298fd71SBarry Smith ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr); 8247c6ae99SBarry Smith 8347c6ae99SBarry Smith /* create interpolation matrix */ 84*ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)dac),&mat);CHKERRQ(ierr); 8547c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr); 8647c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 870298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,2,NULL);CHKERRQ(ierr); 880298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,2,NULL,1,NULL);CHKERRQ(ierr); 8947c6ae99SBarry Smith 9047c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 9185afcc9aSBarry Smith if (!NEWVERSION) { 92b3bd94feSDave May 9347c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 9447c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 9547c6ae99SBarry Smith row = idx_f[dof*(i-i_start_ghost)]/dof; 9647c6ae99SBarry Smith 9747c6ae99SBarry Smith i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 98aa219208SBarry 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\ 9947c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 10047c6ae99SBarry Smith 10147c6ae99SBarry Smith /* 10247c6ae99SBarry Smith Only include those interpolation points that are truly 10347c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 10447c6ae99SBarry Smith in x direction; since they have no right neighbor 10547c6ae99SBarry Smith */ 10647c6ae99SBarry Smith x = ((double)(i - i_c*ratio))/((double)ratio); 10747c6ae99SBarry Smith nc = 0; 10847c6ae99SBarry Smith /* one left and below; or we are right on it */ 10947c6ae99SBarry Smith col = dof*(i_c-i_start_ghost_c); 11047c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 11147c6ae99SBarry Smith v[nc++] = -x + 1.0; 11247c6ae99SBarry Smith /* one right? */ 11347c6ae99SBarry Smith if (i_c*ratio != i) { 11447c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 11547c6ae99SBarry Smith v[nc++] = x; 11647c6ae99SBarry Smith } 11747c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 11847c6ae99SBarry Smith } 119b3bd94feSDave May 120b3bd94feSDave May } else { 121b3bd94feSDave May PetscScalar *xi; 122b3bd94feSDave May PetscInt li,nxi,n; 123b3bd94feSDave May PetscScalar Ni[2]; 124b3bd94feSDave May 125b3bd94feSDave May /* compute local coordinate arrays */ 126b3bd94feSDave May nxi = ratio + 1; 127b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr); 128b3bd94feSDave May for (li=0; li<nxi; li++) { 12952218ef2SJed Brown xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1)); 130b3bd94feSDave May } 131b3bd94feSDave May 132b3bd94feSDave May for (i=i_start; i<i_start+m_f; i++) { 133b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */ 134b3bd94feSDave May row = idx_f[dof*(i-i_start_ghost)]/dof; 135b3bd94feSDave May 136b3bd94feSDave May i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 137b3bd94feSDave 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\ 138b3bd94feSDave May i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 139b3bd94feSDave May 140b3bd94feSDave May /* remainders */ 141b3bd94feSDave May li = i - ratio * (i/ratio); 1428865f1eaSKarl Rupp if (i==mx-1) li = nxi-1; 143b3bd94feSDave May 144b3bd94feSDave May /* corners */ 145b3bd94feSDave May col = dof*(i_c-i_start_ghost_c); 146b3bd94feSDave May cols[0] = idx_c[col]/dof; 147b3bd94feSDave May Ni[0] = 1.0; 148b3bd94feSDave May if ((li==0) || (li==nxi-1)) { 149b3bd94feSDave May ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr); 150b3bd94feSDave May continue; 151b3bd94feSDave May } 152b3bd94feSDave May 153b3bd94feSDave May /* edges + interior */ 154b3bd94feSDave May /* remainders */ 1558865f1eaSKarl Rupp if (i==mx-1) i_c--; 156b3bd94feSDave May 157b3bd94feSDave May col = dof*(i_c-i_start_ghost_c); 158b3bd94feSDave May cols[0] = idx_c[col]/dof; /* one left and below; or we are right on it */ 159b3bd94feSDave May cols[1] = idx_c[col+dof]/dof; 160b3bd94feSDave May 161b3bd94feSDave May Ni[0] = 0.5*(1.0-xi[li]); 162b3bd94feSDave May Ni[1] = 0.5*(1.0+xi[li]); 163b3bd94feSDave May for (n=0; n<2; n++) { 1648865f1eaSKarl Rupp if (PetscAbsScalar(Ni[n])<1.0e-32) cols[n]=-1; 165b3bd94feSDave May } 166b3bd94feSDave May ierr = MatSetValues(mat,1,&row,2,cols,Ni,INSERT_VALUES);CHKERRQ(ierr); 167b3bd94feSDave May } 168b3bd94feSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 169b3bd94feSDave May } 17047c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 17147c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 17247c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 173fcfd50ebSBarry Smith ierr = MatDestroy(&mat);CHKERRQ(ierr); 17447c6ae99SBarry Smith PetscFunctionReturn(0); 17547c6ae99SBarry Smith } 17647c6ae99SBarry Smith 17747c6ae99SBarry Smith #undef __FUNCT__ 178e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA_1D_Q0" 179e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_1D_Q0(DM dac,DM daf,Mat *A) 18047c6ae99SBarry Smith { 18147c6ae99SBarry Smith PetscErrorCode ierr; 18247c6ae99SBarry Smith PetscInt i,i_start,m_f,Mx,*idx_f; 18347c6ae99SBarry Smith PetscInt m_ghost,*idx_c,m_ghost_c; 18447c6ae99SBarry Smith PetscInt row,col,i_start_ghost,mx,m_c,nc,ratio; 18547c6ae99SBarry Smith PetscInt i_c,i_start_c,i_start_ghost_c,cols[2],dof; 18647c6ae99SBarry Smith PetscScalar v[2],x; 18747c6ae99SBarry Smith Mat mat; 1881321219cSEthan Coon DMDABoundaryType bx; 18947c6ae99SBarry Smith 19047c6ae99SBarry Smith PetscFunctionBegin; 1911321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr); 1921321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 1931321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC) { 19447c6ae99SBarry Smith ratio = mx/Mx; 19547c6ae99SBarry 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); 19647c6ae99SBarry Smith } else { 19747c6ae99SBarry Smith ratio = (mx-1)/(Mx-1); 19847c6ae99SBarry 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); 19947c6ae99SBarry Smith } 20047c6ae99SBarry Smith 201aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); 202aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); 2030298fd71SBarry Smith ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr); 20447c6ae99SBarry Smith 205aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr); 206aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr); 2070298fd71SBarry Smith ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr); 20847c6ae99SBarry Smith 20947c6ae99SBarry Smith /* create interpolation matrix */ 210*ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)dac),&mat);CHKERRQ(ierr); 21147c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr); 21247c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 2130298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,2,NULL);CHKERRQ(ierr); 2140298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,2,NULL,0,NULL);CHKERRQ(ierr); 21547c6ae99SBarry Smith 21647c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 21747c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 21847c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 21947c6ae99SBarry Smith row = idx_f[dof*(i-i_start_ghost)]/dof; 22047c6ae99SBarry Smith 22147c6ae99SBarry Smith i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 22247c6ae99SBarry Smith 22347c6ae99SBarry Smith /* 22447c6ae99SBarry Smith Only include those interpolation points that are truly 22547c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 22647c6ae99SBarry Smith in x direction; since they have no right neighbor 22747c6ae99SBarry Smith */ 22847c6ae99SBarry Smith x = ((double)(i - i_c*ratio))/((double)ratio); 22947c6ae99SBarry Smith nc = 0; 23047c6ae99SBarry Smith /* one left and below; or we are right on it */ 23147c6ae99SBarry Smith col = dof*(i_c-i_start_ghost_c); 23247c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 23347c6ae99SBarry Smith v[nc++] = -x + 1.0; 23447c6ae99SBarry Smith /* one right? */ 23547c6ae99SBarry Smith if (i_c*ratio != i) { 23647c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 23747c6ae99SBarry Smith v[nc++] = x; 23847c6ae99SBarry Smith } 23947c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 24047c6ae99SBarry Smith } 24147c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 24247c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 24347c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 244fcfd50ebSBarry Smith ierr = MatDestroy(&mat);CHKERRQ(ierr); 24547c6ae99SBarry Smith ierr = PetscLogFlops(5.0*m_f);CHKERRQ(ierr); 24647c6ae99SBarry Smith PetscFunctionReturn(0); 24747c6ae99SBarry Smith } 24847c6ae99SBarry Smith 24947c6ae99SBarry Smith #undef __FUNCT__ 250e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA_2D_Q1" 251e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_2D_Q1(DM dac,DM daf,Mat *A) 25247c6ae99SBarry Smith { 25347c6ae99SBarry Smith PetscErrorCode ierr; 25447c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 25547c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz; 25647c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj; 25747c6ae99SBarry 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; 25847c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 25947c6ae99SBarry Smith PetscScalar v[4],x,y; 26047c6ae99SBarry Smith Mat mat; 2611321219cSEthan Coon DMDABoundaryType bx,by; 26247c6ae99SBarry Smith 26347c6ae99SBarry Smith PetscFunctionBegin; 2641321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr); 2651321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 2661321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC) { 26747c6ae99SBarry Smith ratioi = mx/Mx; 26847c6ae99SBarry 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); 26947c6ae99SBarry Smith } else { 27047c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 27147c6ae99SBarry 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); 27247c6ae99SBarry Smith } 2731321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC) { 27447c6ae99SBarry Smith ratioj = my/My; 27547c6ae99SBarry 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); 27647c6ae99SBarry Smith } else { 27747c6ae99SBarry Smith ratioj = (my-1)/(My-1); 27847c6ae99SBarry 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); 27947c6ae99SBarry Smith } 28047c6ae99SBarry Smith 28147c6ae99SBarry Smith 282aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 283aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 2840298fd71SBarry Smith ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr); 28547c6ae99SBarry Smith 286aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 287aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 2880298fd71SBarry Smith ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr); 28947c6ae99SBarry Smith 29047c6ae99SBarry Smith /* 291aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 29247c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 29347c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 29447c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 29547c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 29647c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 29747c6ae99SBarry Smith 29847c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 29947c6ae99SBarry Smith */ 300*ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c);CHKERRQ(ierr); 301*ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f);CHKERRQ(ierr); 302*ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f);CHKERRQ(ierr); 30347c6ae99SBarry Smith col_scale = size_f/size_c; 30447c6ae99SBarry Smith col_shift = Mx*My*(rank_f/size_c); 30547c6ae99SBarry Smith 306*ce94432eSBarry Smith ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)daf),m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr); 30747c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 30847c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 30947c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 31047c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 31147c6ae99SBarry Smith 31247c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 31347c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 31447c6ae99SBarry Smith 315aa219208SBarry 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\ 31647c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 317aa219208SBarry 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\ 31847c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 31947c6ae99SBarry Smith 32047c6ae99SBarry Smith /* 32147c6ae99SBarry Smith Only include those interpolation points that are truly 32247c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 32347c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 32447c6ae99SBarry Smith */ 32547c6ae99SBarry Smith nc = 0; 32647c6ae99SBarry Smith /* one left and below; or we are right on it */ 32747c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 32847c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 32947c6ae99SBarry Smith /* one right and below */ 3308865f1eaSKarl Rupp if (i_c*ratioi != i) cols[nc++] = col_shift + idx_c[col+dof]/dof; 33147c6ae99SBarry Smith /* one left and above */ 3328865f1eaSKarl Rupp if (j_c*ratioj != j) cols[nc++] = col_shift + idx_c[col+m_ghost_c*dof]/dof; 33347c6ae99SBarry Smith /* one right and above */ 3348865f1eaSKarl Rupp if (i_c*ratioi != i && j_c*ratioj != j) cols[nc++] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; 33547c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 33647c6ae99SBarry Smith } 33747c6ae99SBarry Smith } 338*ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)daf),&mat);CHKERRQ(ierr); 33947c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr); 34047c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 34147c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 34247c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 34347c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 34447c6ae99SBarry Smith 34547c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 34685afcc9aSBarry Smith if (!NEWVERSION) { 347b3bd94feSDave May 34847c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 34947c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 35047c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 35147c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 35247c6ae99SBarry Smith 35347c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 35447c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 35547c6ae99SBarry Smith 35647c6ae99SBarry Smith /* 35747c6ae99SBarry Smith Only include those interpolation points that are truly 35847c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 35947c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 36047c6ae99SBarry Smith */ 36147c6ae99SBarry Smith x = ((double)(i - i_c*ratioi))/((double)ratioi); 36247c6ae99SBarry Smith y = ((double)(j - j_c*ratioj))/((double)ratioj); 363b3bd94feSDave May 36447c6ae99SBarry Smith nc = 0; 36547c6ae99SBarry Smith /* one left and below; or we are right on it */ 36647c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 36747c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 36847c6ae99SBarry Smith v[nc++] = x*y - x - y + 1.0; 36947c6ae99SBarry Smith /* one right and below */ 37047c6ae99SBarry Smith if (i_c*ratioi != i) { 37147c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+dof]/dof; 37247c6ae99SBarry Smith v[nc++] = -x*y + x; 37347c6ae99SBarry Smith } 37447c6ae99SBarry Smith /* one left and above */ 37547c6ae99SBarry Smith if (j_c*ratioj != j) { 37647c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+m_ghost_c*dof]/dof; 37747c6ae99SBarry Smith v[nc++] = -x*y + y; 37847c6ae99SBarry Smith } 37947c6ae99SBarry Smith /* one right and above */ 38047c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 38147c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; 38247c6ae99SBarry Smith v[nc++] = x*y; 38347c6ae99SBarry Smith } 38447c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 38547c6ae99SBarry Smith } 38647c6ae99SBarry Smith } 387b3bd94feSDave May 388b3bd94feSDave May } else { 389b3bd94feSDave May PetscScalar Ni[4]; 390b3bd94feSDave May PetscScalar *xi,*eta; 391b3bd94feSDave May PetscInt li,nxi,lj,neta; 392b3bd94feSDave May 393b3bd94feSDave May /* compute local coordinate arrays */ 394b3bd94feSDave May nxi = ratioi + 1; 395b3bd94feSDave May neta = ratioj + 1; 396b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr); 397b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*neta,&eta);CHKERRQ(ierr); 398b3bd94feSDave May for (li=0; li<nxi; li++) { 39952218ef2SJed Brown xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1)); 400b3bd94feSDave May } 401b3bd94feSDave May for (lj=0; lj<neta; lj++) { 40252218ef2SJed Brown eta[lj] = -1.0 + (PetscScalar)lj*(2.0/(PetscScalar)(neta-1)); 403b3bd94feSDave May } 404b3bd94feSDave May 405b3bd94feSDave May /* loop over local fine grid nodes setting interpolation for those*/ 406b3bd94feSDave May for (j=j_start; j<j_start+n_f; j++) { 407b3bd94feSDave May for (i=i_start; i<i_start+m_f; i++) { 408b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */ 409b3bd94feSDave May row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 410b3bd94feSDave May 411b3bd94feSDave May i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 412b3bd94feSDave May j_c = (j/ratioj); /* coarse grid node below fine grid node */ 413b3bd94feSDave May 414b3bd94feSDave May /* remainders */ 415b3bd94feSDave May li = i - ratioi * (i/ratioi); 4168865f1eaSKarl Rupp if (i==mx-1) li = nxi-1; 417b3bd94feSDave May lj = j - ratioj * (j/ratioj); 4188865f1eaSKarl Rupp if (j==my-1) lj = neta-1; 419b3bd94feSDave May 420b3bd94feSDave May /* corners */ 421b3bd94feSDave May col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 422b3bd94feSDave May cols[0] = col_shift + idx_c[col]/dof; /* left, below */ 423b3bd94feSDave May Ni[0] = 1.0; 424b3bd94feSDave May if ((li==0) || (li==nxi-1)) { 425b3bd94feSDave May if ((lj==0) || (lj==neta-1)) { 426b3bd94feSDave May ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr); 427b3bd94feSDave May continue; 428b3bd94feSDave May } 429b3bd94feSDave May } 430b3bd94feSDave May 431b3bd94feSDave May /* edges + interior */ 432b3bd94feSDave May /* remainders */ 4338865f1eaSKarl Rupp if (i==mx-1) i_c--; 4348865f1eaSKarl Rupp if (j==my-1) j_c--; 435b3bd94feSDave May 436b3bd94feSDave May col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 437b3bd94feSDave May cols[0] = col_shift + idx_c[col]/dof; /* left, below */ 438b3bd94feSDave May cols[1] = col_shift + idx_c[col+dof]/dof; /* right, below */ 439b3bd94feSDave May cols[2] = col_shift + idx_c[col+m_ghost_c*dof]/dof; /* left, above */ 440b3bd94feSDave May cols[3] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; /* right, above */ 441b3bd94feSDave May 442b3bd94feSDave May Ni[0] = 0.25*(1.0-xi[li])*(1.0-eta[lj]); 443b3bd94feSDave May Ni[1] = 0.25*(1.0+xi[li])*(1.0-eta[lj]); 444b3bd94feSDave May Ni[2] = 0.25*(1.0-xi[li])*(1.0+eta[lj]); 445b3bd94feSDave May Ni[3] = 0.25*(1.0+xi[li])*(1.0+eta[lj]); 446b3bd94feSDave May 447b3bd94feSDave May nc = 0; 4488865f1eaSKarl Rupp if (PetscAbsScalar(Ni[0])<1.0e-32) cols[0]=-1; 4498865f1eaSKarl Rupp if (PetscAbsScalar(Ni[1])<1.0e-32) cols[1]=-1; 4508865f1eaSKarl Rupp if (PetscAbsScalar(Ni[2])<1.0e-32) cols[2]=-1; 4518865f1eaSKarl Rupp if (PetscAbsScalar(Ni[3])<1.0e-32) cols[3]=-1; 452b3bd94feSDave May 453b3bd94feSDave May ierr = MatSetValues(mat,1,&row,4,cols,Ni,INSERT_VALUES);CHKERRQ(ierr); 454b3bd94feSDave May } 455b3bd94feSDave May } 456b3bd94feSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 457b3bd94feSDave May ierr = PetscFree(eta);CHKERRQ(ierr); 458b3bd94feSDave May } 45947c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 46047c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 46147c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 462fcfd50ebSBarry Smith ierr = MatDestroy(&mat);CHKERRQ(ierr); 46347c6ae99SBarry Smith PetscFunctionReturn(0); 46447c6ae99SBarry Smith } 46547c6ae99SBarry Smith 46647c6ae99SBarry Smith /* 46747c6ae99SBarry Smith Contributed by Andrei Draganescu <aidraga@sandia.gov> 46847c6ae99SBarry Smith */ 46947c6ae99SBarry Smith #undef __FUNCT__ 470e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA_2D_Q0" 471e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_2D_Q0(DM dac,DM daf,Mat *A) 47247c6ae99SBarry Smith { 47347c6ae99SBarry Smith PetscErrorCode ierr; 47447c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 47547c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz; 47647c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj; 47747c6ae99SBarry 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; 47847c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 47947c6ae99SBarry Smith PetscScalar v[4]; 48047c6ae99SBarry Smith Mat mat; 4811321219cSEthan Coon DMDABoundaryType bx,by; 48247c6ae99SBarry Smith 48347c6ae99SBarry Smith PetscFunctionBegin; 4841321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr); 4851321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 486*ce94432eSBarry Smith if (bx == DMDA_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x"); 487*ce94432eSBarry Smith if (by == DMDA_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y"); 48847c6ae99SBarry Smith ratioi = mx/Mx; 48947c6ae99SBarry Smith ratioj = my/My; 490*ce94432eSBarry 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"); 491*ce94432eSBarry 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"); 492*ce94432eSBarry Smith if (ratioi != 2) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 2"); 493*ce94432eSBarry Smith if (ratioj != 2) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 2"); 49447c6ae99SBarry Smith 495aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 496aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 4970298fd71SBarry Smith ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr); 49847c6ae99SBarry Smith 499aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 500aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 5010298fd71SBarry Smith ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr); 50247c6ae99SBarry Smith 50347c6ae99SBarry Smith /* 504aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 50547c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 50647c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 50747c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 50847c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 50947c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 51047c6ae99SBarry Smith 51147c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 51247c6ae99SBarry Smith */ 513*ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c);CHKERRQ(ierr); 514*ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f);CHKERRQ(ierr); 515*ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f);CHKERRQ(ierr); 51647c6ae99SBarry Smith col_scale = size_f/size_c; 51747c6ae99SBarry Smith col_shift = Mx*My*(rank_f/size_c); 51847c6ae99SBarry Smith 519*ce94432eSBarry Smith ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)daf),m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr); 52047c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 52147c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 52247c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 52347c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 52447c6ae99SBarry Smith 52547c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 52647c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 52747c6ae99SBarry Smith 528aa219208SBarry 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\ 52947c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 530aa219208SBarry 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\ 53147c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 53247c6ae99SBarry Smith 53347c6ae99SBarry Smith /* 53447c6ae99SBarry Smith Only include those interpolation points that are truly 53547c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 53647c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 53747c6ae99SBarry Smith */ 53847c6ae99SBarry Smith nc = 0; 53947c6ae99SBarry Smith /* one left and below; or we are right on it */ 54047c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 54147c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 54247c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 54347c6ae99SBarry Smith } 54447c6ae99SBarry Smith } 545*ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)daf),&mat);CHKERRQ(ierr); 54647c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr); 54747c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 54847c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 54947c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 55047c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 55147c6ae99SBarry Smith 55247c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 55347c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 55447c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 55547c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 55647c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 55747c6ae99SBarry Smith 55847c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 55947c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 56047c6ae99SBarry Smith nc = 0; 56147c6ae99SBarry Smith /* one left and below; or we are right on it */ 56247c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 56347c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 56447c6ae99SBarry Smith v[nc++] = 1.0; 56547c6ae99SBarry Smith 56647c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 56747c6ae99SBarry Smith } 56847c6ae99SBarry Smith } 56947c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 57047c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 57147c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 572fcfd50ebSBarry Smith ierr = MatDestroy(&mat);CHKERRQ(ierr); 57347c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f);CHKERRQ(ierr); 57447c6ae99SBarry Smith PetscFunctionReturn(0); 57547c6ae99SBarry Smith } 57647c6ae99SBarry Smith 57747c6ae99SBarry Smith /* 57847c6ae99SBarry Smith Contributed by Jianming Yang <jianming-yang@uiowa.edu> 57947c6ae99SBarry Smith */ 58047c6ae99SBarry Smith #undef __FUNCT__ 581e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA_3D_Q0" 582e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_3D_Q0(DM dac,DM daf,Mat *A) 58347c6ae99SBarry Smith { 58447c6ae99SBarry Smith PetscErrorCode ierr; 58547c6ae99SBarry Smith PetscInt i,j,l,i_start,j_start,l_start,m_f,n_f,p_f,Mx,My,Mz,*idx_f,dof; 58647c6ae99SBarry Smith PetscInt m_ghost,n_ghost,p_ghost,*idx_c,m_ghost_c,n_ghost_c,p_ghost_c,nc,*dnz,*onz; 58747c6ae99SBarry 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; 58847c6ae99SBarry 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; 58947c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 59047c6ae99SBarry Smith PetscScalar v[8]; 59147c6ae99SBarry Smith Mat mat; 5921321219cSEthan Coon DMDABoundaryType bx,by,bz; 59347c6ae99SBarry Smith 59447c6ae99SBarry Smith PetscFunctionBegin; 5951321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr); 5961321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 597*ce94432eSBarry Smith if (bx == DMDA_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x"); 598*ce94432eSBarry Smith if (by == DMDA_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y"); 599*ce94432eSBarry Smith if (bz == DMDA_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in z"); 60047c6ae99SBarry Smith ratioi = mx/Mx; 60147c6ae99SBarry Smith ratioj = my/My; 60247c6ae99SBarry Smith ratiol = mz/Mz; 603*ce94432eSBarry 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"); 604*ce94432eSBarry 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"); 605*ce94432eSBarry 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"); 606*ce94432eSBarry Smith if (ratioi != 2 && ratioi != 1) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 1 or 2"); 607*ce94432eSBarry Smith if (ratioj != 2 && ratioj != 1) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 1 or 2"); 608*ce94432eSBarry Smith if (ratiol != 2 && ratiol != 1) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in z must be 1 or 2"); 60947c6ae99SBarry Smith 610aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 611aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 6120298fd71SBarry Smith ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr); 61347c6ae99SBarry Smith 614aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 615aa219208SBarry 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); 6160298fd71SBarry Smith ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr); 61747c6ae99SBarry Smith /* 618aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 61947c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 62047c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 62147c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 62247c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 62347c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 62447c6ae99SBarry Smith 62547c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 62647c6ae99SBarry Smith */ 627*ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c);CHKERRQ(ierr); 628*ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f);CHKERRQ(ierr); 629*ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f);CHKERRQ(ierr); 63047c6ae99SBarry Smith col_scale = size_f/size_c; 63147c6ae99SBarry Smith col_shift = Mx*My*Mz*(rank_f/size_c); 63247c6ae99SBarry Smith 633*ce94432eSBarry Smith ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)daf),m_f*n_f*p_f,col_scale*m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr); 63447c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 63547c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 63647c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 63747c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 63847c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 63947c6ae99SBarry Smith 64047c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 64147c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 64247c6ae99SBarry Smith l_c = (l/ratiol); 64347c6ae99SBarry Smith 644aa219208SBarry 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\ 64547c6ae99SBarry Smith l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c); 646aa219208SBarry 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\ 64747c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 648aa219208SBarry 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\ 64947c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 65047c6ae99SBarry Smith 65147c6ae99SBarry Smith /* 65247c6ae99SBarry Smith Only include those interpolation points that are truly 65347c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 65447c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 65547c6ae99SBarry Smith */ 65647c6ae99SBarry Smith nc = 0; 65747c6ae99SBarry Smith /* one left and below; or we are right on it */ 65847c6ae99SBarry Smith col = dof*(m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c) + m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 65947c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 66047c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 66147c6ae99SBarry Smith } 66247c6ae99SBarry Smith } 66347c6ae99SBarry Smith } 664*ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)daf),&mat);CHKERRQ(ierr); 66547c6ae99SBarry 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); 66647c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 66747c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 66847c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 66947c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 67047c6ae99SBarry Smith 67147c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 67247c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 67347c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 67447c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 67547c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 67647c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 67747c6ae99SBarry Smith 67847c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 67947c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 68047c6ae99SBarry Smith l_c = (l/ratiol); 68147c6ae99SBarry Smith nc = 0; 68247c6ae99SBarry Smith /* one left and below; or we are right on it */ 68347c6ae99SBarry Smith col = dof*(m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c) + m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 68447c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 68547c6ae99SBarry Smith v[nc++] = 1.0; 68647c6ae99SBarry Smith 68747c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 68847c6ae99SBarry Smith } 68947c6ae99SBarry Smith } 69047c6ae99SBarry Smith } 69147c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 69247c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 69347c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 694fcfd50ebSBarry Smith ierr = MatDestroy(&mat);CHKERRQ(ierr); 69547c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f*p_f);CHKERRQ(ierr); 69647c6ae99SBarry Smith PetscFunctionReturn(0); 69747c6ae99SBarry Smith } 69847c6ae99SBarry Smith 69947c6ae99SBarry Smith #undef __FUNCT__ 700e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA_3D_Q1" 701e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA_3D_Q1(DM dac,DM daf,Mat *A) 70247c6ae99SBarry Smith { 70347c6ae99SBarry Smith PetscErrorCode ierr; 70447c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof,l; 70547c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,Mz,mz; 70647c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[8],mx,m_c,my,nc,ratioi,ratioj,ratiok; 70747c6ae99SBarry Smith PetscInt i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c; 70847c6ae99SBarry Smith PetscInt l_start,p_f,l_start_ghost,p_ghost,l_start_c,p_c; 70947c6ae99SBarry Smith PetscInt l_start_ghost_c,p_ghost_c,l_c,*dnz,*onz; 71047c6ae99SBarry Smith PetscScalar v[8],x,y,z; 71147c6ae99SBarry Smith Mat mat; 7121321219cSEthan Coon DMDABoundaryType bx,by,bz; 71347c6ae99SBarry Smith 71447c6ae99SBarry Smith PetscFunctionBegin; 7151321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr); 7161321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 71747c6ae99SBarry Smith if (mx == Mx) { 71847c6ae99SBarry Smith ratioi = 1; 7191321219cSEthan Coon } else if (bx == DMDA_BOUNDARY_PERIODIC) { 72047c6ae99SBarry Smith ratioi = mx/Mx; 72147c6ae99SBarry 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); 72247c6ae99SBarry Smith } else { 72347c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 72447c6ae99SBarry 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); 72547c6ae99SBarry Smith } 72647c6ae99SBarry Smith if (my == My) { 72747c6ae99SBarry Smith ratioj = 1; 7281321219cSEthan Coon } else if (by == DMDA_BOUNDARY_PERIODIC) { 72947c6ae99SBarry Smith ratioj = my/My; 73047c6ae99SBarry 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); 73147c6ae99SBarry Smith } else { 73247c6ae99SBarry Smith ratioj = (my-1)/(My-1); 73347c6ae99SBarry 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); 73447c6ae99SBarry Smith } 73547c6ae99SBarry Smith if (mz == Mz) { 73647c6ae99SBarry Smith ratiok = 1; 7371321219cSEthan Coon } else if (bz == DMDA_BOUNDARY_PERIODIC) { 73847c6ae99SBarry Smith ratiok = mz/Mz; 73947c6ae99SBarry 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); 74047c6ae99SBarry Smith } else { 74147c6ae99SBarry Smith ratiok = (mz-1)/(Mz-1); 74247c6ae99SBarry 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); 74347c6ae99SBarry Smith } 74447c6ae99SBarry Smith 745aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 746aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 7470298fd71SBarry Smith ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr); 74847c6ae99SBarry Smith 749aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 750aa219208SBarry 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); 7510298fd71SBarry Smith ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr); 75247c6ae99SBarry Smith 75347c6ae99SBarry Smith /* create interpolation matrix, determining exact preallocation */ 754*ce94432eSBarry Smith ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)dac),m_f*n_f*p_f,m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr); 75547c6ae99SBarry Smith /* loop over local fine grid nodes counting interpolating points */ 75647c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 75747c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 75847c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 75947c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 76047c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 76147c6ae99SBarry Smith i_c = (i/ratioi); 76247c6ae99SBarry Smith j_c = (j/ratioj); 76347c6ae99SBarry Smith l_c = (l/ratiok); 764aa219208SBarry 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\ 76547c6ae99SBarry Smith l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c); 766aa219208SBarry 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\ 76747c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 768aa219208SBarry 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\ 76947c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 77047c6ae99SBarry Smith 77147c6ae99SBarry Smith /* 77247c6ae99SBarry Smith Only include those interpolation points that are truly 77347c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 77447c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 77547c6ae99SBarry Smith */ 77647c6ae99SBarry Smith nc = 0; 77747c6ae99SBarry Smith col = dof*(m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c) + m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 77847c6ae99SBarry Smith cols[nc++] = idx_c[col]/dof; 77947c6ae99SBarry Smith if (i_c*ratioi != i) { 78047c6ae99SBarry Smith cols[nc++] = idx_c[col+dof]/dof; 78147c6ae99SBarry Smith } 78247c6ae99SBarry Smith if (j_c*ratioj != j) { 78347c6ae99SBarry Smith cols[nc++] = idx_c[col+m_ghost_c*dof]/dof; 78447c6ae99SBarry Smith } 78547c6ae99SBarry Smith if (l_c*ratiok != l) { 78647c6ae99SBarry Smith cols[nc++] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; 78747c6ae99SBarry Smith } 78847c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 78947c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c+1)*dof]/dof; 79047c6ae99SBarry Smith } 79147c6ae99SBarry Smith if (j_c*ratioj != j && l_c*ratiok != l) { 79247c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof; 79347c6ae99SBarry Smith } 79447c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l) { 79547c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; 79647c6ae99SBarry Smith } 79747c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) { 79847c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; 79947c6ae99SBarry Smith } 80047c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 80147c6ae99SBarry Smith } 80247c6ae99SBarry Smith } 80347c6ae99SBarry Smith } 804*ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)dac),&mat);CHKERRQ(ierr); 80547c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f*p_f,m_c*n_c*p_c,mx*my*mz,Mx*My*Mz);CHKERRQ(ierr); 80647c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 80747c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 80847c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 80947c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 81047c6ae99SBarry Smith 81147c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 8122adb9dcfSBarry Smith if (!NEWVERSION) { 813b3bd94feSDave May 81447c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 81547c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 81647c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 81747c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 81847c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 81947c6ae99SBarry Smith 82047c6ae99SBarry Smith i_c = (i/ratioi); 82147c6ae99SBarry Smith j_c = (j/ratioj); 82247c6ae99SBarry Smith l_c = (l/ratiok); 82347c6ae99SBarry Smith 82447c6ae99SBarry Smith /* 82547c6ae99SBarry Smith Only include those interpolation points that are truly 82647c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 82747c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 82847c6ae99SBarry Smith */ 82947c6ae99SBarry Smith x = ((double)(i - i_c*ratioi))/((double)ratioi); 83047c6ae99SBarry Smith y = ((double)(j - j_c*ratioj))/((double)ratioj); 83147c6ae99SBarry Smith z = ((double)(l - l_c*ratiok))/((double)ratiok); 832b3bd94feSDave May 83347c6ae99SBarry Smith nc = 0; 83447c6ae99SBarry Smith /* one left and below; or we are right on it */ 83547c6ae99SBarry Smith col = dof*(m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c)+m_ghost_c*(j_c-j_start_ghost_c)+(i_c-i_start_ghost_c)); 83647c6ae99SBarry Smith 83747c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 83847c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.)); 83947c6ae99SBarry Smith 84047c6ae99SBarry Smith if (i_c*ratioi != i) { 84147c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 84247c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.)); 84347c6ae99SBarry Smith } 84447c6ae99SBarry Smith 84547c6ae99SBarry Smith if (j_c*ratioj != j) { 84647c6ae99SBarry Smith cols[nc] = idx_c[col+m_ghost_c*dof]/dof; 84747c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.)); 84847c6ae99SBarry Smith } 84947c6ae99SBarry Smith 85047c6ae99SBarry Smith if (l_c*ratiok != l) { 85147c6ae99SBarry Smith cols[nc] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; 85247c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.)); 85347c6ae99SBarry Smith } 85447c6ae99SBarry Smith 85547c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 85647c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c+1)*dof]/dof; 85747c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.)); 85847c6ae99SBarry Smith } 85947c6ae99SBarry Smith 86047c6ae99SBarry Smith if (j_c*ratioj != j && l_c*ratiok != l) { 86147c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof; 86247c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.)); 86347c6ae99SBarry Smith } 86447c6ae99SBarry Smith 86547c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l) { 86647c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; 86747c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.)); 86847c6ae99SBarry Smith } 86947c6ae99SBarry Smith 87047c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) { 87147c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; 87247c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.)); 87347c6ae99SBarry Smith } 87447c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 87547c6ae99SBarry Smith } 87647c6ae99SBarry Smith } 87747c6ae99SBarry Smith } 878b3bd94feSDave May 879b3bd94feSDave May } else { 880b3bd94feSDave May PetscScalar *xi,*eta,*zeta; 881b3bd94feSDave May PetscInt li,nxi,lj,neta,lk,nzeta,n; 882b3bd94feSDave May PetscScalar Ni[8]; 883b3bd94feSDave May 884b3bd94feSDave May /* compute local coordinate arrays */ 885b3bd94feSDave May nxi = ratioi + 1; 886b3bd94feSDave May neta = ratioj + 1; 887b3bd94feSDave May nzeta = ratiok + 1; 888b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr); 889b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*neta,&eta);CHKERRQ(ierr); 890b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*nzeta,&zeta);CHKERRQ(ierr); 8918865f1eaSKarl Rupp for (li=0; li<nxi; li++) xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1)); 8928865f1eaSKarl Rupp for (lj=0; lj<neta; lj++) eta[lj] = -1.0 + (PetscScalar)lj*(2.0/(PetscScalar)(neta-1)); 8938865f1eaSKarl Rupp for (lk=0; lk<nzeta; lk++) zeta[lk] = -1.0 + (PetscScalar)lk*(2.0/(PetscScalar)(nzeta-1)); 894b3bd94feSDave May 895b3bd94feSDave May for (l=l_start; l<l_start+p_f; l++) { 896b3bd94feSDave May for (j=j_start; j<j_start+n_f; j++) { 897b3bd94feSDave May for (i=i_start; i<i_start+m_f; i++) { 898b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */ 899b3bd94feSDave May row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 900b3bd94feSDave May 901b3bd94feSDave May i_c = (i/ratioi); 902b3bd94feSDave May j_c = (j/ratioj); 903b3bd94feSDave May l_c = (l/ratiok); 904b3bd94feSDave May 905b3bd94feSDave May /* remainders */ 906b3bd94feSDave May li = i - ratioi * (i/ratioi); 9078865f1eaSKarl Rupp if (i==mx-1) li = nxi-1; 908b3bd94feSDave May lj = j - ratioj * (j/ratioj); 9098865f1eaSKarl Rupp if (j==my-1) lj = neta-1; 910b3bd94feSDave May lk = l - ratiok * (l/ratiok); 9118865f1eaSKarl Rupp if (l==mz-1) lk = nzeta-1; 912b3bd94feSDave May 913b3bd94feSDave May /* corners */ 914b3bd94feSDave May col = dof*(m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c)+m_ghost_c*(j_c-j_start_ghost_c)+(i_c-i_start_ghost_c)); 915b3bd94feSDave May cols[0] = idx_c[col]/dof; 916b3bd94feSDave May Ni[0] = 1.0; 917b3bd94feSDave May if ((li==0) || (li==nxi-1)) { 918b3bd94feSDave May if ((lj==0) || (lj==neta-1)) { 919b3bd94feSDave May if ((lk==0) || (lk==nzeta-1)) { 920b3bd94feSDave May ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr); 921b3bd94feSDave May continue; 922b3bd94feSDave May } 923b3bd94feSDave May } 924b3bd94feSDave May } 925b3bd94feSDave May 926b3bd94feSDave May /* edges + interior */ 927b3bd94feSDave May /* remainders */ 9288865f1eaSKarl Rupp if (i==mx-1) i_c--; 9298865f1eaSKarl Rupp if (j==my-1) j_c--; 9308865f1eaSKarl Rupp if (l==mz-1) l_c--; 931b3bd94feSDave May 932b3bd94feSDave May col = dof*(m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c) + m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 933b3bd94feSDave May cols[0] = idx_c[col]/dof; /* one left and below; or we are right on it */ 934b3bd94feSDave May cols[1] = idx_c[col+dof]/dof; /* one right and below */ 935b3bd94feSDave May cols[2] = idx_c[col+m_ghost_c*dof]/dof; /* one left and above */ 936b3bd94feSDave May cols[3] = idx_c[col+(m_ghost_c+1)*dof]/dof; /* one right and above */ 937b3bd94feSDave May 938b3bd94feSDave May cols[4] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; /* one left and below and front; or we are right on it */ 939b3bd94feSDave May cols[5] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; /* one right and below, and front */ 940b3bd94feSDave May cols[6] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof; /* one left and above and front*/ 941b3bd94feSDave May cols[7] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; /* one right and above and front */ 942b3bd94feSDave May 943b3bd94feSDave May Ni[0] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]); 944b3bd94feSDave May Ni[1] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]); 945b3bd94feSDave May Ni[2] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]); 946b3bd94feSDave May Ni[3] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]); 947b3bd94feSDave May 948b3bd94feSDave May Ni[4] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]); 949b3bd94feSDave May Ni[5] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]); 950b3bd94feSDave May Ni[6] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]); 951b3bd94feSDave May Ni[7] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]); 952b3bd94feSDave May 953b3bd94feSDave May for (n=0; n<8; n++) { 9548865f1eaSKarl Rupp if (PetscAbsScalar(Ni[n])<1.0e-32) cols[n]=-1; 955b3bd94feSDave May } 956b3bd94feSDave May ierr = MatSetValues(mat,1,&row,8,cols,Ni,INSERT_VALUES);CHKERRQ(ierr); 957b3bd94feSDave May 958b3bd94feSDave May } 959b3bd94feSDave May } 960b3bd94feSDave May } 961b3bd94feSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 962b3bd94feSDave May ierr = PetscFree(eta);CHKERRQ(ierr); 963b3bd94feSDave May ierr = PetscFree(zeta);CHKERRQ(ierr); 964b3bd94feSDave May } 965b3bd94feSDave May 96647c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 96747c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 96847c6ae99SBarry Smith 96947c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 970fcfd50ebSBarry Smith ierr = MatDestroy(&mat);CHKERRQ(ierr); 97147c6ae99SBarry Smith PetscFunctionReturn(0); 97247c6ae99SBarry Smith } 97347c6ae99SBarry Smith 97447c6ae99SBarry Smith #undef __FUNCT__ 975e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_DA" 976e727c939SJed Brown PetscErrorCode DMCreateInterpolation_DA(DM dac,DM daf,Mat *A,Vec *scale) 97747c6ae99SBarry Smith { 97847c6ae99SBarry Smith PetscErrorCode ierr; 97947c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 9801321219cSEthan Coon DMDABoundaryType bxc,byc,bzc,bxf,byf,bzf; 981aa219208SBarry Smith DMDAStencilType stc,stf; 98247c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data; 98347c6ae99SBarry Smith 98447c6ae99SBarry Smith PetscFunctionBegin; 98547c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 98647c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 98747c6ae99SBarry Smith PetscValidPointer(A,3); 98847c6ae99SBarry Smith if (scale) PetscValidPointer(scale,4); 98947c6ae99SBarry Smith 9901321219cSEthan Coon ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr); 9911321219cSEthan Coon ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr); 992*ce94432eSBarry Smith if (dimc != dimf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 993*ce94432eSBarry Smith if (dofc != doff) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); 994*ce94432eSBarry Smith if (sc != sf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr); 995*ce94432eSBarry Smith if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr); 996*ce94432eSBarry Smith if (stc != stf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 99747c6ae99SBarry Smith if (Mc < 2 && Mf > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction"); 99847c6ae99SBarry 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"); 99947c6ae99SBarry 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"); 100047c6ae99SBarry Smith 1001aa219208SBarry Smith if (ddc->interptype == DMDA_Q1) { 100247c6ae99SBarry Smith if (dimc == 1) { 1003e727c939SJed Brown ierr = DMCreateInterpolation_DA_1D_Q1(dac,daf,A);CHKERRQ(ierr); 100447c6ae99SBarry Smith } else if (dimc == 2) { 1005e727c939SJed Brown ierr = DMCreateInterpolation_DA_2D_Q1(dac,daf,A);CHKERRQ(ierr); 100647c6ae99SBarry Smith } else if (dimc == 3) { 1007e727c939SJed Brown ierr = DMCreateInterpolation_DA_3D_Q1(dac,daf,A);CHKERRQ(ierr); 1008*ce94432eSBarry Smith } else SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype); 1009aa219208SBarry Smith } else if (ddc->interptype == DMDA_Q0) { 101047c6ae99SBarry Smith if (dimc == 1) { 1011e727c939SJed Brown ierr = DMCreateInterpolation_DA_1D_Q0(dac,daf,A);CHKERRQ(ierr); 101247c6ae99SBarry Smith } else if (dimc == 2) { 1013e727c939SJed Brown ierr = DMCreateInterpolation_DA_2D_Q0(dac,daf,A);CHKERRQ(ierr); 101447c6ae99SBarry Smith } else if (dimc == 3) { 1015e727c939SJed Brown ierr = DMCreateInterpolation_DA_3D_Q0(dac,daf,A);CHKERRQ(ierr); 1016*ce94432eSBarry Smith } else SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype); 101747c6ae99SBarry Smith } 101847c6ae99SBarry Smith if (scale) { 1019e727c939SJed Brown ierr = DMCreateInterpolationScale((DM)dac,(DM)daf,*A,scale);CHKERRQ(ierr); 102047c6ae99SBarry Smith } 102147c6ae99SBarry Smith PetscFunctionReturn(0); 102247c6ae99SBarry Smith } 102347c6ae99SBarry Smith 102447c6ae99SBarry Smith #undef __FUNCT__ 1025e727c939SJed Brown #define __FUNCT__ "DMCreateInjection_DA_1D" 1026e727c939SJed Brown PetscErrorCode DMCreateInjection_DA_1D(DM dac,DM daf,VecScatter *inject) 10270bb2b966SJungho Lee { 10280bb2b966SJungho Lee PetscErrorCode ierr; 10290bb2b966SJungho Lee PetscInt i,i_start,m_f,Mx,*idx_f,dof; 10300bb2b966SJungho Lee PetscInt m_ghost,*idx_c,m_ghost_c; 10310bb2b966SJungho Lee PetscInt row,i_start_ghost,mx,m_c,nc,ratioi; 10320bb2b966SJungho Lee PetscInt i_start_c,i_start_ghost_c; 10330bb2b966SJungho Lee PetscInt *cols; 10340bb2b966SJungho Lee DMDABoundaryType bx; 10350bb2b966SJungho Lee Vec vecf,vecc; 10360bb2b966SJungho Lee IS isf; 10370bb2b966SJungho Lee 10380bb2b966SJungho Lee PetscFunctionBegin; 10390bb2b966SJungho Lee ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr); 10400bb2b966SJungho Lee ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 10410bb2b966SJungho Lee if (bx == DMDA_BOUNDARY_PERIODIC) { 10420bb2b966SJungho Lee ratioi = mx/Mx; 10430bb2b966SJungho 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); 10440bb2b966SJungho Lee } else { 10450bb2b966SJungho Lee ratioi = (mx-1)/(Mx-1); 10460bb2b966SJungho 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); 10470bb2b966SJungho Lee } 10480bb2b966SJungho Lee 10490bb2b966SJungho Lee ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); 10500bb2b966SJungho Lee ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); 10510298fd71SBarry Smith ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr); 10520bb2b966SJungho Lee 10530bb2b966SJungho Lee ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr); 10540bb2b966SJungho Lee ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr); 10550298fd71SBarry Smith ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr); 10560bb2b966SJungho Lee 10570bb2b966SJungho Lee 10580bb2b966SJungho Lee /* loop over local fine grid nodes setting interpolation for those*/ 10590bb2b966SJungho Lee nc = 0; 10600bb2b966SJungho Lee ierr = PetscMalloc(m_f*sizeof(PetscInt),&cols);CHKERRQ(ierr); 10610bb2b966SJungho Lee 10620bb2b966SJungho Lee 10630bb2b966SJungho Lee for (i=i_start_c; i<i_start_c+m_c; i++) { 10640bb2b966SJungho Lee PetscInt i_f = i*ratioi; 10650bb2b966SJungho Lee 10668865f1eaSKarl 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); 10678865f1eaSKarl Rupp 10680bb2b966SJungho Lee row = idx_f[dof*(i_f-i_start_ghost)]; 10690bb2b966SJungho Lee cols[nc++] = row/dof; 10700bb2b966SJungho Lee } 10710bb2b966SJungho Lee 10720bb2b966SJungho Lee 1073*ce94432eSBarry Smith ierr = ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); 10740bb2b966SJungho Lee ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); 10750bb2b966SJungho Lee ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); 10760298fd71SBarry Smith ierr = VecScatterCreate(vecf,isf,vecc,NULL,inject);CHKERRQ(ierr); 10770bb2b966SJungho Lee ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); 10780bb2b966SJungho Lee ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); 10790bb2b966SJungho Lee ierr = ISDestroy(&isf);CHKERRQ(ierr); 10800bb2b966SJungho Lee PetscFunctionReturn(0); 10810bb2b966SJungho Lee } 10820bb2b966SJungho Lee 10830bb2b966SJungho Lee #undef __FUNCT__ 1084e727c939SJed Brown #define __FUNCT__ "DMCreateInjection_DA_2D" 1085e727c939SJed Brown PetscErrorCode DMCreateInjection_DA_2D(DM dac,DM daf,VecScatter *inject) 108647c6ae99SBarry Smith { 108747c6ae99SBarry Smith PetscErrorCode ierr; 108847c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 108947c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c; 109047c6ae99SBarry Smith PetscInt row,i_start_ghost,j_start_ghost,mx,m_c,my,nc,ratioi,ratioj; 1091076202ddSJed Brown PetscInt i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c; 109247c6ae99SBarry Smith PetscInt *cols; 10931321219cSEthan Coon DMDABoundaryType bx,by; 109447c6ae99SBarry Smith Vec vecf,vecc; 109547c6ae99SBarry Smith IS isf; 109647c6ae99SBarry Smith 109747c6ae99SBarry Smith PetscFunctionBegin; 10981321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr); 10991321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 11001321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC) { 110147c6ae99SBarry Smith ratioi = mx/Mx; 110247c6ae99SBarry 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); 110347c6ae99SBarry Smith } else { 110447c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 110547c6ae99SBarry 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); 110647c6ae99SBarry Smith } 11071321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC) { 110847c6ae99SBarry Smith ratioj = my/My; 110947c6ae99SBarry 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); 111047c6ae99SBarry Smith } else { 111147c6ae99SBarry Smith ratioj = (my-1)/(My-1); 111247c6ae99SBarry 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); 111347c6ae99SBarry Smith } 111447c6ae99SBarry Smith 1115aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 1116aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 11170298fd71SBarry Smith ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr); 111847c6ae99SBarry Smith 1119aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 1120aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 11210298fd71SBarry Smith ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr); 112247c6ae99SBarry Smith 112347c6ae99SBarry Smith 112447c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 112547c6ae99SBarry Smith nc = 0; 112647c6ae99SBarry Smith ierr = PetscMalloc(n_f*m_f*sizeof(PetscInt),&cols);CHKERRQ(ierr); 1127076202ddSJed Brown for (j=j_start_c; j<j_start_c+n_c; j++) { 1128076202ddSJed Brown for (i=i_start_c; i<i_start_c+m_c; i++) { 1129076202ddSJed Brown PetscInt i_f = i*ratioi,j_f = j*ratioj; 1130076202ddSJed 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\ 1131076202ddSJed Brown j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost); 1132076202ddSJed 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\ 1133076202ddSJed Brown i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost); 1134076202ddSJed Brown row = idx_f[dof*(m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))]; 113547c6ae99SBarry Smith cols[nc++] = row/dof; 113647c6ae99SBarry Smith } 113747c6ae99SBarry Smith } 113847c6ae99SBarry Smith 1139*ce94432eSBarry Smith ierr = ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); 11409a42bb27SBarry Smith ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); 11419a42bb27SBarry Smith ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); 11420298fd71SBarry Smith ierr = VecScatterCreate(vecf,isf,vecc,NULL,inject);CHKERRQ(ierr); 11439a42bb27SBarry Smith ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); 11449a42bb27SBarry Smith ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); 1145fcfd50ebSBarry Smith ierr = ISDestroy(&isf);CHKERRQ(ierr); 114647c6ae99SBarry Smith PetscFunctionReturn(0); 114747c6ae99SBarry Smith } 114847c6ae99SBarry Smith 114947c6ae99SBarry Smith #undef __FUNCT__ 1150e727c939SJed Brown #define __FUNCT__ "DMCreateInjection_DA_3D" 1151e727c939SJed Brown PetscErrorCode DMCreateInjection_DA_3D(DM dac,DM daf,VecScatter *inject) 1152b1757049SJed Brown { 1153b1757049SJed Brown PetscErrorCode ierr; 1154b1757049SJed Brown PetscInt i,j,k,i_start,j_start,k_start,m_f,n_f,p_f,Mx,My,Mz; 1155b1757049SJed Brown PetscInt m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c; 1156b1757049SJed Brown PetscInt i_start_ghost,j_start_ghost,k_start_ghost; 1157b1757049SJed Brown PetscInt mx,my,mz,ratioi,ratioj,ratiok; 1158b1757049SJed Brown PetscInt i_start_c,j_start_c,k_start_c; 1159b1757049SJed Brown PetscInt m_c,n_c,p_c; 1160b1757049SJed Brown PetscInt i_start_ghost_c,j_start_ghost_c,k_start_ghost_c; 1161b1757049SJed Brown PetscInt row,nc,dof; 1162b1757049SJed Brown PetscInt *idx_c,*idx_f; 1163b1757049SJed Brown PetscInt *cols; 11641321219cSEthan Coon DMDABoundaryType bx,by,bz; 1165b1757049SJed Brown Vec vecf,vecc; 1166b1757049SJed Brown IS isf; 1167b1757049SJed Brown 1168b1757049SJed Brown PetscFunctionBegin; 11691321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr); 11701321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 1171b1757049SJed Brown 11721321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC) { 1173b1757049SJed Brown ratioi = mx/Mx; 1174b1757049SJed 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); 1175b1757049SJed Brown } else { 1176b1757049SJed Brown ratioi = (mx-1)/(Mx-1); 1177b1757049SJed 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); 1178b1757049SJed Brown } 11791321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC) { 1180b1757049SJed Brown ratioj = my/My; 1181b1757049SJed 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); 1182b1757049SJed Brown } else { 1183b1757049SJed Brown ratioj = (my-1)/(My-1); 1184b1757049SJed 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); 1185b1757049SJed Brown } 11861321219cSEthan Coon if (bz == DMDA_BOUNDARY_PERIODIC) { 1187b1757049SJed Brown ratiok = mz/Mz; 1188b1757049SJed 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); 1189b1757049SJed Brown } else { 1190b1757049SJed Brown ratiok = (mz-1)/(Mz-1); 1191b1757049SJed 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); 1192b1757049SJed Brown } 1193b1757049SJed Brown 1194b1757049SJed Brown ierr = DMDAGetCorners(daf,&i_start,&j_start,&k_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 1195b1757049SJed Brown ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&k_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 11960298fd71SBarry Smith ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr); 1197b1757049SJed Brown 1198b1757049SJed Brown ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&k_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 1199b1757049SJed 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); 12000298fd71SBarry Smith ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr); 1201b1757049SJed Brown 1202b1757049SJed Brown 1203b1757049SJed Brown /* loop over local fine grid nodes setting interpolation for those*/ 1204b1757049SJed Brown nc = 0; 1205b1757049SJed Brown ierr = PetscMalloc(n_f*m_f*p_f*sizeof(PetscInt),&cols);CHKERRQ(ierr); 1206b1757049SJed Brown for (k=k_start_c; k<k_start_c+p_c; k++) { 1207b1757049SJed Brown for (j=j_start_c; j<j_start_c+n_c; j++) { 1208b1757049SJed Brown for (i=i_start_c; i<i_start_c+m_c; i++) { 1209b1757049SJed Brown PetscInt i_f = i*ratioi,j_f = j*ratioj,k_f = k*ratiok; 1210b1757049SJed 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 " 1211b1757049SJed Brown "k_c %D k_f %D fine ghost range [%D,%D]",k,k_f,k_start_ghost,k_start_ghost+p_ghost); 1212b1757049SJed 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 " 1213b1757049SJed Brown "j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost); 1214b1757049SJed 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 " 1215b1757049SJed Brown "i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost); 1216b1757049SJed Brown row = idx_f[dof*(m_ghost*n_ghost*(k_f-k_start_ghost) + m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))]; 1217b1757049SJed Brown cols[nc++] = row/dof; 1218b1757049SJed Brown } 1219b1757049SJed Brown } 1220b1757049SJed Brown } 1221b1757049SJed Brown 1222*ce94432eSBarry Smith ierr = ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); 1223b1757049SJed Brown ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); 1224b1757049SJed Brown ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); 12250298fd71SBarry Smith ierr = VecScatterCreate(vecf,isf,vecc,NULL,inject);CHKERRQ(ierr); 1226b1757049SJed Brown ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); 1227b1757049SJed Brown ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); 1228fcfd50ebSBarry Smith ierr = ISDestroy(&isf);CHKERRQ(ierr); 1229b1757049SJed Brown PetscFunctionReturn(0); 1230b1757049SJed Brown } 1231b1757049SJed Brown 1232b1757049SJed Brown #undef __FUNCT__ 1233e727c939SJed Brown #define __FUNCT__ "DMCreateInjection_DA" 1234e727c939SJed Brown PetscErrorCode DMCreateInjection_DA(DM dac,DM daf,VecScatter *inject) 123547c6ae99SBarry Smith { 123647c6ae99SBarry Smith PetscErrorCode ierr; 123747c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 12381321219cSEthan Coon DMDABoundaryType bxc,byc,bzc,bxf,byf,bzf; 1239aa219208SBarry Smith DMDAStencilType stc,stf; 124047c6ae99SBarry Smith 124147c6ae99SBarry Smith PetscFunctionBegin; 124247c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 124347c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 124447c6ae99SBarry Smith PetscValidPointer(inject,3); 124547c6ae99SBarry Smith 12461321219cSEthan Coon ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr); 12471321219cSEthan Coon ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr); 1248*ce94432eSBarry Smith if (dimc != dimf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 1249*ce94432eSBarry Smith if (dofc != doff) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); 1250*ce94432eSBarry Smith if (sc != sf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr); 1251*ce94432eSBarry Smith if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr); 1252*ce94432eSBarry Smith if (stc != stf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 125347c6ae99SBarry Smith if (Mc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction"); 125447c6ae99SBarry Smith if (dimc > 1 && Nc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in y direction"); 125547c6ae99SBarry Smith if (dimc > 2 && Pc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in z direction"); 125647c6ae99SBarry Smith 12570bb2b966SJungho Lee if (dimc == 1) { 1258e727c939SJed Brown ierr = DMCreateInjection_DA_1D(dac,daf,inject);CHKERRQ(ierr); 12590bb2b966SJungho Lee } else if (dimc == 2) { 1260e727c939SJed Brown ierr = DMCreateInjection_DA_2D(dac,daf,inject);CHKERRQ(ierr); 1261b1757049SJed Brown } else if (dimc == 3) { 1262e727c939SJed Brown ierr = DMCreateInjection_DA_3D(dac,daf,inject);CHKERRQ(ierr); 126347c6ae99SBarry Smith } 126447c6ae99SBarry Smith PetscFunctionReturn(0); 126547c6ae99SBarry Smith } 126647c6ae99SBarry Smith 126747c6ae99SBarry Smith #undef __FUNCT__ 1268e727c939SJed Brown #define __FUNCT__ "DMCreateAggregates_DA" 1269e727c939SJed Brown PetscErrorCode DMCreateAggregates_DA(DM dac,DM daf,Mat *rest) 127047c6ae99SBarry Smith { 127147c6ae99SBarry Smith PetscErrorCode ierr; 127247c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc; 127347c6ae99SBarry Smith PetscInt dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 12741321219cSEthan Coon DMDABoundaryType bxc,byc,bzc,bxf,byf,bzf; 1275aa219208SBarry Smith DMDAStencilType stc,stf; 127647c6ae99SBarry Smith PetscInt i,j,l; 127747c6ae99SBarry Smith PetscInt i_start,j_start,l_start, m_f,n_f,p_f; 127847c6ae99SBarry Smith PetscInt i_start_ghost,j_start_ghost,l_start_ghost,m_ghost,n_ghost,p_ghost; 127947c6ae99SBarry Smith PetscInt *idx_f; 128047c6ae99SBarry Smith PetscInt i_c,j_c,l_c; 128147c6ae99SBarry Smith PetscInt i_start_c,j_start_c,l_start_c, m_c,n_c,p_c; 128247c6ae99SBarry Smith PetscInt i_start_ghost_c,j_start_ghost_c,l_start_ghost_c,m_ghost_c,n_ghost_c,p_ghost_c; 128347c6ae99SBarry Smith PetscInt *idx_c; 128447c6ae99SBarry Smith PetscInt d; 128547c6ae99SBarry Smith PetscInt a; 128647c6ae99SBarry Smith PetscInt max_agg_size; 128747c6ae99SBarry Smith PetscInt *fine_nodes; 128847c6ae99SBarry Smith PetscScalar *one_vec; 128947c6ae99SBarry Smith PetscInt fn_idx; 129047c6ae99SBarry Smith 129147c6ae99SBarry Smith PetscFunctionBegin; 129247c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 129347c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 129447c6ae99SBarry Smith PetscValidPointer(rest,3); 129547c6ae99SBarry Smith 12961321219cSEthan Coon ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr); 12971321219cSEthan Coon ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr); 1298*ce94432eSBarry Smith if (dimc != dimf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 1299*ce94432eSBarry Smith if (dofc != doff) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); 1300*ce94432eSBarry Smith if (sc != sf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr); 1301*ce94432eSBarry Smith if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr); 1302*ce94432eSBarry Smith if (stc != stf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 130347c6ae99SBarry Smith 130447c6ae99SBarry 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); 130547c6ae99SBarry 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); 130647c6ae99SBarry 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); 130747c6ae99SBarry Smith 130847c6ae99SBarry Smith if (Pc < 0) Pc = 1; 130947c6ae99SBarry Smith if (Pf < 0) Pf = 1; 131047c6ae99SBarry Smith if (Nc < 0) Nc = 1; 131147c6ae99SBarry Smith if (Nf < 0) Nf = 1; 131247c6ae99SBarry Smith 1313aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 1314aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 13150298fd71SBarry Smith ierr = DMDAGetGlobalIndices(daf,NULL,&idx_f);CHKERRQ(ierr); 131647c6ae99SBarry Smith 1317aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 1318aa219208SBarry 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); 13190298fd71SBarry Smith ierr = DMDAGetGlobalIndices(dac,NULL,&idx_c);CHKERRQ(ierr); 132047c6ae99SBarry Smith 132147c6ae99SBarry Smith /* 132247c6ae99SBarry Smith Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios 132347c6ae99SBarry Smith for dimension 1 and 2 respectively. 132447c6ae99SBarry Smith Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1) 132547c6ae99SBarry Smith and r_y*j and r_y*(j+1) will be grouped into the same coarse grid agregate. 132647c6ae99SBarry Smith Each specific dof on the fine grid is mapped to one dof on the coarse grid. 132747c6ae99SBarry Smith */ 132847c6ae99SBarry Smith 132947c6ae99SBarry Smith max_agg_size = (Mf/Mc+1)*(Nf/Nc+1)*(Pf/Pc+1); 133047c6ae99SBarry Smith 133147c6ae99SBarry Smith /* create the matrix that will contain the restriction operator */ 1332*ce94432eSBarry 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, 13330298fd71SBarry Smith max_agg_size, NULL, max_agg_size, NULL, rest);CHKERRQ(ierr); 133447c6ae99SBarry Smith 133547c6ae99SBarry Smith /* store nodes in the fine grid here */ 133647c6ae99SBarry Smith ierr = PetscMalloc2(max_agg_size,PetscScalar, &one_vec,max_agg_size,PetscInt, &fine_nodes);CHKERRQ(ierr); 133747c6ae99SBarry Smith for (i=0; i<max_agg_size; i++) one_vec[i] = 1.0; 133847c6ae99SBarry Smith 133947c6ae99SBarry Smith /* loop over all coarse nodes */ 134047c6ae99SBarry Smith for (l_c=l_start_c; l_c<l_start_c+p_c; l_c++) { 134147c6ae99SBarry Smith for (j_c=j_start_c; j_c<j_start_c+n_c; j_c++) { 134247c6ae99SBarry Smith for (i_c=i_start_c; i_c<i_start_c+m_c; i_c++) { 134347c6ae99SBarry Smith for (d=0; d<dofc; d++) { 134447c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 134547c6ae99SBarry 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; 134647c6ae99SBarry Smith 134747c6ae99SBarry Smith fn_idx = 0; 134847c6ae99SBarry Smith /* Corresponding fine points are all points (i_f, j_f, l_f) such that 134947c6ae99SBarry Smith i_c*Mf/Mc <= i_f < (i_c+1)*Mf/Mc 135047c6ae99SBarry Smith (same for other dimensions) 135147c6ae99SBarry Smith */ 135247c6ae99SBarry Smith for (l=l_c*Pf/Pc; l<PetscMin((l_c+1)*Pf/Pc,Pf); l++) { 135347c6ae99SBarry Smith for (j=j_c*Nf/Nc; j<PetscMin((j_c+1)*Nf/Nc,Nf); j++) { 135447c6ae99SBarry Smith for (i=i_c*Mf/Mc; i<PetscMin((i_c+1)*Mf/Mc,Mf); i++) { 135547c6ae99SBarry 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; 135647c6ae99SBarry Smith fn_idx++; 135747c6ae99SBarry Smith } 135847c6ae99SBarry Smith } 135947c6ae99SBarry Smith } 136047c6ae99SBarry Smith /* add all these points to one aggregate */ 136147c6ae99SBarry Smith ierr = MatSetValues(*rest, 1, &a, fn_idx, fine_nodes, one_vec, INSERT_VALUES);CHKERRQ(ierr); 136247c6ae99SBarry Smith } 136347c6ae99SBarry Smith } 136447c6ae99SBarry Smith } 136547c6ae99SBarry Smith } 136647c6ae99SBarry Smith ierr = PetscFree2(one_vec,fine_nodes);CHKERRQ(ierr); 136747c6ae99SBarry Smith ierr = MatAssemblyBegin(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 136847c6ae99SBarry Smith ierr = MatAssemblyEnd(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 136947c6ae99SBarry Smith PetscFunctionReturn(0); 137047c6ae99SBarry Smith } 1371