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 20c4762a1bSJed Brown PetscFunctionBegin; 215f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&seqsbaij)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPISBAIJ,&mpisbaij)); 25c4762a1bSJed Brown sbaij = (seqsbaij||mpisbaij) ? PETSC_TRUE : PETSC_FALSE; 26c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 27c4762a1bSJed Brown for (j=(sbaij?i:0); j<N; j++) { 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetValue(A,i,j,&val)); 29c4762a1bSJed Brown v = MakeValue(i,j,M); w = PetscRealPart(val); 302c71b3e2SJacob 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); 31c4762a1bSJed Brown } 32c4762a1bSJed Brown } 33c4762a1bSJed Brown PetscFunctionReturn(0); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown 36c4762a1bSJed Brown int main(int argc,char **args) 37c4762a1bSJed Brown { 38c4762a1bSJed Brown Mat A; 39c4762a1bSJed Brown PetscInt M = 24,N = 24,bs = 3; 40c4762a1bSJed Brown PetscInt rstart,rend,i,j; 41c4762a1bSJed Brown PetscViewer view; 42c4762a1bSJed Brown 43*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,NULL,help)); 44c4762a1bSJed Brown /* 45c4762a1bSJed Brown Create a parallel SBAIJ matrix shared by all processors 46c4762a1bSJed Brown */ 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,M,N,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&A)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* 50c4762a1bSJed Brown Set values into the matrix 51c4762a1bSJed Brown */ 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 54c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 55c4762a1bSJed Brown for (j=0; j<N; j++) { 56c4762a1bSJed Brown PetscReal v = MakeValue(i,j,M); 57c4762a1bSJed Brown if (PetscAbsReal(v) > 0) { 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,i,j,v,INSERT_VALUES)); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown } 61c4762a1bSJed Brown } 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(A,NULL,"-mat_base_view")); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* 67c4762a1bSJed Brown Store the binary matrix to a file 68c4762a1bSJed Brown */ 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view)); 70c4762a1bSJed Brown for (i=0; i<3; i++) { 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,view)); 72c4762a1bSJed Brown } 735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* 77c4762a1bSJed Brown Now reload the matrix and check its values 78c4762a1bSJed Brown */ 795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSBAIJ)); 82c4762a1bSJed Brown for (i=0; i<3; i++) { 835f80ce2aSJacob Faibussowitsch if (i > 0) CHKERRQ(MatZeroEntries(A)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,view)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesAIJ(A)); 86c4762a1bSJed Brown } 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(A,NULL,"-mat_load_view")); 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* 92c4762a1bSJed Brown Reload in SEQSBAIJ matrix and check its values 93c4762a1bSJed Brown */ 945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,"matrix.dat",FILE_MODE_READ,&view)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSEQSBAIJ)); 97c4762a1bSJed Brown for (i=0; i<3; i++) { 985f80ce2aSJacob Faibussowitsch if (i > 0) CHKERRQ(MatZeroEntries(A)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,view)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesAIJ(A)); 101c4762a1bSJed Brown } 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* 106c4762a1bSJed Brown Reload in MPISBAIJ matrix and check its values 107c4762a1bSJed Brown */ 1085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATMPISBAIJ)); 111c4762a1bSJed Brown for (i=0; i<3; i++) { 1125f80ce2aSJacob Faibussowitsch if (i > 0) CHKERRQ(MatZeroEntries(A)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,view)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesAIJ(A)); 115c4762a1bSJed Brown } 1165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 118c4762a1bSJed Brown 119*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 120*b122ec5aSJacob Faibussowitsch return 0; 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown /*TEST 124c4762a1bSJed Brown 125c4762a1bSJed Brown testset: 126c4762a1bSJed Brown args: -viewer_binary_mpiio 0 127c4762a1bSJed Brown output_file: output/ex50.out 128c4762a1bSJed Brown test: 129c4762a1bSJed Brown suffix: stdio_1 130c4762a1bSJed Brown nsize: 1 131c4762a1bSJed Brown test: 132c4762a1bSJed Brown suffix: stdio_2 133c4762a1bSJed Brown nsize: 2 134c4762a1bSJed Brown test: 135c4762a1bSJed Brown suffix: stdio_3 136c4762a1bSJed Brown nsize: 3 137c4762a1bSJed Brown test: 138c4762a1bSJed Brown suffix: stdio_4 139c4762a1bSJed Brown nsize: 4 140c4762a1bSJed Brown test: 141c4762a1bSJed Brown suffix: stdio_5 142c4762a1bSJed Brown nsize: 4 143c4762a1bSJed Brown 144c4762a1bSJed Brown testset: 145c4762a1bSJed Brown requires: mpiio 146c4762a1bSJed Brown args: -viewer_binary_mpiio 1 147c4762a1bSJed Brown output_file: output/ex50.out 148c4762a1bSJed Brown test: 149c4762a1bSJed Brown suffix: mpiio_1 150c4762a1bSJed Brown nsize: 1 151c4762a1bSJed Brown test: 152c4762a1bSJed Brown suffix: mpiio_2 153c4762a1bSJed Brown nsize: 2 154c4762a1bSJed Brown test: 155c4762a1bSJed Brown suffix: mpiio_3 156c4762a1bSJed Brown nsize: 3 157c4762a1bSJed Brown test: 158c4762a1bSJed Brown suffix: mpiio_4 159c4762a1bSJed Brown nsize: 4 160c4762a1bSJed Brown test: 161c4762a1bSJed Brown suffix: mpiio_5 162c4762a1bSJed Brown nsize: 5 163c4762a1bSJed Brown 164c4762a1bSJed Brown TEST*/ 165