147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 3aa219208SBarry Smith Code for interpolating between grids represented by DMDAs 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 6*85afcc9aSBarry Smith #define NEWVERSION 1 7*85afcc9aSBarry Smith 8c6db04a5SJed Brown #include <private/daimpl.h> /*I "petscdmda.h" I*/ 9c6db04a5SJed Brown #include <petscpcmg.h> 1047c6ae99SBarry Smith 1147c6ae99SBarry Smith #undef __FUNCT__ 1247c6ae99SBarry Smith #define __FUNCT__ "DMGetInterpolationScale" 1347c6ae99SBarry Smith /*@ 1447c6ae99SBarry Smith DMGetInterpolationScale - Forms L = R*1/diag(R*1) - L.*v is like a coarse grid average of the 1547c6ae99SBarry Smith nearby fine grid points. 1647c6ae99SBarry Smith 1747c6ae99SBarry Smith Input Parameters: 1847c6ae99SBarry Smith + dac - DM that defines a coarse mesh 1947c6ae99SBarry Smith . daf - DM that defines a fine mesh 2047c6ae99SBarry Smith - mat - the restriction (or interpolation operator) from fine to coarse 2147c6ae99SBarry Smith 2247c6ae99SBarry Smith Output Parameter: 2347c6ae99SBarry Smith . scale - the scaled vector 2447c6ae99SBarry Smith 2547c6ae99SBarry Smith Level: developer 2647c6ae99SBarry Smith 2747c6ae99SBarry Smith .seealso: DMGetInterpolation() 2847c6ae99SBarry Smith 2947c6ae99SBarry Smith @*/ 307087cfbeSBarry Smith PetscErrorCode DMGetInterpolationScale(DM dac,DM daf,Mat mat,Vec *scale) 3147c6ae99SBarry Smith { 3247c6ae99SBarry Smith PetscErrorCode ierr; 3347c6ae99SBarry Smith Vec fine; 3447c6ae99SBarry Smith PetscScalar one = 1.0; 3547c6ae99SBarry Smith 3647c6ae99SBarry Smith PetscFunctionBegin; 3747c6ae99SBarry Smith ierr = DMCreateGlobalVector(daf,&fine);CHKERRQ(ierr); 3847c6ae99SBarry Smith ierr = DMCreateGlobalVector(dac,scale);CHKERRQ(ierr); 3947c6ae99SBarry Smith ierr = VecSet(fine,one);CHKERRQ(ierr); 4047c6ae99SBarry Smith ierr = MatRestrict(mat,fine,*scale);CHKERRQ(ierr); 41fcfd50ebSBarry Smith ierr = VecDestroy(&fine);CHKERRQ(ierr); 4247c6ae99SBarry Smith ierr = VecReciprocal(*scale);CHKERRQ(ierr); 4347c6ae99SBarry Smith PetscFunctionReturn(0); 4447c6ae99SBarry Smith } 4547c6ae99SBarry Smith 4647c6ae99SBarry Smith #undef __FUNCT__ 4794013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_1D_Q1" 4894013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_1D_Q1(DM dac,DM daf,Mat *A) 4947c6ae99SBarry Smith { 5047c6ae99SBarry Smith PetscErrorCode ierr; 5147c6ae99SBarry Smith PetscInt i,i_start,m_f,Mx,*idx_f; 5247c6ae99SBarry Smith PetscInt m_ghost,*idx_c,m_ghost_c; 5347c6ae99SBarry Smith PetscInt row,col,i_start_ghost,mx,m_c,nc,ratio; 5447c6ae99SBarry Smith PetscInt i_c,i_start_c,i_start_ghost_c,cols[2],dof; 55*85afcc9aSBarry Smith PetscScalar v[2],x; 5647c6ae99SBarry Smith Mat mat; 571321219cSEthan Coon DMDABoundaryType bx; 5847c6ae99SBarry Smith 5947c6ae99SBarry Smith PetscFunctionBegin; 601321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr); 611321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 621321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC){ 6347c6ae99SBarry Smith ratio = mx/Mx; 6447c6ae99SBarry Smith if (ratio*Mx != mx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: mx/Mx must be integer: mx %D Mx %D",mx,Mx); 6547c6ae99SBarry Smith } else { 6647c6ae99SBarry Smith ratio = (mx-1)/(Mx-1); 6747c6ae99SBarry Smith if (ratio*(Mx-1) != mx-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %D Mx %D",mx,Mx); 6847c6ae99SBarry Smith } 6947c6ae99SBarry Smith 70aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); 71aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); 72aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 7347c6ae99SBarry Smith 74aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr); 75aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr); 76aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 7747c6ae99SBarry Smith 7847c6ae99SBarry Smith /* create interpolation matrix */ 7947c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr); 8047c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr); 8147c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 8247c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,2,PETSC_NULL);CHKERRQ(ierr); 8347c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,2,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 8447c6ae99SBarry Smith 8547c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 86*85afcc9aSBarry Smith if (!NEWVERSION) { 87b3bd94feSDave May 8847c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 8947c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 9047c6ae99SBarry Smith row = idx_f[dof*(i-i_start_ghost)]/dof; 9147c6ae99SBarry Smith 9247c6ae99SBarry Smith i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 93aa219208SBarry 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\ 9447c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 9547c6ae99SBarry Smith 9647c6ae99SBarry Smith /* 9747c6ae99SBarry Smith Only include those interpolation points that are truly 9847c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 9947c6ae99SBarry Smith in x direction; since they have no right neighbor 10047c6ae99SBarry Smith */ 10147c6ae99SBarry Smith x = ((double)(i - i_c*ratio))/((double)ratio); 10247c6ae99SBarry Smith nc = 0; 10347c6ae99SBarry Smith /* one left and below; or we are right on it */ 10447c6ae99SBarry Smith col = dof*(i_c-i_start_ghost_c); 10547c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 10647c6ae99SBarry Smith v[nc++] = - x + 1.0; 10747c6ae99SBarry Smith /* one right? */ 10847c6ae99SBarry Smith if (i_c*ratio != i) { 10947c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 11047c6ae99SBarry Smith v[nc++] = x; 11147c6ae99SBarry Smith } 11247c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 11347c6ae99SBarry Smith } 114b3bd94feSDave May 115b3bd94feSDave May } else { 116b3bd94feSDave May PetscScalar *xi; 117b3bd94feSDave May PetscInt li,nxi,n; 118b3bd94feSDave May PetscScalar Ni[2]; 119b3bd94feSDave May 120b3bd94feSDave May /* compute local coordinate arrays */ 121b3bd94feSDave May nxi = ratio + 1; 122b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr); 123b3bd94feSDave May for (li=0; li<nxi; li++) { 12452218ef2SJed Brown xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1)); 125b3bd94feSDave May } 126b3bd94feSDave May 127b3bd94feSDave May for (i=i_start; i<i_start+m_f; i++) { 128b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */ 129b3bd94feSDave May row = idx_f[dof*(i-i_start_ghost)]/dof; 130b3bd94feSDave May 131b3bd94feSDave May i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 132b3bd94feSDave 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\ 133b3bd94feSDave May i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 134b3bd94feSDave May 135b3bd94feSDave May /* remainders */ 136b3bd94feSDave May li = i - ratio * (i/ratio); 137b3bd94feSDave May if (i==mx-1){ li = nxi-1; } 138b3bd94feSDave May 139b3bd94feSDave May /* corners */ 140b3bd94feSDave May col = dof*(i_c-i_start_ghost_c); 141b3bd94feSDave May cols[0] = idx_c[col]/dof; 142b3bd94feSDave May Ni[0] = 1.0; 143b3bd94feSDave May if ( (li==0) || (li==nxi-1) ) { 144b3bd94feSDave May ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr); 145b3bd94feSDave May continue; 146b3bd94feSDave May } 147b3bd94feSDave May 148b3bd94feSDave May /* edges + interior */ 149b3bd94feSDave May /* remainders */ 150b3bd94feSDave May if (i==mx-1){ i_c--; } 151b3bd94feSDave May 152b3bd94feSDave May col = dof*(i_c-i_start_ghost_c); 153b3bd94feSDave May cols[0] = idx_c[col]/dof; /* one left and below; or we are right on it */ 154b3bd94feSDave May cols[1] = idx_c[col+dof]/dof; 155b3bd94feSDave May 156b3bd94feSDave May Ni[0] = 0.5*(1.0-xi[li]); 157b3bd94feSDave May Ni[1] = 0.5*(1.0+xi[li]); 158b3bd94feSDave May for (n=0; n<2; n++) { 159b3bd94feSDave May if( PetscAbsScalar(Ni[n])<1.0e-32) { cols[n]=-1; } 160b3bd94feSDave May } 161b3bd94feSDave May ierr = MatSetValues(mat,1,&row,2,cols,Ni,INSERT_VALUES);CHKERRQ(ierr); 162b3bd94feSDave May } 163b3bd94feSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 164b3bd94feSDave May } 16547c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 16647c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 16747c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 168fcfd50ebSBarry Smith ierr = MatDestroy(&mat);CHKERRQ(ierr); 16947c6ae99SBarry Smith PetscFunctionReturn(0); 17047c6ae99SBarry Smith } 17147c6ae99SBarry Smith 17247c6ae99SBarry Smith #undef __FUNCT__ 17394013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_1D_Q0" 17494013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_1D_Q0(DM dac,DM daf,Mat *A) 17547c6ae99SBarry Smith { 17647c6ae99SBarry Smith PetscErrorCode ierr; 17747c6ae99SBarry Smith PetscInt i,i_start,m_f,Mx,*idx_f; 17847c6ae99SBarry Smith PetscInt m_ghost,*idx_c,m_ghost_c; 17947c6ae99SBarry Smith PetscInt row,col,i_start_ghost,mx,m_c,nc,ratio; 18047c6ae99SBarry Smith PetscInt i_c,i_start_c,i_start_ghost_c,cols[2],dof; 18147c6ae99SBarry Smith PetscScalar v[2],x; 18247c6ae99SBarry Smith Mat mat; 1831321219cSEthan Coon DMDABoundaryType bx; 18447c6ae99SBarry Smith 18547c6ae99SBarry Smith PetscFunctionBegin; 1861321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr); 1871321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 1881321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC){ 18947c6ae99SBarry Smith ratio = mx/Mx; 19047c6ae99SBarry 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); 19147c6ae99SBarry Smith } else { 19247c6ae99SBarry Smith ratio = (mx-1)/(Mx-1); 19347c6ae99SBarry 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); 19447c6ae99SBarry Smith } 19547c6ae99SBarry Smith 196aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); 197aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); 198aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 19947c6ae99SBarry Smith 200aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr); 201aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr); 202aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 20347c6ae99SBarry Smith 20447c6ae99SBarry Smith /* create interpolation matrix */ 20547c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr); 20647c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr); 20747c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 20847c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,2,PETSC_NULL);CHKERRQ(ierr); 20947c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,2,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 21047c6ae99SBarry Smith 21147c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 21247c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 21347c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 21447c6ae99SBarry Smith row = idx_f[dof*(i-i_start_ghost)]/dof; 21547c6ae99SBarry Smith 21647c6ae99SBarry Smith i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 21747c6ae99SBarry Smith 21847c6ae99SBarry Smith /* 21947c6ae99SBarry Smith Only include those interpolation points that are truly 22047c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 22147c6ae99SBarry Smith in x direction; since they have no right neighbor 22247c6ae99SBarry Smith */ 22347c6ae99SBarry Smith x = ((double)(i - i_c*ratio))/((double)ratio); 22447c6ae99SBarry Smith nc = 0; 22547c6ae99SBarry Smith /* one left and below; or we are right on it */ 22647c6ae99SBarry Smith col = dof*(i_c-i_start_ghost_c); 22747c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 22847c6ae99SBarry Smith v[nc++] = - x + 1.0; 22947c6ae99SBarry Smith /* one right? */ 23047c6ae99SBarry Smith if (i_c*ratio != i) { 23147c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 23247c6ae99SBarry Smith v[nc++] = x; 23347c6ae99SBarry Smith } 23447c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 23547c6ae99SBarry Smith } 23647c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 23747c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 23847c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 239fcfd50ebSBarry Smith ierr = MatDestroy(&mat);CHKERRQ(ierr); 24047c6ae99SBarry Smith ierr = PetscLogFlops(5.0*m_f);CHKERRQ(ierr); 24147c6ae99SBarry Smith PetscFunctionReturn(0); 24247c6ae99SBarry Smith } 24347c6ae99SBarry Smith 24447c6ae99SBarry Smith #undef __FUNCT__ 24594013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_2D_Q1" 24694013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_2D_Q1(DM dac,DM daf,Mat *A) 24747c6ae99SBarry Smith { 24847c6ae99SBarry Smith PetscErrorCode ierr; 24947c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 25047c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz; 25147c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj; 25247c6ae99SBarry 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; 25347c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 25447c6ae99SBarry Smith PetscScalar v[4],x,y; 25547c6ae99SBarry Smith Mat mat; 2561321219cSEthan Coon DMDABoundaryType bx,by; 25747c6ae99SBarry Smith 25847c6ae99SBarry Smith PetscFunctionBegin; 2591321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr); 2601321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 2611321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC){ 26247c6ae99SBarry Smith ratioi = mx/Mx; 26347c6ae99SBarry 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); 26447c6ae99SBarry Smith } else { 26547c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 26647c6ae99SBarry 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); 26747c6ae99SBarry Smith } 2681321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC){ 26947c6ae99SBarry Smith ratioj = my/My; 27047c6ae99SBarry 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); 27147c6ae99SBarry Smith } else { 27247c6ae99SBarry Smith ratioj = (my-1)/(My-1); 27347c6ae99SBarry 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); 27447c6ae99SBarry Smith } 27547c6ae99SBarry Smith 27647c6ae99SBarry Smith 277aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 278aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 279aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 28047c6ae99SBarry Smith 281aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 282aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 283aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 28447c6ae99SBarry Smith 28547c6ae99SBarry Smith /* 286aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 28747c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 28847c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 28947c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 29047c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 29147c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 29247c6ae99SBarry Smith 29347c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 29447c6ae99SBarry Smith */ 29547c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr); 29647c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr); 29747c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr); 29847c6ae99SBarry Smith col_scale = size_f/size_c; 29947c6ae99SBarry Smith col_shift = Mx*My*(rank_f/size_c); 30047c6ae99SBarry Smith 30147c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr); 30247c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 30347c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 30447c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 30547c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 30647c6ae99SBarry Smith 30747c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 30847c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 30947c6ae99SBarry Smith 310aa219208SBarry 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\ 31147c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 312aa219208SBarry 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\ 31347c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 31447c6ae99SBarry Smith 31547c6ae99SBarry Smith /* 31647c6ae99SBarry Smith Only include those interpolation points that are truly 31747c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 31847c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 31947c6ae99SBarry Smith */ 32047c6ae99SBarry Smith nc = 0; 32147c6ae99SBarry Smith /* one left and below; or we are right on it */ 32247c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 32347c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 32447c6ae99SBarry Smith /* one right and below */ 32547c6ae99SBarry Smith if (i_c*ratioi != i) { 32647c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col+dof]/dof; 32747c6ae99SBarry Smith } 32847c6ae99SBarry Smith /* one left and above */ 32947c6ae99SBarry Smith if (j_c*ratioj != j) { 33047c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col+m_ghost_c*dof]/dof; 33147c6ae99SBarry Smith } 33247c6ae99SBarry Smith /* one right and above */ 3330c82216cSJed Brown if (i_c*ratioi != i && j_c*ratioj != j) { 33447c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; 33547c6ae99SBarry Smith } 33647c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 33747c6ae99SBarry Smith } 33847c6ae99SBarry Smith } 33947c6ae99SBarry Smith ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr); 34047c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr); 34147c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 34247c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 34347c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 34447c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 34547c6ae99SBarry Smith 34647c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 347*85afcc9aSBarry Smith if (!NEWVERSION) { 348b3bd94feSDave May 34947c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 35047c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 35147c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 35247c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 35347c6ae99SBarry Smith 35447c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 35547c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 35647c6ae99SBarry Smith 35747c6ae99SBarry Smith /* 35847c6ae99SBarry Smith Only include those interpolation points that are truly 35947c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 36047c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 36147c6ae99SBarry Smith */ 36247c6ae99SBarry Smith x = ((double)(i - i_c*ratioi))/((double)ratioi); 36347c6ae99SBarry Smith y = ((double)(j - j_c*ratioj))/((double)ratioj); 364b3bd94feSDave May 36547c6ae99SBarry Smith nc = 0; 36647c6ae99SBarry Smith /* one left and below; or we are right on it */ 36747c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 36847c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 36947c6ae99SBarry Smith v[nc++] = x*y - x - y + 1.0; 37047c6ae99SBarry Smith /* one right and below */ 37147c6ae99SBarry Smith if (i_c*ratioi != i) { 37247c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+dof]/dof; 37347c6ae99SBarry Smith v[nc++] = -x*y + x; 37447c6ae99SBarry Smith } 37547c6ae99SBarry Smith /* one left and above */ 37647c6ae99SBarry Smith if (j_c*ratioj != j) { 37747c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+m_ghost_c*dof]/dof; 37847c6ae99SBarry Smith v[nc++] = -x*y + y; 37947c6ae99SBarry Smith } 38047c6ae99SBarry Smith /* one right and above */ 38147c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 38247c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; 38347c6ae99SBarry Smith v[nc++] = x*y; 38447c6ae99SBarry Smith } 38547c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 38647c6ae99SBarry Smith } 38747c6ae99SBarry Smith } 388b3bd94feSDave May 389b3bd94feSDave May } else { 390b3bd94feSDave May PetscScalar Ni[4]; 391b3bd94feSDave May PetscScalar *xi,*eta; 392b3bd94feSDave May PetscInt li,nxi,lj,neta; 393b3bd94feSDave May 394b3bd94feSDave May /* compute local coordinate arrays */ 395b3bd94feSDave May nxi = ratioi + 1; 396b3bd94feSDave May neta = ratioj + 1; 397b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr); 398b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*neta,&eta);CHKERRQ(ierr); 399b3bd94feSDave May for (li=0; li<nxi; li++) { 40052218ef2SJed Brown xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1)); 401b3bd94feSDave May } 402b3bd94feSDave May for (lj=0; lj<neta; lj++) { 40352218ef2SJed Brown eta[lj] = -1.0 + (PetscScalar)lj*(2.0/(PetscScalar)(neta-1)); 404b3bd94feSDave May } 405b3bd94feSDave May 406b3bd94feSDave May /* loop over local fine grid nodes setting interpolation for those*/ 407b3bd94feSDave May for (j=j_start; j<j_start+n_f; j++) { 408b3bd94feSDave May for (i=i_start; i<i_start+m_f; i++) { 409b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */ 410b3bd94feSDave May row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 411b3bd94feSDave May 412b3bd94feSDave May i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 413b3bd94feSDave May j_c = (j/ratioj); /* coarse grid node below fine grid node */ 414b3bd94feSDave May 415b3bd94feSDave May /* remainders */ 416b3bd94feSDave May li = i - ratioi * (i/ratioi); 417b3bd94feSDave May if (i==mx-1){ li = nxi-1; } 418b3bd94feSDave May lj = j - ratioj * (j/ratioj); 419b3bd94feSDave May if (j==my-1){ lj = neta-1; } 420b3bd94feSDave May 421b3bd94feSDave May /* corners */ 422b3bd94feSDave May col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 423b3bd94feSDave May cols[0] = col_shift + idx_c[col]/dof; /* left, below */ 424b3bd94feSDave May Ni[0] = 1.0; 425b3bd94feSDave May if ( (li==0) || (li==nxi-1) ) { 426b3bd94feSDave May if ( (lj==0) || (lj==neta-1) ) { 427b3bd94feSDave May ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr); 428b3bd94feSDave May continue; 429b3bd94feSDave May } 430b3bd94feSDave May } 431b3bd94feSDave May 432b3bd94feSDave May /* edges + interior */ 433b3bd94feSDave May /* remainders */ 434b3bd94feSDave May if (i==mx-1){ i_c--; } 435b3bd94feSDave May if (j==my-1){ j_c--; } 436b3bd94feSDave May 437b3bd94feSDave May col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 438b3bd94feSDave May cols[0] = col_shift + idx_c[col]/dof; /* left, below */ 439b3bd94feSDave May cols[1] = col_shift + idx_c[col+dof]/dof; /* right, below */ 440b3bd94feSDave May cols[2] = col_shift + idx_c[col+m_ghost_c*dof]/dof; /* left, above */ 441b3bd94feSDave May cols[3] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; /* right, above */ 442b3bd94feSDave May 443b3bd94feSDave May Ni[0] = 0.25*(1.0-xi[li])*(1.0-eta[lj]); 444b3bd94feSDave May Ni[1] = 0.25*(1.0+xi[li])*(1.0-eta[lj]); 445b3bd94feSDave May Ni[2] = 0.25*(1.0-xi[li])*(1.0+eta[lj]); 446b3bd94feSDave May Ni[3] = 0.25*(1.0+xi[li])*(1.0+eta[lj]); 447b3bd94feSDave May 448b3bd94feSDave May nc = 0; 449b3bd94feSDave May if( PetscAbsScalar(Ni[0])<1.0e-32) { cols[0]=-1; } 450b3bd94feSDave May if( PetscAbsScalar(Ni[1])<1.0e-32) { cols[1]=-1; } 451b3bd94feSDave May if( PetscAbsScalar(Ni[2])<1.0e-32) { cols[2]=-1; } 452b3bd94feSDave May if( PetscAbsScalar(Ni[3])<1.0e-32) { cols[3]=-1; } 453b3bd94feSDave May 454b3bd94feSDave May ierr = MatSetValues(mat,1,&row,4,cols,Ni,INSERT_VALUES);CHKERRQ(ierr); 455b3bd94feSDave May } 456b3bd94feSDave May } 457b3bd94feSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 458b3bd94feSDave May ierr = PetscFree(eta);CHKERRQ(ierr); 459b3bd94feSDave May } 46047c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 46147c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 46247c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 463fcfd50ebSBarry Smith ierr = MatDestroy(&mat);CHKERRQ(ierr); 46447c6ae99SBarry Smith PetscFunctionReturn(0); 46547c6ae99SBarry Smith } 46647c6ae99SBarry Smith 46747c6ae99SBarry Smith /* 46847c6ae99SBarry Smith Contributed by Andrei Draganescu <aidraga@sandia.gov> 46947c6ae99SBarry Smith */ 47047c6ae99SBarry Smith #undef __FUNCT__ 47194013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_2D_Q0" 47294013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_2D_Q0(DM dac,DM daf,Mat *A) 47347c6ae99SBarry Smith { 47447c6ae99SBarry Smith PetscErrorCode ierr; 47547c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 47647c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz; 47747c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj; 47847c6ae99SBarry 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; 47947c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 48047c6ae99SBarry Smith PetscScalar v[4]; 48147c6ae99SBarry Smith Mat mat; 4821321219cSEthan Coon DMDABoundaryType bx,by; 48347c6ae99SBarry Smith 48447c6ae99SBarry Smith PetscFunctionBegin; 4851321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr); 4861321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 4871321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x"); 4881321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y"); 48947c6ae99SBarry Smith ratioi = mx/Mx; 49047c6ae99SBarry Smith ratioj = my/My; 49147c6ae99SBarry Smith if (ratioi*Mx != mx) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in x"); 49247c6ae99SBarry Smith if (ratioj*My != my) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in y"); 49347c6ae99SBarry Smith if (ratioi != 2) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 2"); 49447c6ae99SBarry Smith if (ratioj != 2) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 2"); 49547c6ae99SBarry Smith 496aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 497aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 498aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 49947c6ae99SBarry Smith 500aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 501aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 502aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 50347c6ae99SBarry Smith 50447c6ae99SBarry Smith /* 505aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 50647c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 50747c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 50847c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 50947c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 51047c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 51147c6ae99SBarry Smith 51247c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 51347c6ae99SBarry Smith */ 51447c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr); 51547c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr); 51647c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr); 51747c6ae99SBarry Smith col_scale = size_f/size_c; 51847c6ae99SBarry Smith col_shift = Mx*My*(rank_f/size_c); 51947c6ae99SBarry Smith 52047c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr); 52147c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 52247c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 52347c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 52447c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 52547c6ae99SBarry Smith 52647c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 52747c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 52847c6ae99SBarry Smith 529aa219208SBarry 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\ 53047c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 531aa219208SBarry 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\ 53247c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 53347c6ae99SBarry Smith 53447c6ae99SBarry Smith /* 53547c6ae99SBarry Smith Only include those interpolation points that are truly 53647c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 53747c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 53847c6ae99SBarry Smith */ 53947c6ae99SBarry Smith nc = 0; 54047c6ae99SBarry Smith /* one left and below; or we are right on it */ 54147c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 54247c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 54347c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 54447c6ae99SBarry Smith } 54547c6ae99SBarry Smith } 54647c6ae99SBarry Smith ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr); 54747c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr); 54847c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 54947c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 55047c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 55147c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 55247c6ae99SBarry Smith 55347c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 55447c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 55547c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 55647c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 55747c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 55847c6ae99SBarry Smith 55947c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 56047c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 56147c6ae99SBarry Smith nc = 0; 56247c6ae99SBarry Smith /* one left and below; or we are right on it */ 56347c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 56447c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 56547c6ae99SBarry Smith v[nc++] = 1.0; 56647c6ae99SBarry Smith 56747c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 56847c6ae99SBarry Smith } 56947c6ae99SBarry Smith } 57047c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 57147c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 57247c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 573fcfd50ebSBarry Smith ierr = MatDestroy(&mat);CHKERRQ(ierr); 57447c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f);CHKERRQ(ierr); 57547c6ae99SBarry Smith PetscFunctionReturn(0); 57647c6ae99SBarry Smith } 57747c6ae99SBarry Smith 57847c6ae99SBarry Smith /* 57947c6ae99SBarry Smith Contributed by Jianming Yang <jianming-yang@uiowa.edu> 58047c6ae99SBarry Smith */ 58147c6ae99SBarry Smith #undef __FUNCT__ 58294013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_3D_Q0" 58394013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_3D_Q0(DM dac,DM daf,Mat *A) 58447c6ae99SBarry Smith { 58547c6ae99SBarry Smith PetscErrorCode ierr; 58647c6ae99SBarry Smith PetscInt i,j,l,i_start,j_start,l_start,m_f,n_f,p_f,Mx,My,Mz,*idx_f,dof; 58747c6ae99SBarry Smith PetscInt m_ghost,n_ghost,p_ghost,*idx_c,m_ghost_c,n_ghost_c,p_ghost_c,nc,*dnz,*onz; 58847c6ae99SBarry 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; 58947c6ae99SBarry 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; 59047c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 59147c6ae99SBarry Smith PetscScalar v[8]; 59247c6ae99SBarry Smith Mat mat; 5931321219cSEthan Coon DMDABoundaryType bx,by,bz; 59447c6ae99SBarry Smith 59547c6ae99SBarry Smith PetscFunctionBegin; 5961321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr); 5971321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 5981321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x"); 5991321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y"); 6001321219cSEthan Coon if (bz == DMDA_BOUNDARY_PERIODIC) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in z"); 60147c6ae99SBarry Smith ratioi = mx/Mx; 60247c6ae99SBarry Smith ratioj = my/My; 60347c6ae99SBarry Smith ratiol = mz/Mz; 60447c6ae99SBarry Smith if (ratioi*Mx != mx) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in x"); 60547c6ae99SBarry Smith if (ratioj*My != my) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in y"); 60647c6ae99SBarry Smith if (ratiol*Mz != mz) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in z"); 60747c6ae99SBarry Smith if (ratioi != 2 && ratioi != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 1 or 2"); 60847c6ae99SBarry Smith if (ratioj != 2 && ratioj != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 1 or 2"); 60947c6ae99SBarry Smith if (ratiol != 2 && ratiol != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in z must be 1 or 2"); 61047c6ae99SBarry Smith 611aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 612aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 613aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 61447c6ae99SBarry Smith 615aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 616aa219208SBarry 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); 617aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 61847c6ae99SBarry Smith /* 619aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 62047c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 62147c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 62247c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 62347c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 62447c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 62547c6ae99SBarry Smith 62647c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 62747c6ae99SBarry Smith */ 62847c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr); 62947c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr); 63047c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr); 63147c6ae99SBarry Smith col_scale = size_f/size_c; 63247c6ae99SBarry Smith col_shift = Mx*My*Mz*(rank_f/size_c); 63347c6ae99SBarry Smith 63447c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f*p_f,col_scale*m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr); 63547c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 63647c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 63747c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 63847c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 63947c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 64047c6ae99SBarry Smith 64147c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 64247c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 64347c6ae99SBarry Smith l_c = (l/ratiol); 64447c6ae99SBarry Smith 645aa219208SBarry 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\ 64647c6ae99SBarry Smith l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c); 647aa219208SBarry 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\ 64847c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 649aa219208SBarry 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\ 65047c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 65147c6ae99SBarry Smith 65247c6ae99SBarry Smith /* 65347c6ae99SBarry Smith Only include those interpolation points that are truly 65447c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 65547c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 65647c6ae99SBarry Smith */ 65747c6ae99SBarry Smith nc = 0; 65847c6ae99SBarry Smith /* one left and below; or we are right on it */ 65947c6ae99SBarry 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)); 66047c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 66147c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 66247c6ae99SBarry Smith } 66347c6ae99SBarry Smith } 66447c6ae99SBarry Smith } 66547c6ae99SBarry Smith ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr); 66647c6ae99SBarry 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); 66747c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 66847c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 66947c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 67047c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 67147c6ae99SBarry Smith 67247c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 67347c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 67447c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 67547c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 67647c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 67747c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 67847c6ae99SBarry Smith 67947c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 68047c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 68147c6ae99SBarry Smith l_c = (l/ratiol); 68247c6ae99SBarry Smith nc = 0; 68347c6ae99SBarry Smith /* one left and below; or we are right on it */ 68447c6ae99SBarry 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)); 68547c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 68647c6ae99SBarry Smith v[nc++] = 1.0; 68747c6ae99SBarry Smith 68847c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 68947c6ae99SBarry Smith } 69047c6ae99SBarry Smith } 69147c6ae99SBarry Smith } 69247c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 69347c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 69447c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 695fcfd50ebSBarry Smith ierr = MatDestroy(&mat);CHKERRQ(ierr); 69647c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f*p_f);CHKERRQ(ierr); 69747c6ae99SBarry Smith PetscFunctionReturn(0); 69847c6ae99SBarry Smith } 69947c6ae99SBarry Smith 70047c6ae99SBarry Smith #undef __FUNCT__ 70194013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_3D_Q1" 70294013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_3D_Q1(DM dac,DM daf,Mat *A) 70347c6ae99SBarry Smith { 70447c6ae99SBarry Smith PetscErrorCode ierr; 70547c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof,l; 70647c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,Mz,mz; 70747c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[8],mx,m_c,my,nc,ratioi,ratioj,ratiok; 70847c6ae99SBarry Smith PetscInt i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c; 70947c6ae99SBarry Smith PetscInt l_start,p_f,l_start_ghost,p_ghost,l_start_c,p_c; 71047c6ae99SBarry Smith PetscInt l_start_ghost_c,p_ghost_c,l_c,*dnz,*onz; 71147c6ae99SBarry Smith PetscScalar v[8],x,y,z; 71247c6ae99SBarry Smith Mat mat; 7131321219cSEthan Coon DMDABoundaryType bx,by,bz; 71447c6ae99SBarry Smith 71547c6ae99SBarry Smith PetscFunctionBegin; 7161321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr); 7171321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 71847c6ae99SBarry Smith if (mx == Mx) { 71947c6ae99SBarry Smith ratioi = 1; 7201321219cSEthan Coon } else if (bx == DMDA_BOUNDARY_PERIODIC) { 72147c6ae99SBarry Smith ratioi = mx/Mx; 72247c6ae99SBarry 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); 72347c6ae99SBarry Smith } else { 72447c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 72547c6ae99SBarry 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); 72647c6ae99SBarry Smith } 72747c6ae99SBarry Smith if (my == My) { 72847c6ae99SBarry Smith ratioj = 1; 7291321219cSEthan Coon } else if (by == DMDA_BOUNDARY_PERIODIC) { 73047c6ae99SBarry Smith ratioj = my/My; 73147c6ae99SBarry 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); 73247c6ae99SBarry Smith } else { 73347c6ae99SBarry Smith ratioj = (my-1)/(My-1); 73447c6ae99SBarry 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); 73547c6ae99SBarry Smith } 73647c6ae99SBarry Smith if (mz == Mz) { 73747c6ae99SBarry Smith ratiok = 1; 7381321219cSEthan Coon } else if (bz == DMDA_BOUNDARY_PERIODIC) { 73947c6ae99SBarry Smith ratiok = mz/Mz; 74047c6ae99SBarry 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); 74147c6ae99SBarry Smith } else { 74247c6ae99SBarry Smith ratiok = (mz-1)/(Mz-1); 74347c6ae99SBarry 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); 74447c6ae99SBarry Smith } 74547c6ae99SBarry Smith 746aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 747aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 748aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 74947c6ae99SBarry Smith 750aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 751aa219208SBarry 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); 752aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 75347c6ae99SBarry Smith 75447c6ae99SBarry Smith /* create interpolation matrix, determining exact preallocation */ 75547c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)dac)->comm,m_f*n_f*p_f,m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr); 75647c6ae99SBarry Smith /* loop over local fine grid nodes counting interpolating points */ 75747c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 75847c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 75947c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 76047c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 76147c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 76247c6ae99SBarry Smith i_c = (i/ratioi); 76347c6ae99SBarry Smith j_c = (j/ratioj); 76447c6ae99SBarry Smith l_c = (l/ratiok); 765aa219208SBarry 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\ 76647c6ae99SBarry Smith l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c); 767aa219208SBarry 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\ 76847c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 769aa219208SBarry 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\ 77047c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 77147c6ae99SBarry Smith 77247c6ae99SBarry Smith /* 77347c6ae99SBarry Smith Only include those interpolation points that are truly 77447c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 77547c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 77647c6ae99SBarry Smith */ 77747c6ae99SBarry Smith nc = 0; 77847c6ae99SBarry 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)); 77947c6ae99SBarry Smith cols[nc++] = idx_c[col]/dof; 78047c6ae99SBarry Smith if (i_c*ratioi != i) { 78147c6ae99SBarry Smith cols[nc++] = idx_c[col+dof]/dof; 78247c6ae99SBarry Smith } 78347c6ae99SBarry Smith if (j_c*ratioj != j) { 78447c6ae99SBarry Smith cols[nc++] = idx_c[col+m_ghost_c*dof]/dof; 78547c6ae99SBarry Smith } 78647c6ae99SBarry Smith if (l_c*ratiok != l) { 78747c6ae99SBarry Smith cols[nc++] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; 78847c6ae99SBarry Smith } 78947c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 79047c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c+1)*dof]/dof; 79147c6ae99SBarry Smith } 79247c6ae99SBarry Smith if (j_c*ratioj != j && l_c*ratiok != l) { 79347c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof; 79447c6ae99SBarry Smith } 79547c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l) { 79647c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; 79747c6ae99SBarry Smith } 79847c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) { 79947c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; 80047c6ae99SBarry Smith } 80147c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 80247c6ae99SBarry Smith } 80347c6ae99SBarry Smith } 80447c6ae99SBarry Smith } 80547c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr); 80647c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f*p_f,m_c*n_c*p_c,mx*my*mz,Mx*My*Mz);CHKERRQ(ierr); 80747c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 80847c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 80947c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 81047c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 81147c6ae99SBarry Smith 81247c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 813*85afcc9aSBarry Smith if (NEWVERSION) { 814b3bd94feSDave May 81547c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 81647c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 81747c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 81847c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 81947c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 82047c6ae99SBarry Smith 82147c6ae99SBarry Smith i_c = (i/ratioi); 82247c6ae99SBarry Smith j_c = (j/ratioj); 82347c6ae99SBarry Smith l_c = (l/ratiok); 82447c6ae99SBarry Smith 82547c6ae99SBarry Smith /* 82647c6ae99SBarry Smith Only include those interpolation points that are truly 82747c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 82847c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 82947c6ae99SBarry Smith */ 83047c6ae99SBarry Smith x = ((double)(i - i_c*ratioi))/((double)ratioi); 83147c6ae99SBarry Smith y = ((double)(j - j_c*ratioj))/((double)ratioj); 83247c6ae99SBarry Smith z = ((double)(l - l_c*ratiok))/((double)ratiok); 833b3bd94feSDave May 83447c6ae99SBarry Smith nc = 0; 83547c6ae99SBarry Smith /* one left and below; or we are right on it */ 83647c6ae99SBarry 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)); 83747c6ae99SBarry Smith 83847c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 83947c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.)); 84047c6ae99SBarry Smith 84147c6ae99SBarry Smith if (i_c*ratioi != i) { 84247c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 84347c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.)); 84447c6ae99SBarry Smith } 84547c6ae99SBarry Smith 84647c6ae99SBarry Smith if (j_c*ratioj != j) { 84747c6ae99SBarry Smith cols[nc] = idx_c[col+m_ghost_c*dof]/dof; 84847c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.)); 84947c6ae99SBarry Smith } 85047c6ae99SBarry Smith 85147c6ae99SBarry Smith if (l_c*ratiok != l) { 85247c6ae99SBarry Smith cols[nc] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; 85347c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.)); 85447c6ae99SBarry Smith } 85547c6ae99SBarry Smith 85647c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 85747c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c+1)*dof]/dof; 85847c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.)); 85947c6ae99SBarry Smith } 86047c6ae99SBarry Smith 86147c6ae99SBarry Smith if (j_c*ratioj != j && l_c*ratiok != l) { 86247c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof; 86347c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.)); 86447c6ae99SBarry Smith } 86547c6ae99SBarry Smith 86647c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l) { 86747c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; 86847c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.)); 86947c6ae99SBarry Smith } 87047c6ae99SBarry Smith 87147c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) { 87247c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; 87347c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.)); 87447c6ae99SBarry Smith } 87547c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 87647c6ae99SBarry Smith } 87747c6ae99SBarry Smith } 87847c6ae99SBarry Smith } 879b3bd94feSDave May 880b3bd94feSDave May } else { 881b3bd94feSDave May PetscScalar *xi,*eta,*zeta; 882b3bd94feSDave May PetscInt li,nxi,lj,neta,lk,nzeta,n; 883b3bd94feSDave May PetscScalar Ni[8]; 884b3bd94feSDave May 885b3bd94feSDave May /* compute local coordinate arrays */ 886b3bd94feSDave May nxi = ratioi + 1; 887b3bd94feSDave May neta = ratioj + 1; 888b3bd94feSDave May nzeta = ratiok + 1; 889b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr); 890b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*neta,&eta);CHKERRQ(ierr); 891b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*nzeta,&zeta);CHKERRQ(ierr); 892b3bd94feSDave May for (li=0; li<nxi; li++) { 89352218ef2SJed Brown xi[li] = -1.0 + (PetscScalar)li*(2.0/(PetscScalar)(nxi-1)); 894b3bd94feSDave May } 895b3bd94feSDave May for (lj=0; lj<neta; lj++) { 89652218ef2SJed Brown eta[lj] = -1.0 + (PetscScalar)lj*(2.0/(PetscScalar)(neta-1)); 897b3bd94feSDave May } 898b3bd94feSDave May for (lk=0; lk<nzeta; lk++) { 89952218ef2SJed Brown zeta[lk] = -1.0 + (PetscScalar)lk*(2.0/(PetscScalar)(nzeta-1)); 900b3bd94feSDave May } 901b3bd94feSDave May 902b3bd94feSDave May for (l=l_start; l<l_start+p_f; l++) { 903b3bd94feSDave May for (j=j_start; j<j_start+n_f; j++) { 904b3bd94feSDave May for (i=i_start; i<i_start+m_f; i++) { 905b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */ 906b3bd94feSDave May row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 907b3bd94feSDave May 908b3bd94feSDave May i_c = (i/ratioi); 909b3bd94feSDave May j_c = (j/ratioj); 910b3bd94feSDave May l_c = (l/ratiok); 911b3bd94feSDave May 912b3bd94feSDave May /* remainders */ 913b3bd94feSDave May li = i - ratioi * (i/ratioi); 914b3bd94feSDave May if (i==mx-1){ li = nxi-1; } 915b3bd94feSDave May lj = j - ratioj * (j/ratioj); 916b3bd94feSDave May if (j==my-1){ lj = neta-1; } 917b3bd94feSDave May lk = l - ratiok * (l/ratiok); 918b3bd94feSDave May if (l==mz-1){ lk = nzeta-1; } 919b3bd94feSDave May 920b3bd94feSDave May /* corners */ 921b3bd94feSDave 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)); 922b3bd94feSDave May cols[0] = idx_c[col]/dof; 923b3bd94feSDave May Ni[0] = 1.0; 924b3bd94feSDave May if ( (li==0) || (li==nxi-1) ) { 925b3bd94feSDave May if ( (lj==0) || (lj==neta-1) ) { 926b3bd94feSDave May if ( (lk==0) || (lk==nzeta-1) ) { 927b3bd94feSDave May ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr); 928b3bd94feSDave May continue; 929b3bd94feSDave May } 930b3bd94feSDave May } 931b3bd94feSDave May } 932b3bd94feSDave May 933b3bd94feSDave May /* edges + interior */ 934b3bd94feSDave May /* remainders */ 935b3bd94feSDave May if (i==mx-1){ i_c--; } 936b3bd94feSDave May if (j==my-1){ j_c--; } 937b3bd94feSDave May if (l==mz-1){ l_c--; } 938b3bd94feSDave May 939b3bd94feSDave 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)); 940b3bd94feSDave May cols[0] = idx_c[col]/dof; /* one left and below; or we are right on it */ 941b3bd94feSDave May cols[1] = idx_c[col+dof]/dof; /* one right and below */ 942b3bd94feSDave May cols[2] = idx_c[col+m_ghost_c*dof]/dof; /* one left and above */ 943b3bd94feSDave May cols[3] = idx_c[col+(m_ghost_c+1)*dof]/dof; /* one right and above */ 944b3bd94feSDave May 945b3bd94feSDave 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 */ 946b3bd94feSDave May cols[5] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; /* one right and below, and front */ 947b3bd94feSDave May cols[6] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof;/* one left and above and front*/ 948b3bd94feSDave May cols[7] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; /* one right and above and front */ 949b3bd94feSDave May 950b3bd94feSDave May Ni[0] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]); 951b3bd94feSDave May Ni[1] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]); 952b3bd94feSDave May Ni[2] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]); 953b3bd94feSDave May Ni[3] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]); 954b3bd94feSDave May 955b3bd94feSDave May Ni[4] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]); 956b3bd94feSDave May Ni[5] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]); 957b3bd94feSDave May Ni[6] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]); 958b3bd94feSDave May Ni[7] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]); 959b3bd94feSDave May 960b3bd94feSDave May for (n=0; n<8; n++) { 961b3bd94feSDave May if( PetscAbsScalar(Ni[n])<1.0e-32) { cols[n]=-1; } 962b3bd94feSDave May } 963b3bd94feSDave May ierr = MatSetValues(mat,1,&row,8,cols,Ni,INSERT_VALUES);CHKERRQ(ierr); 964b3bd94feSDave May 965b3bd94feSDave May } 966b3bd94feSDave May } 967b3bd94feSDave May } 968b3bd94feSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 969b3bd94feSDave May ierr = PetscFree(eta);CHKERRQ(ierr); 970b3bd94feSDave May ierr = PetscFree(zeta);CHKERRQ(ierr); 971b3bd94feSDave May } 972b3bd94feSDave May 97347c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 97447c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 97547c6ae99SBarry Smith 97647c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 977fcfd50ebSBarry Smith ierr = MatDestroy(&mat);CHKERRQ(ierr); 97847c6ae99SBarry Smith PetscFunctionReturn(0); 97947c6ae99SBarry Smith } 98047c6ae99SBarry Smith 98147c6ae99SBarry Smith #undef __FUNCT__ 98294013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA" 9837087cfbeSBarry Smith PetscErrorCode DMGetInterpolation_DA(DM dac,DM daf,Mat *A,Vec *scale) 98447c6ae99SBarry Smith { 98547c6ae99SBarry Smith PetscErrorCode ierr; 98647c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 9871321219cSEthan Coon DMDABoundaryType bxc,byc,bzc,bxf,byf,bzf; 988aa219208SBarry Smith DMDAStencilType stc,stf; 98947c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data; 99047c6ae99SBarry Smith 99147c6ae99SBarry Smith PetscFunctionBegin; 99247c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 99347c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 99447c6ae99SBarry Smith PetscValidPointer(A,3); 99547c6ae99SBarry Smith if (scale) PetscValidPointer(scale,4); 99647c6ae99SBarry Smith 9971321219cSEthan Coon ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr); 9981321219cSEthan Coon ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr); 999aa219208SBarry Smith if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 1000aa219208SBarry Smith if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); 1001aa219208SBarry Smith if (sc != sf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr); 10021321219cSEthan Coon if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr); 1003aa219208SBarry Smith if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 100447c6ae99SBarry Smith if (Mc < 2 && Mf > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction"); 100547c6ae99SBarry 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"); 100647c6ae99SBarry 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"); 100747c6ae99SBarry Smith 1008aa219208SBarry Smith if (ddc->interptype == DMDA_Q1){ 100947c6ae99SBarry Smith if (dimc == 1){ 101094013140SBarry Smith ierr = DMGetInterpolation_DA_1D_Q1(dac,daf,A);CHKERRQ(ierr); 101147c6ae99SBarry Smith } else if (dimc == 2){ 101294013140SBarry Smith ierr = DMGetInterpolation_DA_2D_Q1(dac,daf,A);CHKERRQ(ierr); 101347c6ae99SBarry Smith } else if (dimc == 3){ 101494013140SBarry Smith ierr = DMGetInterpolation_DA_3D_Q1(dac,daf,A);CHKERRQ(ierr); 1015*85afcc9aSBarry Smith } else SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype); 1016aa219208SBarry Smith } else if (ddc->interptype == DMDA_Q0){ 101747c6ae99SBarry Smith if (dimc == 1){ 101894013140SBarry Smith ierr = DMGetInterpolation_DA_1D_Q0(dac,daf,A);CHKERRQ(ierr); 101947c6ae99SBarry Smith } else if (dimc == 2){ 102094013140SBarry Smith ierr = DMGetInterpolation_DA_2D_Q0(dac,daf,A);CHKERRQ(ierr); 102147c6ae99SBarry Smith } else if (dimc == 3){ 102294013140SBarry Smith ierr = DMGetInterpolation_DA_3D_Q0(dac,daf,A);CHKERRQ(ierr); 1023*85afcc9aSBarry Smith } else SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype); 102447c6ae99SBarry Smith } 102547c6ae99SBarry Smith if (scale) { 102647c6ae99SBarry Smith ierr = DMGetInterpolationScale((DM)dac,(DM)daf,*A,scale);CHKERRQ(ierr); 102747c6ae99SBarry Smith } 102847c6ae99SBarry Smith PetscFunctionReturn(0); 102947c6ae99SBarry Smith } 103047c6ae99SBarry Smith 103147c6ae99SBarry Smith #undef __FUNCT__ 10320bb2b966SJungho Lee #define __FUNCT__ "DMGetInjection_DA_1D" 10330bb2b966SJungho Lee PetscErrorCode DMGetInjection_DA_1D(DM dac,DM daf,VecScatter *inject) 10340bb2b966SJungho Lee { 10350bb2b966SJungho Lee PetscErrorCode ierr; 10360bb2b966SJungho Lee PetscInt i,i_start,m_f,Mx,*idx_f,dof; 10370bb2b966SJungho Lee PetscInt m_ghost,*idx_c,m_ghost_c; 10380bb2b966SJungho Lee PetscInt row,i_start_ghost,mx,m_c,nc,ratioi; 10390bb2b966SJungho Lee PetscInt i_start_c,i_start_ghost_c; 10400bb2b966SJungho Lee PetscInt *cols; 10410bb2b966SJungho Lee DMDABoundaryType bx; 10420bb2b966SJungho Lee Vec vecf,vecc; 10430bb2b966SJungho Lee IS isf; 10440bb2b966SJungho Lee 10450bb2b966SJungho Lee PetscFunctionBegin; 10460bb2b966SJungho Lee ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr); 10470bb2b966SJungho Lee ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 10480bb2b966SJungho Lee if (bx == DMDA_BOUNDARY_PERIODIC) { 10490bb2b966SJungho Lee ratioi = mx/Mx; 10500bb2b966SJungho 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); 10510bb2b966SJungho Lee } else { 10520bb2b966SJungho Lee ratioi = (mx-1)/(Mx-1); 10530bb2b966SJungho 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); 10540bb2b966SJungho Lee } 10550bb2b966SJungho Lee 10560bb2b966SJungho Lee ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); 10570bb2b966SJungho Lee ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); 10580bb2b966SJungho Lee ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 10590bb2b966SJungho Lee 10600bb2b966SJungho Lee ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr); 10610bb2b966SJungho Lee ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr); 10620bb2b966SJungho Lee ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 10630bb2b966SJungho Lee 10640bb2b966SJungho Lee 10650bb2b966SJungho Lee /* loop over local fine grid nodes setting interpolation for those*/ 10660bb2b966SJungho Lee nc = 0; 10670bb2b966SJungho Lee ierr = PetscMalloc(m_f*sizeof(PetscInt),&cols);CHKERRQ(ierr); 10680bb2b966SJungho Lee 10690bb2b966SJungho Lee 10700bb2b966SJungho Lee for (i=i_start_c; i<i_start_c+m_c; i++) { 10710bb2b966SJungho Lee PetscInt i_f = i*ratioi; 10720bb2b966SJungho Lee 10730bb2b966SJungho Lee if (i_f < i_start_ghost || i_f >= i_start_ghost+m_ghost) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\ 10740bb2b966SJungho Lee i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost); 10750bb2b966SJungho Lee row = idx_f[dof*(i_f-i_start_ghost)]; 10760bb2b966SJungho Lee cols[nc++] = row/dof; 10770bb2b966SJungho Lee } 10780bb2b966SJungho Lee 10790bb2b966SJungho Lee 10800bb2b966SJungho Lee ierr = ISCreateBlock(((PetscObject)daf)->comm,dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); 10810bb2b966SJungho Lee ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); 10820bb2b966SJungho Lee ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); 10830bb2b966SJungho Lee ierr = VecScatterCreate(vecf,isf,vecc,PETSC_NULL,inject);CHKERRQ(ierr); 10840bb2b966SJungho Lee ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); 10850bb2b966SJungho Lee ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); 10860bb2b966SJungho Lee ierr = ISDestroy(&isf);CHKERRQ(ierr); 10870bb2b966SJungho Lee PetscFunctionReturn(0); 10880bb2b966SJungho Lee } 10890bb2b966SJungho Lee 10900bb2b966SJungho Lee #undef __FUNCT__ 109194013140SBarry Smith #define __FUNCT__ "DMGetInjection_DA_2D" 109294013140SBarry Smith PetscErrorCode DMGetInjection_DA_2D(DM dac,DM daf,VecScatter *inject) 109347c6ae99SBarry Smith { 109447c6ae99SBarry Smith PetscErrorCode ierr; 109547c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 109647c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c; 109747c6ae99SBarry Smith PetscInt row,i_start_ghost,j_start_ghost,mx,m_c,my,nc,ratioi,ratioj; 1098076202ddSJed Brown PetscInt i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c; 109947c6ae99SBarry Smith PetscInt *cols; 11001321219cSEthan Coon DMDABoundaryType bx,by; 110147c6ae99SBarry Smith Vec vecf,vecc; 110247c6ae99SBarry Smith IS isf; 110347c6ae99SBarry Smith 110447c6ae99SBarry Smith PetscFunctionBegin; 11051321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr); 11061321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 11071321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC) { 110847c6ae99SBarry Smith ratioi = mx/Mx; 110947c6ae99SBarry 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); 111047c6ae99SBarry Smith } else { 111147c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 111247c6ae99SBarry 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); 111347c6ae99SBarry Smith } 11141321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC) { 111547c6ae99SBarry Smith ratioj = my/My; 111647c6ae99SBarry 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); 111747c6ae99SBarry Smith } else { 111847c6ae99SBarry Smith ratioj = (my-1)/(My-1); 111947c6ae99SBarry 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); 112047c6ae99SBarry Smith } 112147c6ae99SBarry Smith 1122aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 1123aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 1124aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 112547c6ae99SBarry Smith 1126aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 1127aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 1128aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 112947c6ae99SBarry Smith 113047c6ae99SBarry Smith 113147c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 113247c6ae99SBarry Smith nc = 0; 113347c6ae99SBarry Smith ierr = PetscMalloc(n_f*m_f*sizeof(PetscInt),&cols);CHKERRQ(ierr); 1134076202ddSJed Brown for (j=j_start_c; j<j_start_c+n_c; j++) { 1135076202ddSJed Brown for (i=i_start_c; i<i_start_c+m_c; i++) { 1136076202ddSJed Brown PetscInt i_f = i*ratioi,j_f = j*ratioj; 1137076202ddSJed 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\ 1138076202ddSJed Brown j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost); 1139076202ddSJed 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\ 1140076202ddSJed Brown i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost); 1141076202ddSJed Brown row = idx_f[dof*(m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))]; 114247c6ae99SBarry Smith cols[nc++] = row/dof; 114347c6ae99SBarry Smith } 114447c6ae99SBarry Smith } 114547c6ae99SBarry Smith 114647c6ae99SBarry Smith ierr = ISCreateBlock(((PetscObject)daf)->comm,dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); 11479a42bb27SBarry Smith ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); 11489a42bb27SBarry Smith ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); 114947c6ae99SBarry Smith ierr = VecScatterCreate(vecf,isf,vecc,PETSC_NULL,inject);CHKERRQ(ierr); 11509a42bb27SBarry Smith ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); 11519a42bb27SBarry Smith ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); 1152fcfd50ebSBarry Smith ierr = ISDestroy(&isf);CHKERRQ(ierr); 115347c6ae99SBarry Smith PetscFunctionReturn(0); 115447c6ae99SBarry Smith } 115547c6ae99SBarry Smith 115647c6ae99SBarry Smith #undef __FUNCT__ 1157b1757049SJed Brown #define __FUNCT__ "DMGetInjection_DA_3D" 1158b1757049SJed Brown PetscErrorCode DMGetInjection_DA_3D(DM dac,DM daf,VecScatter *inject) 1159b1757049SJed Brown { 1160b1757049SJed Brown PetscErrorCode ierr; 1161b1757049SJed Brown PetscInt i,j,k,i_start,j_start,k_start,m_f,n_f,p_f,Mx,My,Mz; 1162b1757049SJed Brown PetscInt m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c; 1163b1757049SJed Brown PetscInt i_start_ghost,j_start_ghost,k_start_ghost; 1164b1757049SJed Brown PetscInt mx,my,mz,ratioi,ratioj,ratiok; 1165b1757049SJed Brown PetscInt i_start_c,j_start_c,k_start_c; 1166b1757049SJed Brown PetscInt m_c,n_c,p_c; 1167b1757049SJed Brown PetscInt i_start_ghost_c,j_start_ghost_c,k_start_ghost_c; 1168b1757049SJed Brown PetscInt row,nc,dof; 1169b1757049SJed Brown PetscInt *idx_c,*idx_f; 1170b1757049SJed Brown PetscInt *cols; 11711321219cSEthan Coon DMDABoundaryType bx,by,bz; 1172b1757049SJed Brown Vec vecf,vecc; 1173b1757049SJed Brown IS isf; 1174b1757049SJed Brown 1175b1757049SJed Brown PetscFunctionBegin; 11761321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr); 11771321219cSEthan Coon ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 1178b1757049SJed Brown 11791321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC){ 1180b1757049SJed Brown ratioi = mx/Mx; 1181b1757049SJed 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); 1182b1757049SJed Brown } else { 1183b1757049SJed Brown ratioi = (mx-1)/(Mx-1); 1184b1757049SJed 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); 1185b1757049SJed Brown } 11861321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC){ 1187b1757049SJed Brown ratioj = my/My; 1188b1757049SJed 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); 1189b1757049SJed Brown } else { 1190b1757049SJed Brown ratioj = (my-1)/(My-1); 1191b1757049SJed 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); 1192b1757049SJed Brown } 11931321219cSEthan Coon if (bz == DMDA_BOUNDARY_PERIODIC){ 1194b1757049SJed Brown ratiok = mz/Mz; 1195b1757049SJed 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); 1196b1757049SJed Brown } else { 1197b1757049SJed Brown ratiok = (mz-1)/(Mz-1); 1198b1757049SJed 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); 1199b1757049SJed Brown } 1200b1757049SJed Brown 1201b1757049SJed Brown ierr = DMDAGetCorners(daf,&i_start,&j_start,&k_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 1202b1757049SJed Brown ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&k_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 1203b1757049SJed Brown ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 1204b1757049SJed Brown 1205b1757049SJed Brown ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&k_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 1206b1757049SJed 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); 1207b1757049SJed Brown ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 1208b1757049SJed Brown 1209b1757049SJed Brown 1210b1757049SJed Brown /* loop over local fine grid nodes setting interpolation for those*/ 1211b1757049SJed Brown nc = 0; 1212b1757049SJed Brown ierr = PetscMalloc(n_f*m_f*p_f*sizeof(PetscInt),&cols);CHKERRQ(ierr); 1213b1757049SJed Brown for (k=k_start_c; k<k_start_c+p_c; k++) { 1214b1757049SJed Brown for (j=j_start_c; j<j_start_c+n_c; j++) { 1215b1757049SJed Brown for (i=i_start_c; i<i_start_c+m_c; i++) { 1216b1757049SJed Brown PetscInt i_f = i*ratioi,j_f = j*ratioj,k_f = k*ratiok; 1217b1757049SJed 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 " 1218b1757049SJed Brown "k_c %D k_f %D fine ghost range [%D,%D]",k,k_f,k_start_ghost,k_start_ghost+p_ghost); 1219b1757049SJed 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 " 1220b1757049SJed Brown "j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost); 1221b1757049SJed 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 " 1222b1757049SJed Brown "i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost); 1223b1757049SJed 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))]; 1224b1757049SJed Brown cols[nc++] = row/dof; 1225b1757049SJed Brown } 1226b1757049SJed Brown } 1227b1757049SJed Brown } 1228b1757049SJed Brown 1229b1757049SJed Brown ierr = ISCreateBlock(((PetscObject)daf)->comm,dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); 1230b1757049SJed Brown ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); 1231b1757049SJed Brown ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); 1232b1757049SJed Brown ierr = VecScatterCreate(vecf,isf,vecc,PETSC_NULL,inject);CHKERRQ(ierr); 1233b1757049SJed Brown ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); 1234b1757049SJed Brown ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); 1235fcfd50ebSBarry Smith ierr = ISDestroy(&isf);CHKERRQ(ierr); 1236b1757049SJed Brown PetscFunctionReturn(0); 1237b1757049SJed Brown } 1238b1757049SJed Brown 1239b1757049SJed Brown #undef __FUNCT__ 124094013140SBarry Smith #define __FUNCT__ "DMGetInjection_DA" 12417087cfbeSBarry Smith PetscErrorCode DMGetInjection_DA(DM dac,DM daf,VecScatter *inject) 124247c6ae99SBarry Smith { 124347c6ae99SBarry Smith PetscErrorCode ierr; 124447c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 12451321219cSEthan Coon DMDABoundaryType bxc,byc,bzc,bxf,byf,bzf; 1246aa219208SBarry Smith DMDAStencilType stc,stf; 124747c6ae99SBarry Smith 124847c6ae99SBarry Smith PetscFunctionBegin; 124947c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 125047c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 125147c6ae99SBarry Smith PetscValidPointer(inject,3); 125247c6ae99SBarry Smith 12531321219cSEthan Coon ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr); 12541321219cSEthan Coon ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr); 1255aa219208SBarry Smith if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 1256aa219208SBarry Smith if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); 1257aa219208SBarry Smith if (sc != sf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr); 12581321219cSEthan Coon if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr); 1259aa219208SBarry Smith if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 126047c6ae99SBarry Smith if (Mc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction"); 126147c6ae99SBarry Smith if (dimc > 1 && Nc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in y direction"); 126247c6ae99SBarry Smith if (dimc > 2 && Pc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in z direction"); 126347c6ae99SBarry Smith 12640bb2b966SJungho Lee if (dimc == 1){ 12650bb2b966SJungho Lee ierr = DMGetInjection_DA_1D(dac,daf,inject);CHKERRQ(ierr); 12660bb2b966SJungho Lee } else if (dimc == 2) { 126794013140SBarry Smith ierr = DMGetInjection_DA_2D(dac,daf,inject);CHKERRQ(ierr); 1268b1757049SJed Brown } else if (dimc == 3) { 1269b1757049SJed Brown ierr = DMGetInjection_DA_3D(dac,daf,inject);CHKERRQ(ierr); 127047c6ae99SBarry Smith } 127147c6ae99SBarry Smith PetscFunctionReturn(0); 127247c6ae99SBarry Smith } 127347c6ae99SBarry Smith 127447c6ae99SBarry Smith #undef __FUNCT__ 127594013140SBarry Smith #define __FUNCT__ "DMGetAggregates_DA" 12767087cfbeSBarry Smith PetscErrorCode DMGetAggregates_DA(DM dac,DM daf,Mat *rest) 127747c6ae99SBarry Smith { 127847c6ae99SBarry Smith PetscErrorCode ierr; 127947c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc; 128047c6ae99SBarry Smith PetscInt dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 12811321219cSEthan Coon DMDABoundaryType bxc,byc,bzc,bxf,byf,bzf; 1282aa219208SBarry Smith DMDAStencilType stc,stf; 128347c6ae99SBarry Smith PetscInt i,j,l; 128447c6ae99SBarry Smith PetscInt i_start,j_start,l_start, m_f,n_f,p_f; 128547c6ae99SBarry Smith PetscInt i_start_ghost,j_start_ghost,l_start_ghost,m_ghost,n_ghost,p_ghost; 128647c6ae99SBarry Smith PetscInt *idx_f; 128747c6ae99SBarry Smith PetscInt i_c,j_c,l_c; 128847c6ae99SBarry Smith PetscInt i_start_c,j_start_c,l_start_c, m_c,n_c,p_c; 128947c6ae99SBarry Smith PetscInt i_start_ghost_c,j_start_ghost_c,l_start_ghost_c,m_ghost_c,n_ghost_c,p_ghost_c; 129047c6ae99SBarry Smith PetscInt *idx_c; 129147c6ae99SBarry Smith PetscInt d; 129247c6ae99SBarry Smith PetscInt a; 129347c6ae99SBarry Smith PetscInt max_agg_size; 129447c6ae99SBarry Smith PetscInt *fine_nodes; 129547c6ae99SBarry Smith PetscScalar *one_vec; 129647c6ae99SBarry Smith PetscInt fn_idx; 129747c6ae99SBarry Smith 129847c6ae99SBarry Smith PetscFunctionBegin; 129947c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 130047c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 130147c6ae99SBarry Smith PetscValidPointer(rest,3); 130247c6ae99SBarry Smith 13031321219cSEthan Coon ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr); 13041321219cSEthan Coon ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr); 1305aa219208SBarry Smith if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 1306aa219208SBarry Smith if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); 1307aa219208SBarry Smith if (sc != sf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr); 13081321219cSEthan Coon if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr); 1309aa219208SBarry Smith if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 131047c6ae99SBarry Smith 131147c6ae99SBarry 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); 131247c6ae99SBarry 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); 131347c6ae99SBarry 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); 131447c6ae99SBarry Smith 131547c6ae99SBarry Smith if (Pc < 0) Pc = 1; 131647c6ae99SBarry Smith if (Pf < 0) Pf = 1; 131747c6ae99SBarry Smith if (Nc < 0) Nc = 1; 131847c6ae99SBarry Smith if (Nf < 0) Nf = 1; 131947c6ae99SBarry Smith 1320aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 1321aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 1322aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 132347c6ae99SBarry Smith 1324aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 1325aa219208SBarry 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); 1326aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 132747c6ae99SBarry Smith 132847c6ae99SBarry Smith /* 132947c6ae99SBarry Smith Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios 133047c6ae99SBarry Smith for dimension 1 and 2 respectively. 133147c6ae99SBarry Smith Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1) 133247c6ae99SBarry Smith and r_y*j and r_y*(j+1) will be grouped into the same coarse grid agregate. 133347c6ae99SBarry Smith Each specific dof on the fine grid is mapped to one dof on the coarse grid. 133447c6ae99SBarry Smith */ 133547c6ae99SBarry Smith 133647c6ae99SBarry Smith max_agg_size = (Mf/Mc+1)*(Nf/Nc+1)*(Pf/Pc+1); 133747c6ae99SBarry Smith 133847c6ae99SBarry Smith /* create the matrix that will contain the restriction operator */ 133947c6ae99SBarry Smith ierr = MatCreateMPIAIJ( ((PetscObject)daf)->comm, m_c*n_c*p_c*dofc, m_f*n_f*p_f*doff, Mc*Nc*Pc*dofc, Mf*Nf*Pf*doff, 134047c6ae99SBarry Smith max_agg_size, PETSC_NULL, max_agg_size, PETSC_NULL, rest);CHKERRQ(ierr); 134147c6ae99SBarry Smith 134247c6ae99SBarry Smith /* store nodes in the fine grid here */ 134347c6ae99SBarry Smith ierr = PetscMalloc2(max_agg_size,PetscScalar, &one_vec,max_agg_size,PetscInt, &fine_nodes);CHKERRQ(ierr); 134447c6ae99SBarry Smith for(i=0; i<max_agg_size; i++) one_vec[i] = 1.0; 134547c6ae99SBarry Smith 134647c6ae99SBarry Smith /* loop over all coarse nodes */ 134747c6ae99SBarry Smith for (l_c=l_start_c; l_c<l_start_c+p_c; l_c++) { 134847c6ae99SBarry Smith for (j_c=j_start_c; j_c<j_start_c+n_c; j_c++) { 134947c6ae99SBarry Smith for (i_c=i_start_c; i_c<i_start_c+m_c; i_c++) { 135047c6ae99SBarry Smith for(d=0; d<dofc; d++) { 135147c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 135247c6ae99SBarry 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; 135347c6ae99SBarry Smith 135447c6ae99SBarry Smith fn_idx = 0; 135547c6ae99SBarry Smith /* Corresponding fine points are all points (i_f, j_f, l_f) such that 135647c6ae99SBarry Smith i_c*Mf/Mc <= i_f < (i_c+1)*Mf/Mc 135747c6ae99SBarry Smith (same for other dimensions) 135847c6ae99SBarry Smith */ 135947c6ae99SBarry Smith for (l=l_c*Pf/Pc; l<PetscMin((l_c+1)*Pf/Pc,Pf); l++) { 136047c6ae99SBarry Smith for (j=j_c*Nf/Nc; j<PetscMin((j_c+1)*Nf/Nc,Nf); j++) { 136147c6ae99SBarry Smith for (i=i_c*Mf/Mc; i<PetscMin((i_c+1)*Mf/Mc,Mf); i++) { 136247c6ae99SBarry 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; 136347c6ae99SBarry Smith fn_idx++; 136447c6ae99SBarry Smith } 136547c6ae99SBarry Smith } 136647c6ae99SBarry Smith } 136747c6ae99SBarry Smith /* add all these points to one aggregate */ 136847c6ae99SBarry Smith ierr = MatSetValues(*rest, 1, &a, fn_idx, fine_nodes, one_vec, INSERT_VALUES);CHKERRQ(ierr); 136947c6ae99SBarry Smith } 137047c6ae99SBarry Smith } 137147c6ae99SBarry Smith } 137247c6ae99SBarry Smith } 137347c6ae99SBarry Smith ierr = PetscFree2(one_vec,fine_nodes);CHKERRQ(ierr); 137447c6ae99SBarry Smith ierr = MatAssemblyBegin(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 137547c6ae99SBarry Smith ierr = MatAssemblyEnd(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 137647c6ae99SBarry Smith PetscFunctionReturn(0); 137747c6ae99SBarry Smith } 1378