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