1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests VecView() functionality with DMDA objects when using:"\ 3*c4762a1bSJed Brown "(i) a PetscViewer binary with MPI-IO support; and (ii) when the binary header is skipped.\n\n"; 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown #include <petscdm.h> 6*c4762a1bSJed Brown #include <petscdmda.h> 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown #define DMDA_I 5 10*c4762a1bSJed Brown #define DMDA_J 4 11*c4762a1bSJed Brown #define DMDA_K 6 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown const PetscReal dmda_i_val[] = { 1.10, 2.3006, 2.32444, 3.44006, 66.9009 }; 14*c4762a1bSJed Brown const PetscReal dmda_j_val[] = { 0.0, 0.25, 0.5, 0.75 }; 15*c4762a1bSJed Brown const PetscReal dmda_k_val[] = { 0.0, 1.1, 2.2, 3.3, 4.4, 5.5 }; 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown PetscErrorCode MyVecDump(const char fname[],PetscBool skippheader,PetscBool usempiio,Vec x) 18*c4762a1bSJed Brown { 19*c4762a1bSJed Brown MPI_Comm comm; 20*c4762a1bSJed Brown PetscViewer viewer; 21*c4762a1bSJed Brown PetscBool ismpiio,isskip; 22*c4762a1bSJed Brown PetscErrorCode ierr; 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown PetscFunctionBeginUser; 25*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr); 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown ierr = PetscViewerCreate(comm,&viewer);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = PetscViewerSetType(viewer,PETSCVIEWERBINARY);CHKERRQ(ierr); 29*c4762a1bSJed Brown if (skippheader) { ierr = PetscViewerBinarySetSkipHeader(viewer,PETSC_TRUE);CHKERRQ(ierr); } 30*c4762a1bSJed Brown ierr = PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);CHKERRQ(ierr); 31*c4762a1bSJed Brown if (usempiio) { ierr = PetscViewerBinarySetUseMPIIO(viewer,PETSC_TRUE);CHKERRQ(ierr); } 32*c4762a1bSJed Brown ierr = PetscViewerFileSetName(viewer,fname);CHKERRQ(ierr); 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown ierr = VecView(x,viewer);CHKERRQ(ierr); 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown ierr = PetscViewerBinaryGetUseMPIIO(viewer,&ismpiio);CHKERRQ(ierr); 37*c4762a1bSJed Brown if (ismpiio) { ierr = PetscPrintf(comm,"*** PetscViewer[write] using MPI-IO ***\n");CHKERRQ(ierr); } 38*c4762a1bSJed Brown ierr = PetscViewerBinaryGetSkipHeader(viewer,&isskip);CHKERRQ(ierr); 39*c4762a1bSJed Brown if (isskip) { ierr = PetscPrintf(comm,"*** PetscViewer[write] skipping header ***\n");CHKERRQ(ierr); } 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 42*c4762a1bSJed Brown PetscFunctionReturn(0); 43*c4762a1bSJed Brown } 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown PetscErrorCode MyVecLoad(const char fname[],PetscBool skippheader,PetscBool usempiio,Vec x) 46*c4762a1bSJed Brown { 47*c4762a1bSJed Brown MPI_Comm comm; 48*c4762a1bSJed Brown PetscViewer viewer; 49*c4762a1bSJed Brown PetscBool ismpiio,isskip; 50*c4762a1bSJed Brown PetscErrorCode ierr; 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown PetscFunctionBeginUser; 53*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr); 54*c4762a1bSJed Brown 55*c4762a1bSJed Brown ierr = PetscViewerCreate(comm,&viewer);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = PetscViewerSetType(viewer,PETSCVIEWERBINARY);CHKERRQ(ierr); 57*c4762a1bSJed Brown if (skippheader) { ierr = PetscViewerBinarySetSkipHeader(viewer,PETSC_TRUE);CHKERRQ(ierr); } 58*c4762a1bSJed Brown ierr = PetscViewerFileSetMode(viewer,FILE_MODE_READ);CHKERRQ(ierr); 59*c4762a1bSJed Brown if (usempiio) { ierr = PetscViewerBinarySetUseMPIIO(viewer,PETSC_TRUE);CHKERRQ(ierr); } 60*c4762a1bSJed Brown ierr = PetscViewerFileSetName(viewer,fname);CHKERRQ(ierr); 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown ierr = VecLoad(x,viewer);CHKERRQ(ierr); 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown ierr = PetscViewerBinaryGetSkipHeader(viewer,&isskip);CHKERRQ(ierr); 65*c4762a1bSJed Brown if (isskip) { ierr = PetscPrintf(comm,"*** PetscViewer[load] skipping header ***\n");CHKERRQ(ierr); } 66*c4762a1bSJed Brown ierr = PetscViewerBinaryGetUseMPIIO(viewer,&ismpiio);CHKERRQ(ierr); 67*c4762a1bSJed Brown if (ismpiio) { ierr = PetscPrintf(comm,"*** PetscViewer[load] using MPI-IO ***\n");CHKERRQ(ierr); } 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 70*c4762a1bSJed Brown PetscFunctionReturn(0); 71*c4762a1bSJed Brown } 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown PetscErrorCode DMDAVecGenerateEntries(DM dm,Vec a) 74*c4762a1bSJed Brown { 75*c4762a1bSJed Brown PetscScalar ****LA_v; 76*c4762a1bSJed Brown PetscInt i,j,k,l,si,sj,sk,ni,nj,nk,M,N,dof; 77*c4762a1bSJed Brown PetscErrorCode ierr; 78*c4762a1bSJed Brown 79*c4762a1bSJed Brown PetscFunctionBeginUser; 80*c4762a1bSJed Brown ierr = DMDAGetInfo(dm,NULL,&M,&N,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = DMDAGetCorners(dm,&si,&sj,&sk,&ni,&nj,&nk);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = DMDAVecGetArrayDOF(dm,a,&LA_v);CHKERRQ(ierr); 83*c4762a1bSJed Brown for (k=sk; k<sk+nk; k++) { 84*c4762a1bSJed Brown for (j=sj; j<sj+nj; j++) { 85*c4762a1bSJed Brown for (i=si; i<si+ni; i++) { 86*c4762a1bSJed Brown PetscScalar test_value_s; 87*c4762a1bSJed Brown 88*c4762a1bSJed Brown test_value_s = dmda_i_val[i]*((PetscScalar)i) + dmda_j_val[j]*((PetscScalar)(i+j*M)) + dmda_k_val[k]*((PetscScalar)(i + j*M + k*M*N)); 89*c4762a1bSJed Brown for (l=0; l<dof; l++) { 90*c4762a1bSJed Brown LA_v[k][j][i][l] = (PetscScalar)dof * test_value_s + (PetscScalar)l; 91*c4762a1bSJed Brown } 92*c4762a1bSJed Brown } 93*c4762a1bSJed Brown } 94*c4762a1bSJed Brown } 95*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOF(dm,a,&LA_v);CHKERRQ(ierr); 96*c4762a1bSJed Brown PetscFunctionReturn(0); 97*c4762a1bSJed Brown } 98*c4762a1bSJed Brown 99*c4762a1bSJed Brown PetscErrorCode HeaderlessBinaryReadCheck(DM dm,const char name[]) 100*c4762a1bSJed Brown { 101*c4762a1bSJed Brown PetscErrorCode ierr; 102*c4762a1bSJed Brown int fdes; 103*c4762a1bSJed Brown PetscScalar buffer[DMDA_I*DMDA_J*DMDA_K*10]; 104*c4762a1bSJed Brown PetscInt len,d,i,j,k,M,N,dof; 105*c4762a1bSJed Brown PetscMPIInt rank; 106*c4762a1bSJed Brown PetscBool dataverified = PETSC_TRUE; 107*c4762a1bSJed Brown 108*c4762a1bSJed Brown PetscFunctionBeginUser; 109*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 110*c4762a1bSJed Brown ierr = DMDAGetInfo(dm,NULL,&M,&N,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 111*c4762a1bSJed Brown len = DMDA_I*DMDA_J*DMDA_K*dof; 112*c4762a1bSJed Brown if (!rank) { 113*c4762a1bSJed Brown ierr = PetscBinaryOpen(name,FILE_MODE_READ,&fdes);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = PetscBinaryRead(fdes,buffer,len,NULL,PETSC_SCALAR);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = PetscBinaryClose(fdes);CHKERRQ(ierr); 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown for (k=0; k<DMDA_K; k++) { 118*c4762a1bSJed Brown for (j=0; j<DMDA_J; j++) { 119*c4762a1bSJed Brown for (i=0; i<DMDA_I; i++) { 120*c4762a1bSJed Brown for (d=0; d<dof; d++) { 121*c4762a1bSJed Brown PetscScalar v,test_value_s,test_value; 122*c4762a1bSJed Brown PetscInt index; 123*c4762a1bSJed Brown 124*c4762a1bSJed Brown test_value_s = dmda_i_val[i]*((PetscScalar)i) + dmda_j_val[j]*((PetscScalar)(i+j*M)) + dmda_k_val[k]*((PetscScalar)(i + j*M + k*M*N)); 125*c4762a1bSJed Brown test_value = (PetscScalar)dof * test_value_s + (PetscScalar)d; 126*c4762a1bSJed Brown 127*c4762a1bSJed Brown index = dof*(i + j*M + k*M*N) + d; 128*c4762a1bSJed Brown v = PetscAbsScalar(test_value-buffer[index]); 129*c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 130*c4762a1bSJed Brown if ((PetscRealPart(v) > 1.0e-10) || (PetscImaginaryPart(v) > 1.0e-10)) { 131*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"ERROR: Difference > 1.0e-10 occurred (delta = (%+1.12e,%+1.12e) [loc %D,%D,%D(%D)])\n",(double)PetscRealPart(test_value),(double)PetscImaginaryPart(test_value),i,j,k,d);CHKERRQ(ierr); 132*c4762a1bSJed Brown dataverified = PETSC_FALSE; 133*c4762a1bSJed Brown } 134*c4762a1bSJed Brown #else 135*c4762a1bSJed Brown if (PetscRealPart(v) > 1.0e-10) { 136*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"ERROR: Difference > 1.0e-10 occurred (delta = %+1.12e [loc %D,%D,%D(%D)])\n",(double)PetscRealPart(test_value),i,j,k,d);CHKERRQ(ierr); 137*c4762a1bSJed Brown dataverified = PETSC_FALSE; 138*c4762a1bSJed Brown } 139*c4762a1bSJed Brown #endif 140*c4762a1bSJed Brown } 141*c4762a1bSJed Brown } 142*c4762a1bSJed Brown } 143*c4762a1bSJed Brown } 144*c4762a1bSJed Brown if (dataverified) { 145*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Headerless read of data verified for: %s\n",name);CHKERRQ(ierr); 146*c4762a1bSJed Brown } 147*c4762a1bSJed Brown } 148*c4762a1bSJed Brown PetscFunctionReturn(0); 149*c4762a1bSJed Brown } 150*c4762a1bSJed Brown 151*c4762a1bSJed Brown PetscErrorCode VecCompare(Vec a,Vec b) 152*c4762a1bSJed Brown { 153*c4762a1bSJed Brown PetscInt locmin[2],locmax[2]; 154*c4762a1bSJed Brown PetscReal min[2],max[2]; 155*c4762a1bSJed Brown Vec ref; 156*c4762a1bSJed Brown PetscErrorCode ierr; 157*c4762a1bSJed Brown 158*c4762a1bSJed Brown PetscFunctionBeginUser; 159*c4762a1bSJed Brown ierr = VecMin(a,&locmin[0],&min[0]);CHKERRQ(ierr); 160*c4762a1bSJed Brown ierr = VecMax(a,&locmax[0],&max[0]);CHKERRQ(ierr); 161*c4762a1bSJed Brown 162*c4762a1bSJed Brown ierr = VecMin(b,&locmin[1],&min[1]);CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = VecMax(b,&locmax[1],&max[1]);CHKERRQ(ierr); 164*c4762a1bSJed Brown 165*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"VecCompare\n");CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," min(a) = %+1.2e [loc %D]\n",(double)min[0],locmin[0]);CHKERRQ(ierr); 167*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," max(a) = %+1.2e [loc %D]\n",(double)max[0],locmax[0]);CHKERRQ(ierr); 168*c4762a1bSJed Brown 169*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," min(b) = %+1.2e [loc %D]\n",(double)min[1],locmin[1]);CHKERRQ(ierr); 170*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," max(b) = %+1.2e [loc %D]\n",(double)max[1],locmax[1]);CHKERRQ(ierr); 171*c4762a1bSJed Brown 172*c4762a1bSJed Brown ierr = VecDuplicate(a,&ref);CHKERRQ(ierr); 173*c4762a1bSJed Brown ierr = VecCopy(a,ref);CHKERRQ(ierr); 174*c4762a1bSJed Brown ierr = VecAXPY(ref,-1.0,b);CHKERRQ(ierr); 175*c4762a1bSJed Brown ierr = VecMin(ref,&locmin[0],&min[0]);CHKERRQ(ierr); 176*c4762a1bSJed Brown if (PetscAbsReal(min[0]) > 1.0e-10) { 177*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," ERROR: min(a-b) > 1.0e-10\n");CHKERRQ(ierr); 178*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," min(a-b) = %+1.10e\n",(double)PetscAbsReal(min[0]));CHKERRQ(ierr); 179*c4762a1bSJed Brown } else { 180*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," min(a-b) < 1.0e-10\n");CHKERRQ(ierr); 181*c4762a1bSJed Brown } 182*c4762a1bSJed Brown ierr = VecDestroy(&ref);CHKERRQ(ierr); 183*c4762a1bSJed Brown PetscFunctionReturn(0); 184*c4762a1bSJed Brown } 185*c4762a1bSJed Brown 186*c4762a1bSJed Brown PetscErrorCode TestDMDAVec(PetscBool usempiio) 187*c4762a1bSJed Brown { 188*c4762a1bSJed Brown DM dm; 189*c4762a1bSJed Brown Vec x_ref,x_test; 190*c4762a1bSJed Brown PetscBool skipheader = PETSC_TRUE; 191*c4762a1bSJed Brown PetscErrorCode ierr; 192*c4762a1bSJed Brown 193*c4762a1bSJed Brown PetscFunctionBeginUser; 194*c4762a1bSJed Brown if (!usempiio) { ierr = PetscPrintf(PETSC_COMM_WORLD,"%s\n",PETSC_FUNCTION_NAME);CHKERRQ(ierr); } 195*c4762a1bSJed Brown else { ierr = PetscPrintf(PETSC_COMM_WORLD,"%s [using mpi-io]\n",PETSC_FUNCTION_NAME);CHKERRQ(ierr); } 196*c4762a1bSJed Brown ierr = DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,DMDA_I,DMDA_J,DMDA_K,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE, 197*c4762a1bSJed Brown 3,2,NULL,NULL,NULL,&dm);CHKERRQ(ierr); 198*c4762a1bSJed Brown ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 199*c4762a1bSJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 200*c4762a1bSJed Brown 201*c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm,&x_ref);CHKERRQ(ierr); 202*c4762a1bSJed Brown ierr = DMDAVecGenerateEntries(dm,x_ref);CHKERRQ(ierr); 203*c4762a1bSJed Brown 204*c4762a1bSJed Brown if (!usempiio) { 205*c4762a1bSJed Brown ierr = MyVecDump("dmda.pbvec",skipheader,PETSC_FALSE,x_ref);CHKERRQ(ierr); 206*c4762a1bSJed Brown } else { 207*c4762a1bSJed Brown ierr = MyVecDump("dmda-mpiio.pbvec",skipheader,PETSC_TRUE,x_ref);CHKERRQ(ierr); 208*c4762a1bSJed Brown } 209*c4762a1bSJed Brown 210*c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm,&x_test);CHKERRQ(ierr); 211*c4762a1bSJed Brown 212*c4762a1bSJed Brown if (!usempiio) { 213*c4762a1bSJed Brown ierr = MyVecLoad("dmda.pbvec",skipheader,usempiio,x_test);CHKERRQ(ierr); 214*c4762a1bSJed Brown } else { 215*c4762a1bSJed Brown ierr = MyVecLoad("dmda-mpiio.pbvec",skipheader,usempiio,x_test);CHKERRQ(ierr); 216*c4762a1bSJed Brown } 217*c4762a1bSJed Brown 218*c4762a1bSJed Brown ierr = VecCompare(x_ref,x_test);CHKERRQ(ierr); 219*c4762a1bSJed Brown 220*c4762a1bSJed Brown if (!usempiio) { 221*c4762a1bSJed Brown ierr = HeaderlessBinaryReadCheck(dm,"dmda.pbvec");CHKERRQ(ierr); 222*c4762a1bSJed Brown } else { 223*c4762a1bSJed Brown ierr = HeaderlessBinaryReadCheck(dm,"dmda-mpiio.pbvec");CHKERRQ(ierr); 224*c4762a1bSJed Brown } 225*c4762a1bSJed Brown ierr = VecDestroy(&x_ref);CHKERRQ(ierr); 226*c4762a1bSJed Brown ierr = VecDestroy(&x_test);CHKERRQ(ierr); 227*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 228*c4762a1bSJed Brown PetscFunctionReturn(0); 229*c4762a1bSJed Brown } 230*c4762a1bSJed Brown 231*c4762a1bSJed Brown int main(int argc,char **args) 232*c4762a1bSJed Brown { 233*c4762a1bSJed Brown PetscErrorCode ierr; 234*c4762a1bSJed Brown PetscBool usempiio = PETSC_FALSE; 235*c4762a1bSJed Brown 236*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 237*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-usempiio",&usempiio,NULL);CHKERRQ(ierr); 238*c4762a1bSJed Brown if (!usempiio) { 239*c4762a1bSJed Brown ierr = TestDMDAVec(PETSC_FALSE);CHKERRQ(ierr); 240*c4762a1bSJed Brown } else { 241*c4762a1bSJed Brown #if defined(PETSC_HAVE_MPIIO) 242*c4762a1bSJed Brown ierr = TestDMDAVec(PETSC_TRUE);CHKERRQ(ierr); 243*c4762a1bSJed Brown #else 244*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Executing TestDMDAVec(PETSC_TRUE) requires a working MPI-2 implementation\n");CHKERRQ(ierr); 245*c4762a1bSJed Brown #endif 246*c4762a1bSJed Brown } 247*c4762a1bSJed Brown ierr = PetscFinalize(); 248*c4762a1bSJed Brown return ierr; 249*c4762a1bSJed Brown } 250*c4762a1bSJed Brown 251*c4762a1bSJed Brown 252*c4762a1bSJed Brown /*TEST 253*c4762a1bSJed Brown 254*c4762a1bSJed Brown test: 255*c4762a1bSJed Brown 256*c4762a1bSJed Brown test: 257*c4762a1bSJed Brown suffix: 2 258*c4762a1bSJed Brown nsize: 12 259*c4762a1bSJed Brown 260*c4762a1bSJed Brown test: 261*c4762a1bSJed Brown suffix: 3 262*c4762a1bSJed Brown nsize: 12 263*c4762a1bSJed Brown requires: define(PETSC_HAVE_MPIIO) 264*c4762a1bSJed Brown args: -usempiio 265*c4762a1bSJed Brown 266*c4762a1bSJed Brown TEST*/ 267*c4762a1bSJed Brown 268