147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 3aa219208SBarry Smith Plots vectors obtained with DMDACreate1d() 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 747c6ae99SBarry Smith 847c6ae99SBarry Smith /*@ 9aa219208SBarry Smith DMDASetUniformCoordinates - Sets a DMDA coordinates to be a uniform grid 1047c6ae99SBarry Smith 11d083f849SBarry Smith Collective on da 1247c6ae99SBarry Smith 1347c6ae99SBarry Smith Input Parameters: 1447c6ae99SBarry Smith + da - the distributed array object 1547c6ae99SBarry Smith . xmin,xmax - extremes in the x direction 16362c83fbSPatrick Sanan . ymin,ymax - extremes in the y direction (value ignored for 1 dimensional problems) 17362c83fbSPatrick Sanan - zmin,zmax - extremes in the z direction (value ignored for 1 or 2 dimensional problems) 1847c6ae99SBarry Smith 1947c6ae99SBarry Smith Level: beginner 2047c6ae99SBarry Smith 213bd220d7SPatrick Sanan .seealso: DMSetCoordinates(), DMGetCoordinates(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMStagSetUniformCoordinates() 2247c6ae99SBarry Smith 2347c6ae99SBarry Smith @*/ 247087cfbeSBarry Smith PetscErrorCode DMDASetUniformCoordinates(DM da,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax) 2547c6ae99SBarry Smith { 2647c6ae99SBarry Smith MPI_Comm comm; 279a42bb27SBarry Smith DM cda; 28a5bc1bf3SBarry Smith DM_DA *dd = (DM_DA*)da->data; 29bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by,bz; 3047c6ae99SBarry Smith Vec xcoor; 3147c6ae99SBarry Smith PetscScalar *coors; 3247c6ae99SBarry Smith PetscReal hx,hy,hz_; 3347c6ae99SBarry Smith PetscInt i,j,k,M,N,P,istart,isize,jstart,jsize,kstart,ksize,dim,cnt; 3447c6ae99SBarry Smith 3547c6ae99SBarry Smith PetscFunctionBegin; 36a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 377a8be351SBarry Smith PetscCheck(dd->gtol,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set coordinates until after DMDA has been setup"); 389566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,&M,&N,&P,NULL,NULL,NULL,NULL,NULL,&bx,&by,&bz,NULL)); 39*08401ef6SPierre Jolivet PetscCheck(xmax >= xmin,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"xmax must be larger than xmin %g %g",(double)xmin,(double)xmax); 40*08401ef6SPierre Jolivet PetscCheck(!(dim > 1) || !(ymax < ymin),PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"ymax must be larger than ymin %g %g",(double)ymin,(double)ymax); 41*08401ef6SPierre Jolivet PetscCheck(!(dim > 2) || !(zmax < zmin),PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"zmax must be larger than zmin %g %g",(double)zmin,(double)zmax); 429566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da,&comm)); 439566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&istart,&jstart,&kstart,&isize,&jsize,&ksize)); 449566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda)); 459566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(cda, &xcoor)); 4647c6ae99SBarry Smith if (dim == 1) { 47bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax-xmin)/M; 48ce00eea3SSatish Balay else hx = (xmax-xmin)/(M-1); 499566063dSJacob Faibussowitsch PetscCall(VecGetArray(xcoor,&coors)); 5047c6ae99SBarry Smith for (i=0; i<isize; i++) { 5147c6ae99SBarry Smith coors[i] = xmin + hx*(i+istart); 5247c6ae99SBarry Smith } 539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xcoor,&coors)); 5447c6ae99SBarry Smith } else if (dim == 2) { 55bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax-xmin)/(M); 5647c6ae99SBarry Smith else hx = (xmax-xmin)/(M-1); 57bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) hy = (ymax-ymin)/(N); 5847c6ae99SBarry Smith else hy = (ymax-ymin)/(N-1); 599566063dSJacob Faibussowitsch PetscCall(VecGetArray(xcoor,&coors)); 6047c6ae99SBarry Smith cnt = 0; 6147c6ae99SBarry Smith for (j=0; j<jsize; j++) { 6247c6ae99SBarry Smith for (i=0; i<isize; i++) { 6347c6ae99SBarry Smith coors[cnt++] = xmin + hx*(i+istart); 6447c6ae99SBarry Smith coors[cnt++] = ymin + hy*(j+jstart); 6547c6ae99SBarry Smith } 6647c6ae99SBarry Smith } 679566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xcoor,&coors)); 6847c6ae99SBarry Smith } else if (dim == 3) { 69bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax-xmin)/(M); 7047c6ae99SBarry Smith else hx = (xmax-xmin)/(M-1); 71bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) hy = (ymax-ymin)/(N); 7247c6ae99SBarry Smith else hy = (ymax-ymin)/(N-1); 73bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC) hz_ = (zmax-zmin)/(P); 7447c6ae99SBarry Smith else hz_ = (zmax-zmin)/(P-1); 759566063dSJacob Faibussowitsch PetscCall(VecGetArray(xcoor,&coors)); 7647c6ae99SBarry Smith cnt = 0; 7747c6ae99SBarry Smith for (k=0; k<ksize; k++) { 7847c6ae99SBarry Smith for (j=0; j<jsize; j++) { 7947c6ae99SBarry Smith for (i=0; i<isize; i++) { 8047c6ae99SBarry Smith coors[cnt++] = xmin + hx*(i+istart); 8147c6ae99SBarry Smith coors[cnt++] = ymin + hy*(j+jstart); 8247c6ae99SBarry Smith coors[cnt++] = zmin + hz_*(k+kstart); 8347c6ae99SBarry Smith } 8447c6ae99SBarry Smith } 8547c6ae99SBarry Smith } 869566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xcoor,&coors)); 8798921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot create uniform coordinates for this dimension %D",dim); 889566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(da,xcoor)); 899566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)xcoor)); 909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xcoor)); 9147c6ae99SBarry Smith PetscFunctionReturn(0); 9247c6ae99SBarry Smith } 9347c6ae99SBarry Smith 949ae14b6eSBarry Smith /* 959ae14b6eSBarry Smith Allows a user to select a subset of the fields to be drawn by VecView() when the vector comes from a DMDA 969ae14b6eSBarry Smith */ 974e6118eeSBarry Smith PetscErrorCode DMDASelectFields(DM da,PetscInt *outfields,PetscInt **fields) 984e6118eeSBarry Smith { 994e6118eeSBarry Smith PetscInt step,ndisplayfields,*displayfields,k,j; 1004e6118eeSBarry Smith PetscBool flg; 1014e6118eeSBarry Smith 1024e6118eeSBarry Smith PetscFunctionBegin; 1039566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&step,NULL,NULL,NULL,NULL,NULL)); 1049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(step,&displayfields)); 1054e6118eeSBarry Smith for (k=0; k<step; k++) displayfields[k] = k; 1064e6118eeSBarry Smith ndisplayfields = step; 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL,NULL,"-draw_fields",displayfields,&ndisplayfields,&flg)); 1084e6118eeSBarry Smith if (!ndisplayfields) ndisplayfields = step; 1094e6118eeSBarry Smith if (!flg) { 1104e6118eeSBarry Smith char **fields; 1114e6118eeSBarry Smith const char *fieldname; 1124e6118eeSBarry Smith PetscInt nfields = step; 1139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(step,&fields)); 1149566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(NULL,NULL,"-draw_fields_by_name",fields,&nfields,&flg)); 1154e6118eeSBarry Smith if (flg) { 1164e6118eeSBarry Smith ndisplayfields = 0; 1174e6118eeSBarry Smith for (k=0; k<nfields;k++) { 1184e6118eeSBarry Smith for (j=0; j<step; j++) { 1199566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da,j,&fieldname)); 1209566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fieldname,fields[k],&flg)); 1214e6118eeSBarry Smith if (flg) { 1224e6118eeSBarry Smith goto found; 1234e6118eeSBarry Smith } 1244e6118eeSBarry Smith } 12598921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Unknown fieldname %s",fields[k]); 1264e6118eeSBarry Smith found: displayfields[ndisplayfields++] = j; 1274e6118eeSBarry Smith } 1284e6118eeSBarry Smith } 129e3f3e4b6SBarry Smith for (k=0; k<nfields; k++) { 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(fields[k])); 131e3f3e4b6SBarry Smith } 1329566063dSJacob Faibussowitsch PetscCall(PetscFree(fields)); 1334e6118eeSBarry Smith } 1344e6118eeSBarry Smith *fields = displayfields; 1354e6118eeSBarry Smith *outfields = ndisplayfields; 1364e6118eeSBarry Smith PetscFunctionReturn(0); 1374e6118eeSBarry Smith } 1384e6118eeSBarry Smith 1399804daf3SBarry Smith #include <petscdraw.h> 14073f7a4c5SBarry Smith 14147c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v) 14247c6ae99SBarry Smith { 1439a42bb27SBarry Smith DM da; 144e5ab1681SLisandro Dalcin PetscMPIInt rank,size,tag; 145e5ab1681SLisandro Dalcin PetscInt i,n,N,dof,istart,isize,j,nbounds; 14647c6ae99SBarry Smith MPI_Status status; 147e5ab1681SLisandro Dalcin PetscReal min,max,xmin = 0.0,xmax = 0.0,tmp = 0.0,xgtmp = 0.0; 14847c6ae99SBarry Smith const PetscScalar *array,*xg; 14947c6ae99SBarry Smith PetscDraw draw; 150e5ab1681SLisandro Dalcin PetscBool isnull,useports = PETSC_FALSE,showmarkers = PETSC_FALSE; 15147c6ae99SBarry Smith MPI_Comm comm; 15247c6ae99SBarry Smith PetscDrawAxis axis; 15347c6ae99SBarry Smith Vec xcoor; 154bff4a2f0SMatthew G. Knepley DMBoundaryType bx; 155e5ab1681SLisandro Dalcin const char *tlabel = NULL,*xlabel = NULL; 156a56202b9SBarry Smith const PetscReal *bounds; 15703193ff8SBarry Smith PetscInt *displayfields; 15803193ff8SBarry Smith PetscInt k,ndisplayfields; 1594e6118eeSBarry Smith PetscBool hold; 160e5ab1681SLisandro Dalcin PetscDrawViewPorts *ports = NULL; 161e5ab1681SLisandro Dalcin PetscViewerFormat format; 16247c6ae99SBarry Smith 16347c6ae99SBarry Smith PetscFunctionBegin; 1649566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(v,0,&draw)); 1659566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw,&isnull)); 16645f3bb6eSLisandro Dalcin if (isnull) PetscFunctionReturn(0); 1679566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetBounds(v,&nbounds,&bounds)); 16847c6ae99SBarry Smith 1699566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin,&da)); 1707a8be351SBarry Smith PetscCheck(da,PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 1719566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)xin,&comm)); 1729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 1739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 17447c6ae99SBarry Smith 1759566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-draw_vec_use_markers",&showmarkers,NULL)); 17647c6ae99SBarry Smith 1779566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,NULL,&N,NULL,NULL,NULL,NULL,NULL,&dof,NULL,&bx,NULL,NULL,NULL)); 1789566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&istart,NULL,NULL,&isize,NULL,NULL)); 1799566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xin,&array)); 1809566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xin,&n)); 181e5ab1681SLisandro Dalcin n = n/dof; 18247c6ae99SBarry Smith 183832b7cebSLisandro Dalcin /* Get coordinates of nodes */ 1849566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da,&xcoor)); 18547c6ae99SBarry Smith if (!xcoor) { 1869566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da,0.0,1.0,0.0,0.0,0.0,0.0)); 1879566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da,&xcoor)); 18847c6ae99SBarry Smith } 1899566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xcoor,&xg)); 1909566063dSJacob Faibussowitsch PetscCall(DMDAGetCoordinateName(da,0,&xlabel)); 191832b7cebSLisandro Dalcin 192832b7cebSLisandro Dalcin /* Determine the min and max coordinate in plot */ 193dd400576SPatrick Sanan if (rank == 0) xmin = PetscRealPart(xg[0]); 194e5ab1681SLisandro Dalcin if (rank == size-1) xmax = PetscRealPart(xg[n-1]); 1959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&xmin,1,MPIU_REAL,0,comm)); 1969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&xmax,1,MPIU_REAL,size-1,comm)); 19747c6ae99SBarry Smith 1989566063dSJacob Faibussowitsch PetscCall(DMDASelectFields(da,&ndisplayfields,&displayfields)); 1999566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(v,&format)); 2009566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-draw_ports",&useports,NULL)); 201832b7cebSLisandro Dalcin if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE; 202832b7cebSLisandro Dalcin if (useports) { 2039566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(v,0,&draw)); 2049566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawAxis(v,0,&axis)); 2059566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 2069566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 2079566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsCreate(draw,ndisplayfields,&ports)); 20873f7a4c5SBarry Smith } 209e5ab1681SLisandro Dalcin 210832b7cebSLisandro Dalcin /* Loop over each field; drawing each in a different window */ 21103193ff8SBarry Smith for (k=0; k<ndisplayfields; k++) { 2125f80ce2aSJacob Faibussowitsch PetscErrorCode ierr; 21303193ff8SBarry Smith j = displayfields[k]; 21447c6ae99SBarry Smith 215e5ab1681SLisandro Dalcin /* determine the min and max value in plot */ 2169566063dSJacob Faibussowitsch PetscCall(VecStrideMin(xin,j,NULL,&min)); 2179566063dSJacob Faibussowitsch PetscCall(VecStrideMax(xin,j,NULL,&max)); 218a56202b9SBarry Smith if (j < nbounds) { 219a56202b9SBarry Smith min = PetscMin(min,bounds[2*j]); 220a56202b9SBarry Smith max = PetscMax(max,bounds[2*j+1]); 221a56202b9SBarry Smith } 222e5ab1681SLisandro Dalcin if (min == max) { 223e5ab1681SLisandro Dalcin min -= 1.e-5; 224e5ab1681SLisandro Dalcin max += 1.e-5; 225e5ab1681SLisandro Dalcin } 22647c6ae99SBarry Smith 227832b7cebSLisandro Dalcin if (useports) { 2289566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports,k)); 2299566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da,j,&tlabel)); 230832b7cebSLisandro Dalcin } else { 231832b7cebSLisandro Dalcin const char *title; 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetHold(v,&hold)); 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(v,k,&draw)); 2349566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawAxis(v,k,&axis)); 2359566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da,j,&title)); 2369566063dSJacob Faibussowitsch if (title) PetscCall(PetscDrawSetTitle(draw,title)); 2379566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 2389566063dSJacob Faibussowitsch if (!hold) PetscCall(PetscDrawClear(draw)); 239832b7cebSLisandro Dalcin } 2409566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis,tlabel,xlabel,NULL)); 2419566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLimits(axis,xmin,xmax,min,max)); 2429566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisDraw(axis)); 24347c6ae99SBarry Smith 24447c6ae99SBarry Smith /* draw local part of vector */ 2459566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)xin,&tag)); 24647c6ae99SBarry Smith if (rank < size-1) { /*send value to right */ 2479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void*)&xg[n-1],1,MPIU_REAL,rank+1,tag,comm)); 2489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void*)&array[j+(n-1)*dof],1,MPIU_REAL,rank+1,tag,comm)); 24947c6ae99SBarry Smith } 25047c6ae99SBarry Smith if (rank) { /* receive value from left */ 2519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&xgtmp,1,MPIU_REAL,rank-1,tag,comm,&status)); 2529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag,comm,&status)); 253e5ab1681SLisandro Dalcin } 2549566063dSJacob Faibussowitsch ierr = PetscDrawCollectiveBegin(draw);PetscCall(ierr); 255e5ab1681SLisandro Dalcin if (rank) { 2569566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,xgtmp,tmp,PetscRealPart(xg[0]),PetscRealPart(array[j]),PETSC_DRAW_RED)); 2579566063dSJacob Faibussowitsch if (showmarkers) PetscCall(PetscDrawPoint(draw,xgtmp,tmp,PETSC_DRAW_BLACK)); 25847c6ae99SBarry Smith } 259e5ab1681SLisandro Dalcin for (i=1; i<n; i++) { 2609566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+dof*(i-1)]),PetscRealPart(xg[i]),PetscRealPart(array[j+dof*i]),PETSC_DRAW_RED)); 2619566063dSJacob Faibussowitsch if (showmarkers) PetscCall(PetscDrawMarker(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+dof*(i-1)]),PETSC_DRAW_BLACK)); 26247c6ae99SBarry Smith } 263e5ab1681SLisandro Dalcin if (rank == size-1) { 2649566063dSJacob Faibussowitsch if (showmarkers) PetscCall(PetscDrawMarker(draw,PetscRealPart(xg[n-1]),PetscRealPart(array[j+dof*(n-1)]),PETSC_DRAW_BLACK)); 26547c6ae99SBarry Smith } 2669566063dSJacob Faibussowitsch ierr = PetscDrawCollectiveEnd(draw);PetscCall(ierr); 2679566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 2689566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 2699566063dSJacob Faibussowitsch if (!useports) PetscCall(PetscDrawSave(draw)); 27047c6ae99SBarry Smith } 271832b7cebSLisandro Dalcin if (useports) { 2729566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(v,0,&draw)); 2739566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 274832b7cebSLisandro Dalcin } 275832b7cebSLisandro Dalcin 2769566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsDestroy(ports)); 2779566063dSJacob Faibussowitsch PetscCall(PetscFree(displayfields)); 2789566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xcoor,&xg)); 2799566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xin,&array)); 28047c6ae99SBarry Smith PetscFunctionReturn(0); 28147c6ae99SBarry Smith } 282