1*47c6ae99SBarry Smith #define PETSCDM_DLL 2*47c6ae99SBarry Smith /* 3*47c6ae99SBarry Smith Code for manipulating distributed regular 1d arrays in parallel. 4*47c6ae99SBarry Smith This file was created by Peter Mell 6/30/95 5*47c6ae99SBarry Smith */ 6*47c6ae99SBarry Smith 7*47c6ae99SBarry Smith #include "private/daimpl.h" /*I "petscda.h" I*/ 8*47c6ae99SBarry Smith 9*47c6ae99SBarry Smith const char *DAPeriodicTypes[] = {"NONPERIODIC","XPERIODIC","YPERIODIC","XYPERIODIC", 10*47c6ae99SBarry Smith "XYZPERIODIC","XZPERIODIC","YZPERIODIC","ZPERIODIC","XYZGHOSTED","DAPeriodicType","DA_",0}; 11*47c6ae99SBarry Smith 12*47c6ae99SBarry Smith #undef __FUNCT__ 13*47c6ae99SBarry Smith #define __FUNCT__ "DAView_1d" 14*47c6ae99SBarry Smith PetscErrorCode DAView_1d(DA da,PetscViewer viewer) 15*47c6ae99SBarry Smith { 16*47c6ae99SBarry Smith PetscErrorCode ierr; 17*47c6ae99SBarry Smith PetscMPIInt rank; 18*47c6ae99SBarry Smith PetscBool iascii,isdraw; 19*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 20*47c6ae99SBarry Smith 21*47c6ae99SBarry Smith PetscFunctionBegin; 22*47c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr); 23*47c6ae99SBarry Smith 24*47c6ae99SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 25*47c6ae99SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 26*47c6ae99SBarry Smith if (iascii) { 27*47c6ae99SBarry Smith PetscViewerFormat format; 28*47c6ae99SBarry Smith 29*47c6ae99SBarry Smith ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 30*47c6ae99SBarry Smith if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) { 31*47c6ae99SBarry Smith DALocalInfo info; 32*47c6ae99SBarry Smith ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 33*47c6ae99SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D m %D w %D s %D\n",rank,dd->M,dd->m,dd->w,dd->s);CHKERRQ(ierr); 34*47c6ae99SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D\n",info.xs,info.xs+info.xm);CHKERRQ(ierr); 35*47c6ae99SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 36*47c6ae99SBarry Smith } 37*47c6ae99SBarry Smith } else if (isdraw) { 38*47c6ae99SBarry Smith PetscDraw draw; 39*47c6ae99SBarry Smith double ymin = -1,ymax = 1,xmin = -1,xmax = dd->M,x; 40*47c6ae99SBarry Smith PetscInt base; 41*47c6ae99SBarry Smith char node[10]; 42*47c6ae99SBarry Smith PetscBool isnull; 43*47c6ae99SBarry Smith 44*47c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 45*47c6ae99SBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 46*47c6ae99SBarry Smith 47*47c6ae99SBarry Smith ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 48*47c6ae99SBarry Smith ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr); 49*47c6ae99SBarry Smith 50*47c6ae99SBarry Smith /* first processor draws all node lines */ 51*47c6ae99SBarry Smith if (!rank) { 52*47c6ae99SBarry Smith PetscInt xmin_tmp; 53*47c6ae99SBarry Smith ymin = 0.0; ymax = 0.3; 54*47c6ae99SBarry Smith 55*47c6ae99SBarry Smith /* ADIC doesn't like doubles in a for loop */ 56*47c6ae99SBarry Smith for (xmin_tmp =0; xmin_tmp < dd->M; xmin_tmp++) { 57*47c6ae99SBarry Smith ierr = PetscDrawLine(draw,(double)xmin_tmp,ymin,(double)xmin_tmp,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 58*47c6ae99SBarry Smith } 59*47c6ae99SBarry Smith 60*47c6ae99SBarry Smith xmin = 0.0; xmax = dd->M - 1; 61*47c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 62*47c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 63*47c6ae99SBarry Smith } 64*47c6ae99SBarry Smith 65*47c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 66*47c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 67*47c6ae99SBarry Smith 68*47c6ae99SBarry Smith /* draw my box */ 69*47c6ae99SBarry Smith ymin = 0; ymax = 0.3; xmin = dd->xs / dd->w; xmax = (dd->xe / dd->w) - 1; 70*47c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 71*47c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 72*47c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 73*47c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 74*47c6ae99SBarry Smith 75*47c6ae99SBarry Smith /* Put in index numbers */ 76*47c6ae99SBarry Smith base = dd->base / dd->w; 77*47c6ae99SBarry Smith for (x=xmin; x<=xmax; x++) { 78*47c6ae99SBarry Smith sprintf(node,"%d",(int)base++); 79*47c6ae99SBarry Smith ierr = PetscDrawString(draw,x,ymin,PETSC_DRAW_RED,node);CHKERRQ(ierr); 80*47c6ae99SBarry Smith } 81*47c6ae99SBarry Smith 82*47c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 83*47c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 84*47c6ae99SBarry Smith } else { 85*47c6ae99SBarry Smith SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DA 1d",((PetscObject)viewer)->type_name); 86*47c6ae99SBarry Smith } 87*47c6ae99SBarry Smith PetscFunctionReturn(0); 88*47c6ae99SBarry Smith } 89*47c6ae99SBarry Smith 90*47c6ae99SBarry Smith #undef __FUNCT__ 91*47c6ae99SBarry Smith #define __FUNCT__ "DAView_Private" 92*47c6ae99SBarry Smith /* 93*47c6ae99SBarry Smith Processes command line options to determine if/how a DA 94*47c6ae99SBarry Smith is to be viewed. Called by DACreateXX() 95*47c6ae99SBarry Smith */ 96*47c6ae99SBarry Smith PetscErrorCode DAView_Private(DA da) 97*47c6ae99SBarry Smith { 98*47c6ae99SBarry Smith PetscErrorCode ierr; 99*47c6ae99SBarry Smith PetscBool flg1 = PETSC_FALSE; 100*47c6ae99SBarry Smith PetscViewer view; 101*47c6ae99SBarry Smith 102*47c6ae99SBarry Smith PetscFunctionBegin; 103*47c6ae99SBarry Smith ierr = PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DA viewing options","DA");CHKERRQ(ierr); 104*47c6ae99SBarry Smith ierr = PetscOptionsTruth("-da_view","Print information about the DA's distribution","DAView",PETSC_FALSE,&flg1,PETSC_NULL);CHKERRQ(ierr); 105*47c6ae99SBarry Smith if (flg1) { 106*47c6ae99SBarry Smith ierr = PetscViewerASCIIGetStdout(((PetscObject)da)->comm,&view);CHKERRQ(ierr); 107*47c6ae99SBarry Smith ierr = DAView(da,view);CHKERRQ(ierr); 108*47c6ae99SBarry Smith } 109*47c6ae99SBarry Smith flg1 = PETSC_FALSE; 110*47c6ae99SBarry Smith ierr = PetscOptionsTruth("-da_view_draw","Draw how the DA is distributed","DAView",PETSC_FALSE,&flg1,PETSC_NULL);CHKERRQ(ierr); 111*47c6ae99SBarry Smith if (flg1) {ierr = DAView(da,PETSC_VIEWER_DRAW_(((PetscObject)da)->comm));CHKERRQ(ierr);} 112*47c6ae99SBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 113*47c6ae99SBarry Smith PetscFunctionReturn(0); 114*47c6ae99SBarry Smith } 115*47c6ae99SBarry Smith 116*47c6ae99SBarry Smith EXTERN_C_BEGIN 117*47c6ae99SBarry Smith #undef __FUNCT__ 118*47c6ae99SBarry Smith #define __FUNCT__ "DACreate_1D" 119*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DACreate_1D(DA da) 120*47c6ae99SBarry Smith { 121*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 122*47c6ae99SBarry Smith const PetscInt dim = dd->dim; 123*47c6ae99SBarry Smith const PetscInt M = dd->M; 124*47c6ae99SBarry Smith const PetscInt dof = dd->w; 125*47c6ae99SBarry Smith const PetscInt s = dd->s; 126*47c6ae99SBarry Smith const PetscInt sDist = s*dof; /* absolute stencil distance */ 127*47c6ae99SBarry Smith const PetscInt *lx = dd->lx; 128*47c6ae99SBarry Smith const DAPeriodicType wrap = dd->wrap; 129*47c6ae99SBarry Smith MPI_Comm comm; 130*47c6ae99SBarry Smith Vec local, global; 131*47c6ae99SBarry Smith VecScatter ltog, gtol; 132*47c6ae99SBarry Smith IS to, from; 133*47c6ae99SBarry Smith PetscBool flg1 = PETSC_FALSE, flg2 = PETSC_FALSE; 134*47c6ae99SBarry Smith PetscMPIInt rank, size; 135*47c6ae99SBarry Smith PetscInt i,*idx,nn,left,xs,xe,x,Xs,Xe,start,end,m; 136*47c6ae99SBarry Smith PetscErrorCode ierr; 137*47c6ae99SBarry Smith 138*47c6ae99SBarry Smith PetscFunctionBegin; 139*47c6ae99SBarry Smith if (dim != PETSC_DECIDE && dim != 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Dimension should be 1: %D",dim); 140*47c6ae99SBarry Smith if (dof < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof); 141*47c6ae99SBarry Smith if (s < 0) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s); 142*47c6ae99SBarry Smith 143*47c6ae99SBarry Smith dd->dim = 1; 144*47c6ae99SBarry Smith ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr); 145*47c6ae99SBarry Smith ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr); 146*47c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 147*47c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 148*47c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 149*47c6ae99SBarry Smith 150*47c6ae99SBarry Smith dd->m = size; 151*47c6ae99SBarry Smith m = dd->m; 152*47c6ae99SBarry Smith 153*47c6ae99SBarry Smith if (s > 0) { 154*47c6ae99SBarry Smith /* if not communicating data then should be ok to have nothing on some processes */ 155*47c6ae99SBarry Smith if (M < m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"More processes than data points! %D %D",m,M); 156*47c6ae99SBarry Smith if ((M-1) < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array is too small for stencil! %D %D",M-1,s); 157*47c6ae99SBarry Smith } 158*47c6ae99SBarry Smith 159*47c6ae99SBarry Smith /* 160*47c6ae99SBarry Smith Determine locally owned region 161*47c6ae99SBarry Smith xs is the first local node number, x is the number of local nodes 162*47c6ae99SBarry Smith */ 163*47c6ae99SBarry Smith if (!lx) { 164*47c6ae99SBarry Smith ierr = PetscOptionsGetTruth(PETSC_NULL,"-da_partition_blockcomm",&flg1,PETSC_NULL);CHKERRQ(ierr); 165*47c6ae99SBarry Smith ierr = PetscOptionsGetTruth(PETSC_NULL,"-da_partition_nodes_at_end",&flg2,PETSC_NULL);CHKERRQ(ierr); 166*47c6ae99SBarry Smith if (flg1) { /* Block Comm type Distribution */ 167*47c6ae99SBarry Smith xs = rank*M/m; 168*47c6ae99SBarry Smith x = (rank + 1)*M/m - xs; 169*47c6ae99SBarry Smith } else if (flg2) { /* The odd nodes are evenly distributed across last nodes */ 170*47c6ae99SBarry Smith x = (M + rank)/m; 171*47c6ae99SBarry Smith if (M/m == x) { xs = rank*x; } 172*47c6ae99SBarry Smith else { xs = rank*(x-1) + (M+rank)%(x*m); } 173*47c6ae99SBarry Smith } else { /* The odd nodes are evenly distributed across the first k nodes */ 174*47c6ae99SBarry Smith /* Regular PETSc Distribution */ 175*47c6ae99SBarry Smith x = M/m + ((M % m) > rank); 176*47c6ae99SBarry Smith if (rank >= (M % m)) {xs = (rank * (PetscInt)(M/m) + M % m);} 177*47c6ae99SBarry Smith else {xs = rank * (PetscInt)(M/m) + rank;} 178*47c6ae99SBarry Smith } 179*47c6ae99SBarry Smith } else { 180*47c6ae99SBarry Smith x = lx[rank]; 181*47c6ae99SBarry Smith xs = 0; 182*47c6ae99SBarry Smith for (i=0; i<rank; i++) { 183*47c6ae99SBarry Smith xs += lx[i]; 184*47c6ae99SBarry Smith } 185*47c6ae99SBarry Smith /* verify that data user provided is consistent */ 186*47c6ae99SBarry Smith left = xs; 187*47c6ae99SBarry Smith for (i=rank; i<size; i++) { 188*47c6ae99SBarry Smith left += lx[i]; 189*47c6ae99SBarry Smith } 190*47c6ae99SBarry Smith if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M %D %D",left,M); 191*47c6ae99SBarry Smith } 192*47c6ae99SBarry Smith 193*47c6ae99SBarry Smith /* From now on x,xs,xe,Xs,Xe are the exact location in the array */ 194*47c6ae99SBarry Smith x *= dof; 195*47c6ae99SBarry Smith xs *= dof; 196*47c6ae99SBarry Smith xe = xs + x; 197*47c6ae99SBarry Smith 198*47c6ae99SBarry Smith /* determine ghost region */ 199*47c6ae99SBarry Smith if (wrap == DA_XPERIODIC || wrap == DA_XYZGHOSTED) { 200*47c6ae99SBarry Smith Xs = xs - sDist; 201*47c6ae99SBarry Smith Xe = xe + sDist; 202*47c6ae99SBarry Smith } else { 203*47c6ae99SBarry Smith if ((xs-sDist) >= 0) Xs = xs-sDist; else Xs = 0; 204*47c6ae99SBarry Smith if ((xe+sDist) <= M*dof) Xe = xe+sDist; else Xe = M*dof; 205*47c6ae99SBarry Smith } 206*47c6ae99SBarry Smith 207*47c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 208*47c6ae99SBarry Smith dd->Nlocal = x; 209*47c6ae99SBarry Smith ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr); 210*47c6ae99SBarry Smith ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr); 211*47c6ae99SBarry Smith dd->nlocal = (Xe-Xs); 212*47c6ae99SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr); 213*47c6ae99SBarry Smith ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr); 214*47c6ae99SBarry Smith 215*47c6ae99SBarry Smith /* Create Local to Global Vector Scatter Context */ 216*47c6ae99SBarry Smith /* local to global inserts non-ghost point region into global */ 217*47c6ae99SBarry Smith ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr); 218*47c6ae99SBarry Smith ierr = ISCreateStride(comm,x,start,1,&to);CHKERRQ(ierr); 219*47c6ae99SBarry Smith ierr = ISCreateStride(comm,x,xs-Xs,1,&from);CHKERRQ(ierr); 220*47c6ae99SBarry Smith ierr = VecScatterCreate(local,from,global,to,<og);CHKERRQ(ierr); 221*47c6ae99SBarry Smith ierr = PetscLogObjectParent(da,ltog);CHKERRQ(ierr); 222*47c6ae99SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 223*47c6ae99SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 224*47c6ae99SBarry Smith 225*47c6ae99SBarry Smith /* Create Global to Local Vector Scatter Context */ 226*47c6ae99SBarry Smith /* global to local must retrieve ghost points */ 227*47c6ae99SBarry Smith if (wrap == DA_XYZGHOSTED) { 228*47c6ae99SBarry Smith if (size == 1) { 229*47c6ae99SBarry Smith ierr = ISCreateStride(comm,(xe-xs),sDist,1,&to);CHKERRQ(ierr); 230*47c6ae99SBarry Smith } else if (!rank) { 231*47c6ae99SBarry Smith ierr = ISCreateStride(comm,(Xe-xs),sDist,1,&to);CHKERRQ(ierr); 232*47c6ae99SBarry Smith } else if (rank == size-1) { 233*47c6ae99SBarry Smith ierr = ISCreateStride(comm,(xe-Xs),0,1,&to);CHKERRQ(ierr); 234*47c6ae99SBarry Smith } else { 235*47c6ae99SBarry Smith ierr = ISCreateStride(comm,(Xe-Xs),0,1,&to);CHKERRQ(ierr); 236*47c6ae99SBarry Smith } 237*47c6ae99SBarry Smith } else { 238*47c6ae99SBarry Smith ierr = ISCreateStride(comm,(Xe-Xs),0,1,&to);CHKERRQ(ierr); 239*47c6ae99SBarry Smith } 240*47c6ae99SBarry Smith 241*47c6ae99SBarry Smith ierr = PetscMalloc((x+2*sDist)*sizeof(PetscInt),&idx);CHKERRQ(ierr); 242*47c6ae99SBarry Smith ierr = PetscLogObjectMemory(da,(x+2*sDist)*sizeof(PetscInt));CHKERRQ(ierr); 243*47c6ae99SBarry Smith 244*47c6ae99SBarry Smith nn = 0; 245*47c6ae99SBarry Smith if (wrap == DA_XPERIODIC) { /* Handle all cases with wrap first */ 246*47c6ae99SBarry Smith 247*47c6ae99SBarry Smith for (i=0; i<sDist; i++) { /* Left ghost points */ 248*47c6ae99SBarry Smith if ((xs-sDist+i)>=0) { idx[nn++] = xs-sDist+i;} 249*47c6ae99SBarry Smith else { idx[nn++] = M*dof+(xs-sDist+i);} 250*47c6ae99SBarry Smith } 251*47c6ae99SBarry Smith 252*47c6ae99SBarry Smith for (i=0; i<x; i++) { idx [nn++] = xs + i;} /* Non-ghost points */ 253*47c6ae99SBarry Smith 254*47c6ae99SBarry Smith for (i=0; i<sDist; i++) { /* Right ghost points */ 255*47c6ae99SBarry Smith if ((xe+i)<M*dof) { idx [nn++] = xe+i; } 256*47c6ae99SBarry Smith else { idx [nn++] = (xe+i) - M*dof;} 257*47c6ae99SBarry Smith } 258*47c6ae99SBarry Smith } else if (wrap == DA_XYZGHOSTED) { 259*47c6ae99SBarry Smith 260*47c6ae99SBarry Smith if (sDist <= xs) {for (i=0; i<sDist; i++) {idx[nn++] = xs - sDist + i;}} 261*47c6ae99SBarry Smith 262*47c6ae99SBarry Smith for (i=0; i<x; i++) { idx [nn++] = xs + i;} 263*47c6ae99SBarry Smith 264*47c6ae99SBarry Smith if ((xe+sDist)<=M*dof) {for (i=0; i<sDist; i++) {idx[nn++]=xe+i;}} 265*47c6ae99SBarry Smith 266*47c6ae99SBarry Smith } else { /* Now do all cases with no wrapping */ 267*47c6ae99SBarry Smith 268*47c6ae99SBarry Smith if (sDist <= xs) {for (i=0; i<sDist; i++) {idx[nn++] = xs - sDist + i;}} 269*47c6ae99SBarry Smith else {for (i=0; i<xs; i++) {idx[nn++] = i;}} 270*47c6ae99SBarry Smith 271*47c6ae99SBarry Smith for (i=0; i<x; i++) { idx [nn++] = xs + i;} 272*47c6ae99SBarry Smith 273*47c6ae99SBarry Smith if ((xe+sDist)<=M*dof) {for (i=0; i<sDist; i++) {idx[nn++]=xe+i;}} 274*47c6ae99SBarry Smith else {for (i=xe; i<(M*dof); i++) {idx[nn++]=i;}} 275*47c6ae99SBarry Smith } 276*47c6ae99SBarry Smith 277*47c6ae99SBarry Smith ierr = ISCreateGeneral(comm,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 278*47c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 279*47c6ae99SBarry Smith ierr = PetscLogObjectParent(da,to);CHKERRQ(ierr); 280*47c6ae99SBarry Smith ierr = PetscLogObjectParent(da,from);CHKERRQ(ierr); 281*47c6ae99SBarry Smith ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr); 282*47c6ae99SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 283*47c6ae99SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 284*47c6ae99SBarry Smith ierr = VecDestroy(local);CHKERRQ(ierr); 285*47c6ae99SBarry Smith ierr = VecDestroy(global);CHKERRQ(ierr); 286*47c6ae99SBarry Smith 287*47c6ae99SBarry Smith dd->xs = xs; dd->xe = xe; dd->ys = 0; dd->ye = 1; dd->zs = 0; dd->ze = 1; 288*47c6ae99SBarry Smith dd->Xs = Xs; dd->Xe = Xe; dd->Ys = 0; dd->Ye = 1; dd->Zs = 0; dd->Ze = 1; 289*47c6ae99SBarry Smith 290*47c6ae99SBarry Smith dd->gtol = gtol; 291*47c6ae99SBarry Smith dd->ltog = ltog; 292*47c6ae99SBarry Smith dd->base = xs; 293*47c6ae99SBarry Smith da->ops->view = DAView_1d; 294*47c6ae99SBarry Smith 295*47c6ae99SBarry Smith /* 296*47c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 297*47c6ae99SBarry Smith of VecSetValuesLocal(). 298*47c6ae99SBarry Smith */ 299*47c6ae99SBarry Smith if (wrap == DA_XYZGHOSTED) { 300*47c6ae99SBarry Smith PetscInt *tmpidx; 301*47c6ae99SBarry Smith if (size == 1) { 302*47c6ae99SBarry Smith ierr = PetscMalloc((nn+2*sDist)*sizeof(PetscInt),&tmpidx);CHKERRQ(ierr); 303*47c6ae99SBarry Smith for (i=0; i<sDist; i++) tmpidx[i] = -1; 304*47c6ae99SBarry Smith ierr = PetscMemcpy(tmpidx+sDist,idx,nn*sizeof(PetscInt));CHKERRQ(ierr); 305*47c6ae99SBarry Smith for (i=nn+sDist; i<nn+2*sDist; i++) tmpidx[i] = -1; 306*47c6ae99SBarry Smith ierr = PetscFree(idx);CHKERRQ(ierr); 307*47c6ae99SBarry Smith idx = tmpidx; 308*47c6ae99SBarry Smith nn += 2*sDist; 309*47c6ae99SBarry Smith } else if (!rank) { /* must preprend -1 marker for ghost location that have no global value */ 310*47c6ae99SBarry Smith ierr = PetscMalloc((nn+sDist)*sizeof(PetscInt),&tmpidx);CHKERRQ(ierr); 311*47c6ae99SBarry Smith for (i=0; i<sDist; i++) tmpidx[i] = -1; 312*47c6ae99SBarry Smith ierr = PetscMemcpy(tmpidx+sDist,idx,nn*sizeof(PetscInt));CHKERRQ(ierr); 313*47c6ae99SBarry Smith ierr = PetscFree(idx);CHKERRQ(ierr); 314*47c6ae99SBarry Smith idx = tmpidx; 315*47c6ae99SBarry Smith nn += sDist; 316*47c6ae99SBarry Smith } else if (rank == size-1) { /* must postpend -1 marker for ghost location that have no global value */ 317*47c6ae99SBarry Smith ierr = PetscMalloc((nn+sDist)*sizeof(PetscInt),&tmpidx);CHKERRQ(ierr); 318*47c6ae99SBarry Smith ierr = PetscMemcpy(tmpidx,idx,nn*sizeof(PetscInt));CHKERRQ(ierr); 319*47c6ae99SBarry Smith for (i=nn; i<nn+sDist; i++) tmpidx[i] = -1; 320*47c6ae99SBarry Smith ierr = PetscFree(idx);CHKERRQ(ierr); 321*47c6ae99SBarry Smith idx = tmpidx; 322*47c6ae99SBarry Smith nn += sDist; 323*47c6ae99SBarry Smith } 324*47c6ae99SBarry Smith } 325*47c6ae99SBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,nn,idx,PETSC_OWN_POINTER,&dd->ltogmap);CHKERRQ(ierr); 326*47c6ae99SBarry Smith ierr = ISLocalToGlobalMappingBlock(dd->ltogmap,dd->w,&dd->ltogmapb);CHKERRQ(ierr); 327*47c6ae99SBarry Smith ierr = PetscLogObjectParent(da,dd->ltogmap);CHKERRQ(ierr); 328*47c6ae99SBarry Smith 329*47c6ae99SBarry Smith dd->idx = idx; 330*47c6ae99SBarry Smith dd->Nl = nn; 331*47c6ae99SBarry Smith 332*47c6ae99SBarry Smith PetscFunctionReturn(0); 333*47c6ae99SBarry Smith } 334*47c6ae99SBarry Smith EXTERN_C_END 335*47c6ae99SBarry Smith 336*47c6ae99SBarry Smith #undef __FUNCT__ 337*47c6ae99SBarry Smith #define __FUNCT__ "DACreate1d" 338*47c6ae99SBarry Smith /*@C 339*47c6ae99SBarry Smith DACreate1d - Creates an object that will manage the communication of one-dimensional 340*47c6ae99SBarry Smith regular array data that is distributed across some processors. 341*47c6ae99SBarry Smith 342*47c6ae99SBarry Smith Collective on MPI_Comm 343*47c6ae99SBarry Smith 344*47c6ae99SBarry Smith Input Parameters: 345*47c6ae99SBarry Smith + comm - MPI communicator 346*47c6ae99SBarry Smith . wrap - type of periodicity should the array have, if any. Use 347*47c6ae99SBarry Smith either DA_NONPERIODIC or DA_XPERIODIC 348*47c6ae99SBarry Smith . M - global dimension of the array (use -M to indicate that it may be set to a different value 349*47c6ae99SBarry Smith from the command line with -da_grid_x <M>) 350*47c6ae99SBarry Smith . dof - number of degrees of freedom per node 351*47c6ae99SBarry Smith . s - stencil width 352*47c6ae99SBarry Smith - lx - array containing number of nodes in the X direction on each processor, 353*47c6ae99SBarry Smith or PETSC_NULL. If non-null, must be of length as m. 354*47c6ae99SBarry Smith 355*47c6ae99SBarry Smith Output Parameter: 356*47c6ae99SBarry Smith . da - the resulting distributed array object 357*47c6ae99SBarry Smith 358*47c6ae99SBarry Smith Options Database Key: 359*47c6ae99SBarry Smith + -da_view - Calls DAView() at the conclusion of DACreate1d() 360*47c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction; can set if M < 0 361*47c6ae99SBarry Smith - -da_refine_x - refinement factor 362*47c6ae99SBarry Smith 363*47c6ae99SBarry Smith Level: beginner 364*47c6ae99SBarry Smith 365*47c6ae99SBarry Smith Notes: 366*47c6ae99SBarry Smith The array data itself is NOT stored in the DA, it is stored in Vec objects; 367*47c6ae99SBarry Smith The appropriate vector objects can be obtained with calls to DACreateGlobalVector() 368*47c6ae99SBarry Smith and DACreateLocalVector() and calls to VecDuplicate() if more are needed. 369*47c6ae99SBarry Smith 370*47c6ae99SBarry Smith 371*47c6ae99SBarry Smith .keywords: distributed array, create, one-dimensional 372*47c6ae99SBarry Smith 373*47c6ae99SBarry Smith .seealso: DADestroy(), DAView(), DACreate2d(), DACreate3d(), DAGlobalToLocalBegin(), DASetRefinementFactor(), 374*47c6ae99SBarry Smith DAGlobalToLocalEnd(), DALocalToGlobal(), DALocalToLocalBegin(), DALocalToLocalEnd(), DAGetRefinementFactor(), 375*47c6ae99SBarry Smith DAGetInfo(), DACreateGlobalVector(), DACreateLocalVector(), DACreateNaturalVector(), DALoad(), DAView(), DAGetOwnershipRanges() 376*47c6ae99SBarry Smith 377*47c6ae99SBarry Smith @*/ 378*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DACreate1d(MPI_Comm comm, DAPeriodicType wrap, PetscInt M, PetscInt dof, PetscInt s, const PetscInt lx[], DA *da) 379*47c6ae99SBarry Smith { 380*47c6ae99SBarry Smith PetscErrorCode ierr; 381*47c6ae99SBarry Smith PetscMPIInt size; 382*47c6ae99SBarry Smith 383*47c6ae99SBarry Smith PetscFunctionBegin; 384*47c6ae99SBarry Smith ierr = DACreate(comm, da);CHKERRQ(ierr); 385*47c6ae99SBarry Smith ierr = DASetDim(*da, 1);CHKERRQ(ierr); 386*47c6ae99SBarry Smith ierr = DASetSizes(*da, M, 1, 1);CHKERRQ(ierr); 387*47c6ae99SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 388*47c6ae99SBarry Smith ierr = DASetNumProcs(*da, size, PETSC_DECIDE, PETSC_DECIDE);CHKERRQ(ierr); 389*47c6ae99SBarry Smith ierr = DASetPeriodicity(*da, wrap);CHKERRQ(ierr); 390*47c6ae99SBarry Smith ierr = DASetDof(*da, dof);CHKERRQ(ierr); 391*47c6ae99SBarry Smith ierr = DASetStencilWidth(*da, s);CHKERRQ(ierr); 392*47c6ae99SBarry Smith ierr = DASetOwnershipRanges(*da, lx, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 393*47c6ae99SBarry Smith /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */ 394*47c6ae99SBarry Smith ierr = DASetFromOptions(*da);CHKERRQ(ierr); 395*47c6ae99SBarry Smith ierr = DASetType(*da, DA1D);CHKERRQ(ierr); 396*47c6ae99SBarry Smith PetscFunctionReturn(0); 397*47c6ae99SBarry Smith } 398