1c4762a1bSJed Brown static char help[] = "Tests MatView()/MatLoad() with binary viewers for BAIJ matrices.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown #include <petscviewer.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown #include <petsc/private/hashtable.h> 7c4762a1bSJed Brown static PetscReal MakeValue(PetscInt i,PetscInt j,PetscInt M) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown PetscHash_t h = PetscHashCombine(PetscHashInt(i),PetscHashInt(j)); 10c4762a1bSJed Brown return (PetscReal) ((h % 5 == 0) ? (1 + i + j*M) : 0); 11c4762a1bSJed Brown } 12c4762a1bSJed Brown 13c4762a1bSJed Brown static PetscErrorCode CheckValuesAIJ(Mat A) 14c4762a1bSJed Brown { 15c4762a1bSJed Brown PetscInt M,N,rstart,rend,i,j; 16c4762a1bSJed Brown PetscReal v,w; 17c4762a1bSJed Brown PetscScalar val; 18c4762a1bSJed Brown 19c4762a1bSJed Brown PetscFunctionBegin; 205f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 22c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 23c4762a1bSJed Brown for (j=0; j<N; j++) { 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetValue(A,i,j,&val)); 25c4762a1bSJed Brown v = MakeValue(i,j,M); w = PetscRealPart(val); 262c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsReal(v-w) > 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Matrix entry (%" PetscInt_FMT ",%" PetscInt_FMT ") should be %g, got %g",i,j,(double)v,(double)w); 27c4762a1bSJed Brown } 28c4762a1bSJed Brown } 29c4762a1bSJed Brown PetscFunctionReturn(0); 30c4762a1bSJed Brown } 31c4762a1bSJed Brown 32c4762a1bSJed Brown int main(int argc,char **args) 33c4762a1bSJed Brown { 34c4762a1bSJed Brown Mat A; 35c4762a1bSJed Brown PetscInt M = 24,N = 48,bs = 2; 36c4762a1bSJed Brown PetscInt rstart,rend,i,j; 37c4762a1bSJed Brown PetscErrorCode ierr; 38c4762a1bSJed Brown PetscViewer view; 39c4762a1bSJed Brown 40*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,NULL,help)); 41c4762a1bSJed Brown /* 42c4762a1bSJed Brown Create a parallel BAIJ matrix shared by all processors 43c4762a1bSJed Brown */ 44c4762a1bSJed Brown ierr = MatCreateBAIJ(PETSC_COMM_WORLD, 45c4762a1bSJed Brown bs, 46c4762a1bSJed Brown PETSC_DECIDE,PETSC_DECIDE, 47c4762a1bSJed Brown M,N, 48c4762a1bSJed Brown PETSC_DECIDE,NULL, 49c4762a1bSJed Brown PETSC_DECIDE,NULL, 50c4762a1bSJed Brown &A);CHKERRQ(ierr); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* 53c4762a1bSJed Brown Set values into the matrix 54c4762a1bSJed Brown */ 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 57c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 58c4762a1bSJed Brown for (j=0; j<N; j++) { 59c4762a1bSJed Brown PetscReal v = MakeValue(i,j,M); 60c4762a1bSJed Brown if (PetscAbsReal(v) > 0) { 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,i,j,v,INSERT_VALUES)); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown } 64c4762a1bSJed Brown } 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(A,NULL,"-mat_base_view")); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* 70c4762a1bSJed Brown Store the binary matrix to a file 71c4762a1bSJed Brown */ 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view)); 73c4762a1bSJed Brown for (i=0; i<3; i++) { 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,view)); 75c4762a1bSJed Brown } 765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* 80c4762a1bSJed Brown Now reload the matrix and check its values 81c4762a1bSJed Brown */ 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATBAIJ)); 85c4762a1bSJed Brown for (i=0; i<3; i++) { 865f80ce2aSJacob Faibussowitsch if (i > 0) CHKERRQ(MatZeroEntries(A)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,view)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesAIJ(A)); 89c4762a1bSJed Brown } 905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(A,NULL,"-mat_load_view")); 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* 95c4762a1bSJed Brown Reload in SEQBAIJ matrix and check its values 96c4762a1bSJed Brown */ 975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,"matrix.dat",FILE_MODE_READ,&view)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSEQBAIJ)); 100c4762a1bSJed Brown for (i=0; i<3; i++) { 1015f80ce2aSJacob Faibussowitsch if (i > 0) CHKERRQ(MatZeroEntries(A)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,view)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesAIJ(A)); 104c4762a1bSJed Brown } 1055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* 109c4762a1bSJed Brown Reload in MPIBAIJ matrix and check its values 110c4762a1bSJed Brown */ 1115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATMPIBAIJ)); 114c4762a1bSJed Brown for (i=0; i<3; i++) { 1155f80ce2aSJacob Faibussowitsch if (i > 0) CHKERRQ(MatZeroEntries(A)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,view)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesAIJ(A)); 118c4762a1bSJed Brown } 1195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 121c4762a1bSJed Brown 122*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 123*b122ec5aSJacob Faibussowitsch return 0; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown /*TEST 127c4762a1bSJed Brown 128c4762a1bSJed Brown testset: 129c4762a1bSJed Brown args: -viewer_binary_mpiio 0 130c4762a1bSJed Brown output_file: output/ex45.out 131c4762a1bSJed Brown test: 132c4762a1bSJed Brown suffix: stdio_1 133c4762a1bSJed Brown nsize: 1 134c4762a1bSJed Brown test: 135c4762a1bSJed Brown suffix: stdio_2 136c4762a1bSJed Brown nsize: 2 137c4762a1bSJed Brown test: 138c4762a1bSJed Brown suffix: stdio_3 139c4762a1bSJed Brown nsize: 3 140c4762a1bSJed Brown test: 141c4762a1bSJed Brown suffix: stdio_4 142c4762a1bSJed Brown nsize: 4 143c4762a1bSJed Brown test: 144c4762a1bSJed Brown suffix: stdio_5 145c4762a1bSJed Brown nsize: 4 146c4762a1bSJed Brown 147c4762a1bSJed Brown testset: 148c4762a1bSJed Brown requires: mpiio 149c4762a1bSJed Brown args: -viewer_binary_mpiio 1 150c4762a1bSJed Brown output_file: output/ex45.out 151c4762a1bSJed Brown test: 152c4762a1bSJed Brown suffix: mpiio_1 153c4762a1bSJed Brown nsize: 1 154c4762a1bSJed Brown test: 155c4762a1bSJed Brown suffix: mpiio_2 156c4762a1bSJed Brown nsize: 2 157c4762a1bSJed Brown test: 158c4762a1bSJed Brown suffix: mpiio_3 159c4762a1bSJed Brown nsize: 3 160c4762a1bSJed Brown test: 161c4762a1bSJed Brown suffix: mpiio_4 162c4762a1bSJed Brown nsize: 4 163c4762a1bSJed Brown test: 164c4762a1bSJed Brown suffix: mpiio_5 165c4762a1bSJed Brown nsize: 5 166c4762a1bSJed Brown 167c4762a1bSJed Brown TEST*/ 168