/* Code for interpolating between grids represented by DMDAs */ /* 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 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 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 consider it cleaner, but old version is turned on since it handles periodic case. */ #define NEWVERSION 0 #include /*I "petscdmda.h" I*/ #undef __FUNCT__ #define __FUNCT__ "DMCreateInterpolationScale" /*@ DMCreateInterpolationScale - Forms L = R*1/diag(R*1) - L.*v is like a coarse grid average of the nearby fine grid points. Input Parameters: + dac - DM that defines a coarse mesh . daf - DM that defines a fine mesh - mat - the restriction (or interpolation operator) from fine to coarse Output Parameter: . scale - the scaled vector Level: developer .seealso: DMCreateInterpolation() @*/ PetscErrorCode DMCreateInterpolationScale(DM dac,DM daf,Mat mat,Vec *scale) { PetscErrorCode ierr; Vec fine; PetscScalar one = 1.0; PetscFunctionBegin; ierr = DMCreateGlobalVector(daf,&fine);CHKERRQ(ierr); ierr = DMCreateGlobalVector(dac,scale);CHKERRQ(ierr); ierr = VecSet(fine,one);CHKERRQ(ierr); ierr = MatRestrict(mat,fine,*scale);CHKERRQ(ierr); ierr = VecDestroy(&fine);CHKERRQ(ierr); ierr = VecReciprocal(*scale);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "DMCreateInterpolation_DA_1D_Q1" PetscErrorCode DMCreateInterpolation_DA_1D_Q1(DM dac,DM daf,Mat *A) { PetscErrorCode ierr; PetscInt i,i_start,m_f,Mx; const PetscInt *idx_f,*idx_c; PetscInt m_ghost,m_ghost_c; PetscInt row,col,i_start_ghost,mx,m_c,nc,ratio; PetscInt i_c,i_start_c,i_start_ghost_c,cols[2],dof; PetscScalar v[2],x; Mat mat; DMBoundaryType bx; ISLocalToGlobalMapping ltog_f,ltog_c; PetscFunctionBegin; ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr); ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); if (bx == DM_BOUNDARY_PERIODIC) { ratio = mx/Mx; 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); } else { ratio = (mx-1)/(Mx-1); 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); } ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); ierr = DMGetLocalToGlobalMapping(daf,<og_f);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr); ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr); ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr); ierr = DMGetLocalToGlobalMapping(dac,<og_c);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr); /* create interpolation matrix */ ierr = MatCreate(PetscObjectComm((PetscObject)dac),&mat);CHKERRQ(ierr); ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr); ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); ierr = MatSeqAIJSetPreallocation(mat,2,NULL);CHKERRQ(ierr); ierr = MatMPIAIJSetPreallocation(mat,2,NULL,1,NULL);CHKERRQ(ierr); /* loop over local fine grid nodes setting interpolation for those*/ if (!NEWVERSION) { for (i=i_start; i */ #undef __FUNCT__ #define __FUNCT__ "DMCreateInterpolation_DA_2D_Q0" PetscErrorCode DMCreateInterpolation_DA_2D_Q0(DM dac,DM daf,Mat *A) { PetscErrorCode ierr; PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,dof; const PetscInt *idx_c,*idx_f; ISLocalToGlobalMapping ltog_f,ltog_c; PetscInt m_ghost,n_ghost,m_ghost_c,n_ghost_c,*dnz,*onz; PetscInt row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj; 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; PetscMPIInt size_c,size_f,rank_f; PetscScalar v[4]; Mat mat; DMBoundaryType bx,by; PetscFunctionBegin; ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr); ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); if (bx == DM_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x"); if (by == DM_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y"); ratioi = mx/Mx; ratioj = my/My; if (ratioi*Mx != mx) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in x"); if (ratioj*My != my) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in y"); if (ratioi != 2) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 2"); if (ratioj != 2) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 2"); ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); ierr = DMGetLocalToGlobalMapping(daf,<og_f);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr); ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); ierr = DMGetLocalToGlobalMapping(dac,<og_c);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr); /* Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the processors). It's effective length is hence 4 times its normal length, this is why the col_scale is multiplied by the interpolation matrix column sizes. sol_shift allows each set of 1/4 processors do its own interpolation using ITS copy of the coarse vector. A bit of a hack but you do better. In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 */ ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c);CHKERRQ(ierr); ierr = MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f);CHKERRQ(ierr); ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f);CHKERRQ(ierr); col_scale = size_f/size_c; col_shift = Mx*My*(rank_f/size_c); ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)daf),m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr); for (j=j_start; j */ #undef __FUNCT__ #define __FUNCT__ "DMCreateInterpolation_DA_3D_Q0" PetscErrorCode DMCreateInterpolation_DA_3D_Q0(DM dac,DM daf,Mat *A) { PetscErrorCode ierr; PetscInt i,j,l,i_start,j_start,l_start,m_f,n_f,p_f,Mx,My,Mz,dof; const PetscInt *idx_c,*idx_f; ISLocalToGlobalMapping ltog_f,ltog_c; PetscInt m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c,nc,*dnz,*onz; 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; 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; PetscMPIInt size_c,size_f,rank_f; PetscScalar v[8]; Mat mat; DMBoundaryType bx,by,bz; PetscFunctionBegin; ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr); ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); if (bx == DM_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x"); if (by == DM_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y"); if (bz == DM_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in z"); ratioi = mx/Mx; ratioj = my/My; ratiol = mz/Mz; if (ratioi*Mx != mx) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in x"); if (ratioj*My != my) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in y"); if (ratiol*Mz != mz) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Fine grid points must be multiple of coarse grid points in z"); if (ratioi != 2 && ratioi != 1) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 1 or 2"); if (ratioj != 2 && ratioj != 1) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 1 or 2"); if (ratiol != 2 && ratiol != 1) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_WRONG,"Coarsening factor in z must be 1 or 2"); ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); ierr = DMGetLocalToGlobalMapping(daf,<og_f);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr); ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 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); ierr = DMGetLocalToGlobalMapping(dac,<og_c);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr); /* Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the processors). It's effective length is hence 4 times its normal length, this is why the col_scale is multiplied by the interpolation matrix column sizes. sol_shift allows each set of 1/4 processors do its own interpolation using ITS copy of the coarse vector. A bit of a hack but you do better. In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 */ ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dac),&size_c);CHKERRQ(ierr); ierr = MPI_Comm_size(PetscObjectComm((PetscObject)daf),&size_f);CHKERRQ(ierr); ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)daf),&rank_f);CHKERRQ(ierr); col_scale = size_f/size_c; col_shift = Mx*My*Mz*(rank_f/size_c); ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)daf),m_f*n_f*p_f,col_scale*m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr); for (l=l_start; ldata; PetscFunctionBegin; PetscValidHeaderSpecific(dac,DM_CLASSID,1); PetscValidHeaderSpecific(daf,DM_CLASSID,2); PetscValidPointer(A,3); if (scale) PetscValidPointer(scale,4); ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr); ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr); if (dimc != dimf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); if (dofc != doff) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); if (sc != sf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr); if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr); if (stc != stf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); if (Mc < 2 && Mf > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction"); 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"); 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"); if (ddc->interptype == DMDA_Q1) { if (dimc == 1) { ierr = DMCreateInterpolation_DA_1D_Q1(dac,daf,A);CHKERRQ(ierr); } else if (dimc == 2) { ierr = DMCreateInterpolation_DA_2D_Q1(dac,daf,A);CHKERRQ(ierr); } else if (dimc == 3) { ierr = DMCreateInterpolation_DA_3D_Q1(dac,daf,A);CHKERRQ(ierr); } else SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype); } else if (ddc->interptype == DMDA_Q0) { if (dimc == 1) { ierr = DMCreateInterpolation_DA_1D_Q0(dac,daf,A);CHKERRQ(ierr); } else if (dimc == 2) { ierr = DMCreateInterpolation_DA_2D_Q0(dac,daf,A);CHKERRQ(ierr); } else if (dimc == 3) { ierr = DMCreateInterpolation_DA_3D_Q0(dac,daf,A);CHKERRQ(ierr); } else SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype); } if (scale) { ierr = DMCreateInterpolationScale((DM)dac,(DM)daf,*A,scale);CHKERRQ(ierr); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "DMCreateInjection_DA_1D" PetscErrorCode DMCreateInjection_DA_1D(DM dac,DM daf,VecScatter *inject) { PetscErrorCode ierr; PetscInt i,i_start,m_f,Mx,dof; const PetscInt *idx_f; ISLocalToGlobalMapping ltog_f; PetscInt m_ghost,m_ghost_c; PetscInt row,i_start_ghost,mx,m_c,nc,ratioi; PetscInt i_start_c,i_start_ghost_c; PetscInt *cols; DMBoundaryType bx; Vec vecf,vecc; IS isf; PetscFunctionBegin; ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&bx,0,0,0);CHKERRQ(ierr); ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); if (bx == DM_BOUNDARY_PERIODIC) { ratioi = mx/Mx; 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); } else { ratioi = (mx-1)/(Mx-1); 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); } ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); ierr = DMGetLocalToGlobalMapping(daf,<og_f);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr); ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr); ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr); /* loop over local fine grid nodes setting interpolation for those*/ nc = 0; ierr = PetscMalloc1(m_f,&cols);CHKERRQ(ierr); for (i=i_start_c; i= 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); row = idx_f[(i_f-i_start_ghost)]; cols[nc++] = row; } ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr); ierr = ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); ierr = VecScatterCreate(vecf,isf,vecc,NULL,inject);CHKERRQ(ierr); ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); ierr = ISDestroy(&isf);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "DMCreateInjection_DA_2D" PetscErrorCode DMCreateInjection_DA_2D(DM dac,DM daf,VecScatter *inject) { PetscErrorCode ierr; PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,dof; const PetscInt *idx_c,*idx_f; ISLocalToGlobalMapping ltog_f,ltog_c; PetscInt m_ghost,n_ghost,m_ghost_c,n_ghost_c; PetscInt row,i_start_ghost,j_start_ghost,mx,m_c,my,nc,ratioi,ratioj; PetscInt i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c; PetscInt *cols; DMBoundaryType bx,by; Vec vecf,vecc; IS isf; PetscFunctionBegin; ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&bx,&by,0,0);CHKERRQ(ierr); ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); if (bx == DM_BOUNDARY_PERIODIC) { ratioi = mx/Mx; 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); } else { ratioi = (mx-1)/(Mx-1); 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); } if (by == DM_BOUNDARY_PERIODIC) { ratioj = my/My; 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); } else { ratioj = (my-1)/(My-1); 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); } ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); ierr = DMGetLocalToGlobalMapping(daf,<og_f);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr); ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); ierr = DMGetLocalToGlobalMapping(dac,<og_c);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr); /* loop over local fine grid nodes setting interpolation for those*/ nc = 0; ierr = PetscMalloc1(n_f*m_f,&cols);CHKERRQ(ierr); for (j=j_start_c; j= j_start_ghost+n_ghost) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\ j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost); 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\ i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost); row = idx_f[(m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))]; cols[nc++] = row; } } ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr); ierr = ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); ierr = VecScatterCreate(vecf,isf,vecc,NULL,inject);CHKERRQ(ierr); ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); ierr = ISDestroy(&isf);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "DMCreateInjection_DA_3D" PetscErrorCode DMCreateInjection_DA_3D(DM dac,DM daf,VecScatter *inject) { PetscErrorCode ierr; PetscInt i,j,k,i_start,j_start,k_start,m_f,n_f,p_f,Mx,My,Mz; PetscInt m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c; PetscInt i_start_ghost,j_start_ghost,k_start_ghost; PetscInt mx,my,mz,ratioi,ratioj,ratiok; PetscInt i_start_c,j_start_c,k_start_c; PetscInt m_c,n_c,p_c; PetscInt i_start_ghost_c,j_start_ghost_c,k_start_ghost_c; PetscInt row,nc,dof; const PetscInt *idx_c,*idx_f; ISLocalToGlobalMapping ltog_f,ltog_c; PetscInt *cols; DMBoundaryType bx,by,bz; Vec vecf,vecc; IS isf; PetscFunctionBegin; ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr); ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); if (bx == DM_BOUNDARY_PERIODIC) { ratioi = mx/Mx; 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); } else { ratioi = (mx-1)/(Mx-1); 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); } if (by == DM_BOUNDARY_PERIODIC) { ratioj = my/My; 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); } else { ratioj = (my-1)/(My-1); 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); } if (bz == DM_BOUNDARY_PERIODIC) { ratiok = mz/Mz; 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); } else { ratiok = (mz-1)/(Mz-1); 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); } ierr = DMDAGetCorners(daf,&i_start,&j_start,&k_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&k_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); ierr = DMGetLocalToGlobalMapping(daf,<og_f);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr); ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&k_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 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); ierr = DMGetLocalToGlobalMapping(dac,<og_c);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr); /* loop over local fine grid nodes setting interpolation for those*/ nc = 0; ierr = PetscMalloc1(n_f*m_f*p_f,&cols);CHKERRQ(ierr); for (k=k_start_c; k= k_start_ghost+p_ghost) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA " "k_c %D k_f %D fine ghost range [%D,%D]",k,k_f,k_start_ghost,k_start_ghost+p_ghost); 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 " "j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost); 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 " "i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost); row = idx_f[(m_ghost*n_ghost*(k_f-k_start_ghost) + m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))]; cols[nc++] = row; } } } ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_f,&idx_f);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingRestoreBlockIndices(ltog_c,&idx_c);CHKERRQ(ierr); ierr = ISCreateBlock(PetscObjectComm((PetscObject)daf),dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); ierr = VecScatterCreate(vecf,isf,vecc,NULL,inject);CHKERRQ(ierr); ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); ierr = ISDestroy(&isf);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "DMCreateInjection_DA" PetscErrorCode DMCreateInjection_DA(DM dac,DM daf,Mat *mat) { PetscErrorCode ierr; PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; DMBoundaryType bxc,byc,bzc,bxf,byf,bzf; DMDAStencilType stc,stf; VecScatter inject = NULL; PetscFunctionBegin; PetscValidHeaderSpecific(dac,DM_CLASSID,1); PetscValidHeaderSpecific(daf,DM_CLASSID,2); PetscValidPointer(mat,3); ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr); ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr); if (dimc != dimf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); if (dofc != doff) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); if (sc != sf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr); if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr); if (stc != stf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); if (Mc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction"); if (dimc > 1 && Nc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in y direction"); if (dimc > 2 && Pc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in z direction"); if (dimc == 1) { ierr = DMCreateInjection_DA_1D(dac,daf,&inject);CHKERRQ(ierr); } else if (dimc == 2) { ierr = DMCreateInjection_DA_2D(dac,daf,&inject);CHKERRQ(ierr); } else if (dimc == 3) { ierr = DMCreateInjection_DA_3D(dac,daf,&inject);CHKERRQ(ierr); } ierr = MatCreateScatter(PetscObjectComm((PetscObject)inject), inject, mat);CHKERRQ(ierr); ierr = VecScatterDestroy(&inject);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "DMCreateAggregates_DA" PetscErrorCode DMCreateAggregates_DA(DM dac,DM daf,Mat *rest) { PetscErrorCode ierr; PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc; PetscInt dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; DMBoundaryType bxc,byc,bzc,bxf,byf,bzf; DMDAStencilType stc,stf; PetscInt i,j,l; PetscInt i_start,j_start,l_start, m_f,n_f,p_f; PetscInt i_start_ghost,j_start_ghost,l_start_ghost,m_ghost,n_ghost,p_ghost; const PetscInt *idx_f; PetscInt i_c,j_c,l_c; PetscInt i_start_c,j_start_c,l_start_c, m_c,n_c,p_c; PetscInt i_start_ghost_c,j_start_ghost_c,l_start_ghost_c,m_ghost_c,n_ghost_c,p_ghost_c; const PetscInt *idx_c; PetscInt d; PetscInt a; PetscInt max_agg_size; PetscInt *fine_nodes; PetscScalar *one_vec; PetscInt fn_idx; ISLocalToGlobalMapping ltogmf,ltogmc; PetscFunctionBegin; PetscValidHeaderSpecific(dac,DM_CLASSID,1); PetscValidHeaderSpecific(daf,DM_CLASSID,2); PetscValidPointer(rest,3); ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&bxc,&byc,&bzc,&stc);CHKERRQ(ierr); ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&bxf,&byf,&bzf,&stf);CHKERRQ(ierr); if (dimc != dimf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); if (dofc != doff) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); if (sc != sf) SETERRQ2(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil width of DMDA do not match %D %D",sc,sf);CHKERRQ(ierr); if (bxc != bxf || byc != byf || bzc != bzf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Boundary type different in two DMDAs");CHKERRQ(ierr); if (stc != stf) SETERRQ(PetscObjectComm((PetscObject)daf),PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 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); 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); 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); if (Pc < 0) Pc = 1; if (Pf < 0) Pf = 1; if (Nc < 0) Nc = 1; if (Nf < 0) Nf = 1; ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); ierr = DMGetLocalToGlobalMapping(daf,<ogmf);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetIndices(ltogmf,&idx_f);CHKERRQ(ierr); ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 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); ierr = DMGetLocalToGlobalMapping(dac,<ogmc);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetIndices(ltogmc,&idx_c);CHKERRQ(ierr); /* Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios for dimension 1 and 2 respectively. Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1) and r_y*j and r_y*(j+1) will be grouped into the same coarse grid agregate. Each specific dof on the fine grid is mapped to one dof on the coarse grid. */ max_agg_size = (Mf/Mc+1)*(Nf/Nc+1)*(Pf/Pc+1); /* create the matrix that will contain the restriction operator */ 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, max_agg_size, NULL, max_agg_size, NULL, rest);CHKERRQ(ierr); /* store nodes in the fine grid here */ ierr = PetscMalloc2(max_agg_size, &one_vec,max_agg_size, &fine_nodes);CHKERRQ(ierr); for (i=0; i