1c4762a1bSJed Brown static char help[] = "Tests MatView()/MatLoad() with binary viewers for SBAIJ 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 PetscBool seqsbaij,mpisbaij,sbaij; 19c4762a1bSJed Brown PetscErrorCode ierr; 20c4762a1bSJed Brown 21c4762a1bSJed Brown PetscFunctionBegin; 22c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 23c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 24c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&seqsbaij);CHKERRQ(ierr); 25c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)A,MATMPISBAIJ,&mpisbaij);CHKERRQ(ierr); 26c4762a1bSJed Brown sbaij = (seqsbaij||mpisbaij) ? PETSC_TRUE : PETSC_FALSE; 27c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 28c4762a1bSJed Brown for (j=(sbaij?i:0); j<N; j++) { 29c4762a1bSJed Brown ierr = MatGetValue(A,i,j,&val);CHKERRQ(ierr); 30c4762a1bSJed Brown v = MakeValue(i,j,M); w = PetscRealPart(val); 31*2c71b3e2SJacob 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); 32c4762a1bSJed Brown } 33c4762a1bSJed Brown } 34c4762a1bSJed Brown PetscFunctionReturn(0); 35c4762a1bSJed Brown } 36c4762a1bSJed Brown 37c4762a1bSJed Brown int main(int argc,char **args) 38c4762a1bSJed Brown { 39c4762a1bSJed Brown Mat A; 40c4762a1bSJed Brown PetscInt M = 24,N = 24,bs = 3; 41c4762a1bSJed Brown PetscInt rstart,rend,i,j; 42c4762a1bSJed Brown PetscErrorCode ierr; 43c4762a1bSJed Brown PetscViewer view; 44c4762a1bSJed Brown 45c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr; 46c4762a1bSJed Brown /* 47c4762a1bSJed Brown Create a parallel SBAIJ matrix shared by all processors 48c4762a1bSJed Brown */ 49b5675b0fSBarry Smith ierr = MatCreateSBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,M,N,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&A);CHKERRQ(ierr); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* 52c4762a1bSJed Brown Set values into the matrix 53c4762a1bSJed Brown */ 54c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 56c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 57c4762a1bSJed Brown for (j=0; j<N; j++) { 58c4762a1bSJed Brown PetscReal v = MakeValue(i,j,M); 59c4762a1bSJed Brown if (PetscAbsReal(v) > 0) { 60c4762a1bSJed Brown ierr = MatSetValue(A,i,j,v,INSERT_VALUES);CHKERRQ(ierr); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown } 63c4762a1bSJed Brown } 64c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 65c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = MatViewFromOptions(A,NULL,"-mat_base_view");CHKERRQ(ierr); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* 69c4762a1bSJed Brown Store the binary matrix to a file 70c4762a1bSJed Brown */ 71c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view);CHKERRQ(ierr); 72c4762a1bSJed Brown for (i=0; i<3; i++) { 73c4762a1bSJed Brown ierr = MatView(A,view);CHKERRQ(ierr); 74c4762a1bSJed Brown } 75c4762a1bSJed Brown ierr = PetscViewerDestroy(&view);CHKERRQ(ierr); 76c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* 79c4762a1bSJed Brown Now reload the matrix and check its values 80c4762a1bSJed Brown */ 81c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view);CHKERRQ(ierr); 82c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 83c4762a1bSJed Brown ierr = MatSetType(A,MATSBAIJ);CHKERRQ(ierr); 84c4762a1bSJed Brown for (i=0; i<3; i++) { 85c4762a1bSJed Brown if (i > 0) {ierr = MatZeroEntries(A);CHKERRQ(ierr);} 86c4762a1bSJed Brown ierr = MatLoad(A,view);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = CheckValuesAIJ(A);CHKERRQ(ierr); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown ierr = PetscViewerDestroy(&view);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = MatViewFromOptions(A,NULL,"-mat_load_view");CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* 94c4762a1bSJed Brown Reload in SEQSBAIJ matrix and check its values 95c4762a1bSJed Brown */ 96c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"matrix.dat",FILE_MODE_READ,&view);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&A);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr); 99c4762a1bSJed Brown for (i=0; i<3; i++) { 100c4762a1bSJed Brown if (i > 0) {ierr = MatZeroEntries(A);CHKERRQ(ierr);} 101c4762a1bSJed Brown ierr = MatLoad(A,view);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = CheckValuesAIJ(A);CHKERRQ(ierr); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown ierr = PetscViewerDestroy(&view);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* 108c4762a1bSJed Brown Reload in MPISBAIJ matrix and check its values 109c4762a1bSJed Brown */ 110c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr); 113c4762a1bSJed Brown for (i=0; i<3; i++) { 114c4762a1bSJed Brown if (i > 0) {ierr = MatZeroEntries(A);CHKERRQ(ierr);} 115c4762a1bSJed Brown ierr = MatLoad(A,view);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = CheckValuesAIJ(A);CHKERRQ(ierr); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown ierr = PetscViewerDestroy(&view);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 120c4762a1bSJed Brown 121c4762a1bSJed Brown ierr = PetscFinalize(); 122c4762a1bSJed Brown return ierr; 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125c4762a1bSJed Brown /*TEST 126c4762a1bSJed Brown 127c4762a1bSJed Brown testset: 128c4762a1bSJed Brown args: -viewer_binary_mpiio 0 129c4762a1bSJed Brown output_file: output/ex50.out 130c4762a1bSJed Brown test: 131c4762a1bSJed Brown suffix: stdio_1 132c4762a1bSJed Brown nsize: 1 133c4762a1bSJed Brown test: 134c4762a1bSJed Brown suffix: stdio_2 135c4762a1bSJed Brown nsize: 2 136c4762a1bSJed Brown test: 137c4762a1bSJed Brown suffix: stdio_3 138c4762a1bSJed Brown nsize: 3 139c4762a1bSJed Brown test: 140c4762a1bSJed Brown suffix: stdio_4 141c4762a1bSJed Brown nsize: 4 142c4762a1bSJed Brown test: 143c4762a1bSJed Brown suffix: stdio_5 144c4762a1bSJed Brown nsize: 4 145c4762a1bSJed Brown 146c4762a1bSJed Brown testset: 147c4762a1bSJed Brown requires: mpiio 148c4762a1bSJed Brown args: -viewer_binary_mpiio 1 149c4762a1bSJed Brown output_file: output/ex50.out 150c4762a1bSJed Brown test: 151c4762a1bSJed Brown suffix: mpiio_1 152c4762a1bSJed Brown nsize: 1 153c4762a1bSJed Brown test: 154c4762a1bSJed Brown suffix: mpiio_2 155c4762a1bSJed Brown nsize: 2 156c4762a1bSJed Brown test: 157c4762a1bSJed Brown suffix: mpiio_3 158c4762a1bSJed Brown nsize: 3 159c4762a1bSJed Brown test: 160c4762a1bSJed Brown suffix: mpiio_4 161c4762a1bSJed Brown nsize: 4 162c4762a1bSJed Brown test: 163c4762a1bSJed Brown suffix: mpiio_5 164c4762a1bSJed Brown nsize: 5 165c4762a1bSJed Brown 166c4762a1bSJed Brown TEST*/ 167