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*/ 93*b3bd94feSDave May if (!vcoors) { 94*b3bd94feSDave May 9547c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 9647c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 9747c6ae99SBarry Smith row = idx_f[dof*(i-i_start_ghost)]/dof; 9847c6ae99SBarry Smith 9947c6ae99SBarry Smith i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 100aa219208SBarry 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\ 10147c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 10247c6ae99SBarry Smith 10347c6ae99SBarry Smith /* 10447c6ae99SBarry Smith Only include those interpolation points that are truly 10547c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 10647c6ae99SBarry Smith in x direction; since they have no right neighbor 10747c6ae99SBarry Smith */ 10847c6ae99SBarry Smith x = ((double)(i - i_c*ratio))/((double)ratio); 10947c6ae99SBarry Smith nc = 0; 11047c6ae99SBarry Smith /* one left and below; or we are right on it */ 11147c6ae99SBarry Smith col = dof*(i_c-i_start_ghost_c); 11247c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 11347c6ae99SBarry Smith v[nc++] = - x + 1.0; 11447c6ae99SBarry Smith /* one right? */ 11547c6ae99SBarry Smith if (i_c*ratio != i) { 11647c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 11747c6ae99SBarry Smith v[nc++] = x; 11847c6ae99SBarry Smith } 11947c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 12047c6ae99SBarry Smith } 121*b3bd94feSDave May 122*b3bd94feSDave May } else { 123*b3bd94feSDave May PetscScalar *xi; 124*b3bd94feSDave May PetscInt li,nxi,n; 125*b3bd94feSDave May PetscScalar Ni[2]; 126*b3bd94feSDave May 127*b3bd94feSDave May /* compute local coordinate arrays */ 128*b3bd94feSDave May nxi = ratio + 1; 129*b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr); 130*b3bd94feSDave May for (li=0; li<nxi; li++) { 131*b3bd94feSDave May xi[li] = -1.0 + li*(2.0/(PetscScalar)(nxi-1)); 132*b3bd94feSDave May } 133*b3bd94feSDave May 134*b3bd94feSDave May for (i=i_start; i<i_start+m_f; i++) { 135*b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */ 136*b3bd94feSDave May row = idx_f[dof*(i-i_start_ghost)]/dof; 137*b3bd94feSDave May 138*b3bd94feSDave May i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 139*b3bd94feSDave May if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DMDA must lie over fine DMDA\n\ 140*b3bd94feSDave May i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 141*b3bd94feSDave May 142*b3bd94feSDave May /* remainders */ 143*b3bd94feSDave May li = i - ratio * (i/ratio); 144*b3bd94feSDave May if (i==mx-1){ li = nxi-1; } 145*b3bd94feSDave May 146*b3bd94feSDave May /* corners */ 147*b3bd94feSDave May col = dof*(i_c-i_start_ghost_c); 148*b3bd94feSDave May cols[0] = idx_c[col]/dof; 149*b3bd94feSDave May Ni[0] = 1.0; 150*b3bd94feSDave May if ( (li==0) || (li==nxi-1) ) { 151*b3bd94feSDave May ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr); 152*b3bd94feSDave May continue; 153*b3bd94feSDave May } 154*b3bd94feSDave May 155*b3bd94feSDave May /* edges + interior */ 156*b3bd94feSDave May /* remainders */ 157*b3bd94feSDave May if (i==mx-1){ i_c--; } 158*b3bd94feSDave May 159*b3bd94feSDave May col = dof*(i_c-i_start_ghost_c); 160*b3bd94feSDave May cols[0] = idx_c[col]/dof; /* one left and below; or we are right on it */ 161*b3bd94feSDave May cols[1] = idx_c[col+dof]/dof; 162*b3bd94feSDave May 163*b3bd94feSDave May Ni[0] = 0.5*(1.0-xi[li]); 164*b3bd94feSDave May Ni[1] = 0.5*(1.0+xi[li]); 165*b3bd94feSDave May for (n=0; n<2; n++) { 166*b3bd94feSDave May if( PetscAbsScalar(Ni[n])<1.0e-32) { cols[n]=-1; } 167*b3bd94feSDave May } 168*b3bd94feSDave May ierr = MatSetValues(mat,1,&row,2,cols,Ni,INSERT_VALUES);CHKERRQ(ierr); 169*b3bd94feSDave May } 170*b3bd94feSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 171*b3bd94feSDave May } 17247c6ae99SBarry Smith if (vcoors) { 173aa219208SBarry Smith ierr = DMDAVecRestoreArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 174aa219208SBarry Smith ierr = DMDAVecRestoreArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 17547c6ae99SBarry Smith } 17647c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 17747c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 17847c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 17947c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 18047c6ae99SBarry Smith PetscFunctionReturn(0); 18147c6ae99SBarry Smith } 18247c6ae99SBarry Smith 18347c6ae99SBarry Smith #undef __FUNCT__ 18494013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_1D_Q0" 18594013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_1D_Q0(DM dac,DM daf,Mat *A) 18647c6ae99SBarry Smith { 18747c6ae99SBarry Smith PetscErrorCode ierr; 18847c6ae99SBarry Smith PetscInt i,i_start,m_f,Mx,*idx_f; 18947c6ae99SBarry Smith PetscInt m_ghost,*idx_c,m_ghost_c; 19047c6ae99SBarry Smith PetscInt row,col,i_start_ghost,mx,m_c,nc,ratio; 19147c6ae99SBarry Smith PetscInt i_c,i_start_c,i_start_ghost_c,cols[2],dof; 19247c6ae99SBarry Smith PetscScalar v[2],x; 19347c6ae99SBarry Smith Mat mat; 194aa219208SBarry Smith DMDAPeriodicType pt; 19547c6ae99SBarry Smith 19647c6ae99SBarry Smith PetscFunctionBegin; 197aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 198aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 199aa219208SBarry Smith if (pt == DMDA_XPERIODIC) { 20047c6ae99SBarry Smith ratio = mx/Mx; 20147c6ae99SBarry 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); 20247c6ae99SBarry Smith } else { 20347c6ae99SBarry Smith ratio = (mx-1)/(Mx-1); 20447c6ae99SBarry 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); 20547c6ae99SBarry Smith } 20647c6ae99SBarry Smith 207aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); 208aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); 209aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 21047c6ae99SBarry Smith 211aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr); 212aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr); 213aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 21447c6ae99SBarry Smith 21547c6ae99SBarry Smith /* create interpolation matrix */ 21647c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr); 21747c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr); 21847c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 21947c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,2,PETSC_NULL);CHKERRQ(ierr); 22047c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,2,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 22147c6ae99SBarry Smith 22247c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 22347c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 22447c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 22547c6ae99SBarry Smith row = idx_f[dof*(i-i_start_ghost)]/dof; 22647c6ae99SBarry Smith 22747c6ae99SBarry Smith i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 22847c6ae99SBarry Smith 22947c6ae99SBarry Smith /* 23047c6ae99SBarry Smith Only include those interpolation points that are truly 23147c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 23247c6ae99SBarry Smith in x direction; since they have no right neighbor 23347c6ae99SBarry Smith */ 23447c6ae99SBarry Smith x = ((double)(i - i_c*ratio))/((double)ratio); 23547c6ae99SBarry Smith nc = 0; 23647c6ae99SBarry Smith /* one left and below; or we are right on it */ 23747c6ae99SBarry Smith col = dof*(i_c-i_start_ghost_c); 23847c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 23947c6ae99SBarry Smith v[nc++] = - x + 1.0; 24047c6ae99SBarry Smith /* one right? */ 24147c6ae99SBarry Smith if (i_c*ratio != i) { 24247c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 24347c6ae99SBarry Smith v[nc++] = x; 24447c6ae99SBarry Smith } 24547c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 24647c6ae99SBarry Smith } 24747c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 24847c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 24947c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 25047c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 25147c6ae99SBarry Smith ierr = PetscLogFlops(5.0*m_f);CHKERRQ(ierr); 25247c6ae99SBarry Smith PetscFunctionReturn(0); 25347c6ae99SBarry Smith } 25447c6ae99SBarry Smith 25547c6ae99SBarry Smith #undef __FUNCT__ 25694013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_2D_Q1" 25794013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_2D_Q1(DM dac,DM daf,Mat *A) 25847c6ae99SBarry Smith { 25947c6ae99SBarry Smith PetscErrorCode ierr; 26047c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 26147c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz; 26247c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj; 26347c6ae99SBarry 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; 26447c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 26547c6ae99SBarry Smith PetscScalar v[4],x,y; 26647c6ae99SBarry Smith Mat mat; 267aa219208SBarry Smith DMDAPeriodicType pt; 268aa219208SBarry Smith DMDACoor2d **coors = 0,**ccoors; 26947c6ae99SBarry Smith Vec vcoors,cvcoors; 27047c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data, *ddf = (DM_DA*)daf->data; 27147c6ae99SBarry Smith 27247c6ae99SBarry Smith PetscFunctionBegin; 273aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 274aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 275aa219208SBarry Smith if (DMDAXPeriodic(pt)){ 27647c6ae99SBarry Smith ratioi = mx/Mx; 27747c6ae99SBarry 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); 27847c6ae99SBarry Smith } else { 27947c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 28047c6ae99SBarry 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); 28147c6ae99SBarry Smith } 282aa219208SBarry Smith if (DMDAYPeriodic(pt)){ 28347c6ae99SBarry Smith ratioj = my/My; 28447c6ae99SBarry 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); 28547c6ae99SBarry Smith } else { 28647c6ae99SBarry Smith ratioj = (my-1)/(My-1); 28747c6ae99SBarry 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); 28847c6ae99SBarry Smith } 28947c6ae99SBarry Smith 29047c6ae99SBarry Smith 291aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 292aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 293aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 29447c6ae99SBarry Smith 295aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 296aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 297aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 29847c6ae99SBarry Smith 29947c6ae99SBarry Smith /* 300aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 30147c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 30247c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 30347c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 30447c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 30547c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 30647c6ae99SBarry Smith 30747c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 30847c6ae99SBarry Smith */ 30947c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr); 31047c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr); 31147c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr); 31247c6ae99SBarry Smith col_scale = size_f/size_c; 31347c6ae99SBarry Smith col_shift = Mx*My*(rank_f/size_c); 31447c6ae99SBarry Smith 31547c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr); 31647c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 31747c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 31847c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 31947c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 32047c6ae99SBarry Smith 32147c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 32247c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 32347c6ae99SBarry Smith 324aa219208SBarry 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\ 32547c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 326aa219208SBarry 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\ 32747c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 32847c6ae99SBarry Smith 32947c6ae99SBarry Smith /* 33047c6ae99SBarry Smith Only include those interpolation points that are truly 33147c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 33247c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 33347c6ae99SBarry Smith */ 33447c6ae99SBarry Smith nc = 0; 33547c6ae99SBarry Smith /* one left and below; or we are right on it */ 33647c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 33747c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 33847c6ae99SBarry Smith /* one right and below */ 33947c6ae99SBarry Smith if (i_c*ratioi != i) { 34047c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col+dof]/dof; 34147c6ae99SBarry Smith } 34247c6ae99SBarry Smith /* one left and above */ 34347c6ae99SBarry Smith if (j_c*ratioj != j) { 34447c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col+m_ghost_c*dof]/dof; 34547c6ae99SBarry Smith } 34647c6ae99SBarry Smith /* one right and above */ 34747c6ae99SBarry Smith if (j_c*ratioi != j && i_c*ratioj != i) { 34847c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; 34947c6ae99SBarry Smith } 35047c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 35147c6ae99SBarry Smith } 35247c6ae99SBarry Smith } 35347c6ae99SBarry Smith ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr); 35447c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr); 35547c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 35647c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 35747c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 35847c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 35947c6ae99SBarry Smith 360aa219208SBarry Smith ierr = DMDAGetCoordinates(daf,&vcoors);CHKERRQ(ierr); 36147c6ae99SBarry Smith if (vcoors) { 362aa219208SBarry Smith ierr = DMDAGetGhostedCoordinates(dac,&cvcoors);CHKERRQ(ierr); 363aa219208SBarry Smith ierr = DMDAVecGetArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 364aa219208SBarry Smith ierr = DMDAVecGetArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 36547c6ae99SBarry Smith } 36647c6ae99SBarry Smith 36747c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 368*b3bd94feSDave May if (!vcoors) { 369*b3bd94feSDave May 37047c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 37147c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 37247c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 37347c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 37447c6ae99SBarry Smith 37547c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 37647c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 37747c6ae99SBarry Smith 37847c6ae99SBarry Smith /* 37947c6ae99SBarry Smith Only include those interpolation points that are truly 38047c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 38147c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 38247c6ae99SBarry Smith */ 38347c6ae99SBarry Smith x = ((double)(i - i_c*ratioi))/((double)ratioi); 38447c6ae99SBarry Smith y = ((double)(j - j_c*ratioj))/((double)ratioj); 385*b3bd94feSDave May 38647c6ae99SBarry Smith nc = 0; 38747c6ae99SBarry Smith /* one left and below; or we are right on it */ 38847c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 38947c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 39047c6ae99SBarry Smith v[nc++] = x*y - x - y + 1.0; 39147c6ae99SBarry Smith /* one right and below */ 39247c6ae99SBarry Smith if (i_c*ratioi != i) { 39347c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+dof]/dof; 39447c6ae99SBarry Smith v[nc++] = -x*y + x; 39547c6ae99SBarry Smith } 39647c6ae99SBarry Smith /* one left and above */ 39747c6ae99SBarry Smith if (j_c*ratioj != j) { 39847c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+m_ghost_c*dof]/dof; 39947c6ae99SBarry Smith v[nc++] = -x*y + y; 40047c6ae99SBarry Smith } 40147c6ae99SBarry Smith /* one right and above */ 40247c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 40347c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; 40447c6ae99SBarry Smith v[nc++] = x*y; 40547c6ae99SBarry Smith } 40647c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 40747c6ae99SBarry Smith } 40847c6ae99SBarry Smith } 409*b3bd94feSDave May 410*b3bd94feSDave May } else { 411*b3bd94feSDave May PetscScalar Ni[4]; 412*b3bd94feSDave May PetscScalar *xi,*eta; 413*b3bd94feSDave May PetscInt li,nxi,lj,neta; 414*b3bd94feSDave May 415*b3bd94feSDave May /* compute local coordinate arrays */ 416*b3bd94feSDave May nxi = ratioi + 1; 417*b3bd94feSDave May neta = ratioj + 1; 418*b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr); 419*b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*neta,&eta);CHKERRQ(ierr); 420*b3bd94feSDave May for (li=0; li<nxi; li++) { 421*b3bd94feSDave May xi[li] = -1.0 + li*(2.0/(PetscScalar)(nxi-1)); 422*b3bd94feSDave May } 423*b3bd94feSDave May for (lj=0; lj<neta; lj++) { 424*b3bd94feSDave May eta[lj] = -1.0 + lj*(2.0/(PetscScalar)(neta-1)); 425*b3bd94feSDave May } 426*b3bd94feSDave May 427*b3bd94feSDave May /* loop over local fine grid nodes setting interpolation for those*/ 428*b3bd94feSDave May for (j=j_start; j<j_start+n_f; j++) { 429*b3bd94feSDave May for (i=i_start; i<i_start+m_f; i++) { 430*b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */ 431*b3bd94feSDave May row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 432*b3bd94feSDave May 433*b3bd94feSDave May i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 434*b3bd94feSDave May j_c = (j/ratioj); /* coarse grid node below fine grid node */ 435*b3bd94feSDave May 436*b3bd94feSDave May /* remainders */ 437*b3bd94feSDave May li = i - ratioi * (i/ratioi); 438*b3bd94feSDave May if (i==mx-1){ li = nxi-1; } 439*b3bd94feSDave May lj = j - ratioj * (j/ratioj); 440*b3bd94feSDave May if (j==my-1){ lj = neta-1; } 441*b3bd94feSDave May 442*b3bd94feSDave May /* corners */ 443*b3bd94feSDave May col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 444*b3bd94feSDave May cols[0] = col_shift + idx_c[col]/dof; /* left, below */ 445*b3bd94feSDave May Ni[0] = 1.0; 446*b3bd94feSDave May if ( (li==0) || (li==nxi-1) ) { 447*b3bd94feSDave May if ( (lj==0) || (lj==neta-1) ) { 448*b3bd94feSDave May ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr); 449*b3bd94feSDave May continue; 450*b3bd94feSDave May } 451*b3bd94feSDave May } 452*b3bd94feSDave May 453*b3bd94feSDave May /* edges + interior */ 454*b3bd94feSDave May /* remainders */ 455*b3bd94feSDave May if (i==mx-1){ i_c--; } 456*b3bd94feSDave May if (j==my-1){ j_c--; } 457*b3bd94feSDave May 458*b3bd94feSDave May col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 459*b3bd94feSDave May cols[0] = col_shift + idx_c[col]/dof; /* left, below */ 460*b3bd94feSDave May cols[1] = col_shift + idx_c[col+dof]/dof; /* right, below */ 461*b3bd94feSDave May cols[2] = col_shift + idx_c[col+m_ghost_c*dof]/dof; /* left, above */ 462*b3bd94feSDave May cols[3] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; /* right, above */ 463*b3bd94feSDave May 464*b3bd94feSDave May Ni[0] = 0.25*(1.0-xi[li])*(1.0-eta[lj]); 465*b3bd94feSDave May Ni[1] = 0.25*(1.0+xi[li])*(1.0-eta[lj]); 466*b3bd94feSDave May Ni[2] = 0.25*(1.0-xi[li])*(1.0+eta[lj]); 467*b3bd94feSDave May Ni[3] = 0.25*(1.0+xi[li])*(1.0+eta[lj]); 468*b3bd94feSDave May 469*b3bd94feSDave May nc = 0; 470*b3bd94feSDave May if( PetscAbsScalar(Ni[0])<1.0e-32) { cols[0]=-1; } 471*b3bd94feSDave May if( PetscAbsScalar(Ni[1])<1.0e-32) { cols[1]=-1; } 472*b3bd94feSDave May if( PetscAbsScalar(Ni[2])<1.0e-32) { cols[2]=-1; } 473*b3bd94feSDave May if( PetscAbsScalar(Ni[3])<1.0e-32) { cols[3]=-1; } 474*b3bd94feSDave May 475*b3bd94feSDave May ierr = MatSetValues(mat,1,&row,4,cols,Ni,INSERT_VALUES);CHKERRQ(ierr); 476*b3bd94feSDave May } 477*b3bd94feSDave May } 478*b3bd94feSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 479*b3bd94feSDave May ierr = PetscFree(eta);CHKERRQ(ierr); 480*b3bd94feSDave May } 48147c6ae99SBarry Smith if (vcoors) { 482aa219208SBarry Smith ierr = DMDAVecRestoreArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 483aa219208SBarry Smith ierr = DMDAVecRestoreArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 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 PetscFunctionReturn(0); 49047c6ae99SBarry Smith } 49147c6ae99SBarry Smith 49247c6ae99SBarry Smith /* 49347c6ae99SBarry Smith Contributed by Andrei Draganescu <aidraga@sandia.gov> 49447c6ae99SBarry Smith */ 49547c6ae99SBarry Smith #undef __FUNCT__ 49694013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_2D_Q0" 49794013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_2D_Q0(DM dac,DM daf,Mat *A) 49847c6ae99SBarry Smith { 49947c6ae99SBarry Smith PetscErrorCode ierr; 50047c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 50147c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz; 50247c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj; 50347c6ae99SBarry 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; 50447c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 50547c6ae99SBarry Smith PetscScalar v[4]; 50647c6ae99SBarry Smith Mat mat; 507aa219208SBarry Smith DMDAPeriodicType pt; 50847c6ae99SBarry Smith 50947c6ae99SBarry Smith PetscFunctionBegin; 510aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 511aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 512aa219208SBarry Smith if (DMDAXPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x"); 513aa219208SBarry Smith if (DMDAYPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y"); 51447c6ae99SBarry Smith ratioi = mx/Mx; 51547c6ae99SBarry Smith ratioj = my/My; 51647c6ae99SBarry 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"); 51747c6ae99SBarry 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"); 51847c6ae99SBarry Smith if (ratioi != 2) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 2"); 51947c6ae99SBarry Smith if (ratioj != 2) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 2"); 52047c6ae99SBarry Smith 521aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 522aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 523aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 52447c6ae99SBarry Smith 525aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 526aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 527aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 52847c6ae99SBarry Smith 52947c6ae99SBarry Smith /* 530aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 53147c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 53247c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 53347c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 53447c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 53547c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 53647c6ae99SBarry Smith 53747c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 53847c6ae99SBarry Smith */ 53947c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr); 54047c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr); 54147c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr); 54247c6ae99SBarry Smith col_scale = size_f/size_c; 54347c6ae99SBarry Smith col_shift = Mx*My*(rank_f/size_c); 54447c6ae99SBarry Smith 54547c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr); 54647c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 54747c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 54847c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 54947c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 55047c6ae99SBarry Smith 55147c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 55247c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 55347c6ae99SBarry Smith 554aa219208SBarry 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\ 55547c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 556aa219208SBarry 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\ 55747c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 55847c6ae99SBarry Smith 55947c6ae99SBarry Smith /* 56047c6ae99SBarry Smith Only include those interpolation points that are truly 56147c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 56247c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 56347c6ae99SBarry Smith */ 56447c6ae99SBarry Smith nc = 0; 56547c6ae99SBarry Smith /* one left and below; or we are right on it */ 56647c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 56747c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 56847c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 56947c6ae99SBarry Smith } 57047c6ae99SBarry Smith } 57147c6ae99SBarry Smith ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr); 57247c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr); 57347c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 57447c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 57547c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 57647c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 57747c6ae99SBarry Smith 57847c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 57947c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 58047c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 58147c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 58247c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 58347c6ae99SBarry Smith 58447c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 58547c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 58647c6ae99SBarry Smith nc = 0; 58747c6ae99SBarry Smith /* one left and below; or we are right on it */ 58847c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 58947c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 59047c6ae99SBarry Smith v[nc++] = 1.0; 59147c6ae99SBarry Smith 59247c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 59347c6ae99SBarry Smith } 59447c6ae99SBarry Smith } 59547c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 59647c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 59747c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 59847c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 59947c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f);CHKERRQ(ierr); 60047c6ae99SBarry Smith PetscFunctionReturn(0); 60147c6ae99SBarry Smith } 60247c6ae99SBarry Smith 60347c6ae99SBarry Smith /* 60447c6ae99SBarry Smith Contributed by Jianming Yang <jianming-yang@uiowa.edu> 60547c6ae99SBarry Smith */ 60647c6ae99SBarry Smith #undef __FUNCT__ 60794013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_3D_Q0" 60894013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_3D_Q0(DM dac,DM daf,Mat *A) 60947c6ae99SBarry Smith { 61047c6ae99SBarry Smith PetscErrorCode ierr; 61147c6ae99SBarry Smith PetscInt i,j,l,i_start,j_start,l_start,m_f,n_f,p_f,Mx,My,Mz,*idx_f,dof; 61247c6ae99SBarry Smith PetscInt m_ghost,n_ghost,p_ghost,*idx_c,m_ghost_c,n_ghost_c,p_ghost_c,nc,*dnz,*onz; 61347c6ae99SBarry 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; 61447c6ae99SBarry 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; 61547c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 61647c6ae99SBarry Smith PetscScalar v[8]; 61747c6ae99SBarry Smith Mat mat; 618aa219208SBarry Smith DMDAPeriodicType pt; 61947c6ae99SBarry Smith 62047c6ae99SBarry Smith PetscFunctionBegin; 621aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 622aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 623aa219208SBarry Smith if (DMDAXPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x"); 624aa219208SBarry Smith if (DMDAYPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y"); 625aa219208SBarry Smith if (DMDAZPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in z"); 62647c6ae99SBarry Smith ratioi = mx/Mx; 62747c6ae99SBarry Smith ratioj = my/My; 62847c6ae99SBarry Smith ratiol = mz/Mz; 62947c6ae99SBarry 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"); 63047c6ae99SBarry 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"); 63147c6ae99SBarry 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"); 63247c6ae99SBarry Smith if (ratioi != 2 && ratioi != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 1 or 2"); 63347c6ae99SBarry Smith if (ratioj != 2 && ratioj != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 1 or 2"); 63447c6ae99SBarry Smith if (ratiol != 2 && ratiol != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in z must be 1 or 2"); 63547c6ae99SBarry Smith 636aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 637aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 638aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 63947c6ae99SBarry Smith 640aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 641aa219208SBarry 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); 642aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 64347c6ae99SBarry Smith /* 644aa219208SBarry Smith Used for handling a coarse DMDA that lives on 1/4 the processors of the fine DMDA. 64547c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 64647c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 64747c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 64847c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 64947c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 65047c6ae99SBarry Smith 65147c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 65247c6ae99SBarry Smith */ 65347c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr); 65447c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr); 65547c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr); 65647c6ae99SBarry Smith col_scale = size_f/size_c; 65747c6ae99SBarry Smith col_shift = Mx*My*Mz*(rank_f/size_c); 65847c6ae99SBarry Smith 65947c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f*p_f,col_scale*m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr); 66047c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 66147c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 66247c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 66347c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 66447c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 66547c6ae99SBarry Smith 66647c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 66747c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 66847c6ae99SBarry Smith l_c = (l/ratiol); 66947c6ae99SBarry Smith 670aa219208SBarry 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\ 67147c6ae99SBarry Smith l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c); 672aa219208SBarry 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\ 67347c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 674aa219208SBarry 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\ 67547c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 67647c6ae99SBarry Smith 67747c6ae99SBarry Smith /* 67847c6ae99SBarry Smith Only include those interpolation points that are truly 67947c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 68047c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 68147c6ae99SBarry Smith */ 68247c6ae99SBarry Smith nc = 0; 68347c6ae99SBarry Smith /* one left and below; or we are right on it */ 68447c6ae99SBarry Smith col = dof*(m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c) + m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 68547c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 68647c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 68747c6ae99SBarry Smith } 68847c6ae99SBarry Smith } 68947c6ae99SBarry Smith } 69047c6ae99SBarry Smith ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr); 69147c6ae99SBarry 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); 69247c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 69347c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 69447c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 69547c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 69647c6ae99SBarry Smith 69747c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 69847c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 69947c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 70047c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 70147c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 70247c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 70347c6ae99SBarry Smith 70447c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 70547c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 70647c6ae99SBarry Smith l_c = (l/ratiol); 70747c6ae99SBarry Smith nc = 0; 70847c6ae99SBarry Smith /* one left and below; or we are right on it */ 70947c6ae99SBarry 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)); 71047c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 71147c6ae99SBarry Smith v[nc++] = 1.0; 71247c6ae99SBarry Smith 71347c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 71447c6ae99SBarry Smith } 71547c6ae99SBarry Smith } 71647c6ae99SBarry Smith } 71747c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 71847c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 71947c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 72047c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 72147c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f*p_f);CHKERRQ(ierr); 72247c6ae99SBarry Smith PetscFunctionReturn(0); 72347c6ae99SBarry Smith } 72447c6ae99SBarry Smith 72547c6ae99SBarry Smith #undef __FUNCT__ 72694013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA_3D_Q1" 72794013140SBarry Smith PetscErrorCode DMGetInterpolation_DA_3D_Q1(DM dac,DM daf,Mat *A) 72847c6ae99SBarry Smith { 72947c6ae99SBarry Smith PetscErrorCode ierr; 73047c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof,l; 73147c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,Mz,mz; 73247c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[8],mx,m_c,my,nc,ratioi,ratioj,ratiok; 73347c6ae99SBarry Smith PetscInt i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c; 73447c6ae99SBarry Smith PetscInt l_start,p_f,l_start_ghost,p_ghost,l_start_c,p_c; 73547c6ae99SBarry Smith PetscInt l_start_ghost_c,p_ghost_c,l_c,*dnz,*onz; 73647c6ae99SBarry Smith PetscScalar v[8],x,y,z; 73747c6ae99SBarry Smith Mat mat; 738aa219208SBarry Smith DMDAPeriodicType pt; 739aa219208SBarry Smith DMDACoor3d ***coors = 0,***ccoors; 74047c6ae99SBarry Smith Vec vcoors,cvcoors; 74147c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data, *ddf = (DM_DA*)daf->data; 74247c6ae99SBarry Smith 74347c6ae99SBarry Smith PetscFunctionBegin; 744aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 745aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 74647c6ae99SBarry Smith if (mx == Mx) { 74747c6ae99SBarry Smith ratioi = 1; 748aa219208SBarry Smith } else if (DMDAXPeriodic(pt)){ 74947c6ae99SBarry Smith ratioi = mx/Mx; 75047c6ae99SBarry 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); 75147c6ae99SBarry Smith } else { 75247c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 75347c6ae99SBarry 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); 75447c6ae99SBarry Smith } 75547c6ae99SBarry Smith if (my == My) { 75647c6ae99SBarry Smith ratioj = 1; 757aa219208SBarry Smith } else if (DMDAYPeriodic(pt)){ 75847c6ae99SBarry Smith ratioj = my/My; 75947c6ae99SBarry 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); 76047c6ae99SBarry Smith } else { 76147c6ae99SBarry Smith ratioj = (my-1)/(My-1); 76247c6ae99SBarry Smith if (ratioj*(My-1) != my-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Ratio between levels: (my - 1)/(My - 1) must be integer: my %D My %D",my,My); 76347c6ae99SBarry Smith } 76447c6ae99SBarry Smith if (mz == Mz) { 76547c6ae99SBarry Smith ratiok = 1; 766aa219208SBarry Smith } else if (DMDAZPeriodic(pt)){ 76747c6ae99SBarry Smith ratiok = mz/Mz; 76847c6ae99SBarry 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); 76947c6ae99SBarry Smith } else { 77047c6ae99SBarry Smith ratiok = (mz-1)/(Mz-1); 77147c6ae99SBarry 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); 77247c6ae99SBarry Smith } 77347c6ae99SBarry Smith 774aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 775aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 776aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 77747c6ae99SBarry Smith 778aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 779aa219208SBarry 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); 780aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 78147c6ae99SBarry Smith 78247c6ae99SBarry Smith /* create interpolation matrix, determining exact preallocation */ 78347c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)dac)->comm,m_f*n_f*p_f,m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr); 78447c6ae99SBarry Smith /* loop over local fine grid nodes counting interpolating points */ 78547c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 78647c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 78747c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 78847c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 78947c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 79047c6ae99SBarry Smith i_c = (i/ratioi); 79147c6ae99SBarry Smith j_c = (j/ratioj); 79247c6ae99SBarry Smith l_c = (l/ratiok); 793aa219208SBarry 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\ 79447c6ae99SBarry Smith l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c); 795aa219208SBarry 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\ 79647c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 797aa219208SBarry 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\ 79847c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 79947c6ae99SBarry Smith 80047c6ae99SBarry Smith /* 80147c6ae99SBarry Smith Only include those interpolation points that are truly 80247c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 80347c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 80447c6ae99SBarry Smith */ 80547c6ae99SBarry Smith nc = 0; 80647c6ae99SBarry 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)); 80747c6ae99SBarry Smith cols[nc++] = idx_c[col]/dof; 80847c6ae99SBarry Smith if (i_c*ratioi != i) { 80947c6ae99SBarry Smith cols[nc++] = idx_c[col+dof]/dof; 81047c6ae99SBarry Smith } 81147c6ae99SBarry Smith if (j_c*ratioj != j) { 81247c6ae99SBarry Smith cols[nc++] = idx_c[col+m_ghost_c*dof]/dof; 81347c6ae99SBarry Smith } 81447c6ae99SBarry Smith if (l_c*ratiok != l) { 81547c6ae99SBarry Smith cols[nc++] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; 81647c6ae99SBarry Smith } 81747c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 81847c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c+1)*dof]/dof; 81947c6ae99SBarry Smith } 82047c6ae99SBarry Smith if (j_c*ratioj != j && l_c*ratiok != l) { 82147c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof; 82247c6ae99SBarry Smith } 82347c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l) { 82447c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; 82547c6ae99SBarry Smith } 82647c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) { 82747c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; 82847c6ae99SBarry Smith } 82947c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 83047c6ae99SBarry Smith } 83147c6ae99SBarry Smith } 83247c6ae99SBarry Smith } 83347c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr); 83447c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f*p_f,m_c*n_c*p_c,mx*my*mz,Mx*My*Mz);CHKERRQ(ierr); 83547c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 83647c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 83747c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 83847c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 83947c6ae99SBarry Smith 840aa219208SBarry Smith ierr = DMDAGetCoordinates(daf,&vcoors);CHKERRQ(ierr); 84147c6ae99SBarry Smith if (vcoors) { 842aa219208SBarry Smith ierr = DMDAGetGhostedCoordinates(dac,&cvcoors);CHKERRQ(ierr); 843aa219208SBarry Smith ierr = DMDAVecGetArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 844aa219208SBarry Smith ierr = DMDAVecGetArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 84547c6ae99SBarry Smith } 84647c6ae99SBarry Smith 84747c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 848*b3bd94feSDave May if (!vcoors) { 849*b3bd94feSDave May 85047c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 85147c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 85247c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 85347c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 85447c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 85547c6ae99SBarry Smith 85647c6ae99SBarry Smith i_c = (i/ratioi); 85747c6ae99SBarry Smith j_c = (j/ratioj); 85847c6ae99SBarry Smith l_c = (l/ratiok); 85947c6ae99SBarry Smith 86047c6ae99SBarry Smith /* 86147c6ae99SBarry Smith Only include those interpolation points that are truly 86247c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 86347c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 86447c6ae99SBarry Smith */ 86547c6ae99SBarry Smith x = ((double)(i - i_c*ratioi))/((double)ratioi); 86647c6ae99SBarry Smith y = ((double)(j - j_c*ratioj))/((double)ratioj); 86747c6ae99SBarry Smith z = ((double)(l - l_c*ratiok))/((double)ratiok); 868*b3bd94feSDave May 86947c6ae99SBarry Smith nc = 0; 87047c6ae99SBarry Smith /* one left and below; or we are right on it */ 87147c6ae99SBarry 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)); 87247c6ae99SBarry Smith 87347c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 87447c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.)); 87547c6ae99SBarry Smith 87647c6ae99SBarry Smith if (i_c*ratioi != i) { 87747c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 87847c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.)); 87947c6ae99SBarry Smith } 88047c6ae99SBarry Smith 88147c6ae99SBarry Smith if (j_c*ratioj != j) { 88247c6ae99SBarry Smith cols[nc] = idx_c[col+m_ghost_c*dof]/dof; 88347c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.)); 88447c6ae99SBarry Smith } 88547c6ae99SBarry Smith 88647c6ae99SBarry Smith if (l_c*ratiok != l) { 88747c6ae99SBarry Smith cols[nc] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; 88847c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.)); 88947c6ae99SBarry Smith } 89047c6ae99SBarry Smith 89147c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 89247c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c+1)*dof]/dof; 89347c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.)); 89447c6ae99SBarry Smith } 89547c6ae99SBarry Smith 89647c6ae99SBarry Smith if (j_c*ratioj != j && l_c*ratiok != l) { 89747c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof; 89847c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.)); 89947c6ae99SBarry Smith } 90047c6ae99SBarry Smith 90147c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l) { 90247c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; 90347c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.)); 90447c6ae99SBarry Smith } 90547c6ae99SBarry Smith 90647c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) { 90747c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; 90847c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.)); 90947c6ae99SBarry Smith } 91047c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 91147c6ae99SBarry Smith } 91247c6ae99SBarry Smith } 91347c6ae99SBarry Smith } 914*b3bd94feSDave May 915*b3bd94feSDave May } else { 916*b3bd94feSDave May PetscScalar *xi,*eta,*zeta; 917*b3bd94feSDave May PetscInt li,nxi,lj,neta,lk,nzeta,n; 918*b3bd94feSDave May PetscScalar Ni[8]; 919*b3bd94feSDave May 920*b3bd94feSDave May /* compute local coordinate arrays */ 921*b3bd94feSDave May nxi = ratioi + 1; 922*b3bd94feSDave May neta = ratioj + 1; 923*b3bd94feSDave May nzeta = ratiok + 1; 924*b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*nxi,&xi);CHKERRQ(ierr); 925*b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*neta,&eta);CHKERRQ(ierr); 926*b3bd94feSDave May ierr = PetscMalloc(sizeof(PetscScalar)*nzeta,&zeta);CHKERRQ(ierr); 927*b3bd94feSDave May for (li=0; li<nxi; li++) { 928*b3bd94feSDave May xi[li] = -1.0 + li*(2.0/(PetscScalar)(nxi-1)); 929*b3bd94feSDave May } 930*b3bd94feSDave May for (lj=0; lj<neta; lj++) { 931*b3bd94feSDave May eta[lj] = -1.0 + lj*(2.0/(PetscScalar)(neta-1)); 932*b3bd94feSDave May } 933*b3bd94feSDave May for (lk=0; lk<nzeta; lk++) { 934*b3bd94feSDave May zeta[lk] = -1.0 + lk*(2.0/(PetscScalar)(nzeta-1)); 935*b3bd94feSDave May } 936*b3bd94feSDave May 937*b3bd94feSDave May for (l=l_start; l<l_start+p_f; l++) { 938*b3bd94feSDave May for (j=j_start; j<j_start+n_f; j++) { 939*b3bd94feSDave May for (i=i_start; i<i_start+m_f; i++) { 940*b3bd94feSDave May /* convert to local "natural" numbering and then to PETSc global numbering */ 941*b3bd94feSDave May row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 942*b3bd94feSDave May 943*b3bd94feSDave May i_c = (i/ratioi); 944*b3bd94feSDave May j_c = (j/ratioj); 945*b3bd94feSDave May l_c = (l/ratiok); 946*b3bd94feSDave May 947*b3bd94feSDave May /* remainders */ 948*b3bd94feSDave May li = i - ratioi * (i/ratioi); 949*b3bd94feSDave May if (i==mx-1){ li = nxi-1; } 950*b3bd94feSDave May lj = j - ratioj * (j/ratioj); 951*b3bd94feSDave May if (j==my-1){ lj = neta-1; } 952*b3bd94feSDave May lk = l - ratiok * (l/ratiok); 953*b3bd94feSDave May if (l==mz-1){ lk = nzeta-1; } 954*b3bd94feSDave May 955*b3bd94feSDave May /* corners */ 956*b3bd94feSDave May col = dof*(m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c)+m_ghost_c*(j_c-j_start_ghost_c)+(i_c-i_start_ghost_c)); 957*b3bd94feSDave May cols[0] = idx_c[col]/dof; 958*b3bd94feSDave May Ni[0] = 1.0; 959*b3bd94feSDave May if ( (li==0) || (li==nxi-1) ) { 960*b3bd94feSDave May if ( (lj==0) || (lj==neta-1) ) { 961*b3bd94feSDave May if ( (lk==0) || (lk==nzeta-1) ) { 962*b3bd94feSDave May ierr = MatSetValue(mat,row,cols[0],Ni[0],INSERT_VALUES);CHKERRQ(ierr); 963*b3bd94feSDave May continue; 964*b3bd94feSDave May } 965*b3bd94feSDave May } 966*b3bd94feSDave May } 967*b3bd94feSDave May 968*b3bd94feSDave May /* edges + interior */ 969*b3bd94feSDave May /* remainders */ 970*b3bd94feSDave May if (i==mx-1){ i_c--; } 971*b3bd94feSDave May if (j==my-1){ j_c--; } 972*b3bd94feSDave May if (l==mz-1){ l_c--; } 973*b3bd94feSDave May 974*b3bd94feSDave May col = dof*(m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c) + m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 975*b3bd94feSDave May cols[0] = idx_c[col]/dof; /* one left and below; or we are right on it */ 976*b3bd94feSDave May cols[1] = idx_c[col+dof]/dof; /* one right and below */ 977*b3bd94feSDave May cols[2] = idx_c[col+m_ghost_c*dof]/dof; /* one left and above */ 978*b3bd94feSDave May cols[3] = idx_c[col+(m_ghost_c+1)*dof]/dof; /* one right and above */ 979*b3bd94feSDave May 980*b3bd94feSDave May cols[4] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; /* one left and below and front; or we are right on it */ 981*b3bd94feSDave May cols[5] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; /* one right and below, and front */ 982*b3bd94feSDave May cols[6] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof;/* one left and above and front*/ 983*b3bd94feSDave May cols[7] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; /* one right and above and front */ 984*b3bd94feSDave May 985*b3bd94feSDave May Ni[0] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]); 986*b3bd94feSDave May Ni[1] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0-zeta[lk]); 987*b3bd94feSDave May Ni[2] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]); 988*b3bd94feSDave May Ni[3] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0-zeta[lk]); 989*b3bd94feSDave May 990*b3bd94feSDave May Ni[4] = 0.125*(1.0-xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]); 991*b3bd94feSDave May Ni[5] = 0.125*(1.0+xi[li])*(1.0-eta[lj])*(1.0+zeta[lk]); 992*b3bd94feSDave May Ni[6] = 0.125*(1.0-xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]); 993*b3bd94feSDave May Ni[7] = 0.125*(1.0+xi[li])*(1.0+eta[lj])*(1.0+zeta[lk]); 994*b3bd94feSDave May 995*b3bd94feSDave May for (n=0; n<8; n++) { 996*b3bd94feSDave May if( PetscAbsScalar(Ni[n])<1.0e-32) { cols[n]=-1; } 997*b3bd94feSDave May } 998*b3bd94feSDave May ierr = MatSetValues(mat,1,&row,8,cols,Ni,INSERT_VALUES);CHKERRQ(ierr); 999*b3bd94feSDave May 1000*b3bd94feSDave May } 1001*b3bd94feSDave May } 1002*b3bd94feSDave May } 1003*b3bd94feSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 1004*b3bd94feSDave May ierr = PetscFree(eta);CHKERRQ(ierr); 1005*b3bd94feSDave May ierr = PetscFree(zeta);CHKERRQ(ierr); 1006*b3bd94feSDave May } 1007*b3bd94feSDave May 100847c6ae99SBarry Smith if (vcoors) { 1009aa219208SBarry Smith ierr = DMDAVecRestoreArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 1010aa219208SBarry Smith ierr = DMDAVecRestoreArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 101147c6ae99SBarry Smith } 101247c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 101347c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 101447c6ae99SBarry Smith 101547c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 101647c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 101747c6ae99SBarry Smith PetscFunctionReturn(0); 101847c6ae99SBarry Smith } 101947c6ae99SBarry Smith 102047c6ae99SBarry Smith #undef __FUNCT__ 102194013140SBarry Smith #define __FUNCT__ "DMGetInterpolation_DA" 102294013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolation_DA(DM dac,DM daf,Mat *A,Vec *scale) 102347c6ae99SBarry Smith { 102447c6ae99SBarry Smith PetscErrorCode ierr; 102547c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 1026aa219208SBarry Smith DMDAPeriodicType wrapc,wrapf; 1027aa219208SBarry Smith DMDAStencilType stc,stf; 102847c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data; 102947c6ae99SBarry Smith 103047c6ae99SBarry Smith PetscFunctionBegin; 103147c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 103247c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 103347c6ae99SBarry Smith PetscValidPointer(A,3); 103447c6ae99SBarry Smith if (scale) PetscValidPointer(scale,4); 103547c6ae99SBarry Smith 1036aa219208SBarry Smith ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&wrapc,&stc);CHKERRQ(ierr); 1037aa219208SBarry Smith ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&wrapf,&stf);CHKERRQ(ierr); 1038aa219208SBarry Smith if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 1039aa219208SBarry Smith if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); 1040aa219208SBarry 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); 1041aa219208SBarry Smith if (wrapc != wrapf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Periodic type different in two DMDAs");CHKERRQ(ierr); 1042aa219208SBarry Smith if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 104347c6ae99SBarry Smith if (Mc < 2 && Mf > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction"); 104447c6ae99SBarry 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"); 104547c6ae99SBarry 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"); 104647c6ae99SBarry Smith 1047aa219208SBarry Smith if (ddc->interptype == DMDA_Q1){ 104847c6ae99SBarry Smith if (dimc == 1){ 104994013140SBarry Smith ierr = DMGetInterpolation_DA_1D_Q1(dac,daf,A);CHKERRQ(ierr); 105047c6ae99SBarry Smith } else if (dimc == 2){ 105194013140SBarry Smith ierr = DMGetInterpolation_DA_2D_Q1(dac,daf,A);CHKERRQ(ierr); 105247c6ae99SBarry Smith } else if (dimc == 3){ 105394013140SBarry Smith ierr = DMGetInterpolation_DA_3D_Q1(dac,daf,A);CHKERRQ(ierr); 105447c6ae99SBarry Smith } else { 1055aa219208SBarry Smith SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype); 105647c6ae99SBarry Smith } 1057aa219208SBarry Smith } else if (ddc->interptype == DMDA_Q0){ 105847c6ae99SBarry Smith if (dimc == 1){ 105994013140SBarry Smith ierr = DMGetInterpolation_DA_1D_Q0(dac,daf,A);CHKERRQ(ierr); 106047c6ae99SBarry Smith } else if (dimc == 2){ 106194013140SBarry Smith ierr = DMGetInterpolation_DA_2D_Q0(dac,daf,A);CHKERRQ(ierr); 106247c6ae99SBarry Smith } else if (dimc == 3){ 106394013140SBarry Smith ierr = DMGetInterpolation_DA_3D_Q0(dac,daf,A);CHKERRQ(ierr); 106447c6ae99SBarry Smith } else { 1065aa219208SBarry Smith SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DMDA dimension %D for interpolation type %d",dimc,(int)ddc->interptype); 106647c6ae99SBarry Smith } 106747c6ae99SBarry Smith } 106847c6ae99SBarry Smith if (scale) { 106947c6ae99SBarry Smith ierr = DMGetInterpolationScale((DM)dac,(DM)daf,*A,scale);CHKERRQ(ierr); 107047c6ae99SBarry Smith } 107147c6ae99SBarry Smith PetscFunctionReturn(0); 107247c6ae99SBarry Smith } 107347c6ae99SBarry Smith 107447c6ae99SBarry Smith #undef __FUNCT__ 107594013140SBarry Smith #define __FUNCT__ "DMGetInjection_DA_2D" 107694013140SBarry Smith PetscErrorCode DMGetInjection_DA_2D(DM dac,DM daf,VecScatter *inject) 107747c6ae99SBarry Smith { 107847c6ae99SBarry Smith PetscErrorCode ierr; 107947c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 108047c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c; 108147c6ae99SBarry Smith PetscInt row,i_start_ghost,j_start_ghost,mx,m_c,my,nc,ratioi,ratioj; 1082076202ddSJed Brown PetscInt i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c; 108347c6ae99SBarry Smith PetscInt *cols; 1084aa219208SBarry Smith DMDAPeriodicType pt; 108547c6ae99SBarry Smith Vec vecf,vecc; 108647c6ae99SBarry Smith IS isf; 108747c6ae99SBarry Smith 108847c6ae99SBarry Smith PetscFunctionBegin; 1089aa219208SBarry Smith ierr = DMDAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 1090aa219208SBarry Smith ierr = DMDAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 1091aa219208SBarry Smith if (DMDAXPeriodic(pt)){ 109247c6ae99SBarry Smith ratioi = mx/Mx; 109347c6ae99SBarry 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); 109447c6ae99SBarry Smith } else { 109547c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 109647c6ae99SBarry 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); 109747c6ae99SBarry Smith } 1098aa219208SBarry Smith if (DMDAYPeriodic(pt)){ 109947c6ae99SBarry Smith ratioj = my/My; 110047c6ae99SBarry 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); 110147c6ae99SBarry Smith } else { 110247c6ae99SBarry Smith ratioj = (my-1)/(My-1); 110347c6ae99SBarry 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); 110447c6ae99SBarry Smith } 110547c6ae99SBarry Smith 1106aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 1107aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 1108aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 110947c6ae99SBarry Smith 1110aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 1111aa219208SBarry Smith ierr = DMDAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 1112aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 111347c6ae99SBarry Smith 111447c6ae99SBarry Smith 111547c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 111647c6ae99SBarry Smith nc = 0; 111747c6ae99SBarry Smith ierr = PetscMalloc(n_f*m_f*sizeof(PetscInt),&cols);CHKERRQ(ierr); 1118076202ddSJed Brown for (j=j_start_c; j<j_start_c+n_c; j++) { 1119076202ddSJed Brown for (i=i_start_c; i<i_start_c+m_c; i++) { 1120076202ddSJed Brown PetscInt i_f = i*ratioi,j_f = j*ratioj; 1121076202ddSJed 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\ 1122076202ddSJed Brown j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost); 1123076202ddSJed 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\ 1124076202ddSJed Brown i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost); 1125076202ddSJed Brown row = idx_f[dof*(m_ghost*(j_f-j_start_ghost) + (i_f-i_start_ghost))]; 112647c6ae99SBarry Smith cols[nc++] = row/dof; 112747c6ae99SBarry Smith } 112847c6ae99SBarry Smith } 112947c6ae99SBarry Smith 113047c6ae99SBarry Smith ierr = ISCreateBlock(((PetscObject)daf)->comm,dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); 11319a42bb27SBarry Smith ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); 11329a42bb27SBarry Smith ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); 113347c6ae99SBarry Smith ierr = VecScatterCreate(vecf,isf,vecc,PETSC_NULL,inject);CHKERRQ(ierr); 11349a42bb27SBarry Smith ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); 11359a42bb27SBarry Smith ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); 113647c6ae99SBarry Smith ierr = ISDestroy(isf);CHKERRQ(ierr); 113747c6ae99SBarry Smith PetscFunctionReturn(0); 113847c6ae99SBarry Smith } 113947c6ae99SBarry Smith 114047c6ae99SBarry Smith #undef __FUNCT__ 1141b1757049SJed Brown #define __FUNCT__ "DMGetInjection_DA_3D" 1142b1757049SJed Brown PetscErrorCode DMGetInjection_DA_3D(DM dac,DM daf,VecScatter *inject) 1143b1757049SJed Brown { 1144b1757049SJed Brown PetscErrorCode ierr; 1145b1757049SJed Brown PetscInt i,j,k,i_start,j_start,k_start,m_f,n_f,p_f,Mx,My,Mz; 1146b1757049SJed Brown PetscInt m_ghost,n_ghost,p_ghost,m_ghost_c,n_ghost_c,p_ghost_c; 1147b1757049SJed Brown PetscInt i_start_ghost,j_start_ghost,k_start_ghost; 1148b1757049SJed Brown PetscInt mx,my,mz,ratioi,ratioj,ratiok; 1149b1757049SJed Brown PetscInt i_start_c,j_start_c,k_start_c; 1150b1757049SJed Brown PetscInt m_c,n_c,p_c; 1151b1757049SJed Brown PetscInt i_start_ghost_c,j_start_ghost_c,k_start_ghost_c; 1152b1757049SJed Brown PetscInt row,nc,dof; 1153b1757049SJed Brown PetscInt *idx_c,*idx_f; 1154b1757049SJed Brown PetscInt *cols; 1155b1757049SJed Brown DMDAPeriodicType pt; 1156b1757049SJed Brown Vec vecf,vecc; 1157b1757049SJed Brown IS isf; 1158b1757049SJed Brown 1159b1757049SJed Brown PetscFunctionBegin; 1160b1757049SJed Brown ierr = DMDAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 1161b1757049SJed Brown ierr = DMDAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 1162b1757049SJed Brown 1163b1757049SJed Brown if (DMDAXPeriodic(pt)){ 1164b1757049SJed Brown ratioi = mx/Mx; 1165b1757049SJed 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); 1166b1757049SJed Brown } else { 1167b1757049SJed Brown ratioi = (mx-1)/(Mx-1); 1168b1757049SJed 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); 1169b1757049SJed Brown } 1170b1757049SJed Brown if (DMDAYPeriodic(pt)){ 1171b1757049SJed Brown ratioj = my/My; 1172b1757049SJed 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); 1173b1757049SJed Brown } else { 1174b1757049SJed Brown ratioj = (my-1)/(My-1); 1175b1757049SJed 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); 1176b1757049SJed Brown } 1177b1757049SJed Brown if (DMDAZPeriodic(pt)){ 1178b1757049SJed Brown ratiok = mz/Mz; 1179b1757049SJed 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); 1180b1757049SJed Brown } else { 1181b1757049SJed Brown ratiok = (mz-1)/(Mz-1); 1182b1757049SJed 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); 1183b1757049SJed Brown } 1184b1757049SJed Brown 1185b1757049SJed Brown ierr = DMDAGetCorners(daf,&i_start,&j_start,&k_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 1186b1757049SJed Brown ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&k_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 1187b1757049SJed Brown ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 1188b1757049SJed Brown 1189b1757049SJed Brown ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&k_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 1190b1757049SJed 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); 1191b1757049SJed Brown ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 1192b1757049SJed Brown 1193b1757049SJed Brown 1194b1757049SJed Brown /* loop over local fine grid nodes setting interpolation for those*/ 1195b1757049SJed Brown nc = 0; 1196b1757049SJed Brown ierr = PetscMalloc(n_f*m_f*p_f*sizeof(PetscInt),&cols);CHKERRQ(ierr); 1197b1757049SJed Brown for (k=k_start_c; k<k_start_c+p_c; k++) { 1198b1757049SJed Brown for (j=j_start_c; j<j_start_c+n_c; j++) { 1199b1757049SJed Brown for (i=i_start_c; i<i_start_c+m_c; i++) { 1200b1757049SJed Brown PetscInt i_f = i*ratioi,j_f = j*ratioj,k_f = k*ratiok; 1201b1757049SJed 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 " 1202b1757049SJed Brown "k_c %D k_f %D fine ghost range [%D,%D]",k,k_f,k_start_ghost,k_start_ghost+p_ghost); 1203b1757049SJed 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 " 1204b1757049SJed Brown "j_c %D j_f %D fine ghost range [%D,%D]",j,j_f,j_start_ghost,j_start_ghost+n_ghost); 1205b1757049SJed 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 " 1206b1757049SJed Brown "i_c %D i_f %D fine ghost range [%D,%D]",i,i_f,i_start_ghost,i_start_ghost+m_ghost); 1207b1757049SJed 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))]; 1208b1757049SJed Brown cols[nc++] = row/dof; 1209b1757049SJed Brown } 1210b1757049SJed Brown } 1211b1757049SJed Brown } 1212b1757049SJed Brown 1213b1757049SJed Brown ierr = ISCreateBlock(((PetscObject)daf)->comm,dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); 1214b1757049SJed Brown ierr = DMGetGlobalVector(dac,&vecc);CHKERRQ(ierr); 1215b1757049SJed Brown ierr = DMGetGlobalVector(daf,&vecf);CHKERRQ(ierr); 1216b1757049SJed Brown ierr = VecScatterCreate(vecf,isf,vecc,PETSC_NULL,inject);CHKERRQ(ierr); 1217b1757049SJed Brown ierr = DMRestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); 1218b1757049SJed Brown ierr = DMRestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); 1219b1757049SJed Brown ierr = ISDestroy(isf);CHKERRQ(ierr); 1220b1757049SJed Brown PetscFunctionReturn(0); 1221b1757049SJed Brown } 1222b1757049SJed Brown 1223b1757049SJed Brown #undef __FUNCT__ 122494013140SBarry Smith #define __FUNCT__ "DMGetInjection_DA" 122594013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetInjection_DA(DM dac,DM daf,VecScatter *inject) 122647c6ae99SBarry Smith { 122747c6ae99SBarry Smith PetscErrorCode ierr; 122847c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 1229aa219208SBarry Smith DMDAPeriodicType wrapc,wrapf; 1230aa219208SBarry Smith DMDAStencilType stc,stf; 123147c6ae99SBarry Smith 123247c6ae99SBarry Smith PetscFunctionBegin; 123347c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 123447c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 123547c6ae99SBarry Smith PetscValidPointer(inject,3); 123647c6ae99SBarry Smith 1237aa219208SBarry Smith ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&wrapc,&stc);CHKERRQ(ierr); 1238aa219208SBarry Smith ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&wrapf,&stf);CHKERRQ(ierr); 1239aa219208SBarry Smith if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 1240aa219208SBarry Smith if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); 1241aa219208SBarry 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); 1242aa219208SBarry Smith if (wrapc != wrapf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Periodic type different in two DMDAs");CHKERRQ(ierr); 1243aa219208SBarry Smith if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 124447c6ae99SBarry Smith if (Mc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction"); 124547c6ae99SBarry Smith if (dimc > 1 && Nc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in y direction"); 124647c6ae99SBarry Smith if (dimc > 2 && Pc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in z direction"); 124747c6ae99SBarry Smith 124847c6ae99SBarry Smith if (dimc == 2){ 124994013140SBarry Smith ierr = DMGetInjection_DA_2D(dac,daf,inject);CHKERRQ(ierr); 1250b1757049SJed Brown } else if (dimc == 3) { 1251b1757049SJed Brown ierr = DMGetInjection_DA_3D(dac,daf,inject);CHKERRQ(ierr); 125247c6ae99SBarry Smith } else { 1253aa219208SBarry Smith SETERRQ1(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DMDA dimension %D",dimc); 125447c6ae99SBarry Smith } 125547c6ae99SBarry Smith PetscFunctionReturn(0); 125647c6ae99SBarry Smith } 125747c6ae99SBarry Smith 125847c6ae99SBarry Smith #undef __FUNCT__ 125994013140SBarry Smith #define __FUNCT__ "DMGetAggregates_DA" 126094013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetAggregates_DA(DM dac,DM daf,Mat *rest) 126147c6ae99SBarry Smith { 126247c6ae99SBarry Smith PetscErrorCode ierr; 126347c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc; 126447c6ae99SBarry Smith PetscInt dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 1265aa219208SBarry Smith DMDAPeriodicType wrapc,wrapf; 1266aa219208SBarry Smith DMDAStencilType stc,stf; 126747c6ae99SBarry Smith PetscInt i,j,l; 126847c6ae99SBarry Smith PetscInt i_start,j_start,l_start, m_f,n_f,p_f; 126947c6ae99SBarry Smith PetscInt i_start_ghost,j_start_ghost,l_start_ghost,m_ghost,n_ghost,p_ghost; 127047c6ae99SBarry Smith PetscInt *idx_f; 127147c6ae99SBarry Smith PetscInt i_c,j_c,l_c; 127247c6ae99SBarry Smith PetscInt i_start_c,j_start_c,l_start_c, m_c,n_c,p_c; 127347c6ae99SBarry Smith PetscInt i_start_ghost_c,j_start_ghost_c,l_start_ghost_c,m_ghost_c,n_ghost_c,p_ghost_c; 127447c6ae99SBarry Smith PetscInt *idx_c; 127547c6ae99SBarry Smith PetscInt d; 127647c6ae99SBarry Smith PetscInt a; 127747c6ae99SBarry Smith PetscInt max_agg_size; 127847c6ae99SBarry Smith PetscInt *fine_nodes; 127947c6ae99SBarry Smith PetscScalar *one_vec; 128047c6ae99SBarry Smith PetscInt fn_idx; 128147c6ae99SBarry Smith 128247c6ae99SBarry Smith PetscFunctionBegin; 128347c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 128447c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 128547c6ae99SBarry Smith PetscValidPointer(rest,3); 128647c6ae99SBarry Smith 1287aa219208SBarry Smith ierr = DMDAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&wrapc,&stc);CHKERRQ(ierr); 1288aa219208SBarry Smith ierr = DMDAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&wrapf,&stf);CHKERRQ(ierr); 1289aa219208SBarry Smith if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DMDA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 1290aa219208SBarry Smith if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DMDA do not match %D %D",dofc,doff);CHKERRQ(ierr); 1291aa219208SBarry 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); 1292aa219208SBarry Smith if (wrapc != wrapf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Periodic type different in two DMDAs");CHKERRQ(ierr); 1293aa219208SBarry Smith if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DMDAs");CHKERRQ(ierr); 129447c6ae99SBarry Smith 129547c6ae99SBarry 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); 129647c6ae99SBarry 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); 129747c6ae99SBarry 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); 129847c6ae99SBarry Smith 129947c6ae99SBarry Smith if (Pc < 0) Pc = 1; 130047c6ae99SBarry Smith if (Pf < 0) Pf = 1; 130147c6ae99SBarry Smith if (Nc < 0) Nc = 1; 130247c6ae99SBarry Smith if (Nf < 0) Nf = 1; 130347c6ae99SBarry Smith 1304aa219208SBarry Smith ierr = DMDAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 1305aa219208SBarry Smith ierr = DMDAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 1306aa219208SBarry Smith ierr = DMDAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 130747c6ae99SBarry Smith 1308aa219208SBarry Smith ierr = DMDAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 1309aa219208SBarry 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); 1310aa219208SBarry Smith ierr = DMDAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 131147c6ae99SBarry Smith 131247c6ae99SBarry Smith /* 131347c6ae99SBarry Smith Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios 131447c6ae99SBarry Smith for dimension 1 and 2 respectively. 131547c6ae99SBarry Smith Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1) 131647c6ae99SBarry Smith and r_y*j and r_y*(j+1) will be grouped into the same coarse grid agregate. 131747c6ae99SBarry Smith Each specific dof on the fine grid is mapped to one dof on the coarse grid. 131847c6ae99SBarry Smith */ 131947c6ae99SBarry Smith 132047c6ae99SBarry Smith max_agg_size = (Mf/Mc+1)*(Nf/Nc+1)*(Pf/Pc+1); 132147c6ae99SBarry Smith 132247c6ae99SBarry Smith /* create the matrix that will contain the restriction operator */ 132347c6ae99SBarry 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, 132447c6ae99SBarry Smith max_agg_size, PETSC_NULL, max_agg_size, PETSC_NULL, rest);CHKERRQ(ierr); 132547c6ae99SBarry Smith 132647c6ae99SBarry Smith /* store nodes in the fine grid here */ 132747c6ae99SBarry Smith ierr = PetscMalloc2(max_agg_size,PetscScalar, &one_vec,max_agg_size,PetscInt, &fine_nodes);CHKERRQ(ierr); 132847c6ae99SBarry Smith for(i=0; i<max_agg_size; i++) one_vec[i] = 1.0; 132947c6ae99SBarry Smith 133047c6ae99SBarry Smith /* loop over all coarse nodes */ 133147c6ae99SBarry Smith for (l_c=l_start_c; l_c<l_start_c+p_c; l_c++) { 133247c6ae99SBarry Smith for (j_c=j_start_c; j_c<j_start_c+n_c; j_c++) { 133347c6ae99SBarry Smith for (i_c=i_start_c; i_c<i_start_c+m_c; i_c++) { 133447c6ae99SBarry Smith for(d=0; d<dofc; d++) { 133547c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 133647c6ae99SBarry 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; 133747c6ae99SBarry Smith 133847c6ae99SBarry Smith fn_idx = 0; 133947c6ae99SBarry Smith /* Corresponding fine points are all points (i_f, j_f, l_f) such that 134047c6ae99SBarry Smith i_c*Mf/Mc <= i_f < (i_c+1)*Mf/Mc 134147c6ae99SBarry Smith (same for other dimensions) 134247c6ae99SBarry Smith */ 134347c6ae99SBarry Smith for (l=l_c*Pf/Pc; l<PetscMin((l_c+1)*Pf/Pc,Pf); l++) { 134447c6ae99SBarry Smith for (j=j_c*Nf/Nc; j<PetscMin((j_c+1)*Nf/Nc,Nf); j++) { 134547c6ae99SBarry Smith for (i=i_c*Mf/Mc; i<PetscMin((i_c+1)*Mf/Mc,Mf); i++) { 134647c6ae99SBarry 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; 134747c6ae99SBarry Smith fn_idx++; 134847c6ae99SBarry Smith } 134947c6ae99SBarry Smith } 135047c6ae99SBarry Smith } 135147c6ae99SBarry Smith /* add all these points to one aggregate */ 135247c6ae99SBarry Smith ierr = MatSetValues(*rest, 1, &a, fn_idx, fine_nodes, one_vec, INSERT_VALUES);CHKERRQ(ierr); 135347c6ae99SBarry Smith } 135447c6ae99SBarry Smith } 135547c6ae99SBarry Smith } 135647c6ae99SBarry Smith } 135747c6ae99SBarry Smith ierr = PetscFree2(one_vec,fine_nodes);CHKERRQ(ierr); 135847c6ae99SBarry Smith ierr = MatAssemblyBegin(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 135947c6ae99SBarry Smith ierr = MatAssemblyEnd(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 136047c6ae99SBarry Smith PetscFunctionReturn(0); 136147c6ae99SBarry Smith } 1362