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