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; 209566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&M,&N)); 219566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 22c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 23c4762a1bSJed Brown for (j=0; j<N; j++) { 249566063dSJacob Faibussowitsch PetscCall(MatGetValue(A,i,j,&val)); 25c4762a1bSJed Brown v = MakeValue(i,j,M); w = PetscRealPart(val); 2608401ef6SPierre Jolivet PetscCheck(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 PetscViewer view; 38c4762a1bSJed Brown 399566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,NULL,help)); 40c4762a1bSJed Brown /* 41c4762a1bSJed Brown Create a parallel BAIJ matrix shared by all processors 42c4762a1bSJed Brown */ 43*d0609cedSBarry Smith PetscCall(MatCreateBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,M,N,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&A)); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* 46c4762a1bSJed Brown Set values into the matrix 47c4762a1bSJed Brown */ 489566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&M,&N)); 499566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 50c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 51c4762a1bSJed Brown for (j=0; j<N; j++) { 52c4762a1bSJed Brown PetscReal v = MakeValue(i,j,M); 53c4762a1bSJed Brown if (PetscAbsReal(v) > 0) { 549566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,j,v,INSERT_VALUES)); 55c4762a1bSJed Brown } 56c4762a1bSJed Brown } 57c4762a1bSJed Brown } 589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 609566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-mat_base_view")); 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* 63c4762a1bSJed Brown Store the binary matrix to a file 64c4762a1bSJed Brown */ 659566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view)); 66c4762a1bSJed Brown for (i=0; i<3; i++) { 679566063dSJacob Faibussowitsch PetscCall(MatView(A,view)); 68c4762a1bSJed Brown } 699566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* 73c4762a1bSJed Brown Now reload the matrix and check its values 74c4762a1bSJed Brown */ 759566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view)); 769566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 779566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATBAIJ)); 78c4762a1bSJed Brown for (i=0; i<3; i++) { 799566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 809566063dSJacob Faibussowitsch PetscCall(MatLoad(A,view)); 819566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A)); 82c4762a1bSJed Brown } 839566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 849566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-mat_load_view")); 859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* 88c4762a1bSJed Brown Reload in SEQBAIJ matrix and check its values 89c4762a1bSJed Brown */ 909566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF,"matrix.dat",FILE_MODE_READ,&view)); 919566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&A)); 929566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSEQBAIJ)); 93c4762a1bSJed Brown for (i=0; i<3; i++) { 949566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 959566063dSJacob Faibussowitsch PetscCall(MatLoad(A,view)); 969566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A)); 97c4762a1bSJed Brown } 989566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* 102c4762a1bSJed Brown Reload in MPIBAIJ matrix and check its values 103c4762a1bSJed Brown */ 1049566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view)); 1059566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 1069566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATMPIBAIJ)); 107c4762a1bSJed Brown for (i=0; i<3; i++) { 1089566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 1099566063dSJacob Faibussowitsch PetscCall(MatLoad(A,view)); 1109566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A)); 111c4762a1bSJed Brown } 1129566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 114c4762a1bSJed Brown 1159566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 116b122ec5aSJacob Faibussowitsch return 0; 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown /*TEST 120c4762a1bSJed Brown 121c4762a1bSJed Brown testset: 122c4762a1bSJed Brown args: -viewer_binary_mpiio 0 123c4762a1bSJed Brown output_file: output/ex45.out 124c4762a1bSJed Brown test: 125c4762a1bSJed Brown suffix: stdio_1 126c4762a1bSJed Brown nsize: 1 127c4762a1bSJed Brown test: 128c4762a1bSJed Brown suffix: stdio_2 129c4762a1bSJed Brown nsize: 2 130c4762a1bSJed Brown test: 131c4762a1bSJed Brown suffix: stdio_3 132c4762a1bSJed Brown nsize: 3 133c4762a1bSJed Brown test: 134c4762a1bSJed Brown suffix: stdio_4 135c4762a1bSJed Brown nsize: 4 136c4762a1bSJed Brown test: 137c4762a1bSJed Brown suffix: stdio_5 138c4762a1bSJed Brown nsize: 4 139c4762a1bSJed Brown 140c4762a1bSJed Brown testset: 141c4762a1bSJed Brown requires: mpiio 142c4762a1bSJed Brown args: -viewer_binary_mpiio 1 143c4762a1bSJed Brown output_file: output/ex45.out 144c4762a1bSJed Brown test: 145c4762a1bSJed Brown suffix: mpiio_1 146c4762a1bSJed Brown nsize: 1 147c4762a1bSJed Brown test: 148c4762a1bSJed Brown suffix: mpiio_2 149c4762a1bSJed Brown nsize: 2 150c4762a1bSJed Brown test: 151c4762a1bSJed Brown suffix: mpiio_3 152c4762a1bSJed Brown nsize: 3 153c4762a1bSJed Brown test: 154c4762a1bSJed Brown suffix: mpiio_4 155c4762a1bSJed Brown nsize: 4 156c4762a1bSJed Brown test: 157c4762a1bSJed Brown suffix: mpiio_5 158c4762a1bSJed Brown nsize: 5 159c4762a1bSJed Brown 160c4762a1bSJed Brown TEST*/ 161