1b45d2f2cSJed Brown #include <petsc-private/daimpl.h> /*I "petscdmda.h" I*/ 247c6ae99SBarry Smith 347c6ae99SBarry Smith 447c6ae99SBarry Smith #undef __FUNCT__ 5aa219208SBarry Smith #define __FUNCT__ "DMDASetDim" 647c6ae99SBarry Smith /*@ 7aa219208SBarry Smith DMDASetDim - Sets the dimension 847c6ae99SBarry Smith 9aa219208SBarry Smith Collective on DMDA 1047c6ae99SBarry Smith 1147c6ae99SBarry Smith Input Parameters: 12aa219208SBarry Smith + da - the DMDA 1347c6ae99SBarry Smith - dim - the dimension (or PETSC_DECIDE) 1447c6ae99SBarry Smith 1547c6ae99SBarry Smith Level: intermediate 1647c6ae99SBarry Smith 17aa219208SBarry Smith .seealso: DaGetDim(), DMDASetSizes() 1847c6ae99SBarry Smith @*/ 197087cfbeSBarry Smith PetscErrorCode DMDASetDim(DM da, PetscInt dim) 2047c6ae99SBarry Smith { 2147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 2247c6ae99SBarry Smith 2347c6ae99SBarry Smith PetscFunctionBegin; 2447c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 25aa219208SBarry Smith if (dd->dim > 0 && dim != dd->dim) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change DMDA dim from %D after it was set to %D",dd->dim,dim); 2647c6ae99SBarry Smith dd->dim = dim; 2747c6ae99SBarry Smith PetscFunctionReturn(0); 2847c6ae99SBarry Smith } 2947c6ae99SBarry Smith 3047c6ae99SBarry Smith #undef __FUNCT__ 31aa219208SBarry Smith #define __FUNCT__ "DMDASetSizes" 3247c6ae99SBarry Smith /*@ 33aa219208SBarry Smith DMDASetSizes - Sets the global sizes 3447c6ae99SBarry Smith 35aa219208SBarry Smith Logically Collective on DMDA 3647c6ae99SBarry Smith 3747c6ae99SBarry Smith Input Parameters: 38aa219208SBarry Smith + da - the DMDA 3947c6ae99SBarry Smith . M - the global X size (or PETSC_DECIDE) 4047c6ae99SBarry Smith . N - the global Y size (or PETSC_DECIDE) 4147c6ae99SBarry Smith - P - the global Z size (or PETSC_DECIDE) 4247c6ae99SBarry Smith 4347c6ae99SBarry Smith Level: intermediate 4447c6ae99SBarry Smith 45aa219208SBarry Smith .seealso: DMDAGetSize(), PetscSplitOwnership() 4647c6ae99SBarry Smith @*/ 477087cfbeSBarry Smith PetscErrorCode DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P) 4847c6ae99SBarry Smith { 4947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5047c6ae99SBarry Smith 5147c6ae99SBarry Smith PetscFunctionBegin; 5247c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 5347c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,M,2); 5447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,N,3); 5547c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,P,4); 56cb630486SJed Brown if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 5747c6ae99SBarry Smith 5847c6ae99SBarry Smith dd->M = M; 5947c6ae99SBarry Smith dd->N = N; 6047c6ae99SBarry Smith dd->P = P; 6147c6ae99SBarry Smith PetscFunctionReturn(0); 6247c6ae99SBarry Smith } 6347c6ae99SBarry Smith 6447c6ae99SBarry Smith #undef __FUNCT__ 65aa219208SBarry Smith #define __FUNCT__ "DMDASetNumProcs" 6647c6ae99SBarry Smith /*@ 67aa219208SBarry Smith DMDASetNumProcs - Sets the number of processes in each dimension 6847c6ae99SBarry Smith 69aa219208SBarry Smith Logically Collective on DMDA 7047c6ae99SBarry Smith 7147c6ae99SBarry Smith Input Parameters: 72aa219208SBarry Smith + da - the DMDA 7347c6ae99SBarry Smith . m - the number of X procs (or PETSC_DECIDE) 7447c6ae99SBarry Smith . n - the number of Y procs (or PETSC_DECIDE) 7547c6ae99SBarry Smith - p - the number of Z procs (or PETSC_DECIDE) 7647c6ae99SBarry Smith 7747c6ae99SBarry Smith Level: intermediate 7847c6ae99SBarry Smith 79aa219208SBarry Smith .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership() 8047c6ae99SBarry Smith @*/ 817087cfbeSBarry Smith PetscErrorCode DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p) 8247c6ae99SBarry Smith { 8347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 8447c6ae99SBarry Smith 8547c6ae99SBarry Smith PetscFunctionBegin; 8647c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 8747c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,m,2); 8847c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,n,3); 8947c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,p,4); 90cb630486SJed Brown if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 9147c6ae99SBarry Smith dd->m = m; 9247c6ae99SBarry Smith dd->n = n; 9347c6ae99SBarry Smith dd->p = p; 9447c6ae99SBarry Smith PetscFunctionReturn(0); 9547c6ae99SBarry Smith } 9647c6ae99SBarry Smith 9747c6ae99SBarry Smith #undef __FUNCT__ 983e24977eSEthan Coon #define __FUNCT__ "DMDASetBoundaryType" 9947c6ae99SBarry Smith /*@ 1001321219cSEthan Coon DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries. 10147c6ae99SBarry Smith 10247c6ae99SBarry Smith Not collective 10347c6ae99SBarry Smith 10447c6ae99SBarry Smith Input Parameter: 105aa219208SBarry Smith + da - The DMDA 1061321219cSEthan Coon - bx,by,bz - One of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC 10747c6ae99SBarry Smith 10847c6ae99SBarry Smith Level: intermediate 10947c6ae99SBarry Smith 11047c6ae99SBarry Smith .keywords: distributed array, periodicity 11196e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA, DMDABoundaryType 11247c6ae99SBarry Smith @*/ 1131321219cSEthan Coon PetscErrorCode DMDASetBoundaryType(DM da,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz) 11447c6ae99SBarry Smith { 11547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 11647c6ae99SBarry Smith 11747c6ae99SBarry Smith PetscFunctionBegin; 11847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1195a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,bx,2); 1205a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,by,3); 1215a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,bz,4); 122cb630486SJed Brown if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 1231321219cSEthan Coon dd->bx = bx; 1241321219cSEthan Coon dd->by = by; 1251321219cSEthan Coon dd->bz = bz; 12647c6ae99SBarry Smith PetscFunctionReturn(0); 12747c6ae99SBarry Smith } 12847c6ae99SBarry Smith 12947c6ae99SBarry Smith #undef __FUNCT__ 130aa219208SBarry Smith #define __FUNCT__ "DMDASetDof" 13147c6ae99SBarry Smith /*@ 132aa219208SBarry Smith DMDASetDof - Sets the number of degrees of freedom per vertex 13347c6ae99SBarry Smith 13447c6ae99SBarry Smith Not collective 13547c6ae99SBarry Smith 13647c6ae99SBarry Smith Input Parameter: 137aa219208SBarry Smith + da - The DMDA 13847c6ae99SBarry Smith - dof - Number of degrees of freedom 13947c6ae99SBarry Smith 14047c6ae99SBarry Smith Level: intermediate 14147c6ae99SBarry Smith 14247c6ae99SBarry Smith .keywords: distributed array, degrees of freedom 14396e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 14447c6ae99SBarry Smith @*/ 14554cfb0beSLisandro Dalcin PetscErrorCode DMDASetDof(DM da, PetscInt dof) 14647c6ae99SBarry Smith { 14747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 14847c6ae99SBarry Smith 14947c6ae99SBarry Smith PetscFunctionBegin; 15047c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 15154cfb0beSLisandro Dalcin PetscValidLogicalCollectiveInt(da,dof,2); 152cb630486SJed Brown if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 15347c6ae99SBarry Smith dd->w = dof; 1541411c6eeSJed Brown da->bs = dof; 15547c6ae99SBarry Smith PetscFunctionReturn(0); 15647c6ae99SBarry Smith } 15747c6ae99SBarry Smith 15847c6ae99SBarry Smith #undef __FUNCT__ 159aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType" 16047c6ae99SBarry Smith /*@ 161aa219208SBarry Smith DMDASetStencilType - Sets the type of the communication stencil 16247c6ae99SBarry Smith 163aa219208SBarry Smith Logically Collective on DMDA 16447c6ae99SBarry Smith 16547c6ae99SBarry Smith Input Parameter: 166aa219208SBarry Smith + da - The DMDA 167aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 16847c6ae99SBarry Smith 16947c6ae99SBarry Smith Level: intermediate 17047c6ae99SBarry Smith 17147c6ae99SBarry Smith .keywords: distributed array, stencil 17296e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 17347c6ae99SBarry Smith @*/ 1747087cfbeSBarry Smith PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) 17547c6ae99SBarry Smith { 17647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 17747c6ae99SBarry Smith 17847c6ae99SBarry Smith PetscFunctionBegin; 17947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 18047c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,stype,2); 181cb630486SJed Brown if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 18247c6ae99SBarry Smith dd->stencil_type = stype; 18347c6ae99SBarry Smith PetscFunctionReturn(0); 18447c6ae99SBarry Smith } 18547c6ae99SBarry Smith 18647c6ae99SBarry Smith #undef __FUNCT__ 187aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth" 18847c6ae99SBarry Smith /*@ 189aa219208SBarry Smith DMDASetStencilWidth - Sets the width of the communication stencil 19047c6ae99SBarry Smith 191aa219208SBarry Smith Logically Collective on DMDA 19247c6ae99SBarry Smith 19347c6ae99SBarry Smith Input Parameter: 194aa219208SBarry Smith + da - The DMDA 19547c6ae99SBarry Smith - width - The stencil width 19647c6ae99SBarry Smith 19747c6ae99SBarry Smith Level: intermediate 19847c6ae99SBarry Smith 19947c6ae99SBarry Smith .keywords: distributed array, stencil 20096e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 20147c6ae99SBarry Smith @*/ 2027087cfbeSBarry Smith PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) 20347c6ae99SBarry Smith { 20447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 20547c6ae99SBarry Smith 20647c6ae99SBarry Smith PetscFunctionBegin; 20747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 20847c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,width,2); 209cb630486SJed Brown if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 21047c6ae99SBarry Smith dd->s = width; 21147c6ae99SBarry Smith PetscFunctionReturn(0); 21247c6ae99SBarry Smith } 21347c6ae99SBarry Smith 21447c6ae99SBarry Smith #undef __FUNCT__ 215aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private" 216aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[]) 21747c6ae99SBarry Smith { 21847c6ae99SBarry Smith PetscInt i,sum; 21947c6ae99SBarry Smith 22047c6ae99SBarry Smith PetscFunctionBegin; 22147c6ae99SBarry Smith if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set"); 22247c6ae99SBarry Smith for (i=sum=0; i<m; i++) sum += lx[i]; 22347c6ae99SBarry Smith if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M); 22447c6ae99SBarry Smith PetscFunctionReturn(0); 22547c6ae99SBarry Smith } 22647c6ae99SBarry Smith 22747c6ae99SBarry Smith #undef __FUNCT__ 228aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges" 22947c6ae99SBarry Smith /*@ 230aa219208SBarry Smith DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 23147c6ae99SBarry Smith 232aa219208SBarry Smith Logically Collective on DMDA 23347c6ae99SBarry Smith 23447c6ae99SBarry Smith Input Parameter: 235aa219208SBarry Smith + da - The DMDA 23647c6ae99SBarry Smith . lx - array containing number of nodes in the X direction on each process, or PETSC_NULL. If non-null, must be of length da->m 23747c6ae99SBarry Smith . ly - array containing number of nodes in the Y direction on each process, or PETSC_NULL. If non-null, must be of length da->n 23847c6ae99SBarry Smith - lz - array containing number of nodes in the Z direction on each process, or PETSC_NULL. If non-null, must be of length da->p. 23947c6ae99SBarry Smith 24047c6ae99SBarry Smith Level: intermediate 24147c6ae99SBarry Smith 24247c6ae99SBarry Smith .keywords: distributed array 24396e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 24447c6ae99SBarry Smith @*/ 2457087cfbeSBarry Smith PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 24647c6ae99SBarry Smith { 24747c6ae99SBarry Smith PetscErrorCode ierr; 24847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 24947c6ae99SBarry Smith 25047c6ae99SBarry Smith PetscFunctionBegin; 25147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 252cb630486SJed Brown if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 25347c6ae99SBarry Smith if (lx) { 25447c6ae99SBarry Smith if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 255aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr); 25647c6ae99SBarry Smith if (!dd->lx) { 25747c6ae99SBarry Smith ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 25847c6ae99SBarry Smith } 25947c6ae99SBarry Smith ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr); 26047c6ae99SBarry Smith } 26147c6ae99SBarry Smith if (ly) { 26247c6ae99SBarry Smith if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 263aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr); 26447c6ae99SBarry Smith if (!dd->ly) { 26547c6ae99SBarry Smith ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 26647c6ae99SBarry Smith } 26747c6ae99SBarry Smith ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr); 26847c6ae99SBarry Smith } 26947c6ae99SBarry Smith if (lz) { 27047c6ae99SBarry Smith if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 271aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr); 27247c6ae99SBarry Smith if (!dd->lz) { 27347c6ae99SBarry Smith ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr); 27447c6ae99SBarry Smith } 27547c6ae99SBarry Smith ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr); 27647c6ae99SBarry Smith } 27747c6ae99SBarry Smith PetscFunctionReturn(0); 27847c6ae99SBarry Smith } 27947c6ae99SBarry Smith 28047c6ae99SBarry Smith #undef __FUNCT__ 281aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType" 28247c6ae99SBarry Smith /*@ 283aa219208SBarry Smith DMDASetInterpolationType - Sets the type of interpolation that will be 284e727c939SJed Brown returned by DMCreateInterpolation() 28547c6ae99SBarry Smith 286aa219208SBarry Smith Logically Collective on DMDA 28747c6ae99SBarry Smith 28847c6ae99SBarry Smith Input Parameter: 28947c6ae99SBarry Smith + da - initial distributed array 290aa219208SBarry Smith . ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms 29147c6ae99SBarry Smith 29247c6ae99SBarry Smith Level: intermediate 29347c6ae99SBarry Smith 294e727c939SJed Brown Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation() 29547c6ae99SBarry Smith 29647c6ae99SBarry Smith .keywords: distributed array, interpolation 29747c6ae99SBarry Smith 29896e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType 29947c6ae99SBarry Smith @*/ 3007087cfbeSBarry Smith PetscErrorCode DMDASetInterpolationType(DM da,DMDAInterpolationType ctype) 30147c6ae99SBarry Smith { 30247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 30347c6ae99SBarry Smith 30447c6ae99SBarry Smith PetscFunctionBegin; 30547c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 30647c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,ctype,2); 30747c6ae99SBarry Smith dd->interptype = ctype; 30847c6ae99SBarry Smith PetscFunctionReturn(0); 30947c6ae99SBarry Smith } 31047c6ae99SBarry Smith 3112dde6fd4SLisandro Dalcin #undef __FUNCT__ 3122dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetInterpolationType" 3132dde6fd4SLisandro Dalcin /*@ 3142dde6fd4SLisandro Dalcin DMDAGetInterpolationType - Gets the type of interpolation that will be 315e727c939SJed Brown used by DMCreateInterpolation() 3162dde6fd4SLisandro Dalcin 3172dde6fd4SLisandro Dalcin Not Collective 3182dde6fd4SLisandro Dalcin 3192dde6fd4SLisandro Dalcin Input Parameter: 3202dde6fd4SLisandro Dalcin . da - distributed array 3212dde6fd4SLisandro Dalcin 3222dde6fd4SLisandro Dalcin Output Parameter: 3232dde6fd4SLisandro Dalcin . ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms) 3242dde6fd4SLisandro Dalcin 3252dde6fd4SLisandro Dalcin Level: intermediate 3262dde6fd4SLisandro Dalcin 3272dde6fd4SLisandro Dalcin .keywords: distributed array, interpolation 3282dde6fd4SLisandro Dalcin 329e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation() 3302dde6fd4SLisandro Dalcin @*/ 3312dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype) 3322dde6fd4SLisandro Dalcin { 3332dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 3342dde6fd4SLisandro Dalcin 3352dde6fd4SLisandro Dalcin PetscFunctionBegin; 3362dde6fd4SLisandro Dalcin PetscValidHeaderSpecific(da,DM_CLASSID,1); 3372dde6fd4SLisandro Dalcin PetscValidPointer(ctype,2); 3382dde6fd4SLisandro Dalcin *ctype = dd->interptype; 3392dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 3402dde6fd4SLisandro Dalcin } 34147c6ae99SBarry Smith 34247c6ae99SBarry Smith #undef __FUNCT__ 343aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors" 34447c6ae99SBarry Smith /*@C 345aa219208SBarry Smith DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 34647c6ae99SBarry Smith processes neighbors. 34747c6ae99SBarry Smith 34847c6ae99SBarry Smith Not Collective 34947c6ae99SBarry Smith 35047c6ae99SBarry Smith Input Parameter: 351aa219208SBarry Smith . da - the DMDA object 35247c6ae99SBarry Smith 35347c6ae99SBarry Smith Output Parameters: 35447c6ae99SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. 35547c6ae99SBarry Smith this process itself is in the list 35647c6ae99SBarry Smith 35747c6ae99SBarry Smith Notes: In 2d the array is of length 9, in 3d of length 27 35847c6ae99SBarry Smith Not supported in 1d 359aa219208SBarry Smith Do not free the array, it is freed when the DMDA is destroyed. 36047c6ae99SBarry Smith 36147c6ae99SBarry Smith Fortran Notes: In fortran you must pass in an array of the appropriate length. 36247c6ae99SBarry Smith 36347c6ae99SBarry Smith Level: intermediate 36447c6ae99SBarry Smith 36547c6ae99SBarry Smith @*/ 3667087cfbeSBarry Smith PetscErrorCode DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[]) 36747c6ae99SBarry Smith { 36847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 36947c6ae99SBarry Smith PetscFunctionBegin; 37047c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 37147c6ae99SBarry Smith *ranks = dd->neighbors; 37247c6ae99SBarry Smith PetscFunctionReturn(0); 37347c6ae99SBarry Smith } 37447c6ae99SBarry Smith 37547c6ae99SBarry Smith #undef __FUNCT__ 376aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges" 37747c6ae99SBarry Smith /*@C 378aa219208SBarry Smith DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 37947c6ae99SBarry Smith 38047c6ae99SBarry Smith Not Collective 38147c6ae99SBarry Smith 38247c6ae99SBarry Smith Input Parameter: 383aa219208SBarry Smith . da - the DMDA object 38447c6ae99SBarry Smith 38547c6ae99SBarry Smith Output Parameter: 38647c6ae99SBarry Smith + lx - ownership along x direction (optional) 38747c6ae99SBarry Smith . ly - ownership along y direction (optional) 38847c6ae99SBarry Smith - lz - ownership along z direction (optional) 38947c6ae99SBarry Smith 39047c6ae99SBarry Smith Level: intermediate 39147c6ae99SBarry Smith 392aa219208SBarry Smith Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d() 39347c6ae99SBarry Smith 39447c6ae99SBarry Smith In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and 395aa219208SBarry Smith eighth arguments from DMDAGetInfo() 39647c6ae99SBarry Smith 39747c6ae99SBarry Smith In C you should not free these arrays, nor change the values in them. They will only have valid values while the 398aa219208SBarry Smith DMDA they came from still exists (has not been destroyed). 39947c6ae99SBarry Smith 400aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges() 40147c6ae99SBarry Smith @*/ 4027087cfbeSBarry Smith PetscErrorCode DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[]) 40347c6ae99SBarry Smith { 40447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 40547c6ae99SBarry Smith 40647c6ae99SBarry Smith PetscFunctionBegin; 40747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 40847c6ae99SBarry Smith if (lx) *lx = dd->lx; 40947c6ae99SBarry Smith if (ly) *ly = dd->ly; 41047c6ae99SBarry Smith if (lz) *lz = dd->lz; 41147c6ae99SBarry Smith PetscFunctionReturn(0); 41247c6ae99SBarry Smith } 41347c6ae99SBarry Smith 41447c6ae99SBarry Smith #undef __FUNCT__ 415aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor" 41647c6ae99SBarry Smith /*@ 417aa219208SBarry Smith DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined 41847c6ae99SBarry Smith 419aa219208SBarry Smith Logically Collective on DMDA 42047c6ae99SBarry Smith 42147c6ae99SBarry Smith Input Parameters: 422aa219208SBarry Smith + da - the DMDA object 42347c6ae99SBarry Smith . refine_x - ratio of fine grid to coarse in x direction (2 by default) 42447c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 42547c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 42647c6ae99SBarry Smith 42747c6ae99SBarry Smith Options Database: 42847c6ae99SBarry Smith + -da_refine_x - refinement ratio in x direction 42947c6ae99SBarry Smith . -da_refine_y - refinement ratio in y direction 43047c6ae99SBarry Smith - -da_refine_z - refinement ratio in z direction 43147c6ae99SBarry Smith 43247c6ae99SBarry Smith Level: intermediate 43347c6ae99SBarry Smith 43447c6ae99SBarry Smith Notes: Pass PETSC_IGNORE to leave a value unchanged 43547c6ae99SBarry Smith 436aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor() 43747c6ae99SBarry Smith @*/ 4387087cfbeSBarry Smith PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z) 43947c6ae99SBarry Smith { 44047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 44147c6ae99SBarry Smith 44247c6ae99SBarry Smith PetscFunctionBegin; 44347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 44447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_x,2); 44547c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_y,3); 44647c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_z,4); 44747c6ae99SBarry Smith 44847c6ae99SBarry Smith if (refine_x > 0) dd->refine_x = refine_x; 44947c6ae99SBarry Smith if (refine_y > 0) dd->refine_y = refine_y; 45047c6ae99SBarry Smith if (refine_z > 0) dd->refine_z = refine_z; 45147c6ae99SBarry Smith PetscFunctionReturn(0); 45247c6ae99SBarry Smith } 45347c6ae99SBarry Smith 45447c6ae99SBarry Smith #undef __FUNCT__ 455aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor" 45647c6ae99SBarry Smith /*@C 457aa219208SBarry Smith DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined 45847c6ae99SBarry Smith 45947c6ae99SBarry Smith Not Collective 46047c6ae99SBarry Smith 46147c6ae99SBarry Smith Input Parameter: 462aa219208SBarry Smith . da - the DMDA object 46347c6ae99SBarry Smith 46447c6ae99SBarry Smith Output Parameters: 46547c6ae99SBarry Smith + refine_x - ratio of fine grid to coarse in x direction (2 by default) 46647c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 46747c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 46847c6ae99SBarry Smith 46947c6ae99SBarry Smith Level: intermediate 47047c6ae99SBarry Smith 47147c6ae99SBarry Smith Notes: Pass PETSC_NULL for values you do not need 47247c6ae99SBarry Smith 473aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor() 47447c6ae99SBarry Smith @*/ 4757087cfbeSBarry Smith PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z) 47647c6ae99SBarry Smith { 47747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 47847c6ae99SBarry Smith 47947c6ae99SBarry Smith PetscFunctionBegin; 48047c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 48147c6ae99SBarry Smith if (refine_x) *refine_x = dd->refine_x; 48247c6ae99SBarry Smith if (refine_y) *refine_y = dd->refine_y; 48347c6ae99SBarry Smith if (refine_z) *refine_z = dd->refine_z; 48447c6ae99SBarry Smith PetscFunctionReturn(0); 48547c6ae99SBarry Smith } 48647c6ae99SBarry Smith 48747c6ae99SBarry Smith #undef __FUNCT__ 488aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix" 48947c6ae99SBarry Smith /*@C 490aa219208SBarry Smith DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix. 49147c6ae99SBarry Smith 492aa219208SBarry Smith Logically Collective on DMDA 49347c6ae99SBarry Smith 49447c6ae99SBarry Smith Input Parameters: 495aa219208SBarry Smith + da - the DMDA object 496aa219208SBarry Smith - f - the function that allocates the matrix for that specific DMDA 49747c6ae99SBarry Smith 49847c6ae99SBarry Smith Level: developer 49947c6ae99SBarry Smith 500aa219208SBarry Smith Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for 50147c6ae99SBarry Smith the diagonal and off-diagonal blocks of the matrix 50247c6ae99SBarry Smith 503950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills() 50447c6ae99SBarry Smith @*/ 5057087cfbeSBarry Smith PetscErrorCode DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, const MatType,Mat*)) 50647c6ae99SBarry Smith { 50747c6ae99SBarry Smith PetscFunctionBegin; 50847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 50925296bd5SBarry Smith da->ops->creatematrix = f; 51047c6ae99SBarry Smith PetscFunctionReturn(0); 51147c6ae99SBarry Smith } 51247c6ae99SBarry Smith 51347c6ae99SBarry Smith #undef __FUNCT__ 51494013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges" 51547c6ae99SBarry Smith /* 51647c6ae99SBarry Smith Creates "balanced" ownership ranges after refinement, constrained by the need for the 51747c6ae99SBarry Smith fine grid boundaries to fall within one stencil width of the coarse partition. 51847c6ae99SBarry Smith 51947c6ae99SBarry Smith Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 52047c6ae99SBarry Smith */ 52194013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[]) 52247c6ae99SBarry Smith { 52347c6ae99SBarry Smith PetscInt i,totalc = 0,remaining,startc = 0,startf = 0; 52447c6ae99SBarry Smith PetscErrorCode ierr; 52547c6ae99SBarry Smith 52647c6ae99SBarry Smith PetscFunctionBegin; 52747c6ae99SBarry Smith if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 52847c6ae99SBarry Smith if (ratio == 1) { 52947c6ae99SBarry Smith ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr); 53047c6ae99SBarry Smith PetscFunctionReturn(0); 53147c6ae99SBarry Smith } 53247c6ae99SBarry Smith for (i=0; i<m; i++) totalc += lc[i]; 53347c6ae99SBarry Smith remaining = (!periodic) + ratio * (totalc - (!periodic)); 53447c6ae99SBarry Smith for (i=0; i<m; i++) { 53547c6ae99SBarry Smith PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 53647c6ae99SBarry Smith if (i == m-1) lf[i] = want; 53747c6ae99SBarry Smith else { 5387aca7175SJed Brown const PetscInt nextc = startc + lc[i]; 5397aca7175SJed Brown /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one 5407aca7175SJed Brown * coarse stencil width of the first coarse node in the next subdomain. */ 5417aca7175SJed Brown while ((startf+want)/ratio < nextc - stencil_width) want++; 5427aca7175SJed Brown /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one 5437aca7175SJed Brown * coarse stencil width of the last coarse node in the current subdomain. */ 5447aca7175SJed Brown while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--; 5457aca7175SJed Brown /* Make sure all constraints are satisfied */ 546*30729d88SBarry Smith if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width) 547*30729d88SBarry Smith || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range"); 54847c6ae99SBarry Smith } 54947c6ae99SBarry Smith lf[i] = want; 55047c6ae99SBarry Smith startc += lc[i]; 55147c6ae99SBarry Smith startf += lf[i]; 55247c6ae99SBarry Smith remaining -= lf[i]; 55347c6ae99SBarry Smith } 55447c6ae99SBarry Smith PetscFunctionReturn(0); 55547c6ae99SBarry Smith } 55647c6ae99SBarry Smith 55747c6ae99SBarry Smith #undef __FUNCT__ 5582be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges" 5592be375d4SJed Brown /* 5602be375d4SJed Brown Creates "balanced" ownership ranges after coarsening, constrained by the need for the 5612be375d4SJed Brown fine grid boundaries to fall within one stencil width of the coarse partition. 5622be375d4SJed Brown 5632be375d4SJed Brown Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 5642be375d4SJed Brown */ 5652be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[]) 5662be375d4SJed Brown { 5672be375d4SJed Brown PetscInt i,totalf,remaining,startc,startf; 5682be375d4SJed Brown PetscErrorCode ierr; 5692be375d4SJed Brown 5702be375d4SJed Brown PetscFunctionBegin; 5712be375d4SJed Brown if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 5722be375d4SJed Brown if (ratio == 1) { 5732be375d4SJed Brown ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr); 5742be375d4SJed Brown PetscFunctionReturn(0); 5752be375d4SJed Brown } 5762be375d4SJed Brown for (i=0,totalf=0; i<m; i++) totalf += lf[i]; 5772be375d4SJed Brown remaining = (!periodic) + (totalf - (!periodic)) / ratio; 5782be375d4SJed Brown for (i=0,startc=0,startf=0; i<m; i++) { 5792be375d4SJed Brown PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 5802be375d4SJed Brown if (i == m-1) lc[i] = want; 5812be375d4SJed Brown else { 5822be375d4SJed Brown const PetscInt nextf = startf+lf[i]; 5832be375d4SJed Brown /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 5842be375d4SJed Brown * node is within one stencil width. */ 5852be375d4SJed Brown while (nextf/ratio < startc+want-stencil_width) want--; 5862be375d4SJed Brown /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 5872be375d4SJed Brown * fine node is within one stencil width. */ 5882be375d4SJed Brown while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++; 5892be375d4SJed Brown if (want < 0 || want > remaining 590*30729d88SBarry Smith || (nextf/ratio < startc+want-stencil_width) || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range"); 5912be375d4SJed Brown } 5922be375d4SJed Brown lc[i] = want; 5932be375d4SJed Brown startc += lc[i]; 5942be375d4SJed Brown startf += lf[i]; 5952be375d4SJed Brown remaining -= lc[i]; 5962be375d4SJed Brown } 5972be375d4SJed Brown PetscFunctionReturn(0); 5982be375d4SJed Brown } 5992be375d4SJed Brown 6002be375d4SJed Brown #undef __FUNCT__ 60194013140SBarry Smith #define __FUNCT__ "DMRefine_DA" 6027087cfbeSBarry Smith PetscErrorCode DMRefine_DA(DM da,MPI_Comm comm,DM *daref) 60347c6ae99SBarry Smith { 60447c6ae99SBarry Smith PetscErrorCode ierr; 605f3141302SJed Brown PetscInt M,N,P,i; 6069a42bb27SBarry Smith DM da2; 60747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 60847c6ae99SBarry Smith 60947c6ae99SBarry Smith PetscFunctionBegin; 61047c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 61147c6ae99SBarry Smith PetscValidPointer(daref,3); 61247c6ae99SBarry Smith 6131321219cSEthan Coon if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 61447c6ae99SBarry Smith M = dd->refine_x*dd->M; 61547c6ae99SBarry Smith } else { 61647c6ae99SBarry Smith M = 1 + dd->refine_x*(dd->M - 1); 61747c6ae99SBarry Smith } 6181321219cSEthan Coon if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 61947c6ae99SBarry Smith N = dd->refine_y*dd->N; 62047c6ae99SBarry Smith } else { 62147c6ae99SBarry Smith N = 1 + dd->refine_y*(dd->N - 1); 62247c6ae99SBarry Smith } 6231321219cSEthan Coon if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 62447c6ae99SBarry Smith P = dd->refine_z*dd->P; 62547c6ae99SBarry Smith } else { 62647c6ae99SBarry Smith P = 1 + dd->refine_z*(dd->P - 1); 62747c6ae99SBarry Smith } 62847c6ae99SBarry Smith if (dd->dim == 3) { 62947c6ae99SBarry Smith PetscInt *lx,*ly,*lz; 63047c6ae99SBarry Smith ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr); 6311321219cSEthan Coon ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 6321321219cSEthan Coon ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 6331321219cSEthan Coon ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 6341321219cSEthan Coon ierr = DMDACreate3d(((PetscObject)da)->comm,dd->bx,dd->by,dd->bz,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr); 63547c6ae99SBarry Smith ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 63647c6ae99SBarry Smith } else if (dd->dim == 2) { 63747c6ae99SBarry Smith PetscInt *lx,*ly; 63847c6ae99SBarry Smith ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr); 6391321219cSEthan Coon ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 6401321219cSEthan Coon ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 6411321219cSEthan Coon ierr = DMDACreate2d(((PetscObject)da)->comm,dd->bx,dd->by,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr); 64247c6ae99SBarry Smith ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 64347c6ae99SBarry Smith } else if (dd->dim == 1) { 64447c6ae99SBarry Smith PetscInt *lx; 64547c6ae99SBarry Smith ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); 6461321219cSEthan Coon ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 6471321219cSEthan Coon ierr = DMDACreate1d(((PetscObject)da)->comm,dd->bx,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr); 64847c6ae99SBarry Smith ierr = PetscFree(lx);CHKERRQ(ierr); 64947c6ae99SBarry Smith } 65047c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 65147c6ae99SBarry Smith 65247c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 65325296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 65425296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */ 65547c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 65647c6ae99SBarry Smith dd2->interptype = dd->interptype; 65747c6ae99SBarry Smith 65847c6ae99SBarry Smith /* copy fill information if given */ 65947c6ae99SBarry Smith if (dd->dfill) { 66047c6ae99SBarry Smith ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); 66147c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 66247c6ae99SBarry Smith } 66347c6ae99SBarry Smith if (dd->ofill) { 66447c6ae99SBarry Smith ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); 66547c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 66647c6ae99SBarry Smith } 66747c6ae99SBarry Smith /* copy the refine information */ 66847c6ae99SBarry Smith dd2->refine_x = dd->refine_x; 66947c6ae99SBarry Smith dd2->refine_y = dd->refine_y; 67047c6ae99SBarry Smith dd2->refine_z = dd->refine_z; 67147c6ae99SBarry Smith 67247c6ae99SBarry Smith /* copy vector type information */ 67347c6ae99SBarry Smith ierr = PetscFree(da2->vectype);CHKERRQ(ierr); 67447c6ae99SBarry Smith ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr); 675ddcf8b74SDave May 676efd51863SBarry Smith dd2->lf = dd->lf; 677efd51863SBarry Smith dd2->lj = dd->lj; 678efd51863SBarry Smith 679ddcf8b74SDave May /* interpolate coordinates if they are set on the coarse grid */ 680ddcf8b74SDave May if (dd->coordinates) { 681ddcf8b74SDave May DM cdaf,cdac; 682ddcf8b74SDave May Vec coordsc,coordsf; 683ddcf8b74SDave May Mat II; 684ddcf8b74SDave May 685ddcf8b74SDave May ierr = DMDAGetCoordinateDA(da,&cdac);CHKERRQ(ierr); 686ddcf8b74SDave May ierr = DMDAGetCoordinates(da,&coordsc);CHKERRQ(ierr); 687b61d3410SDave May ierr = DMDAGetCoordinateDA(da2,&cdaf);CHKERRQ(ierr); 688b61d3410SDave May /* force creation of the coordinate vector */ 689b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 690ddcf8b74SDave May ierr = DMDAGetCoordinates(da2,&coordsf);CHKERRQ(ierr); 691e727c939SJed Brown ierr = DMCreateInterpolation(cdac,cdaf,&II,PETSC_NULL);CHKERRQ(ierr); 692ddcf8b74SDave May ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr); 693fcfd50ebSBarry Smith ierr = MatDestroy(&II);CHKERRQ(ierr); 694ddcf8b74SDave May } 695f3141302SJed Brown for (i=0; i<da->bs; i++) { 696f3141302SJed Brown const char *fieldname; 697f3141302SJed Brown ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 698f3141302SJed Brown ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 699f3141302SJed Brown } 70047c6ae99SBarry Smith *daref = da2; 70147c6ae99SBarry Smith PetscFunctionReturn(0); 70247c6ae99SBarry Smith } 70347c6ae99SBarry Smith 70447c6ae99SBarry Smith #undef __FUNCT__ 70594013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA" 7067087cfbeSBarry Smith PetscErrorCode DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref) 70747c6ae99SBarry Smith { 70847c6ae99SBarry Smith PetscErrorCode ierr; 70947c6ae99SBarry Smith PetscInt M,N,P; 7109a42bb27SBarry Smith DM da2; 71147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 71247c6ae99SBarry Smith 71347c6ae99SBarry Smith PetscFunctionBegin; 71447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 71547c6ae99SBarry Smith PetscValidPointer(daref,3); 71647c6ae99SBarry Smith 7171321219cSEthan Coon if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 71847c6ae99SBarry Smith M = dd->M / dd->refine_x; 71947c6ae99SBarry Smith } else { 72047c6ae99SBarry Smith M = 1 + (dd->M - 1) / dd->refine_x; 72147c6ae99SBarry Smith } 7221321219cSEthan Coon if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 72347c6ae99SBarry Smith N = dd->N / dd->refine_y; 72447c6ae99SBarry Smith } else { 72547c6ae99SBarry Smith N = 1 + (dd->N - 1) / dd->refine_y; 72647c6ae99SBarry Smith } 7271321219cSEthan Coon if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 72847c6ae99SBarry Smith P = dd->P / dd->refine_z; 72947c6ae99SBarry Smith } else { 73047c6ae99SBarry Smith P = 1 + (dd->P - 1) / dd->refine_z; 73147c6ae99SBarry Smith } 73247c6ae99SBarry Smith if (dd->dim == 3) { 7332be375d4SJed Brown PetscInt *lx,*ly,*lz; 7342be375d4SJed Brown ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr); 7351321219cSEthan Coon ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 7361321219cSEthan Coon ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 7371321219cSEthan Coon ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 7381321219cSEthan Coon ierr = DMDACreate3d(((PetscObject)da)->comm,dd->bx,dd->by,dd->bz,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr); 7392be375d4SJed Brown ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 74047c6ae99SBarry Smith } else if (dd->dim == 2) { 7412be375d4SJed Brown PetscInt *lx,*ly; 7422be375d4SJed Brown ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr); 7431321219cSEthan Coon ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 7441321219cSEthan Coon ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 7451321219cSEthan Coon ierr = DMDACreate2d(((PetscObject)da)->comm,dd->bx,dd->by,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr); 7462be375d4SJed Brown ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 74747c6ae99SBarry Smith } else if (dd->dim == 1) { 7482be375d4SJed Brown PetscInt *lx; 7492be375d4SJed Brown ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); 7501321219cSEthan Coon ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 7511321219cSEthan Coon ierr = DMDACreate1d(((PetscObject)da)->comm,dd->bx,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr); 7522be375d4SJed Brown ierr = PetscFree(lx);CHKERRQ(ierr); 75347c6ae99SBarry Smith } 75447c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 75547c6ae99SBarry Smith 7564dcab191SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */ 75725296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */ 75825296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 75947c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 76047c6ae99SBarry Smith dd2->interptype = dd->interptype; 76147c6ae99SBarry Smith 76247c6ae99SBarry Smith /* copy fill information if given */ 76347c6ae99SBarry Smith if (dd->dfill) { 76447c6ae99SBarry Smith ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); 76547c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 76647c6ae99SBarry Smith } 76747c6ae99SBarry Smith if (dd->ofill) { 76847c6ae99SBarry Smith ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); 76947c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 77047c6ae99SBarry Smith } 77147c6ae99SBarry Smith /* copy the refine information */ 77247c6ae99SBarry Smith dd2->refine_x = dd->refine_x; 77347c6ae99SBarry Smith dd2->refine_y = dd->refine_y; 77447c6ae99SBarry Smith dd2->refine_z = dd->refine_z; 77547c6ae99SBarry Smith 77647c6ae99SBarry Smith /* copy vector type information */ 77747c6ae99SBarry Smith ierr = PetscFree(da2->vectype);CHKERRQ(ierr); 77847c6ae99SBarry Smith ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr); 77947c6ae99SBarry Smith 780644e2e5bSBarry Smith dd2->lf = dd->lf; 781644e2e5bSBarry Smith dd2->lj = dd->lj; 782644e2e5bSBarry Smith 783ddcf8b74SDave May /* inject coordinates if they are set on the fine grid */ 784ddcf8b74SDave May if (dd->coordinates) { 785ddcf8b74SDave May DM cdaf,cdac; 786ddcf8b74SDave May Vec coordsc,coordsf; 787ddcf8b74SDave May VecScatter inject; 788ddcf8b74SDave May 789ddcf8b74SDave May ierr = DMDAGetCoordinateDA(da,&cdaf);CHKERRQ(ierr); 790ddcf8b74SDave May ierr = DMDAGetCoordinates(da,&coordsf);CHKERRQ(ierr); 791b61d3410SDave May ierr = DMDAGetCoordinateDA(da2,&cdac);CHKERRQ(ierr); 792b61d3410SDave May /* force creation of the coordinate vector */ 793b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 794ddcf8b74SDave May ierr = DMDAGetCoordinates(da2,&coordsc);CHKERRQ(ierr); 795ddcf8b74SDave May 796e727c939SJed Brown ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr); 797ddcf8b74SDave May ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 798ddcf8b74SDave May ierr = VecScatterEnd(inject ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 799fcfd50ebSBarry Smith ierr = VecScatterDestroy(&inject);CHKERRQ(ierr); 800ddcf8b74SDave May } 80147c6ae99SBarry Smith *daref = da2; 80247c6ae99SBarry Smith PetscFunctionReturn(0); 80347c6ae99SBarry Smith } 80447c6ae99SBarry Smith 80547c6ae99SBarry Smith #undef __FUNCT__ 80694013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA" 8077087cfbeSBarry Smith PetscErrorCode DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[]) 80847c6ae99SBarry Smith { 80947c6ae99SBarry Smith PetscErrorCode ierr; 81047c6ae99SBarry Smith PetscInt i,n,*refx,*refy,*refz; 81147c6ae99SBarry Smith 81247c6ae99SBarry Smith PetscFunctionBegin; 81347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 81447c6ae99SBarry Smith if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 81547c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 81647c6ae99SBarry Smith PetscValidPointer(daf,3); 81747c6ae99SBarry Smith 818aa219208SBarry Smith /* Get refinement factors, defaults taken from the coarse DMDA */ 81947c6ae99SBarry Smith ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr); 82047c6ae99SBarry Smith for (i=0; i<nlevels; i++) { 821aa219208SBarry Smith ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr); 82247c6ae99SBarry Smith } 82347c6ae99SBarry Smith n = nlevels; 82447c6ae99SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr); 82547c6ae99SBarry Smith n = nlevels; 82647c6ae99SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr); 82747c6ae99SBarry Smith n = nlevels; 82847c6ae99SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr); 82947c6ae99SBarry Smith 830aa219208SBarry Smith ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr); 83194013140SBarry Smith ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr); 83247c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 833aa219208SBarry Smith ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr); 83494013140SBarry Smith ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr); 83547c6ae99SBarry Smith } 83647c6ae99SBarry Smith ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr); 83747c6ae99SBarry Smith PetscFunctionReturn(0); 83847c6ae99SBarry Smith } 83947c6ae99SBarry Smith 84047c6ae99SBarry Smith #undef __FUNCT__ 84194013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA" 8427087cfbeSBarry Smith PetscErrorCode DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[]) 84347c6ae99SBarry Smith { 84447c6ae99SBarry Smith PetscErrorCode ierr; 84547c6ae99SBarry Smith PetscInt i; 84647c6ae99SBarry Smith 84747c6ae99SBarry Smith PetscFunctionBegin; 84847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 84947c6ae99SBarry Smith if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 85047c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 85147c6ae99SBarry Smith PetscValidPointer(dac,3); 85294013140SBarry Smith ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr); 85347c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 85494013140SBarry Smith ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr); 85547c6ae99SBarry Smith } 85647c6ae99SBarry Smith PetscFunctionReturn(0); 85747c6ae99SBarry Smith } 858