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