147c6ae99SBarry Smith #define PETSCDM_DLL 247c6ae99SBarry Smith 347c6ae99SBarry Smith /* 4aa219208SBarry Smith Code for interpolating between grids represented by DMDAs 547c6ae99SBarry Smith */ 647c6ae99SBarry Smith 7e1589f56SBarry Smith #include "private/daimpl.h" /*I "petscdm.h" I*/ 847c6ae99SBarry Smith #include "petscmg.h" 947c6ae99SBarry Smith 1047c6ae99SBarry Smith #undef __FUNCT__ 1147c6ae99SBarry Smith #define __FUNCT__ "DMGetInterpolationScale" 1247c6ae99SBarry Smith /*@ 1347c6ae99SBarry Smith DMGetInterpolationScale - Forms L = R*1/diag(R*1) - L.*v is like a coarse grid average of the 1447c6ae99SBarry Smith nearby fine grid points. 1547c6ae99SBarry Smith 1647c6ae99SBarry Smith Input Parameters: 1747c6ae99SBarry Smith + dac - DM that defines a coarse mesh 1847c6ae99SBarry Smith . daf - DM that defines a fine mesh 1947c6ae99SBarry Smith - mat - the restriction (or interpolation operator) from fine to coarse 2047c6ae99SBarry Smith 2147c6ae99SBarry Smith Output Parameter: 2247c6ae99SBarry Smith . scale - the scaled vector 2347c6ae99SBarry Smith 2447c6ae99SBarry Smith Level: developer 2547c6ae99SBarry Smith 2647c6ae99SBarry Smith .seealso: DMGetInterpolation() 2747c6ae99SBarry Smith 2847c6ae99SBarry Smith @*/ 2947c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolationScale(DM dac,DM daf,Mat mat,Vec *scale) 3047c6ae99SBarry Smith { 3147c6ae99SBarry Smith PetscErrorCode ierr; 3247c6ae99SBarry Smith Vec fine; 3347c6ae99SBarry Smith PetscScalar one = 1.0; 3447c6ae99SBarry Smith 3547c6ae99SBarry Smith PetscFunctionBegin; 3647c6ae99SBarry Smith ierr = DMCreateGlobalVector(daf,&fine);CHKERRQ(ierr); 3747c6ae99SBarry Smith ierr = DMCreateGlobalVector(dac,scale);CHKERRQ(ierr); 3847c6ae99SBarry Smith ierr = VecSet(fine,one);CHKERRQ(ierr); 3947c6ae99SBarry Smith ierr = MatRestrict(mat,fine,*scale);CHKERRQ(ierr); 4047c6ae99SBarry Smith ierr = VecDestroy(fine);CHKERRQ(ierr); 4147c6ae99SBarry Smith ierr = VecReciprocal(*scale);CHKERRQ(ierr); 4247c6ae99SBarry Smith PetscFunctionReturn(0); 4347c6ae99SBarry Smith } 4447c6ae99SBarry Smith 4547c6ae99SBarry Smith #undef __FUNCT__ 4694013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_1D_Q1" 4794013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_1D_Q1(DM dac,DM daf,Mat *A) 4847c6ae99SBarry Smith { 4947c6ae99SBarry Smith PetscErrorCode ierr; 5047c6ae99SBarry Smith PetscInt i,i_start,m_f,Mx,*idx_f; 5147c6ae99SBarry Smith PetscInt m_ghost,*idx_c,m_ghost_c; 5247c6ae99SBarry Smith PetscInt row,col,i_start_ghost,mx,m_c,nc,ratio; 5347c6ae99SBarry Smith PetscInt i_c,i_start_c,i_start_ghost_c,cols[2],dof; 5447c6ae99SBarry Smith PetscScalar v[2],x,*coors = 0,*ccoors; 5547c6ae99SBarry Smith Mat mat; 56aa219208SBarry Smith DMDAPeriodicType pt; 5747c6ae99SBarry Smith Vec vcoors,cvcoors; 5847c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data, *ddf = (DM_DA*)daf->data; 5947c6ae99SBarry Smith 6047c6ae99SBarry Smith PetscFunctionBegin; 61aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 62aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 63aa219208SBarry Smith if (pt == DMDA_XPERIODIC) { 6447c6ae99SBarry Smith ratio = mx/Mx; 6547c6ae99SBarry 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); 6647c6ae99SBarry Smith } else { 6747c6ae99SBarry Smith ratio = (mx-1)/(Mx-1); 6847c6ae99SBarry 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); 6947c6ae99SBarry Smith } 7047c6ae99SBarry Smith 71aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); 72aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); 73aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 7447c6ae99SBarry Smith 75aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr); 76aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr); 77aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 7847c6ae99SBarry Smith 7947c6ae99SBarry Smith /* create interpolation matrix */ 8047c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr); 8147c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr); 8247c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 8347c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,2,PETSC_NULL);CHKERRQ(ierr); 8447c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,2,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 8547c6ae99SBarry Smith 86aa219208SBarry Smith ierr = DMDAGetCoordinates(daf,&vcoors);CHKERRQ(ierr); 8747c6ae99SBarry Smith if (vcoors) { 88aa219208SBarry Smith ierr = DMDAGetGhostedCoordinates(dac,&cvcoors);CHKERRQ(ierr); 89aa219208SBarry Smith ierr = DMDAVecGetArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 90aa219208SBarry Smith ierr = DMDAVecGetArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 9147c6ae99SBarry Smith } 9247c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 9347c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 9447c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 9547c6ae99SBarry Smith row = idx_f[dof*(i-i_start_ghost)]/dof; 9647c6ae99SBarry Smith 9747c6ae99SBarry Smith i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 98aa219208SBarry Smith if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\ 9947c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 10047c6ae99SBarry Smith 10147c6ae99SBarry Smith /* 10247c6ae99SBarry Smith Only include those interpolation points that are truly 10347c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 10447c6ae99SBarry Smith in x direction; since they have no right neighbor 10547c6ae99SBarry Smith */ 10647c6ae99SBarry Smith if (coors) { 10747c6ae99SBarry Smith x = (coors[i] - ccoors[i_c]); 10847c6ae99SBarry Smith /* only access the next coors point if we know there is one */ 10947c6ae99SBarry Smith /* note this is dangerous because x may not exactly equal ZERO */ 11047c6ae99SBarry Smith if (PetscAbsScalar(x) != 0.0) x = x/(ccoors[i_c+1] - ccoors[i_c]); 11147c6ae99SBarry Smith } else { 11247c6ae99SBarry Smith x = ((double)(i - i_c*ratio))/((double)ratio); 11347c6ae99SBarry Smith } 11447c6ae99SBarry Smith nc = 0; 11547c6ae99SBarry Smith /* one left and below; or we are right on it */ 11647c6ae99SBarry Smith col = dof*(i_c-i_start_ghost_c); 11747c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 11847c6ae99SBarry Smith v[nc++] = - x + 1.0; 11947c6ae99SBarry Smith /* one right? */ 12047c6ae99SBarry Smith if (i_c*ratio != i) { 12147c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 12247c6ae99SBarry Smith v[nc++] = x; 12347c6ae99SBarry Smith } 12447c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 12547c6ae99SBarry Smith } 12647c6ae99SBarry Smith if (vcoors) { 127aa219208SBarry Smith ierr = DMDAVecRestoreArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 128aa219208SBarry Smith ierr = DMDAVecRestoreArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 12947c6ae99SBarry Smith } 13047c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13147c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13247c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 13347c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 13447c6ae99SBarry Smith ierr = PetscLogFlops(5.0*m_f);CHKERRQ(ierr); 13547c6ae99SBarry Smith PetscFunctionReturn(0); 13647c6ae99SBarry Smith } 13747c6ae99SBarry Smith 13847c6ae99SBarry Smith #undef __FUNCT__ 13994013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_1D_Q0" 14094013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_1D_Q0(DM dac,DM daf,Mat *A) 14147c6ae99SBarry Smith { 14247c6ae99SBarry Smith PetscErrorCode ierr; 14347c6ae99SBarry Smith PetscInt i,i_start,m_f,Mx,*idx_f; 14447c6ae99SBarry Smith PetscInt m_ghost,*idx_c,m_ghost_c; 14547c6ae99SBarry Smith PetscInt row,col,i_start_ghost,mx,m_c,nc,ratio; 14647c6ae99SBarry Smith PetscInt i_c,i_start_c,i_start_ghost_c,cols[2],dof; 14747c6ae99SBarry Smith PetscScalar v[2],x; 14847c6ae99SBarry Smith Mat mat; 149aa219208SBarry Smith DMDAPeriodicType pt; 15047c6ae99SBarry Smith 15147c6ae99SBarry Smith PetscFunctionBegin; 152aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 153aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 154aa219208SBarry Smith if (pt == DMDA_XPERIODIC) { 15547c6ae99SBarry Smith ratio = mx/Mx; 15647c6ae99SBarry 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); 15747c6ae99SBarry Smith } else { 15847c6ae99SBarry Smith ratio = (mx-1)/(Mx-1); 15947c6ae99SBarry 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); 16047c6ae99SBarry Smith } 16147c6ae99SBarry Smith 162aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); 163aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); 164aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 16547c6ae99SBarry Smith 166aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr); 167aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr); 168aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 16947c6ae99SBarry Smith 17047c6ae99SBarry Smith /* create interpolation matrix */ 17147c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr); 17247c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr); 17347c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 17447c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,2,PETSC_NULL);CHKERRQ(ierr); 17547c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,2,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 17647c6ae99SBarry Smith 17747c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 17847c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 17947c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 18047c6ae99SBarry Smith row = idx_f[dof*(i-i_start_ghost)]/dof; 18147c6ae99SBarry Smith 18247c6ae99SBarry Smith i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 18347c6ae99SBarry Smith 18447c6ae99SBarry Smith /* 18547c6ae99SBarry Smith Only include those interpolation points that are truly 18647c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 18747c6ae99SBarry Smith in x direction; since they have no right neighbor 18847c6ae99SBarry Smith */ 18947c6ae99SBarry Smith x = ((double)(i - i_c*ratio))/((double)ratio); 19047c6ae99SBarry Smith nc = 0; 19147c6ae99SBarry Smith /* one left and below; or we are right on it */ 19247c6ae99SBarry Smith col = dof*(i_c-i_start_ghost_c); 19347c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 19447c6ae99SBarry Smith v[nc++] = - x + 1.0; 19547c6ae99SBarry Smith /* one right? */ 19647c6ae99SBarry Smith if (i_c*ratio != i) { 19747c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 19847c6ae99SBarry Smith v[nc++] = x; 19947c6ae99SBarry Smith } 20047c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 20147c6ae99SBarry Smith } 20247c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20347c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20447c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 20547c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 20647c6ae99SBarry Smith ierr = PetscLogFlops(5.0*m_f);CHKERRQ(ierr); 20747c6ae99SBarry Smith PetscFunctionReturn(0); 20847c6ae99SBarry Smith } 20947c6ae99SBarry Smith 21047c6ae99SBarry Smith #undef __FUNCT__ 21194013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_2D_Q1" 21294013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_2D_Q1(DM dac,DM daf,Mat *A) 21347c6ae99SBarry Smith { 21447c6ae99SBarry Smith PetscErrorCode ierr; 21547c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 21647c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz; 21747c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj; 21847c6ae99SBarry 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; 21947c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 22047c6ae99SBarry Smith PetscScalar v[4],x,y; 22147c6ae99SBarry Smith Mat mat; 222aa219208SBarry Smith DMDAPeriodicType pt; 223aa219208SBarry Smith DMDACoor2d **coors = 0,**ccoors; 22447c6ae99SBarry Smith Vec vcoors,cvcoors; 22547c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data, *ddf = (DM_DA*)daf->data; 22647c6ae99SBarry Smith 22747c6ae99SBarry Smith PetscFunctionBegin; 228aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 229aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 230aa219208SBarry Smith if (DMDAXPeriodic(pt)){ 23147c6ae99SBarry Smith ratioi = mx/Mx; 23247c6ae99SBarry 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); 23347c6ae99SBarry Smith } else { 23447c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 23547c6ae99SBarry 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); 23647c6ae99SBarry Smith } 237aa219208SBarry Smith if (DMDAYPeriodic(pt)){ 23847c6ae99SBarry Smith ratioj = my/My; 23947c6ae99SBarry 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); 24047c6ae99SBarry Smith } else { 24147c6ae99SBarry Smith ratioj = (my-1)/(My-1); 24247c6ae99SBarry 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); 24347c6ae99SBarry Smith } 24447c6ae99SBarry Smith 24547c6ae99SBarry Smith 246aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 247aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 248aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 24947c6ae99SBarry Smith 250aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 251aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 252aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 25347c6ae99SBarry Smith 25447c6ae99SBarry Smith /* 255aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 25647c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 25747c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 25847c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 25947c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 26047c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 26147c6ae99SBarry Smith 26247c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 26347c6ae99SBarry Smith */ 26447c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr); 26547c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr); 26647c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr); 26747c6ae99SBarry Smith col_scale = size_f/size_c; 26847c6ae99SBarry Smith col_shift = Mx*My*(rank_f/size_c); 26947c6ae99SBarry Smith 27047c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr); 27147c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 27247c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 27347c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 27447c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 27547c6ae99SBarry Smith 27647c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 27747c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 27847c6ae99SBarry Smith 279aa219208SBarry 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\ 28047c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 281aa219208SBarry 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\ 28247c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 28347c6ae99SBarry Smith 28447c6ae99SBarry Smith /* 28547c6ae99SBarry Smith Only include those interpolation points that are truly 28647c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 28747c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 28847c6ae99SBarry Smith */ 28947c6ae99SBarry Smith nc = 0; 29047c6ae99SBarry Smith /* one left and below; or we are right on it */ 29147c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 29247c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 29347c6ae99SBarry Smith /* one right and below */ 29447c6ae99SBarry Smith if (i_c*ratioi != i) { 29547c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col+dof]/dof; 29647c6ae99SBarry Smith } 29747c6ae99SBarry Smith /* one left and above */ 29847c6ae99SBarry Smith if (j_c*ratioj != j) { 29947c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col+m_ghost_c*dof]/dof; 30047c6ae99SBarry Smith } 30147c6ae99SBarry Smith /* one right and above */ 30247c6ae99SBarry Smith if (j_c*ratioi != j && i_c*ratioj != i) { 30347c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; 30447c6ae99SBarry Smith } 30547c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 30647c6ae99SBarry Smith } 30747c6ae99SBarry Smith } 30847c6ae99SBarry Smith ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr); 30947c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr); 31047c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 31147c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 31247c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 31347c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 31447c6ae99SBarry Smith 315aa219208SBarry Smith ierr = DMDAGetCoordinates(daf,&vcoors);CHKERRQ(ierr); 31647c6ae99SBarry Smith if (vcoors) { 317aa219208SBarry Smith ierr = DMDAGetGhostedCoordinates(dac,&cvcoors);CHKERRQ(ierr); 318aa219208SBarry Smith ierr = DMDAVecGetArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 319aa219208SBarry Smith ierr = DMDAVecGetArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 32047c6ae99SBarry Smith } 32147c6ae99SBarry Smith 32247c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 32347c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 32447c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 32547c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 32647c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 32747c6ae99SBarry Smith 32847c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 32947c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 33047c6ae99SBarry Smith 33147c6ae99SBarry Smith /* 33247c6ae99SBarry Smith Only include those interpolation points that are truly 33347c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 33447c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 33547c6ae99SBarry Smith */ 33647c6ae99SBarry Smith if (coors) { 33747c6ae99SBarry Smith /* only access the next coors point if we know there is one */ 33847c6ae99SBarry Smith /* note this is dangerous because x may not exactly equal ZERO */ 33947c6ae99SBarry Smith x = (coors[j][i].x - ccoors[j_c][i_c].x); 34047c6ae99SBarry Smith if (PetscAbsScalar(x) != 0.0) x = x/(ccoors[j_c][i_c+1].x - ccoors[j_c][i_c].x); 34147c6ae99SBarry Smith y = (coors[j][i].y - ccoors[j_c][i_c].y); 34247c6ae99SBarry Smith if (PetscAbsScalar(y) != 0.0) y = y/(ccoors[j_c+1][i_c].y - ccoors[j_c][i_c].y); 34347c6ae99SBarry Smith } else { 34447c6ae99SBarry Smith x = ((double)(i - i_c*ratioi))/((double)ratioi); 34547c6ae99SBarry Smith y = ((double)(j - j_c*ratioj))/((double)ratioj); 34647c6ae99SBarry Smith } 34747c6ae99SBarry Smith nc = 0; 34847c6ae99SBarry Smith /* one left and below; or we are right on it */ 34947c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 35047c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 35147c6ae99SBarry Smith v[nc++] = x*y - x - y + 1.0; 35247c6ae99SBarry Smith /* one right and below */ 35347c6ae99SBarry Smith if (i_c*ratioi != i) { 35447c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+dof]/dof; 35547c6ae99SBarry Smith v[nc++] = -x*y + x; 35647c6ae99SBarry Smith } 35747c6ae99SBarry Smith /* one left and above */ 35847c6ae99SBarry Smith if (j_c*ratioj != j) { 35947c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+m_ghost_c*dof]/dof; 36047c6ae99SBarry Smith v[nc++] = -x*y + y; 36147c6ae99SBarry Smith } 36247c6ae99SBarry Smith /* one right and above */ 36347c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 36447c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; 36547c6ae99SBarry Smith v[nc++] = x*y; 36647c6ae99SBarry Smith } 36747c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 36847c6ae99SBarry Smith } 36947c6ae99SBarry Smith } 37047c6ae99SBarry Smith if (vcoors) { 371aa219208SBarry Smith ierr = DMDAVecRestoreArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 372aa219208SBarry Smith ierr = DMDAVecRestoreArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 37347c6ae99SBarry Smith } 37447c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 37547c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 37647c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 37747c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 37847c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f);CHKERRQ(ierr); 37947c6ae99SBarry Smith PetscFunctionReturn(0); 38047c6ae99SBarry Smith } 38147c6ae99SBarry Smith 38247c6ae99SBarry Smith /* 38347c6ae99SBarry Smith Contributed by Andrei Draganescu <aidraga@sandia.gov> 38447c6ae99SBarry Smith */ 38547c6ae99SBarry Smith #undef __FUNCT__ 38694013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_2D_Q0" 38794013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_2D_Q0(DM dac,DM daf,Mat *A) 38847c6ae99SBarry Smith { 38947c6ae99SBarry Smith PetscErrorCode ierr; 39047c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 39147c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz; 39247c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj; 39347c6ae99SBarry 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; 39447c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 39547c6ae99SBarry Smith PetscScalar v[4]; 39647c6ae99SBarry Smith Mat mat; 397aa219208SBarry Smith DMDAPeriodicType pt; 39847c6ae99SBarry Smith 39947c6ae99SBarry Smith PetscFunctionBegin; 400aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 401aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 402aa219208SBarry Smith if (DMDAXPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x"); 403aa219208SBarry Smith if (DMDAYPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y"); 40447c6ae99SBarry Smith ratioi = mx/Mx; 40547c6ae99SBarry Smith ratioj = my/My; 40647c6ae99SBarry 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"); 40747c6ae99SBarry 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"); 40847c6ae99SBarry Smith if (ratioi != 2) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 2"); 40947c6ae99SBarry Smith if (ratioj != 2) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 2"); 41047c6ae99SBarry Smith 411aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 412aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 413aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 41447c6ae99SBarry Smith 415aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 416aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 417aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 41847c6ae99SBarry Smith 41947c6ae99SBarry Smith /* 420aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 42147c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 42247c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 42347c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 42447c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 42547c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 42647c6ae99SBarry Smith 42747c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 42847c6ae99SBarry Smith */ 42947c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr); 43047c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr); 43147c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr); 43247c6ae99SBarry Smith col_scale = size_f/size_c; 43347c6ae99SBarry Smith col_shift = Mx*My*(rank_f/size_c); 43447c6ae99SBarry Smith 43547c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr); 43647c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 43747c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 43847c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 43947c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 44047c6ae99SBarry Smith 44147c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 44247c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 44347c6ae99SBarry Smith 444aa219208SBarry 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\ 44547c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 446aa219208SBarry 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\ 44747c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 44847c6ae99SBarry Smith 44947c6ae99SBarry Smith /* 45047c6ae99SBarry Smith Only include those interpolation points that are truly 45147c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 45247c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 45347c6ae99SBarry Smith */ 45447c6ae99SBarry Smith nc = 0; 45547c6ae99SBarry Smith /* one left and below; or we are right on it */ 45647c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 45747c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 45847c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 45947c6ae99SBarry Smith } 46047c6ae99SBarry Smith } 46147c6ae99SBarry Smith ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr); 46247c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr); 46347c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 46447c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 46547c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 46647c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 46747c6ae99SBarry Smith 46847c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 46947c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 47047c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 47147c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 47247c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 47347c6ae99SBarry Smith 47447c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 47547c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 47647c6ae99SBarry Smith nc = 0; 47747c6ae99SBarry Smith /* one left and below; or we are right on it */ 47847c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 47947c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 48047c6ae99SBarry Smith v[nc++] = 1.0; 48147c6ae99SBarry Smith 48247c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 48347c6ae99SBarry Smith } 48447c6ae99SBarry Smith } 48547c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 48647c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 48747c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 48847c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 48947c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f);CHKERRQ(ierr); 49047c6ae99SBarry Smith PetscFunctionReturn(0); 49147c6ae99SBarry Smith } 49247c6ae99SBarry Smith 49347c6ae99SBarry Smith /* 49447c6ae99SBarry Smith Contributed by Jianming Yang <jianming-yang@uiowa.edu> 49547c6ae99SBarry Smith */ 49647c6ae99SBarry Smith #undef __FUNCT__ 49794013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_3D_Q0" 49894013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_3D_Q0(DM dac,DM daf,Mat *A) 49947c6ae99SBarry Smith { 50047c6ae99SBarry Smith PetscErrorCode ierr; 50147c6ae99SBarry Smith PetscInt i,j,l,i_start,j_start,l_start,m_f,n_f,p_f,Mx,My,Mz,*idx_f,dof; 50247c6ae99SBarry Smith PetscInt m_ghost,n_ghost,p_ghost,*idx_c,m_ghost_c,n_ghost_c,p_ghost_c,nc,*dnz,*onz; 50347c6ae99SBarry 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; 50447c6ae99SBarry 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; 50547c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 50647c6ae99SBarry Smith PetscScalar v[8]; 50747c6ae99SBarry Smith Mat mat; 508aa219208SBarry Smith DMDAPeriodicType pt; 50947c6ae99SBarry Smith 51047c6ae99SBarry Smith PetscFunctionBegin; 511aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 512aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 513aa219208SBarry Smith if (DMDAXPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x"); 514aa219208SBarry Smith if (DMDAYPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y"); 515aa219208SBarry Smith if (DMDAZPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in z"); 51647c6ae99SBarry Smith ratioi = mx/Mx; 51747c6ae99SBarry Smith ratioj = my/My; 51847c6ae99SBarry Smith ratiol = mz/Mz; 51947c6ae99SBarry 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"); 52047c6ae99SBarry 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"); 52147c6ae99SBarry 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"); 52247c6ae99SBarry Smith if (ratioi != 2 && ratioi != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 1 or 2"); 52347c6ae99SBarry Smith if (ratioj != 2 && ratioj != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 1 or 2"); 52447c6ae99SBarry Smith if (ratiol != 2 && ratiol != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in z must be 1 or 2"); 52547c6ae99SBarry Smith 526aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 527aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 528aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 52947c6ae99SBarry Smith 530aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 531aa219208SBarry 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); 532aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 53347c6ae99SBarry Smith /* 534aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 53547c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 53647c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 53747c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 53847c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 53947c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 54047c6ae99SBarry Smith 54147c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 54247c6ae99SBarry Smith */ 54347c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr); 54447c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr); 54547c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr); 54647c6ae99SBarry Smith col_scale = size_f/size_c; 54747c6ae99SBarry Smith col_shift = Mx*My*Mz*(rank_f/size_c); 54847c6ae99SBarry Smith 54947c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f*p_f,col_scale*m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr); 55047c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 55147c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 55247c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 55347c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 55447c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 55547c6ae99SBarry Smith 55647c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 55747c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 55847c6ae99SBarry Smith l_c = (l/ratiol); 55947c6ae99SBarry Smith 560aa219208SBarry 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\ 56147c6ae99SBarry Smith l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c); 562aa219208SBarry 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\ 56347c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 564aa219208SBarry 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\ 56547c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 56647c6ae99SBarry Smith 56747c6ae99SBarry Smith /* 56847c6ae99SBarry Smith Only include those interpolation points that are truly 56947c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 57047c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 57147c6ae99SBarry Smith */ 57247c6ae99SBarry Smith nc = 0; 57347c6ae99SBarry Smith /* one left and below; or we are right on it */ 57447c6ae99SBarry 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)); 57547c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 57647c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 57747c6ae99SBarry Smith } 57847c6ae99SBarry Smith } 57947c6ae99SBarry Smith } 58047c6ae99SBarry Smith ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr); 58147c6ae99SBarry 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); 58247c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 58347c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 58447c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 58547c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 58647c6ae99SBarry Smith 58747c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 58847c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 58947c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 59047c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 59147c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 59247c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 59347c6ae99SBarry Smith 59447c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 59547c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 59647c6ae99SBarry Smith l_c = (l/ratiol); 59747c6ae99SBarry Smith nc = 0; 59847c6ae99SBarry Smith /* one left and below; or we are right on it */ 59947c6ae99SBarry 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)); 60047c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 60147c6ae99SBarry Smith v[nc++] = 1.0; 60247c6ae99SBarry Smith 60347c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 60447c6ae99SBarry Smith } 60547c6ae99SBarry Smith } 60647c6ae99SBarry Smith } 60747c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60847c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60947c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 61047c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 61147c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f*p_f);CHKERRQ(ierr); 61247c6ae99SBarry Smith PetscFunctionReturn(0); 61347c6ae99SBarry Smith } 61447c6ae99SBarry Smith 61547c6ae99SBarry Smith #undef __FUNCT__ 61694013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_3D_Q1" 61794013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_3D_Q1(DM dac,DM daf,Mat *A) 61847c6ae99SBarry Smith { 61947c6ae99SBarry Smith PetscErrorCode ierr; 62047c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof,l; 62147c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,Mz,mz; 62247c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[8],mx,m_c,my,nc,ratioi,ratioj,ratiok; 62347c6ae99SBarry Smith PetscInt i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c; 62447c6ae99SBarry Smith PetscInt l_start,p_f,l_start_ghost,p_ghost,l_start_c,p_c; 62547c6ae99SBarry Smith PetscInt l_start_ghost_c,p_ghost_c,l_c,*dnz,*onz; 62647c6ae99SBarry Smith PetscScalar v[8],x,y,z; 62747c6ae99SBarry Smith Mat mat; 628aa219208SBarry Smith DMDAPeriodicType pt; 629aa219208SBarry Smith DMDACoor3d ***coors = 0,***ccoors; 63047c6ae99SBarry Smith Vec vcoors,cvcoors; 63147c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data, *ddf = (DM_DA*)daf->data; 63247c6ae99SBarry Smith 63347c6ae99SBarry Smith PetscFunctionBegin; 634aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 635aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 63647c6ae99SBarry Smith if (mx == Mx) { 63747c6ae99SBarry Smith ratioi = 1; 638aa219208SBarry Smith } else if (DMDAXPeriodic(pt)){ 63947c6ae99SBarry Smith ratioi = mx/Mx; 64047c6ae99SBarry 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); 64147c6ae99SBarry Smith } else { 64247c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 64347c6ae99SBarry 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); 64447c6ae99SBarry Smith } 64547c6ae99SBarry Smith if (my == My) { 64647c6ae99SBarry Smith ratioj = 1; 647aa219208SBarry Smith } else if (DMDAYPeriodic(pt)){ 64847c6ae99SBarry Smith ratioj = my/My; 64947c6ae99SBarry 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); 65047c6ae99SBarry Smith } else { 65147c6ae99SBarry Smith ratioj = (my-1)/(My-1); 65247c6ae99SBarry 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); 65347c6ae99SBarry Smith } 65447c6ae99SBarry Smith if (mz == Mz) { 65547c6ae99SBarry Smith ratiok = 1; 656aa219208SBarry Smith } else if (DMDAZPeriodic(pt)){ 65747c6ae99SBarry Smith ratiok = mz/Mz; 65847c6ae99SBarry 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); 65947c6ae99SBarry Smith } else { 66047c6ae99SBarry Smith ratiok = (mz-1)/(Mz-1); 66147c6ae99SBarry 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); 66247c6ae99SBarry Smith } 66347c6ae99SBarry Smith 664aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 665aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 666aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 66747c6ae99SBarry Smith 668aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 669aa219208SBarry 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); 670aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 67147c6ae99SBarry Smith 67247c6ae99SBarry Smith /* create interpolation matrix, determining exact preallocation */ 67347c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)dac)->comm,m_f*n_f*p_f,m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr); 67447c6ae99SBarry Smith /* loop over local fine grid nodes counting interpolating points */ 67547c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 67647c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 67747c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 67847c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 67947c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 68047c6ae99SBarry Smith i_c = (i/ratioi); 68147c6ae99SBarry Smith j_c = (j/ratioj); 68247c6ae99SBarry Smith l_c = (l/ratiok); 683aa219208SBarry 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\ 68447c6ae99SBarry Smith l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c); 685aa219208SBarry 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\ 68647c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 687aa219208SBarry 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\ 68847c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 68947c6ae99SBarry Smith 69047c6ae99SBarry Smith /* 69147c6ae99SBarry Smith Only include those interpolation points that are truly 69247c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 69347c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 69447c6ae99SBarry Smith */ 69547c6ae99SBarry Smith nc = 0; 69647c6ae99SBarry 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)); 69747c6ae99SBarry Smith cols[nc++] = idx_c[col]/dof; 69847c6ae99SBarry Smith if (i_c*ratioi != i) { 69947c6ae99SBarry Smith cols[nc++] = idx_c[col+dof]/dof; 70047c6ae99SBarry Smith } 70147c6ae99SBarry Smith if (j_c*ratioj != j) { 70247c6ae99SBarry Smith cols[nc++] = idx_c[col+m_ghost_c*dof]/dof; 70347c6ae99SBarry Smith } 70447c6ae99SBarry Smith if (l_c*ratiok != l) { 70547c6ae99SBarry Smith cols[nc++] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; 70647c6ae99SBarry Smith } 70747c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 70847c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c+1)*dof]/dof; 70947c6ae99SBarry Smith } 71047c6ae99SBarry Smith if (j_c*ratioj != j && l_c*ratiok != l) { 71147c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof; 71247c6ae99SBarry Smith } 71347c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l) { 71447c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; 71547c6ae99SBarry Smith } 71647c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) { 71747c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; 71847c6ae99SBarry Smith } 71947c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 72047c6ae99SBarry Smith } 72147c6ae99SBarry Smith } 72247c6ae99SBarry Smith } 72347c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr); 72447c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f*p_f,m_c*n_c*p_c,mx*my*mz,Mx*My*Mz);CHKERRQ(ierr); 72547c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 72647c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 72747c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 72847c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 72947c6ae99SBarry Smith 730aa219208SBarry Smith ierr = DMDAGetCoordinates(daf,&vcoors);CHKERRQ(ierr); 73147c6ae99SBarry Smith if (vcoors) { 732aa219208SBarry Smith ierr = DMDAGetGhostedCoordinates(dac,&cvcoors);CHKERRQ(ierr); 733aa219208SBarry Smith ierr = DMDAVecGetArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 734aa219208SBarry Smith ierr = DMDAVecGetArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 73547c6ae99SBarry Smith } 73647c6ae99SBarry Smith 73747c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 73847c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 73947c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 74047c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 74147c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 74247c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 74347c6ae99SBarry Smith 74447c6ae99SBarry Smith i_c = (i/ratioi); 74547c6ae99SBarry Smith j_c = (j/ratioj); 74647c6ae99SBarry Smith l_c = (l/ratiok); 74747c6ae99SBarry Smith 74847c6ae99SBarry Smith /* 74947c6ae99SBarry Smith Only include those interpolation points that are truly 75047c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 75147c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 75247c6ae99SBarry Smith */ 75347c6ae99SBarry Smith if (coors) { 75447c6ae99SBarry Smith /* only access the next coors point if we know there is one */ 75547c6ae99SBarry Smith /* note this is dangerous because x may not exactly equal ZERO */ 75647c6ae99SBarry Smith x = (coors[l][j][i].x - ccoors[l_c][j_c][i_c].x); 75747c6ae99SBarry Smith if (PetscAbsScalar(x) != 0.0) x = x/(ccoors[l_c][j_c][i_c+1].x - ccoors[l_c][j_c][i_c].x); 75847c6ae99SBarry Smith y = (coors[l][j][i].y - ccoors[l_c][j_c][i_c].y); 75947c6ae99SBarry Smith if (PetscAbsScalar(y) != 0.0) y = y/(ccoors[l_c][j_c+1][i_c].y - ccoors[l_c][j_c][i_c].y); 76047c6ae99SBarry Smith z = (coors[l][j][i].z - ccoors[l_c][j_c][i_c].z); 76147c6ae99SBarry Smith if (PetscAbsScalar(z) != 0.0) z = z/(ccoors[l_c+1][j_c][i_c].z - ccoors[l_c][j_c][i_c].z); 76247c6ae99SBarry Smith } else { 76347c6ae99SBarry Smith x = ((double)(i - i_c*ratioi))/((double)ratioi); 76447c6ae99SBarry Smith y = ((double)(j - j_c*ratioj))/((double)ratioj); 76547c6ae99SBarry Smith z = ((double)(l - l_c*ratiok))/((double)ratiok); 76647c6ae99SBarry Smith } 76747c6ae99SBarry Smith nc = 0; 76847c6ae99SBarry Smith /* one left and below; or we are right on it */ 76947c6ae99SBarry 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)); 77047c6ae99SBarry Smith 77147c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 77247c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.)); 77347c6ae99SBarry Smith 77447c6ae99SBarry Smith if (i_c*ratioi != i) { 77547c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 77647c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.)); 77747c6ae99SBarry Smith } 77847c6ae99SBarry Smith 77947c6ae99SBarry Smith if (j_c*ratioj != j) { 78047c6ae99SBarry Smith cols[nc] = idx_c[col+m_ghost_c*dof]/dof; 78147c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.)); 78247c6ae99SBarry Smith } 78347c6ae99SBarry Smith 78447c6ae99SBarry Smith if (l_c*ratiok != l) { 78547c6ae99SBarry Smith cols[nc] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; 78647c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.)); 78747c6ae99SBarry Smith } 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 v[nc++] = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.)); 79247c6ae99SBarry Smith } 79347c6ae99SBarry Smith 79447c6ae99SBarry Smith if (j_c*ratioj != j && l_c*ratiok != l) { 79547c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof; 79647c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.)); 79747c6ae99SBarry Smith } 79847c6ae99SBarry Smith 79947c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l) { 80047c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; 80147c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.)); 80247c6ae99SBarry Smith } 80347c6ae99SBarry Smith 80447c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) { 80547c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; 80647c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.)); 80747c6ae99SBarry Smith } 80847c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 80947c6ae99SBarry Smith } 81047c6ae99SBarry Smith } 81147c6ae99SBarry Smith } 81247c6ae99SBarry Smith if (vcoors) { 813aa219208SBarry Smith ierr = DMDAVecRestoreArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 814aa219208SBarry Smith ierr = DMDAVecRestoreArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 81547c6ae99SBarry Smith } 81647c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 81747c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 81847c6ae99SBarry Smith 81947c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 82047c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 82147c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f);CHKERRQ(ierr); 82247c6ae99SBarry Smith PetscFunctionReturn(0); 82347c6ae99SBarry Smith } 82447c6ae99SBarry Smith 82547c6ae99SBarry Smith #undef __FUNCT__ 82694013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA" 82794013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolation_DA(DM dac,DM daf,Mat *A,Vec *scale) 82847c6ae99SBarry Smith { 82947c6ae99SBarry Smith PetscErrorCode ierr; 83047c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 831aa219208SBarry Smith DMDAPeriodicType wrapc,wrapf; 832aa219208SBarry Smith DMDAStencilType stc,stf; 83347c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data; 83447c6ae99SBarry Smith 83547c6ae99SBarry Smith PetscFunctionBegin; 83647c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 83747c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 83847c6ae99SBarry Smith PetscValidPointer(A,3); 83947c6ae99SBarry Smith if (scale) PetscValidPointer(scale,4); 84047c6ae99SBarry Smith 841aa219208SBarry Smith ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&wrapc,&stc);CHKERRQ(ierr); 842aa219208SBarry Smith ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&wrapf,&stf);CHKERRQ(ierr); 843aa219208SBarry Smith if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 844aa219208SBarry Smith if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); 845aa219208SBarry 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); 846aa219208SBarry Smith if (wrapc != wrapf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Periodic type different in two DMDAs");CHKERRQ(ierr); 847aa219208SBarry Smith if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 84847c6ae99SBarry Smith if (Mc < 2 && Mf > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction"); 84947c6ae99SBarry 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"); 85047c6ae99SBarry 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"); 85147c6ae99SBarry Smith 852aa219208SBarry Smith if (ddc->interptype == DMDA_Q1){ 85347c6ae99SBarry Smith if (dimc == 1){ 85494013140SBarry Smith ierr = DMGetInterpolation_DA_1D_Q1(dac,daf,A);CHKERRQ(ierr); 85547c6ae99SBarry Smith } else if (dimc == 2){ 85694013140SBarry Smith ierr = DMGetInterpolation_DA_2D_Q1(dac,daf,A);CHKERRQ(ierr); 85747c6ae99SBarry Smith } else if (dimc == 3){ 85894013140SBarry Smith ierr = DMGetInterpolation_DA_3D_Q1(dac,daf,A);CHKERRQ(ierr); 85947c6ae99SBarry Smith } else { 860aa219208SBarry Smith SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype); 86147c6ae99SBarry Smith } 862aa219208SBarry Smith } else if (ddc->interptype == DMDA_Q0){ 86347c6ae99SBarry Smith if (dimc == 1){ 86494013140SBarry Smith ierr = DMGetInterpolation_DA_1D_Q0(dac,daf,A);CHKERRQ(ierr); 86547c6ae99SBarry Smith } else if (dimc == 2){ 86694013140SBarry Smith ierr = DMGetInterpolation_DA_2D_Q0(dac,daf,A);CHKERRQ(ierr); 86747c6ae99SBarry Smith } else if (dimc == 3){ 86894013140SBarry Smith ierr = DMGetInterpolation_DA_3D_Q0(dac,daf,A);CHKERRQ(ierr); 86947c6ae99SBarry Smith } else { 870aa219208SBarry Smith SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype); 87147c6ae99SBarry Smith } 87247c6ae99SBarry Smith } 87347c6ae99SBarry Smith if (scale) { 87447c6ae99SBarry Smith ierr = DMGetInterpolationScale((DM)dac,(DM)daf,*A,scale);CHKERRQ(ierr); 87547c6ae99SBarry Smith } 87647c6ae99SBarry Smith PetscFunctionReturn(0); 87747c6ae99SBarry Smith } 87847c6ae99SBarry Smith 87947c6ae99SBarry Smith #undef __FUNCT__ 88094013140SBarry Smith #define __FUNCT__ "DMGetInjection_DA_2D" 88194013140SBarry Smith PetscErrorCode DMGetInjection_DA_2D(DM dac,DM daf,VecScatter *inject) 88247c6ae99SBarry Smith { 88347c6ae99SBarry Smith PetscErrorCode ierr; 88447c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 88547c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c; 88647c6ae99SBarry Smith PetscInt row,i_start_ghost,j_start_ghost,mx,m_c,my,nc,ratioi,ratioj; 887076202ddSJed Brown PetscInt i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c; 88847c6ae99SBarry Smith PetscInt *cols; 889aa219208SBarry Smith DMDAPeriodicType pt; 89047c6ae99SBarry Smith Vec vecf,vecc; 89147c6ae99SBarry Smith IS isf; 89247c6ae99SBarry Smith 89347c6ae99SBarry Smith PetscFunctionBegin; 89447c6ae99SBarry Smith 895aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 896aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 897aa219208SBarry Smith if (DMDAXPeriodic(pt)){ 89847c6ae99SBarry Smith ratioi = mx/Mx; 89947c6ae99SBarry 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); 90047c6ae99SBarry Smith } else { 90147c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 90247c6ae99SBarry 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); 90347c6ae99SBarry Smith } 904aa219208SBarry Smith if (DMDAYPeriodic(pt)){ 90547c6ae99SBarry Smith ratioj = my/My; 90647c6ae99SBarry 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); 90747c6ae99SBarry Smith } else { 90847c6ae99SBarry Smith ratioj = (my-1)/(My-1); 90947c6ae99SBarry 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); 91047c6ae99SBarry Smith } 91147c6ae99SBarry Smith 912aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 913aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 914aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 91547c6ae99SBarry Smith 916aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 917aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 918aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 91947c6ae99SBarry Smith 92047c6ae99SBarry Smith 92147c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 92247c6ae99SBarry Smith nc = 0; 92347c6ae99SBarry Smith ierr = PetscMalloc(n_f*m_f*sizeof(PetscInt),&cols);CHKERRQ(ierr); 924076202ddSJed Brown for (j=j_start_c; j<j_start_c+n_c; j++) { 925076202ddSJed Brown for (i=i_start_c; i<i_start_c+m_c; i++) { 926076202ddSJed Brown PetscInt i_f = i*ratioi,j_f = j*ratioj; 927076202ddSJed 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\ 928076202ddSJed Brown j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost); 929076202ddSJed 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\ 930076202ddSJed Brown i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost); 931076202ddSJed Brown row = idx_f[dof*(m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))]; 93247c6ae99SBarry Smith cols[nc++] = row/dof; 93347c6ae99SBarry Smith } 93447c6ae99SBarry Smith } 93547c6ae99SBarry Smith 93647c6ae99SBarry Smith ierr = ISCreateBlock(((PetscObject)daf)->comm,dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); 9379a42bb27SBarry Smith ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); 9389a42bb27SBarry Smith ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); 93947c6ae99SBarry Smith ierr = VecScatterCreate(vecf,isf,vecc,PETSC_NULL,inject);CHKERRQ(ierr); 9409a42bb27SBarry Smith ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); 9419a42bb27SBarry Smith ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); 94247c6ae99SBarry Smith ierr = ISDestroy(isf);CHKERRQ(ierr); 94347c6ae99SBarry Smith PetscFunctionReturn(0); 94447c6ae99SBarry Smith } 94547c6ae99SBarry Smith 94647c6ae99SBarry Smith #undef __FUNCT__ 947*b1757049SJed Brown #define __FUNCT__ "DMGetInjection_DA_3D" 948*b1757049SJed Brown PetscErrorCode DMGetInjection_DA_3D(DM dac,DM daf,VecScatter *inject) 949*b1757049SJed Brown { 950*b1757049SJed Brown PetscErrorCode ierr; 951*b1757049SJed Brown PetscInt i,j,k,i_start,j_start,k_start,m_f,n_f,p_f,Mx,My,Mz; 952*b1757049SJed Brown PetscInt m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c; 953*b1757049SJed Brown PetscInt i_start_ghost,j_start_ghost,k_start_ghost; 954*b1757049SJed Brown PetscInt mx,my,mz,ratioi,ratioj,ratiok; 955*b1757049SJed Brown PetscInt i_start_c,j_start_c,k_start_c; 956*b1757049SJed Brown PetscInt m_c,n_c,p_c; 957*b1757049SJed Brown PetscInt i_start_ghost_c,j_start_ghost_c,k_start_ghost_c; 958*b1757049SJed Brown PetscInt row,nc,dof; 959*b1757049SJed Brown PetscInt *idx_c,*idx_f; 960*b1757049SJed Brown PetscInt *cols; 961*b1757049SJed Brown DMDAPeriodicType pt; 962*b1757049SJed Brown Vec vecf,vecc; 963*b1757049SJed Brown IS isf; 964*b1757049SJed Brown 965*b1757049SJed Brown PetscFunctionBegin; 966*b1757049SJed Brown ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 967*b1757049SJed Brown ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 968*b1757049SJed Brown 969*b1757049SJed Brown if (DMDAXPeriodic(pt)){ 970*b1757049SJed Brown ratioi = mx/Mx; 971*b1757049SJed 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); 972*b1757049SJed Brown } else { 973*b1757049SJed Brown ratioi = (mx-1)/(Mx-1); 974*b1757049SJed 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); 975*b1757049SJed Brown } 976*b1757049SJed Brown if (DMDAYPeriodic(pt)){ 977*b1757049SJed Brown ratioj = my/My; 978*b1757049SJed 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); 979*b1757049SJed Brown } else { 980*b1757049SJed Brown ratioj = (my-1)/(My-1); 981*b1757049SJed 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); 982*b1757049SJed Brown } 983*b1757049SJed Brown if (DMDAZPeriodic(pt)){ 984*b1757049SJed Brown ratiok = mz/Mz; 985*b1757049SJed 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); 986*b1757049SJed Brown } else { 987*b1757049SJed Brown ratiok = (mz-1)/(Mz-1); 988*b1757049SJed 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); 989*b1757049SJed Brown } 990*b1757049SJed Brown 991*b1757049SJed Brown ierr = DMDAGetCorners(daf,&i_start,&j_start,&k_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 992*b1757049SJed Brown ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&k_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 993*b1757049SJed Brown ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 994*b1757049SJed Brown 995*b1757049SJed Brown ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&k_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 996*b1757049SJed 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); 997*b1757049SJed Brown ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 998*b1757049SJed Brown 999*b1757049SJed Brown 1000*b1757049SJed Brown /* loop over local fine grid nodes setting interpolation for those*/ 1001*b1757049SJed Brown nc = 0; 1002*b1757049SJed Brown ierr = PetscMalloc(n_f*m_f*p_f*sizeof(PetscInt),&cols);CHKERRQ(ierr); 1003*b1757049SJed Brown for (k=k_start_c; k<k_start_c+p_c; k++) { 1004*b1757049SJed Brown for (j=j_start_c; j<j_start_c+n_c; j++) { 1005*b1757049SJed Brown for (i=i_start_c; i<i_start_c+m_c; i++) { 1006*b1757049SJed Brown PetscInt i_f = i*ratioi,j_f = j*ratioj,k_f = k*ratiok; 1007*b1757049SJed 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 " 1008*b1757049SJed Brown "k_c %D k_f %D fine ghost range [%D,%D]",k,k_f,k_start_ghost,k_start_ghost+p_ghost); 1009*b1757049SJed 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 " 1010*b1757049SJed Brown "j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost); 1011*b1757049SJed 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 " 1012*b1757049SJed Brown "i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost); 1013*b1757049SJed 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))]; 1014*b1757049SJed Brown cols[nc++] = row/dof; 1015*b1757049SJed Brown } 1016*b1757049SJed Brown } 1017*b1757049SJed Brown } 1018*b1757049SJed Brown 1019*b1757049SJed Brown ierr = ISCreateBlock(((PetscObject)daf)->comm,dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); 1020*b1757049SJed Brown ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); 1021*b1757049SJed Brown ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); 1022*b1757049SJed Brown ierr = VecScatterCreate(vecf,isf,vecc,PETSC_NULL,inject);CHKERRQ(ierr); 1023*b1757049SJed Brown ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); 1024*b1757049SJed Brown ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); 1025*b1757049SJed Brown ierr = ISDestroy(isf);CHKERRQ(ierr); 1026*b1757049SJed Brown PetscFunctionReturn(0); 1027*b1757049SJed Brown } 1028*b1757049SJed Brown 1029*b1757049SJed Brown #undef __FUNCT__ 103094013140SBarry Smith #define __FUNCT__ "DMGetInjection_DA" 103194013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetInjection_DA(DM dac,DM daf,VecScatter *inject) 103247c6ae99SBarry Smith { 103347c6ae99SBarry Smith PetscErrorCode ierr; 103447c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 1035aa219208SBarry Smith DMDAPeriodicType wrapc,wrapf; 1036aa219208SBarry Smith DMDAStencilType stc,stf; 103747c6ae99SBarry Smith 103847c6ae99SBarry Smith PetscFunctionBegin; 103947c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 104047c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 104147c6ae99SBarry Smith PetscValidPointer(inject,3); 104247c6ae99SBarry Smith 1043aa219208SBarry Smith ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&wrapc,&stc);CHKERRQ(ierr); 1044aa219208SBarry Smith ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&wrapf,&stf);CHKERRQ(ierr); 1045aa219208SBarry Smith if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 1046aa219208SBarry Smith if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); 1047aa219208SBarry 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); 1048aa219208SBarry Smith if (wrapc != wrapf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Periodic type different in two DMDAs");CHKERRQ(ierr); 1049aa219208SBarry Smith if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 105047c6ae99SBarry Smith if (Mc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction"); 105147c6ae99SBarry Smith if (dimc > 1 && Nc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in y direction"); 105247c6ae99SBarry Smith if (dimc > 2 && Pc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in z direction"); 105347c6ae99SBarry Smith 105447c6ae99SBarry Smith if (dimc == 2){ 105594013140SBarry Smith ierr = DMGetInjection_DA_2D(dac,daf,inject);CHKERRQ(ierr); 1056*b1757049SJed Brown } else if (dimc == 3) { 1057*b1757049SJed Brown ierr = DMGetInjection_DA_3D(dac,daf,inject);CHKERRQ(ierr); 105847c6ae99SBarry Smith } else { 1059aa219208SBarry Smith SETERRQ1(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DMDA dimension %D",dimc); 106047c6ae99SBarry Smith } 106147c6ae99SBarry Smith PetscFunctionReturn(0); 106247c6ae99SBarry Smith } 106347c6ae99SBarry Smith 106447c6ae99SBarry Smith #undef __FUNCT__ 106594013140SBarry Smith #define __FUNCT__ "DMGetAggregates_DA" 106694013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetAggregates_DA(DM dac,DM daf,Mat *rest) 106747c6ae99SBarry Smith { 106847c6ae99SBarry Smith PetscErrorCode ierr; 106947c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc; 107047c6ae99SBarry Smith PetscInt dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 1071aa219208SBarry Smith DMDAPeriodicType wrapc,wrapf; 1072aa219208SBarry Smith DMDAStencilType stc,stf; 107347c6ae99SBarry Smith /* PetscReal r_x, r_y, r_z; */ 107447c6ae99SBarry Smith 107547c6ae99SBarry Smith PetscInt i,j,l; 107647c6ae99SBarry Smith PetscInt i_start,j_start,l_start, m_f,n_f,p_f; 107747c6ae99SBarry Smith PetscInt i_start_ghost,j_start_ghost,l_start_ghost,m_ghost,n_ghost,p_ghost; 107847c6ae99SBarry Smith PetscInt *idx_f; 107947c6ae99SBarry Smith PetscInt i_c,j_c,l_c; 108047c6ae99SBarry Smith PetscInt i_start_c,j_start_c,l_start_c, m_c,n_c,p_c; 108147c6ae99SBarry Smith PetscInt i_start_ghost_c,j_start_ghost_c,l_start_ghost_c,m_ghost_c,n_ghost_c,p_ghost_c; 108247c6ae99SBarry Smith PetscInt *idx_c; 108347c6ae99SBarry Smith PetscInt d; 108447c6ae99SBarry Smith PetscInt a; 108547c6ae99SBarry Smith PetscInt max_agg_size; 108647c6ae99SBarry Smith PetscInt *fine_nodes; 108747c6ae99SBarry Smith PetscScalar *one_vec; 108847c6ae99SBarry Smith PetscInt fn_idx; 108947c6ae99SBarry Smith 109047c6ae99SBarry Smith PetscFunctionBegin; 109147c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 109247c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 109347c6ae99SBarry Smith PetscValidPointer(rest,3); 109447c6ae99SBarry Smith 1095aa219208SBarry Smith ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&wrapc,&stc);CHKERRQ(ierr); 1096aa219208SBarry Smith ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&wrapf,&stf);CHKERRQ(ierr); 1097aa219208SBarry Smith if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 1098aa219208SBarry Smith if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); 1099aa219208SBarry 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); 1100aa219208SBarry Smith if (wrapc != wrapf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Periodic type different in two DMDAs");CHKERRQ(ierr); 1101aa219208SBarry Smith if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 110247c6ae99SBarry Smith 110347c6ae99SBarry 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); 110447c6ae99SBarry 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); 110547c6ae99SBarry 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); 110647c6ae99SBarry Smith 110747c6ae99SBarry Smith if (Pc < 0) Pc = 1; 110847c6ae99SBarry Smith if (Pf < 0) Pf = 1; 110947c6ae99SBarry Smith if (Nc < 0) Nc = 1; 111047c6ae99SBarry Smith if (Nf < 0) Nf = 1; 111147c6ae99SBarry Smith 1112aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 1113aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 1114aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 111547c6ae99SBarry Smith 1116aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 1117aa219208SBarry 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); 1118aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 111947c6ae99SBarry Smith 112047c6ae99SBarry Smith /* 112147c6ae99SBarry Smith Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios 112247c6ae99SBarry Smith for dimension 1 and 2 respectively. 112347c6ae99SBarry Smith Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1) 112447c6ae99SBarry Smith and r_y*j and r_y*(j+1) will be grouped into the same coarse grid agregate. 112547c6ae99SBarry Smith Each specific dof on the fine grid is mapped to one dof on the coarse grid. 112647c6ae99SBarry Smith */ 112747c6ae99SBarry Smith 112847c6ae99SBarry Smith max_agg_size = (Mf/Mc+1)*(Nf/Nc+1)*(Pf/Pc+1); 112947c6ae99SBarry Smith 113047c6ae99SBarry Smith /* create the matrix that will contain the restriction operator */ 113147c6ae99SBarry 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, 113247c6ae99SBarry Smith max_agg_size, PETSC_NULL, max_agg_size, PETSC_NULL, rest);CHKERRQ(ierr); 113347c6ae99SBarry Smith 113447c6ae99SBarry Smith /* store nodes in the fine grid here */ 113547c6ae99SBarry Smith ierr = PetscMalloc2(max_agg_size,PetscScalar, &one_vec,max_agg_size,PetscInt, &fine_nodes);CHKERRQ(ierr); 113647c6ae99SBarry Smith for(i=0; i<max_agg_size; i++) one_vec[i] = 1.0; 113747c6ae99SBarry Smith 113847c6ae99SBarry Smith /* loop over all coarse nodes */ 113947c6ae99SBarry Smith for (l_c=l_start_c; l_c<l_start_c+p_c; l_c++) { 114047c6ae99SBarry Smith for (j_c=j_start_c; j_c<j_start_c+n_c; j_c++) { 114147c6ae99SBarry Smith for (i_c=i_start_c; i_c<i_start_c+m_c; i_c++) { 114247c6ae99SBarry Smith for(d=0; d<dofc; d++) { 114347c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 114447c6ae99SBarry 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; 114547c6ae99SBarry Smith 114647c6ae99SBarry Smith fn_idx = 0; 114747c6ae99SBarry Smith /* Corresponding fine points are all points (i_f, j_f, l_f) such that 114847c6ae99SBarry Smith i_c*Mf/Mc <= i_f < (i_c+1)*Mf/Mc 114947c6ae99SBarry Smith (same for other dimensions) 115047c6ae99SBarry Smith */ 115147c6ae99SBarry Smith for (l=l_c*Pf/Pc; l<PetscMin((l_c+1)*Pf/Pc,Pf); l++) { 115247c6ae99SBarry Smith for (j=j_c*Nf/Nc; j<PetscMin((j_c+1)*Nf/Nc,Nf); j++) { 115347c6ae99SBarry Smith for (i=i_c*Mf/Mc; i<PetscMin((i_c+1)*Mf/Mc,Mf); i++) { 115447c6ae99SBarry 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; 115547c6ae99SBarry Smith fn_idx++; 115647c6ae99SBarry Smith } 115747c6ae99SBarry Smith } 115847c6ae99SBarry Smith } 115947c6ae99SBarry Smith /* add all these points to one aggregate */ 116047c6ae99SBarry Smith ierr = MatSetValues(*rest, 1, &a, fn_idx, fine_nodes, one_vec, INSERT_VALUES);CHKERRQ(ierr); 116147c6ae99SBarry Smith } 116247c6ae99SBarry Smith } 116347c6ae99SBarry Smith } 116447c6ae99SBarry Smith } 116547c6ae99SBarry Smith ierr = PetscFree2(one_vec,fine_nodes);CHKERRQ(ierr); 116647c6ae99SBarry Smith ierr = MatAssemblyBegin(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 116747c6ae99SBarry Smith ierr = MatAssemblyEnd(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 116847c6ae99SBarry Smith PetscFunctionReturn(0); 116947c6ae99SBarry Smith } 1170