1*47c6ae99SBarry Smith #define PETSCDM_DLL 2*47c6ae99SBarry Smith 3*47c6ae99SBarry Smith /* 4*47c6ae99SBarry Smith Code for interpolating between grids represented by DAs 5*47c6ae99SBarry Smith */ 6*47c6ae99SBarry Smith 7*47c6ae99SBarry Smith #include "private/daimpl.h" /*I "petscda.h" I*/ 8*47c6ae99SBarry Smith #include "petscmg.h" 9*47c6ae99SBarry Smith 10*47c6ae99SBarry Smith #undef __FUNCT__ 11*47c6ae99SBarry Smith #define __FUNCT__ "DMGetInterpolationScale" 12*47c6ae99SBarry Smith /*@ 13*47c6ae99SBarry Smith DMGetInterpolationScale - Forms L = R*1/diag(R*1) - L.*v is like a coarse grid average of the 14*47c6ae99SBarry Smith nearby fine grid points. 15*47c6ae99SBarry Smith 16*47c6ae99SBarry Smith Input Parameters: 17*47c6ae99SBarry Smith + dac - DM that defines a coarse mesh 18*47c6ae99SBarry Smith . daf - DM that defines a fine mesh 19*47c6ae99SBarry Smith - mat - the restriction (or interpolation operator) from fine to coarse 20*47c6ae99SBarry Smith 21*47c6ae99SBarry Smith Output Parameter: 22*47c6ae99SBarry Smith . scale - the scaled vector 23*47c6ae99SBarry Smith 24*47c6ae99SBarry Smith Level: developer 25*47c6ae99SBarry Smith 26*47c6ae99SBarry Smith .seealso: DMGetInterpolation() 27*47c6ae99SBarry Smith 28*47c6ae99SBarry Smith @*/ 29*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolationScale(DM dac,DM daf,Mat mat,Vec *scale) 30*47c6ae99SBarry Smith { 31*47c6ae99SBarry Smith PetscErrorCode ierr; 32*47c6ae99SBarry Smith Vec fine; 33*47c6ae99SBarry Smith PetscScalar one = 1.0; 34*47c6ae99SBarry Smith 35*47c6ae99SBarry Smith PetscFunctionBegin; 36*47c6ae99SBarry Smith ierr = DMCreateGlobalVector(daf,&fine);CHKERRQ(ierr); 37*47c6ae99SBarry Smith ierr = DMCreateGlobalVector(dac,scale);CHKERRQ(ierr); 38*47c6ae99SBarry Smith ierr = VecSet(fine,one);CHKERRQ(ierr); 39*47c6ae99SBarry Smith ierr = MatRestrict(mat,fine,*scale);CHKERRQ(ierr); 40*47c6ae99SBarry Smith ierr = VecDestroy(fine);CHKERRQ(ierr); 41*47c6ae99SBarry Smith ierr = VecReciprocal(*scale);CHKERRQ(ierr); 42*47c6ae99SBarry Smith PetscFunctionReturn(0); 43*47c6ae99SBarry Smith } 44*47c6ae99SBarry Smith 45*47c6ae99SBarry Smith #undef __FUNCT__ 46*47c6ae99SBarry Smith #define __FUNCT__ "DAGetInterpolation_1D_Q1" 47*47c6ae99SBarry Smith PetscErrorCode DAGetInterpolation_1D_Q1(DA dac,DA daf,Mat *A) 48*47c6ae99SBarry Smith { 49*47c6ae99SBarry Smith PetscErrorCode ierr; 50*47c6ae99SBarry Smith PetscInt i,i_start,m_f,Mx,*idx_f; 51*47c6ae99SBarry Smith PetscInt m_ghost,*idx_c,m_ghost_c; 52*47c6ae99SBarry Smith PetscInt row,col,i_start_ghost,mx,m_c,nc,ratio; 53*47c6ae99SBarry Smith PetscInt i_c,i_start_c,i_start_ghost_c,cols[2],dof; 54*47c6ae99SBarry Smith PetscScalar v[2],x,*coors = 0,*ccoors; 55*47c6ae99SBarry Smith Mat mat; 56*47c6ae99SBarry Smith DAPeriodicType pt; 57*47c6ae99SBarry Smith Vec vcoors,cvcoors; 58*47c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data, *ddf = (DM_DA*)daf->data; 59*47c6ae99SBarry Smith 60*47c6ae99SBarry Smith PetscFunctionBegin; 61*47c6ae99SBarry Smith ierr = DAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 62*47c6ae99SBarry Smith ierr = DAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 63*47c6ae99SBarry Smith if (pt == DA_XPERIODIC) { 64*47c6ae99SBarry Smith ratio = mx/Mx; 65*47c6ae99SBarry 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); 66*47c6ae99SBarry Smith } else { 67*47c6ae99SBarry Smith ratio = (mx-1)/(Mx-1); 68*47c6ae99SBarry 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); 69*47c6ae99SBarry Smith } 70*47c6ae99SBarry Smith 71*47c6ae99SBarry Smith ierr = DAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); 72*47c6ae99SBarry Smith ierr = DAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); 73*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 74*47c6ae99SBarry Smith 75*47c6ae99SBarry Smith ierr = DAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr); 76*47c6ae99SBarry Smith ierr = DAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr); 77*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 78*47c6ae99SBarry Smith 79*47c6ae99SBarry Smith /* create interpolation matrix */ 80*47c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr); 81*47c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr); 82*47c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 83*47c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,2,PETSC_NULL);CHKERRQ(ierr); 84*47c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,2,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 85*47c6ae99SBarry Smith 86*47c6ae99SBarry Smith ierr = DAGetCoordinates(daf,&vcoors);CHKERRQ(ierr); 87*47c6ae99SBarry Smith if (vcoors) { 88*47c6ae99SBarry Smith ierr = DAGetGhostedCoordinates(dac,&cvcoors);CHKERRQ(ierr); 89*47c6ae99SBarry Smith ierr = DAVecGetArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 90*47c6ae99SBarry Smith ierr = DAVecGetArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 91*47c6ae99SBarry Smith } 92*47c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 93*47c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 94*47c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 95*47c6ae99SBarry Smith row = idx_f[dof*(i-i_start_ghost)]/dof; 96*47c6ae99SBarry Smith 97*47c6ae99SBarry Smith i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 98*47c6ae99SBarry Smith if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 99*47c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 100*47c6ae99SBarry Smith 101*47c6ae99SBarry Smith /* 102*47c6ae99SBarry Smith Only include those interpolation points that are truly 103*47c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 104*47c6ae99SBarry Smith in x direction; since they have no right neighbor 105*47c6ae99SBarry Smith */ 106*47c6ae99SBarry Smith if (coors) { 107*47c6ae99SBarry Smith x = (coors[i] - ccoors[i_c]); 108*47c6ae99SBarry Smith /* only access the next coors point if we know there is one */ 109*47c6ae99SBarry Smith /* note this is dangerous because x may not exactly equal ZERO */ 110*47c6ae99SBarry Smith if (PetscAbsScalar(x) != 0.0) x = x/(ccoors[i_c+1] - ccoors[i_c]); 111*47c6ae99SBarry Smith } else { 112*47c6ae99SBarry Smith x = ((double)(i - i_c*ratio))/((double)ratio); 113*47c6ae99SBarry Smith } 114*47c6ae99SBarry Smith nc = 0; 115*47c6ae99SBarry Smith /* one left and below; or we are right on it */ 116*47c6ae99SBarry Smith col = dof*(i_c-i_start_ghost_c); 117*47c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 118*47c6ae99SBarry Smith v[nc++] = - x + 1.0; 119*47c6ae99SBarry Smith /* one right? */ 120*47c6ae99SBarry Smith if (i_c*ratio != i) { 121*47c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 122*47c6ae99SBarry Smith v[nc++] = x; 123*47c6ae99SBarry Smith } 124*47c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 125*47c6ae99SBarry Smith } 126*47c6ae99SBarry Smith if (vcoors) { 127*47c6ae99SBarry Smith ierr = DAVecRestoreArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 128*47c6ae99SBarry Smith ierr = DAVecRestoreArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 129*47c6ae99SBarry Smith } 130*47c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 131*47c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 132*47c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 133*47c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 134*47c6ae99SBarry Smith ierr = PetscLogFlops(5.0*m_f);CHKERRQ(ierr); 135*47c6ae99SBarry Smith PetscFunctionReturn(0); 136*47c6ae99SBarry Smith } 137*47c6ae99SBarry Smith 138*47c6ae99SBarry Smith #undef __FUNCT__ 139*47c6ae99SBarry Smith #define __FUNCT__ "DAGetInterpolation_1D_Q0" 140*47c6ae99SBarry Smith PetscErrorCode DAGetInterpolation_1D_Q0(DA dac,DA daf,Mat *A) 141*47c6ae99SBarry Smith { 142*47c6ae99SBarry Smith PetscErrorCode ierr; 143*47c6ae99SBarry Smith PetscInt i,i_start,m_f,Mx,*idx_f; 144*47c6ae99SBarry Smith PetscInt m_ghost,*idx_c,m_ghost_c; 145*47c6ae99SBarry Smith PetscInt row,col,i_start_ghost,mx,m_c,nc,ratio; 146*47c6ae99SBarry Smith PetscInt i_c,i_start_c,i_start_ghost_c,cols[2],dof; 147*47c6ae99SBarry Smith PetscScalar v[2],x; 148*47c6ae99SBarry Smith Mat mat; 149*47c6ae99SBarry Smith DAPeriodicType pt; 150*47c6ae99SBarry Smith 151*47c6ae99SBarry Smith PetscFunctionBegin; 152*47c6ae99SBarry Smith ierr = DAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 153*47c6ae99SBarry Smith ierr = DAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 154*47c6ae99SBarry Smith if (pt == DA_XPERIODIC) { 155*47c6ae99SBarry Smith ratio = mx/Mx; 156*47c6ae99SBarry 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); 157*47c6ae99SBarry Smith } else { 158*47c6ae99SBarry Smith ratio = (mx-1)/(Mx-1); 159*47c6ae99SBarry 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); 160*47c6ae99SBarry Smith } 161*47c6ae99SBarry Smith 162*47c6ae99SBarry Smith ierr = DAGetCorners(daf,&i_start,0,0,&m_f,0,0);CHKERRQ(ierr); 163*47c6ae99SBarry Smith ierr = DAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);CHKERRQ(ierr); 164*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 165*47c6ae99SBarry Smith 166*47c6ae99SBarry Smith ierr = DAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);CHKERRQ(ierr); 167*47c6ae99SBarry Smith ierr = DAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);CHKERRQ(ierr); 168*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 169*47c6ae99SBarry Smith 170*47c6ae99SBarry Smith /* create interpolation matrix */ 171*47c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr); 172*47c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f,m_c,mx,Mx);CHKERRQ(ierr); 173*47c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 174*47c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,2,PETSC_NULL);CHKERRQ(ierr); 175*47c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,2,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 176*47c6ae99SBarry Smith 177*47c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 178*47c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 179*47c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 180*47c6ae99SBarry Smith row = idx_f[dof*(i-i_start_ghost)]/dof; 181*47c6ae99SBarry Smith 182*47c6ae99SBarry Smith i_c = (i/ratio); /* coarse grid node to left of fine grid node */ 183*47c6ae99SBarry Smith 184*47c6ae99SBarry Smith /* 185*47c6ae99SBarry Smith Only include those interpolation points that are truly 186*47c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 187*47c6ae99SBarry Smith in x direction; since they have no right neighbor 188*47c6ae99SBarry Smith */ 189*47c6ae99SBarry Smith x = ((double)(i - i_c*ratio))/((double)ratio); 190*47c6ae99SBarry Smith nc = 0; 191*47c6ae99SBarry Smith /* one left and below; or we are right on it */ 192*47c6ae99SBarry Smith col = dof*(i_c-i_start_ghost_c); 193*47c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 194*47c6ae99SBarry Smith v[nc++] = - x + 1.0; 195*47c6ae99SBarry Smith /* one right? */ 196*47c6ae99SBarry Smith if (i_c*ratio != i) { 197*47c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 198*47c6ae99SBarry Smith v[nc++] = x; 199*47c6ae99SBarry Smith } 200*47c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 201*47c6ae99SBarry Smith } 202*47c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 203*47c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 204*47c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 205*47c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 206*47c6ae99SBarry Smith ierr = PetscLogFlops(5.0*m_f);CHKERRQ(ierr); 207*47c6ae99SBarry Smith PetscFunctionReturn(0); 208*47c6ae99SBarry Smith } 209*47c6ae99SBarry Smith 210*47c6ae99SBarry Smith #undef __FUNCT__ 211*47c6ae99SBarry Smith #define __FUNCT__ "DAGetInterpolation_2D_Q1" 212*47c6ae99SBarry Smith PetscErrorCode DAGetInterpolation_2D_Q1(DA dac,DA daf,Mat *A) 213*47c6ae99SBarry Smith { 214*47c6ae99SBarry Smith PetscErrorCode ierr; 215*47c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 216*47c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz; 217*47c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj; 218*47c6ae99SBarry 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; 219*47c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 220*47c6ae99SBarry Smith PetscScalar v[4],x,y; 221*47c6ae99SBarry Smith Mat mat; 222*47c6ae99SBarry Smith DAPeriodicType pt; 223*47c6ae99SBarry Smith DACoor2d **coors = 0,**ccoors; 224*47c6ae99SBarry Smith Vec vcoors,cvcoors; 225*47c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data, *ddf = (DM_DA*)daf->data; 226*47c6ae99SBarry Smith 227*47c6ae99SBarry Smith PetscFunctionBegin; 228*47c6ae99SBarry Smith ierr = DAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 229*47c6ae99SBarry Smith ierr = DAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 230*47c6ae99SBarry Smith if (DAXPeriodic(pt)){ 231*47c6ae99SBarry Smith ratioi = mx/Mx; 232*47c6ae99SBarry 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); 233*47c6ae99SBarry Smith } else { 234*47c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 235*47c6ae99SBarry 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); 236*47c6ae99SBarry Smith } 237*47c6ae99SBarry Smith if (DAYPeriodic(pt)){ 238*47c6ae99SBarry Smith ratioj = my/My; 239*47c6ae99SBarry 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); 240*47c6ae99SBarry Smith } else { 241*47c6ae99SBarry Smith ratioj = (my-1)/(My-1); 242*47c6ae99SBarry 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); 243*47c6ae99SBarry Smith } 244*47c6ae99SBarry Smith 245*47c6ae99SBarry Smith 246*47c6ae99SBarry Smith ierr = DAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 247*47c6ae99SBarry Smith ierr = DAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 248*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 249*47c6ae99SBarry Smith 250*47c6ae99SBarry Smith ierr = DAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 251*47c6ae99SBarry Smith ierr = DAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 252*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 253*47c6ae99SBarry Smith 254*47c6ae99SBarry Smith /* 255*47c6ae99SBarry Smith Used for handling a coarse DA that lives on 1/4 the processors of the fine DA. 256*47c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 257*47c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 258*47c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 259*47c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 260*47c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 261*47c6ae99SBarry Smith 262*47c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 263*47c6ae99SBarry Smith */ 264*47c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr); 265*47c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr); 266*47c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr); 267*47c6ae99SBarry Smith col_scale = size_f/size_c; 268*47c6ae99SBarry Smith col_shift = Mx*My*(rank_f/size_c); 269*47c6ae99SBarry Smith 270*47c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr); 271*47c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 272*47c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 273*47c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 274*47c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 275*47c6ae99SBarry Smith 276*47c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 277*47c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 278*47c6ae99SBarry Smith 279*47c6ae99SBarry Smith if (j_c < j_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 280*47c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 281*47c6ae99SBarry Smith if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 282*47c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 283*47c6ae99SBarry Smith 284*47c6ae99SBarry Smith /* 285*47c6ae99SBarry Smith Only include those interpolation points that are truly 286*47c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 287*47c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 288*47c6ae99SBarry Smith */ 289*47c6ae99SBarry Smith nc = 0; 290*47c6ae99SBarry Smith /* one left and below; or we are right on it */ 291*47c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 292*47c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 293*47c6ae99SBarry Smith /* one right and below */ 294*47c6ae99SBarry Smith if (i_c*ratioi != i) { 295*47c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col+dof]/dof; 296*47c6ae99SBarry Smith } 297*47c6ae99SBarry Smith /* one left and above */ 298*47c6ae99SBarry Smith if (j_c*ratioj != j) { 299*47c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col+m_ghost_c*dof]/dof; 300*47c6ae99SBarry Smith } 301*47c6ae99SBarry Smith /* one right and above */ 302*47c6ae99SBarry Smith if (j_c*ratioi != j && i_c*ratioj != i) { 303*47c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; 304*47c6ae99SBarry Smith } 305*47c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 306*47c6ae99SBarry Smith } 307*47c6ae99SBarry Smith } 308*47c6ae99SBarry Smith ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr); 309*47c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr); 310*47c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 311*47c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 312*47c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 313*47c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 314*47c6ae99SBarry Smith 315*47c6ae99SBarry Smith ierr = DAGetCoordinates(daf,&vcoors);CHKERRQ(ierr); 316*47c6ae99SBarry Smith if (vcoors) { 317*47c6ae99SBarry Smith ierr = DAGetGhostedCoordinates(dac,&cvcoors);CHKERRQ(ierr); 318*47c6ae99SBarry Smith ierr = DAVecGetArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 319*47c6ae99SBarry Smith ierr = DAVecGetArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 320*47c6ae99SBarry Smith } 321*47c6ae99SBarry Smith 322*47c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 323*47c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 324*47c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 325*47c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 326*47c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 327*47c6ae99SBarry Smith 328*47c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 329*47c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 330*47c6ae99SBarry Smith 331*47c6ae99SBarry Smith /* 332*47c6ae99SBarry Smith Only include those interpolation points that are truly 333*47c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 334*47c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 335*47c6ae99SBarry Smith */ 336*47c6ae99SBarry Smith if (coors) { 337*47c6ae99SBarry Smith /* only access the next coors point if we know there is one */ 338*47c6ae99SBarry Smith /* note this is dangerous because x may not exactly equal ZERO */ 339*47c6ae99SBarry Smith x = (coors[j][i].x - ccoors[j_c][i_c].x); 340*47c6ae99SBarry Smith if (PetscAbsScalar(x) != 0.0) x = x/(ccoors[j_c][i_c+1].x - ccoors[j_c][i_c].x); 341*47c6ae99SBarry Smith y = (coors[j][i].y - ccoors[j_c][i_c].y); 342*47c6ae99SBarry Smith if (PetscAbsScalar(y) != 0.0) y = y/(ccoors[j_c+1][i_c].y - ccoors[j_c][i_c].y); 343*47c6ae99SBarry Smith } else { 344*47c6ae99SBarry Smith x = ((double)(i - i_c*ratioi))/((double)ratioi); 345*47c6ae99SBarry Smith y = ((double)(j - j_c*ratioj))/((double)ratioj); 346*47c6ae99SBarry Smith } 347*47c6ae99SBarry Smith nc = 0; 348*47c6ae99SBarry Smith /* one left and below; or we are right on it */ 349*47c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 350*47c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 351*47c6ae99SBarry Smith v[nc++] = x*y - x - y + 1.0; 352*47c6ae99SBarry Smith /* one right and below */ 353*47c6ae99SBarry Smith if (i_c*ratioi != i) { 354*47c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+dof]/dof; 355*47c6ae99SBarry Smith v[nc++] = -x*y + x; 356*47c6ae99SBarry Smith } 357*47c6ae99SBarry Smith /* one left and above */ 358*47c6ae99SBarry Smith if (j_c*ratioj != j) { 359*47c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+m_ghost_c*dof]/dof; 360*47c6ae99SBarry Smith v[nc++] = -x*y + y; 361*47c6ae99SBarry Smith } 362*47c6ae99SBarry Smith /* one right and above */ 363*47c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 364*47c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof; 365*47c6ae99SBarry Smith v[nc++] = x*y; 366*47c6ae99SBarry Smith } 367*47c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 368*47c6ae99SBarry Smith } 369*47c6ae99SBarry Smith } 370*47c6ae99SBarry Smith if (vcoors) { 371*47c6ae99SBarry Smith ierr = DAVecRestoreArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 372*47c6ae99SBarry Smith ierr = DAVecRestoreArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 373*47c6ae99SBarry Smith } 374*47c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 375*47c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 376*47c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 377*47c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 378*47c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f);CHKERRQ(ierr); 379*47c6ae99SBarry Smith PetscFunctionReturn(0); 380*47c6ae99SBarry Smith } 381*47c6ae99SBarry Smith 382*47c6ae99SBarry Smith /* 383*47c6ae99SBarry Smith Contributed by Andrei Draganescu <aidraga@sandia.gov> 384*47c6ae99SBarry Smith */ 385*47c6ae99SBarry Smith #undef __FUNCT__ 386*47c6ae99SBarry Smith #define __FUNCT__ "DAGetInterpolation_2D_Q0" 387*47c6ae99SBarry Smith PetscErrorCode DAGetInterpolation_2D_Q0(DA dac,DA daf,Mat *A) 388*47c6ae99SBarry Smith { 389*47c6ae99SBarry Smith PetscErrorCode ierr; 390*47c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 391*47c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz; 392*47c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj; 393*47c6ae99SBarry 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; 394*47c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 395*47c6ae99SBarry Smith PetscScalar v[4]; 396*47c6ae99SBarry Smith Mat mat; 397*47c6ae99SBarry Smith DAPeriodicType pt; 398*47c6ae99SBarry Smith 399*47c6ae99SBarry Smith PetscFunctionBegin; 400*47c6ae99SBarry Smith ierr = DAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 401*47c6ae99SBarry Smith ierr = DAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 402*47c6ae99SBarry Smith if (DAXPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x"); 403*47c6ae99SBarry Smith if (DAYPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y"); 404*47c6ae99SBarry Smith ratioi = mx/Mx; 405*47c6ae99SBarry Smith ratioj = my/My; 406*47c6ae99SBarry 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"); 407*47c6ae99SBarry 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"); 408*47c6ae99SBarry Smith if (ratioi != 2) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 2"); 409*47c6ae99SBarry Smith if (ratioj != 2) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 2"); 410*47c6ae99SBarry Smith 411*47c6ae99SBarry Smith ierr = DAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 412*47c6ae99SBarry Smith ierr = DAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 413*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 414*47c6ae99SBarry Smith 415*47c6ae99SBarry Smith ierr = DAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 416*47c6ae99SBarry Smith ierr = DAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 417*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 418*47c6ae99SBarry Smith 419*47c6ae99SBarry Smith /* 420*47c6ae99SBarry Smith Used for handling a coarse DA that lives on 1/4 the processors of the fine DA. 421*47c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 422*47c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 423*47c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 424*47c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 425*47c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 426*47c6ae99SBarry Smith 427*47c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 428*47c6ae99SBarry Smith */ 429*47c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr); 430*47c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr); 431*47c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr); 432*47c6ae99SBarry Smith col_scale = size_f/size_c; 433*47c6ae99SBarry Smith col_shift = Mx*My*(rank_f/size_c); 434*47c6ae99SBarry Smith 435*47c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f,col_scale*m_c*n_c,dnz,onz);CHKERRQ(ierr); 436*47c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 437*47c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 438*47c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 439*47c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 440*47c6ae99SBarry Smith 441*47c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 442*47c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 443*47c6ae99SBarry Smith 444*47c6ae99SBarry Smith if (j_c < j_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 445*47c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 446*47c6ae99SBarry Smith if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 447*47c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 448*47c6ae99SBarry Smith 449*47c6ae99SBarry Smith /* 450*47c6ae99SBarry Smith Only include those interpolation points that are truly 451*47c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 452*47c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 453*47c6ae99SBarry Smith */ 454*47c6ae99SBarry Smith nc = 0; 455*47c6ae99SBarry Smith /* one left and below; or we are right on it */ 456*47c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 457*47c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 458*47c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 459*47c6ae99SBarry Smith } 460*47c6ae99SBarry Smith } 461*47c6ae99SBarry Smith ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr); 462*47c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My);CHKERRQ(ierr); 463*47c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 464*47c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 465*47c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 466*47c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 467*47c6ae99SBarry Smith 468*47c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 469*47c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 470*47c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 471*47c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 472*47c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 473*47c6ae99SBarry Smith 474*47c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 475*47c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 476*47c6ae99SBarry Smith nc = 0; 477*47c6ae99SBarry Smith /* one left and below; or we are right on it */ 478*47c6ae99SBarry Smith col = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c)); 479*47c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 480*47c6ae99SBarry Smith v[nc++] = 1.0; 481*47c6ae99SBarry Smith 482*47c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 483*47c6ae99SBarry Smith } 484*47c6ae99SBarry Smith } 485*47c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 486*47c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 487*47c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 488*47c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 489*47c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f);CHKERRQ(ierr); 490*47c6ae99SBarry Smith PetscFunctionReturn(0); 491*47c6ae99SBarry Smith } 492*47c6ae99SBarry Smith 493*47c6ae99SBarry Smith /* 494*47c6ae99SBarry Smith Contributed by Jianming Yang <jianming-yang@uiowa.edu> 495*47c6ae99SBarry Smith */ 496*47c6ae99SBarry Smith #undef __FUNCT__ 497*47c6ae99SBarry Smith #define __FUNCT__ "DAGetInterpolation_3D_Q0" 498*47c6ae99SBarry Smith PetscErrorCode DAGetInterpolation_3D_Q0(DA dac,DA daf,Mat *A) 499*47c6ae99SBarry Smith { 500*47c6ae99SBarry Smith PetscErrorCode ierr; 501*47c6ae99SBarry Smith PetscInt i,j,l,i_start,j_start,l_start,m_f,n_f,p_f,Mx,My,Mz,*idx_f,dof; 502*47c6ae99SBarry Smith PetscInt m_ghost,n_ghost,p_ghost,*idx_c,m_ghost_c,n_ghost_c,p_ghost_c,nc,*dnz,*onz; 503*47c6ae99SBarry 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; 504*47c6ae99SBarry 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; 505*47c6ae99SBarry Smith PetscMPIInt size_c,size_f,rank_f; 506*47c6ae99SBarry Smith PetscScalar v[8]; 507*47c6ae99SBarry Smith Mat mat; 508*47c6ae99SBarry Smith DAPeriodicType pt; 509*47c6ae99SBarry Smith 510*47c6ae99SBarry Smith PetscFunctionBegin; 511*47c6ae99SBarry Smith ierr = DAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 512*47c6ae99SBarry Smith ierr = DAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 513*47c6ae99SBarry Smith if (DAXPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in x"); 514*47c6ae99SBarry Smith if (DAYPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in y"); 515*47c6ae99SBarry Smith if (DAZPeriodic(pt)) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Cannot handle periodic grid in z"); 516*47c6ae99SBarry Smith ratioi = mx/Mx; 517*47c6ae99SBarry Smith ratioj = my/My; 518*47c6ae99SBarry Smith ratiol = mz/Mz; 519*47c6ae99SBarry 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"); 520*47c6ae99SBarry 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"); 521*47c6ae99SBarry 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"); 522*47c6ae99SBarry Smith if (ratioi != 2 && ratioi != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in x must be 1 or 2"); 523*47c6ae99SBarry Smith if (ratioj != 2 && ratioj != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in y must be 1 or 2"); 524*47c6ae99SBarry Smith if (ratiol != 2 && ratiol != 1) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_WRONG,"Coarsening factor in z must be 1 or 2"); 525*47c6ae99SBarry Smith 526*47c6ae99SBarry Smith ierr = DAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 527*47c6ae99SBarry Smith ierr = DAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 528*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 529*47c6ae99SBarry Smith 530*47c6ae99SBarry Smith ierr = DAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 531*47c6ae99SBarry Smith ierr = DAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,&l_start_ghost_c,&m_ghost_c,&n_ghost_c,&p_ghost_c);CHKERRQ(ierr); 532*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 533*47c6ae99SBarry Smith /* 534*47c6ae99SBarry Smith Used for handling a coarse DA that lives on 1/4 the processors of the fine DA. 535*47c6ae99SBarry Smith The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 536*47c6ae99SBarry Smith processors). It's effective length is hence 4 times its normal length, this is 537*47c6ae99SBarry Smith why the col_scale is multiplied by the interpolation matrix column sizes. 538*47c6ae99SBarry Smith sol_shift allows each set of 1/4 processors do its own interpolation using ITS 539*47c6ae99SBarry Smith copy of the coarse vector. A bit of a hack but you do better. 540*47c6ae99SBarry Smith 541*47c6ae99SBarry Smith In the standard case when size_f == size_c col_scale == 1 and col_shift == 0 542*47c6ae99SBarry Smith */ 543*47c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dac)->comm,&size_c);CHKERRQ(ierr); 544*47c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)daf)->comm,&size_f);CHKERRQ(ierr); 545*47c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)daf)->comm,&rank_f);CHKERRQ(ierr); 546*47c6ae99SBarry Smith col_scale = size_f/size_c; 547*47c6ae99SBarry Smith col_shift = Mx*My*Mz*(rank_f/size_c); 548*47c6ae99SBarry Smith 549*47c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)daf)->comm,m_f*n_f*p_f,col_scale*m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr); 550*47c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 551*47c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 552*47c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 553*47c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 554*47c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 555*47c6ae99SBarry Smith 556*47c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 557*47c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 558*47c6ae99SBarry Smith l_c = (l/ratiol); 559*47c6ae99SBarry Smith 560*47c6ae99SBarry Smith if (l_c < l_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 561*47c6ae99SBarry Smith l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c); 562*47c6ae99SBarry Smith if (j_c < j_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 563*47c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 564*47c6ae99SBarry Smith if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 565*47c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 566*47c6ae99SBarry Smith 567*47c6ae99SBarry Smith /* 568*47c6ae99SBarry Smith Only include those interpolation points that are truly 569*47c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 570*47c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 571*47c6ae99SBarry Smith */ 572*47c6ae99SBarry Smith nc = 0; 573*47c6ae99SBarry Smith /* one left and below; or we are right on it */ 574*47c6ae99SBarry 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)); 575*47c6ae99SBarry Smith cols[nc++] = col_shift + idx_c[col]/dof; 576*47c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 577*47c6ae99SBarry Smith } 578*47c6ae99SBarry Smith } 579*47c6ae99SBarry Smith } 580*47c6ae99SBarry Smith ierr = MatCreate(((PetscObject)daf)->comm,&mat);CHKERRQ(ierr); 581*47c6ae99SBarry 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); 582*47c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 583*47c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 584*47c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 585*47c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 586*47c6ae99SBarry Smith 587*47c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 588*47c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 589*47c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 590*47c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 591*47c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 592*47c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 593*47c6ae99SBarry Smith 594*47c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 595*47c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 596*47c6ae99SBarry Smith l_c = (l/ratiol); 597*47c6ae99SBarry Smith nc = 0; 598*47c6ae99SBarry Smith /* one left and below; or we are right on it */ 599*47c6ae99SBarry 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)); 600*47c6ae99SBarry Smith cols[nc] = col_shift + idx_c[col]/dof; 601*47c6ae99SBarry Smith v[nc++] = 1.0; 602*47c6ae99SBarry Smith 603*47c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 604*47c6ae99SBarry Smith } 605*47c6ae99SBarry Smith } 606*47c6ae99SBarry Smith } 607*47c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 608*47c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 609*47c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 610*47c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 611*47c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f*p_f);CHKERRQ(ierr); 612*47c6ae99SBarry Smith PetscFunctionReturn(0); 613*47c6ae99SBarry Smith } 614*47c6ae99SBarry Smith 615*47c6ae99SBarry Smith #undef __FUNCT__ 616*47c6ae99SBarry Smith #define __FUNCT__ "DAGetInterpolation_3D_Q1" 617*47c6ae99SBarry Smith PetscErrorCode DAGetInterpolation_3D_Q1(DA dac,DA daf,Mat *A) 618*47c6ae99SBarry Smith { 619*47c6ae99SBarry Smith PetscErrorCode ierr; 620*47c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof,l; 621*47c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,Mz,mz; 622*47c6ae99SBarry Smith PetscInt row,col,i_start_ghost,j_start_ghost,cols[8],mx,m_c,my,nc,ratioi,ratioj,ratiok; 623*47c6ae99SBarry Smith PetscInt i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c; 624*47c6ae99SBarry Smith PetscInt l_start,p_f,l_start_ghost,p_ghost,l_start_c,p_c; 625*47c6ae99SBarry Smith PetscInt l_start_ghost_c,p_ghost_c,l_c,*dnz,*onz; 626*47c6ae99SBarry Smith PetscScalar v[8],x,y,z; 627*47c6ae99SBarry Smith Mat mat; 628*47c6ae99SBarry Smith DAPeriodicType pt; 629*47c6ae99SBarry Smith DACoor3d ***coors = 0,***ccoors; 630*47c6ae99SBarry Smith Vec vcoors,cvcoors; 631*47c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data, *ddf = (DM_DA*)daf->data; 632*47c6ae99SBarry Smith 633*47c6ae99SBarry Smith PetscFunctionBegin; 634*47c6ae99SBarry Smith ierr = DAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 635*47c6ae99SBarry Smith ierr = DAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 636*47c6ae99SBarry Smith if (mx == Mx) { 637*47c6ae99SBarry Smith ratioi = 1; 638*47c6ae99SBarry Smith } else if (DAXPeriodic(pt)){ 639*47c6ae99SBarry Smith ratioi = mx/Mx; 640*47c6ae99SBarry 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); 641*47c6ae99SBarry Smith } else { 642*47c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 643*47c6ae99SBarry 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); 644*47c6ae99SBarry Smith } 645*47c6ae99SBarry Smith if (my == My) { 646*47c6ae99SBarry Smith ratioj = 1; 647*47c6ae99SBarry Smith } else if (DAYPeriodic(pt)){ 648*47c6ae99SBarry Smith ratioj = my/My; 649*47c6ae99SBarry 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); 650*47c6ae99SBarry Smith } else { 651*47c6ae99SBarry Smith ratioj = (my-1)/(My-1); 652*47c6ae99SBarry 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); 653*47c6ae99SBarry Smith } 654*47c6ae99SBarry Smith if (mz == Mz) { 655*47c6ae99SBarry Smith ratiok = 1; 656*47c6ae99SBarry Smith } else if (DAZPeriodic(pt)){ 657*47c6ae99SBarry Smith ratiok = mz/Mz; 658*47c6ae99SBarry 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); 659*47c6ae99SBarry Smith } else { 660*47c6ae99SBarry Smith ratiok = (mz-1)/(Mz-1); 661*47c6ae99SBarry 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); 662*47c6ae99SBarry Smith } 663*47c6ae99SBarry Smith 664*47c6ae99SBarry Smith ierr = DAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 665*47c6ae99SBarry Smith ierr = DAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 666*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 667*47c6ae99SBarry Smith 668*47c6ae99SBarry Smith ierr = DAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 669*47c6ae99SBarry Smith ierr = DAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,&l_start_ghost_c,&m_ghost_c,&n_ghost_c,&p_ghost_c);CHKERRQ(ierr); 670*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 671*47c6ae99SBarry Smith 672*47c6ae99SBarry Smith /* create interpolation matrix, determining exact preallocation */ 673*47c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)dac)->comm,m_f*n_f*p_f,m_c*n_c*p_c,dnz,onz);CHKERRQ(ierr); 674*47c6ae99SBarry Smith /* loop over local fine grid nodes counting interpolating points */ 675*47c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 676*47c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 677*47c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 678*47c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 679*47c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 680*47c6ae99SBarry Smith i_c = (i/ratioi); 681*47c6ae99SBarry Smith j_c = (j/ratioj); 682*47c6ae99SBarry Smith l_c = (l/ratiok); 683*47c6ae99SBarry Smith if (l_c < l_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 684*47c6ae99SBarry Smith l_start %D l_c %D l_start_ghost_c %D",l_start,l_c,l_start_ghost_c); 685*47c6ae99SBarry Smith if (j_c < j_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 686*47c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 687*47c6ae99SBarry Smith if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 688*47c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 689*47c6ae99SBarry Smith 690*47c6ae99SBarry Smith /* 691*47c6ae99SBarry Smith Only include those interpolation points that are truly 692*47c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 693*47c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 694*47c6ae99SBarry Smith */ 695*47c6ae99SBarry Smith nc = 0; 696*47c6ae99SBarry 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)); 697*47c6ae99SBarry Smith cols[nc++] = idx_c[col]/dof; 698*47c6ae99SBarry Smith if (i_c*ratioi != i) { 699*47c6ae99SBarry Smith cols[nc++] = idx_c[col+dof]/dof; 700*47c6ae99SBarry Smith } 701*47c6ae99SBarry Smith if (j_c*ratioj != j) { 702*47c6ae99SBarry Smith cols[nc++] = idx_c[col+m_ghost_c*dof]/dof; 703*47c6ae99SBarry Smith } 704*47c6ae99SBarry Smith if (l_c*ratiok != l) { 705*47c6ae99SBarry Smith cols[nc++] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; 706*47c6ae99SBarry Smith } 707*47c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 708*47c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c+1)*dof]/dof; 709*47c6ae99SBarry Smith } 710*47c6ae99SBarry Smith if (j_c*ratioj != j && l_c*ratiok != l) { 711*47c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof; 712*47c6ae99SBarry Smith } 713*47c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l) { 714*47c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; 715*47c6ae99SBarry Smith } 716*47c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) { 717*47c6ae99SBarry Smith cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; 718*47c6ae99SBarry Smith } 719*47c6ae99SBarry Smith ierr = MatPreallocateSet(row,nc,cols,dnz,onz);CHKERRQ(ierr); 720*47c6ae99SBarry Smith } 721*47c6ae99SBarry Smith } 722*47c6ae99SBarry Smith } 723*47c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dac)->comm,&mat);CHKERRQ(ierr); 724*47c6ae99SBarry Smith ierr = MatSetSizes(mat,m_f*n_f*p_f,m_c*n_c*p_c,mx*my*mz,Mx*My*Mz);CHKERRQ(ierr); 725*47c6ae99SBarry Smith ierr = MatSetType(mat,MATAIJ);CHKERRQ(ierr); 726*47c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(mat,0,dnz);CHKERRQ(ierr); 727*47c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(mat,0,dnz,0,onz);CHKERRQ(ierr); 728*47c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 729*47c6ae99SBarry Smith 730*47c6ae99SBarry Smith ierr = DAGetCoordinates(daf,&vcoors);CHKERRQ(ierr); 731*47c6ae99SBarry Smith if (vcoors) { 732*47c6ae99SBarry Smith ierr = DAGetGhostedCoordinates(dac,&cvcoors);CHKERRQ(ierr); 733*47c6ae99SBarry Smith ierr = DAVecGetArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 734*47c6ae99SBarry Smith ierr = DAVecGetArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 735*47c6ae99SBarry Smith } 736*47c6ae99SBarry Smith 737*47c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 738*47c6ae99SBarry Smith for (l=l_start; l<l_start+p_f; l++) { 739*47c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 740*47c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 741*47c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 742*47c6ae99SBarry Smith row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof; 743*47c6ae99SBarry Smith 744*47c6ae99SBarry Smith i_c = (i/ratioi); 745*47c6ae99SBarry Smith j_c = (j/ratioj); 746*47c6ae99SBarry Smith l_c = (l/ratiok); 747*47c6ae99SBarry Smith 748*47c6ae99SBarry Smith /* 749*47c6ae99SBarry Smith Only include those interpolation points that are truly 750*47c6ae99SBarry Smith nonzero. Note this is very important for final grid lines 751*47c6ae99SBarry Smith in x and y directions; since they have no right/top neighbors 752*47c6ae99SBarry Smith */ 753*47c6ae99SBarry Smith if (coors) { 754*47c6ae99SBarry Smith /* only access the next coors point if we know there is one */ 755*47c6ae99SBarry Smith /* note this is dangerous because x may not exactly equal ZERO */ 756*47c6ae99SBarry Smith x = (coors[l][j][i].x - ccoors[l_c][j_c][i_c].x); 757*47c6ae99SBarry 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); 758*47c6ae99SBarry Smith y = (coors[l][j][i].y - ccoors[l_c][j_c][i_c].y); 759*47c6ae99SBarry 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); 760*47c6ae99SBarry Smith z = (coors[l][j][i].z - ccoors[l_c][j_c][i_c].z); 761*47c6ae99SBarry 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); 762*47c6ae99SBarry Smith } else { 763*47c6ae99SBarry Smith x = ((double)(i - i_c*ratioi))/((double)ratioi); 764*47c6ae99SBarry Smith y = ((double)(j - j_c*ratioj))/((double)ratioj); 765*47c6ae99SBarry Smith z = ((double)(l - l_c*ratiok))/((double)ratiok); 766*47c6ae99SBarry Smith } 767*47c6ae99SBarry Smith nc = 0; 768*47c6ae99SBarry Smith /* one left and below; or we are right on it */ 769*47c6ae99SBarry 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)); 770*47c6ae99SBarry Smith 771*47c6ae99SBarry Smith cols[nc] = idx_c[col]/dof; 772*47c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.)); 773*47c6ae99SBarry Smith 774*47c6ae99SBarry Smith if (i_c*ratioi != i) { 775*47c6ae99SBarry Smith cols[nc] = idx_c[col+dof]/dof; 776*47c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.)); 777*47c6ae99SBarry Smith } 778*47c6ae99SBarry Smith 779*47c6ae99SBarry Smith if (j_c*ratioj != j) { 780*47c6ae99SBarry Smith cols[nc] = idx_c[col+m_ghost_c*dof]/dof; 781*47c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.)); 782*47c6ae99SBarry Smith } 783*47c6ae99SBarry Smith 784*47c6ae99SBarry Smith if (l_c*ratiok != l) { 785*47c6ae99SBarry Smith cols[nc] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof; 786*47c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.)); 787*47c6ae99SBarry Smith } 788*47c6ae99SBarry Smith 789*47c6ae99SBarry Smith if (j_c*ratioj != j && i_c*ratioi != i) { 790*47c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c+1)*dof]/dof; 791*47c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.)); 792*47c6ae99SBarry Smith } 793*47c6ae99SBarry Smith 794*47c6ae99SBarry Smith if (j_c*ratioj != j && l_c*ratiok != l) { 795*47c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof; 796*47c6ae99SBarry Smith v[nc++] = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.)); 797*47c6ae99SBarry Smith } 798*47c6ae99SBarry Smith 799*47c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l) { 800*47c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof; 801*47c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.)); 802*47c6ae99SBarry Smith } 803*47c6ae99SBarry Smith 804*47c6ae99SBarry Smith if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) { 805*47c6ae99SBarry Smith cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof; 806*47c6ae99SBarry Smith v[nc++] = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.)); 807*47c6ae99SBarry Smith } 808*47c6ae99SBarry Smith ierr = MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);CHKERRQ(ierr); 809*47c6ae99SBarry Smith } 810*47c6ae99SBarry Smith } 811*47c6ae99SBarry Smith } 812*47c6ae99SBarry Smith if (vcoors) { 813*47c6ae99SBarry Smith ierr = DAVecRestoreArray(ddf->da_coordinates,vcoors,&coors);CHKERRQ(ierr); 814*47c6ae99SBarry Smith ierr = DAVecRestoreArray(ddc->da_coordinates,cvcoors,&ccoors);CHKERRQ(ierr); 815*47c6ae99SBarry Smith } 816*47c6ae99SBarry Smith ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 817*47c6ae99SBarry Smith ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 818*47c6ae99SBarry Smith 819*47c6ae99SBarry Smith ierr = MatCreateMAIJ(mat,dof,A);CHKERRQ(ierr); 820*47c6ae99SBarry Smith ierr = MatDestroy(mat);CHKERRQ(ierr); 821*47c6ae99SBarry Smith ierr = PetscLogFlops(13.0*m_f*n_f);CHKERRQ(ierr); 822*47c6ae99SBarry Smith PetscFunctionReturn(0); 823*47c6ae99SBarry Smith } 824*47c6ae99SBarry Smith 825*47c6ae99SBarry Smith #undef __FUNCT__ 826*47c6ae99SBarry Smith #define __FUNCT__ "DAGetInterpolation" 827*47c6ae99SBarry Smith /*@C 828*47c6ae99SBarry Smith DAGetInterpolation - Gets an interpolation matrix that maps between 829*47c6ae99SBarry Smith grids associated with two DAs. 830*47c6ae99SBarry Smith 831*47c6ae99SBarry Smith Collective on DA 832*47c6ae99SBarry Smith 833*47c6ae99SBarry Smith Input Parameters: 834*47c6ae99SBarry Smith + dac - the coarse grid DA 835*47c6ae99SBarry Smith - daf - the fine grid DA 836*47c6ae99SBarry Smith 837*47c6ae99SBarry Smith Output Parameters: 838*47c6ae99SBarry Smith + A - the interpolation matrix 839*47c6ae99SBarry Smith - scale - a scaling vector used to scale the coarse grid restricted vector before applying the 840*47c6ae99SBarry Smith grid function or grid Jacobian to it. 841*47c6ae99SBarry Smith 842*47c6ae99SBarry Smith Level: intermediate 843*47c6ae99SBarry Smith 844*47c6ae99SBarry Smith .keywords: interpolation, restriction, multigrid 845*47c6ae99SBarry Smith 846*47c6ae99SBarry Smith .seealso: DARefine(), DAGetInjection() 847*47c6ae99SBarry Smith @*/ 848*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetInterpolation(DA dac,DA daf,Mat *A,Vec *scale) 849*47c6ae99SBarry Smith { 850*47c6ae99SBarry Smith PetscErrorCode ierr; 851*47c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 852*47c6ae99SBarry Smith DAPeriodicType wrapc,wrapf; 853*47c6ae99SBarry Smith DAStencilType stc,stf; 854*47c6ae99SBarry Smith DM_DA *ddc = (DM_DA*)dac->data; 855*47c6ae99SBarry Smith 856*47c6ae99SBarry Smith PetscFunctionBegin; 857*47c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 858*47c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 859*47c6ae99SBarry Smith PetscValidPointer(A,3); 860*47c6ae99SBarry Smith if (scale) PetscValidPointer(scale,4); 861*47c6ae99SBarry Smith 862*47c6ae99SBarry Smith ierr = DAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&wrapc,&stc);CHKERRQ(ierr); 863*47c6ae99SBarry Smith ierr = DAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&wrapf,&stf);CHKERRQ(ierr); 864*47c6ae99SBarry Smith if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 865*47c6ae99SBarry Smith if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DA do not match %D %D",dofc,doff);CHKERRQ(ierr); 866*47c6ae99SBarry Smith if (sc != sf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil width of DA do not match %D %D",sc,sf);CHKERRQ(ierr); 867*47c6ae99SBarry Smith if (wrapc != wrapf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Periodic type different in two DAs");CHKERRQ(ierr); 868*47c6ae99SBarry Smith if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DAs");CHKERRQ(ierr); 869*47c6ae99SBarry Smith if (Mc < 2 && Mf > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction"); 870*47c6ae99SBarry 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"); 871*47c6ae99SBarry 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"); 872*47c6ae99SBarry Smith 873*47c6ae99SBarry Smith if (ddc->interptype == DA_Q1){ 874*47c6ae99SBarry Smith if (dimc == 1){ 875*47c6ae99SBarry Smith ierr = DAGetInterpolation_1D_Q1(dac,daf,A);CHKERRQ(ierr); 876*47c6ae99SBarry Smith } else if (dimc == 2){ 877*47c6ae99SBarry Smith ierr = DAGetInterpolation_2D_Q1(dac,daf,A);CHKERRQ(ierr); 878*47c6ae99SBarry Smith } else if (dimc == 3){ 879*47c6ae99SBarry Smith ierr = DAGetInterpolation_3D_Q1(dac,daf,A);CHKERRQ(ierr); 880*47c6ae99SBarry Smith } else { 881*47c6ae99SBarry Smith SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DA dimension %D for interpolation type %d",dimc,(int)ddc->interptype); 882*47c6ae99SBarry Smith } 883*47c6ae99SBarry Smith } else if (ddc->interptype == DA_Q0){ 884*47c6ae99SBarry Smith if (dimc == 1){ 885*47c6ae99SBarry Smith ierr = DAGetInterpolation_1D_Q0(dac,daf,A);CHKERRQ(ierr); 886*47c6ae99SBarry Smith } else if (dimc == 2){ 887*47c6ae99SBarry Smith ierr = DAGetInterpolation_2D_Q0(dac,daf,A);CHKERRQ(ierr); 888*47c6ae99SBarry Smith } else if (dimc == 3){ 889*47c6ae99SBarry Smith ierr = DAGetInterpolation_3D_Q0(dac,daf,A);CHKERRQ(ierr); 890*47c6ae99SBarry Smith } else { 891*47c6ae99SBarry Smith SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DA dimension %D for interpolation type %d",dimc,(int)ddc->interptype); 892*47c6ae99SBarry Smith } 893*47c6ae99SBarry Smith } 894*47c6ae99SBarry Smith if (scale) { 895*47c6ae99SBarry Smith ierr = DMGetInterpolationScale((DM)dac,(DM)daf,*A,scale);CHKERRQ(ierr); 896*47c6ae99SBarry Smith } 897*47c6ae99SBarry Smith PetscFunctionReturn(0); 898*47c6ae99SBarry Smith } 899*47c6ae99SBarry Smith 900*47c6ae99SBarry Smith #undef __FUNCT__ 901*47c6ae99SBarry Smith #define __FUNCT__ "DAGetInjection_2D" 902*47c6ae99SBarry Smith PetscErrorCode DAGetInjection_2D(DA dac,DA daf,VecScatter *inject) 903*47c6ae99SBarry Smith { 904*47c6ae99SBarry Smith PetscErrorCode ierr; 905*47c6ae99SBarry Smith PetscInt i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof; 906*47c6ae99SBarry Smith PetscInt m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c; 907*47c6ae99SBarry Smith PetscInt row,i_start_ghost,j_start_ghost,mx,m_c,my,nc,ratioi,ratioj; 908*47c6ae99SBarry Smith PetscInt i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c; 909*47c6ae99SBarry Smith PetscInt *cols; 910*47c6ae99SBarry Smith DAPeriodicType pt; 911*47c6ae99SBarry Smith Vec vecf,vecc; 912*47c6ae99SBarry Smith IS isf; 913*47c6ae99SBarry Smith 914*47c6ae99SBarry Smith PetscFunctionBegin; 915*47c6ae99SBarry Smith 916*47c6ae99SBarry Smith ierr = DAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&pt,0);CHKERRQ(ierr); 917*47c6ae99SBarry Smith ierr = DAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 918*47c6ae99SBarry Smith if (DAXPeriodic(pt)){ 919*47c6ae99SBarry Smith ratioi = mx/Mx; 920*47c6ae99SBarry 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); 921*47c6ae99SBarry Smith } else { 922*47c6ae99SBarry Smith ratioi = (mx-1)/(Mx-1); 923*47c6ae99SBarry 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); 924*47c6ae99SBarry Smith } 925*47c6ae99SBarry Smith if (DAYPeriodic(pt)){ 926*47c6ae99SBarry Smith ratioj = my/My; 927*47c6ae99SBarry 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); 928*47c6ae99SBarry Smith } else { 929*47c6ae99SBarry Smith ratioj = (my-1)/(My-1); 930*47c6ae99SBarry 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); 931*47c6ae99SBarry Smith } 932*47c6ae99SBarry Smith 933*47c6ae99SBarry Smith 934*47c6ae99SBarry Smith ierr = DAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);CHKERRQ(ierr); 935*47c6ae99SBarry Smith ierr = DAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);CHKERRQ(ierr); 936*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 937*47c6ae99SBarry Smith 938*47c6ae99SBarry Smith ierr = DAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);CHKERRQ(ierr); 939*47c6ae99SBarry Smith ierr = DAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);CHKERRQ(ierr); 940*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 941*47c6ae99SBarry Smith 942*47c6ae99SBarry Smith 943*47c6ae99SBarry Smith /* loop over local fine grid nodes setting interpolation for those*/ 944*47c6ae99SBarry Smith nc = 0; 945*47c6ae99SBarry Smith ierr = PetscMalloc(n_f*m_f*sizeof(PetscInt),&cols);CHKERRQ(ierr); 946*47c6ae99SBarry Smith for (j=j_start; j<j_start+n_f; j++) { 947*47c6ae99SBarry Smith for (i=i_start; i<i_start+m_f; i++) { 948*47c6ae99SBarry Smith 949*47c6ae99SBarry Smith i_c = (i/ratioi); /* coarse grid node to left of fine grid node */ 950*47c6ae99SBarry Smith j_c = (j/ratioj); /* coarse grid node below fine grid node */ 951*47c6ae99SBarry Smith 952*47c6ae99SBarry Smith if (j_c < j_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 953*47c6ae99SBarry Smith j_start %D j_c %D j_start_ghost_c %D",j_start,j_c,j_start_ghost_c); 954*47c6ae99SBarry Smith if (i_c < i_start_ghost_c) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Processor's coarse DA must lie over fine DA\n\ 955*47c6ae99SBarry Smith i_start %D i_c %D i_start_ghost_c %D",i_start,i_c,i_start_ghost_c); 956*47c6ae99SBarry Smith 957*47c6ae99SBarry Smith if (i_c*ratioi == i && j_c*ratioj == j) { 958*47c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 959*47c6ae99SBarry Smith row = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]; 960*47c6ae99SBarry Smith cols[nc++] = row/dof; 961*47c6ae99SBarry Smith } 962*47c6ae99SBarry Smith } 963*47c6ae99SBarry Smith } 964*47c6ae99SBarry Smith 965*47c6ae99SBarry Smith ierr = ISCreateBlock(((PetscObject)daf)->comm,dof,nc,cols,PETSC_OWN_POINTER,&isf);CHKERRQ(ierr); 966*47c6ae99SBarry Smith ierr = DAGetGlobalVector(dac,&vecc);CHKERRQ(ierr); 967*47c6ae99SBarry Smith ierr = DAGetGlobalVector(daf,&vecf);CHKERRQ(ierr); 968*47c6ae99SBarry Smith ierr = VecScatterCreate(vecf,isf,vecc,PETSC_NULL,inject);CHKERRQ(ierr); 969*47c6ae99SBarry Smith ierr = DARestoreGlobalVector(dac,&vecc);CHKERRQ(ierr); 970*47c6ae99SBarry Smith ierr = DARestoreGlobalVector(daf,&vecf);CHKERRQ(ierr); 971*47c6ae99SBarry Smith ierr = ISDestroy(isf);CHKERRQ(ierr); 972*47c6ae99SBarry Smith PetscFunctionReturn(0); 973*47c6ae99SBarry Smith } 974*47c6ae99SBarry Smith 975*47c6ae99SBarry Smith #undef __FUNCT__ 976*47c6ae99SBarry Smith #define __FUNCT__ "DAGetInjection" 977*47c6ae99SBarry Smith /*@ 978*47c6ae99SBarry Smith DAGetInjection - Gets an injection matrix that maps between 979*47c6ae99SBarry Smith grids associated with two DAs. 980*47c6ae99SBarry Smith 981*47c6ae99SBarry Smith Collective on DA 982*47c6ae99SBarry Smith 983*47c6ae99SBarry Smith Input Parameters: 984*47c6ae99SBarry Smith + dac - the coarse grid DA 985*47c6ae99SBarry Smith - daf - the fine grid DA 986*47c6ae99SBarry Smith 987*47c6ae99SBarry Smith Output Parameters: 988*47c6ae99SBarry Smith . inject - the injection scatter 989*47c6ae99SBarry Smith 990*47c6ae99SBarry Smith Level: intermediate 991*47c6ae99SBarry Smith 992*47c6ae99SBarry Smith .keywords: interpolation, restriction, multigrid, injection 993*47c6ae99SBarry Smith 994*47c6ae99SBarry Smith .seealso: DARefine(), DAGetInterpolation() 995*47c6ae99SBarry Smith @*/ 996*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetInjection(DA dac,DA daf,VecScatter *inject) 997*47c6ae99SBarry Smith { 998*47c6ae99SBarry Smith PetscErrorCode ierr; 999*47c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 1000*47c6ae99SBarry Smith DAPeriodicType wrapc,wrapf; 1001*47c6ae99SBarry Smith DAStencilType stc,stf; 1002*47c6ae99SBarry Smith 1003*47c6ae99SBarry Smith PetscFunctionBegin; 1004*47c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 1005*47c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 1006*47c6ae99SBarry Smith PetscValidPointer(inject,3); 1007*47c6ae99SBarry Smith 1008*47c6ae99SBarry Smith ierr = DAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&wrapc,&stc);CHKERRQ(ierr); 1009*47c6ae99SBarry Smith ierr = DAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&wrapf,&stf);CHKERRQ(ierr); 1010*47c6ae99SBarry Smith if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 1011*47c6ae99SBarry Smith if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DA do not match %D %D",dofc,doff);CHKERRQ(ierr); 1012*47c6ae99SBarry Smith if (sc != sf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil width of DA do not match %D %D",sc,sf);CHKERRQ(ierr); 1013*47c6ae99SBarry Smith if (wrapc != wrapf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Periodic type different in two DAs");CHKERRQ(ierr); 1014*47c6ae99SBarry Smith if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DAs");CHKERRQ(ierr); 1015*47c6ae99SBarry Smith if (Mc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in x direction"); 1016*47c6ae99SBarry Smith if (dimc > 1 && Nc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in y direction"); 1017*47c6ae99SBarry Smith if (dimc > 2 && Pc < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coarse grid requires at least 2 points in z direction"); 1018*47c6ae99SBarry Smith 1019*47c6ae99SBarry Smith if (dimc == 2){ 1020*47c6ae99SBarry Smith ierr = DAGetInjection_2D(dac,daf,inject);CHKERRQ(ierr); 1021*47c6ae99SBarry Smith } else { 1022*47c6ae99SBarry Smith SETERRQ1(((PetscObject)daf)->comm,PETSC_ERR_SUP,"No support for this DA dimension %D",dimc); 1023*47c6ae99SBarry Smith } 1024*47c6ae99SBarry Smith PetscFunctionReturn(0); 1025*47c6ae99SBarry Smith } 1026*47c6ae99SBarry Smith 1027*47c6ae99SBarry Smith #undef __FUNCT__ 1028*47c6ae99SBarry Smith #define __FUNCT__ "DAGetAggregates" 1029*47c6ae99SBarry Smith /*@ 1030*47c6ae99SBarry Smith DAGetAggregates - Gets the aggregates that map between 1031*47c6ae99SBarry Smith grids associated with two DAs. 1032*47c6ae99SBarry Smith 1033*47c6ae99SBarry Smith Collective on DA 1034*47c6ae99SBarry Smith 1035*47c6ae99SBarry Smith Input Parameters: 1036*47c6ae99SBarry Smith + dac - the coarse grid DA 1037*47c6ae99SBarry Smith - daf - the fine grid DA 1038*47c6ae99SBarry Smith 1039*47c6ae99SBarry Smith Output Parameters: 1040*47c6ae99SBarry Smith . rest - the restriction matrix (transpose of the projection matrix) 1041*47c6ae99SBarry Smith 1042*47c6ae99SBarry Smith Level: intermediate 1043*47c6ae99SBarry Smith 1044*47c6ae99SBarry Smith .keywords: interpolation, restriction, multigrid 1045*47c6ae99SBarry Smith 1046*47c6ae99SBarry Smith .seealso: DARefine(), DAGetInjection(), DAGetInterpolation() 1047*47c6ae99SBarry Smith @*/ 1048*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetAggregates(DA dac,DA daf,Mat *rest) 1049*47c6ae99SBarry Smith { 1050*47c6ae99SBarry Smith PetscErrorCode ierr; 1051*47c6ae99SBarry Smith PetscInt dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc; 1052*47c6ae99SBarry Smith PetscInt dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf; 1053*47c6ae99SBarry Smith DAPeriodicType wrapc,wrapf; 1054*47c6ae99SBarry Smith DAStencilType stc,stf; 1055*47c6ae99SBarry Smith /* PetscReal r_x, r_y, r_z; */ 1056*47c6ae99SBarry Smith 1057*47c6ae99SBarry Smith PetscInt i,j,l; 1058*47c6ae99SBarry Smith PetscInt i_start,j_start,l_start, m_f,n_f,p_f; 1059*47c6ae99SBarry Smith PetscInt i_start_ghost,j_start_ghost,l_start_ghost,m_ghost,n_ghost,p_ghost; 1060*47c6ae99SBarry Smith PetscInt *idx_f; 1061*47c6ae99SBarry Smith PetscInt i_c,j_c,l_c; 1062*47c6ae99SBarry Smith PetscInt i_start_c,j_start_c,l_start_c, m_c,n_c,p_c; 1063*47c6ae99SBarry Smith PetscInt i_start_ghost_c,j_start_ghost_c,l_start_ghost_c,m_ghost_c,n_ghost_c,p_ghost_c; 1064*47c6ae99SBarry Smith PetscInt *idx_c; 1065*47c6ae99SBarry Smith PetscInt d; 1066*47c6ae99SBarry Smith PetscInt a; 1067*47c6ae99SBarry Smith PetscInt max_agg_size; 1068*47c6ae99SBarry Smith PetscInt *fine_nodes; 1069*47c6ae99SBarry Smith PetscScalar *one_vec; 1070*47c6ae99SBarry Smith PetscInt fn_idx; 1071*47c6ae99SBarry Smith 1072*47c6ae99SBarry Smith PetscFunctionBegin; 1073*47c6ae99SBarry Smith PetscValidHeaderSpecific(dac,DM_CLASSID,1); 1074*47c6ae99SBarry Smith PetscValidHeaderSpecific(daf,DM_CLASSID,2); 1075*47c6ae99SBarry Smith PetscValidPointer(rest,3); 1076*47c6ae99SBarry Smith 1077*47c6ae99SBarry Smith ierr = DAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&wrapc,&stc);CHKERRQ(ierr); 1078*47c6ae99SBarry Smith ierr = DAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&wrapf,&stf);CHKERRQ(ierr); 1079*47c6ae99SBarry Smith if (dimc != dimf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of DA do not match %D %D",dimc,dimf);CHKERRQ(ierr); 1080*47c6ae99SBarry Smith if (dofc != doff) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"DOF of DA do not match %D %D",dofc,doff);CHKERRQ(ierr); 1081*47c6ae99SBarry Smith if (sc != sf) SETERRQ2(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil width of DA do not match %D %D",sc,sf);CHKERRQ(ierr); 1082*47c6ae99SBarry Smith if (wrapc != wrapf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Periodic type different in two DAs");CHKERRQ(ierr); 1083*47c6ae99SBarry Smith if (stc != stf) SETERRQ(((PetscObject)daf)->comm,PETSC_ERR_ARG_INCOMP,"Stencil type different in two DAs");CHKERRQ(ierr); 1084*47c6ae99SBarry Smith 1085*47c6ae99SBarry 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); 1086*47c6ae99SBarry 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); 1087*47c6ae99SBarry 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); 1088*47c6ae99SBarry Smith 1089*47c6ae99SBarry Smith if (Pc < 0) Pc = 1; 1090*47c6ae99SBarry Smith if (Pf < 0) Pf = 1; 1091*47c6ae99SBarry Smith if (Nc < 0) Nc = 1; 1092*47c6ae99SBarry Smith if (Nf < 0) Nf = 1; 1093*47c6ae99SBarry Smith 1094*47c6ae99SBarry Smith ierr = DAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);CHKERRQ(ierr); 1095*47c6ae99SBarry Smith ierr = DAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);CHKERRQ(ierr); 1096*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(daf,PETSC_NULL,&idx_f);CHKERRQ(ierr); 1097*47c6ae99SBarry Smith 1098*47c6ae99SBarry Smith ierr = DAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);CHKERRQ(ierr); 1099*47c6ae99SBarry Smith ierr = DAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,&l_start_ghost_c,&m_ghost_c,&n_ghost_c,&p_ghost_c);CHKERRQ(ierr); 1100*47c6ae99SBarry Smith ierr = DAGetGlobalIndices(dac,PETSC_NULL,&idx_c);CHKERRQ(ierr); 1101*47c6ae99SBarry Smith 1102*47c6ae99SBarry Smith /* 1103*47c6ae99SBarry Smith Basic idea is as follows. Here's a 2D example, suppose r_x, r_y are the ratios 1104*47c6ae99SBarry Smith for dimension 1 and 2 respectively. 1105*47c6ae99SBarry Smith Let (i,j) be a coarse grid node. All the fine grid nodes between r_x*i and r_x*(i+1) 1106*47c6ae99SBarry Smith and r_y*j and r_y*(j+1) will be grouped into the same coarse grid agregate. 1107*47c6ae99SBarry Smith Each specific dof on the fine grid is mapped to one dof on the coarse grid. 1108*47c6ae99SBarry Smith */ 1109*47c6ae99SBarry Smith 1110*47c6ae99SBarry Smith max_agg_size = (Mf/Mc+1)*(Nf/Nc+1)*(Pf/Pc+1); 1111*47c6ae99SBarry Smith 1112*47c6ae99SBarry Smith /* create the matrix that will contain the restriction operator */ 1113*47c6ae99SBarry 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, 1114*47c6ae99SBarry Smith max_agg_size, PETSC_NULL, max_agg_size, PETSC_NULL, rest);CHKERRQ(ierr); 1115*47c6ae99SBarry Smith 1116*47c6ae99SBarry Smith /* store nodes in the fine grid here */ 1117*47c6ae99SBarry Smith ierr = PetscMalloc2(max_agg_size,PetscScalar, &one_vec,max_agg_size,PetscInt, &fine_nodes);CHKERRQ(ierr); 1118*47c6ae99SBarry Smith for(i=0; i<max_agg_size; i++) one_vec[i] = 1.0; 1119*47c6ae99SBarry Smith 1120*47c6ae99SBarry Smith /* loop over all coarse nodes */ 1121*47c6ae99SBarry Smith for (l_c=l_start_c; l_c<l_start_c+p_c; l_c++) { 1122*47c6ae99SBarry Smith for (j_c=j_start_c; j_c<j_start_c+n_c; j_c++) { 1123*47c6ae99SBarry Smith for (i_c=i_start_c; i_c<i_start_c+m_c; i_c++) { 1124*47c6ae99SBarry Smith for(d=0; d<dofc; d++) { 1125*47c6ae99SBarry Smith /* convert to local "natural" numbering and then to PETSc global numbering */ 1126*47c6ae99SBarry 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; 1127*47c6ae99SBarry Smith 1128*47c6ae99SBarry Smith fn_idx = 0; 1129*47c6ae99SBarry Smith /* Corresponding fine points are all points (i_f, j_f, l_f) such that 1130*47c6ae99SBarry Smith i_c*Mf/Mc <= i_f < (i_c+1)*Mf/Mc 1131*47c6ae99SBarry Smith (same for other dimensions) 1132*47c6ae99SBarry Smith */ 1133*47c6ae99SBarry Smith for (l=l_c*Pf/Pc; l<PetscMin((l_c+1)*Pf/Pc,Pf); l++) { 1134*47c6ae99SBarry Smith for (j=j_c*Nf/Nc; j<PetscMin((j_c+1)*Nf/Nc,Nf); j++) { 1135*47c6ae99SBarry Smith for (i=i_c*Mf/Mc; i<PetscMin((i_c+1)*Mf/Mc,Mf); i++) { 1136*47c6ae99SBarry 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; 1137*47c6ae99SBarry Smith fn_idx++; 1138*47c6ae99SBarry Smith } 1139*47c6ae99SBarry Smith } 1140*47c6ae99SBarry Smith } 1141*47c6ae99SBarry Smith /* add all these points to one aggregate */ 1142*47c6ae99SBarry Smith ierr = MatSetValues(*rest, 1, &a, fn_idx, fine_nodes, one_vec, INSERT_VALUES);CHKERRQ(ierr); 1143*47c6ae99SBarry Smith } 1144*47c6ae99SBarry Smith } 1145*47c6ae99SBarry Smith } 1146*47c6ae99SBarry Smith } 1147*47c6ae99SBarry Smith ierr = PetscFree2(one_vec,fine_nodes);CHKERRQ(ierr); 1148*47c6ae99SBarry Smith ierr = MatAssemblyBegin(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1149*47c6ae99SBarry Smith ierr = MatAssemblyEnd(*rest, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1150*47c6ae99SBarry Smith PetscFunctionReturn(0); 1151*47c6ae99SBarry Smith } 1152