1 2 static char help[] = "Tests various routines in MatKAIJ format.\n"; 3 4 #include <petscmat.h> 5 #define IMAX 15 6 7 int main(int argc,char **args) 8 { 9 Mat A,B,TA; 10 PetscScalar *S,*T; 11 PetscViewer fd; 12 char file[PETSC_MAX_PATH_LEN]; 13 PetscInt m,n,M,N,p=1,q=1,i,j; 14 PetscMPIInt rank,size; 15 PetscErrorCode ierr; 16 PetscBool flg; 17 18 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 19 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 20 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 21 22 /* Load AIJ matrix A */ 23 CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg)); 24 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 25 CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 26 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 27 CHKERRQ(MatLoad(A,fd)); 28 CHKERRQ(PetscViewerDestroy(&fd)); 29 30 /* Get dof, then create S and T */ 31 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL)); 32 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-q",&q,NULL)); 33 CHKERRQ(PetscMalloc2(p*q,&S,p*q,&T)); 34 for (i=0; i<p*q; i++) S[i] = 0; 35 36 for (i=0; i<p; i++) { 37 for (j=0; j<q; j++) { 38 /* Set some random non-zero values */ 39 S[i+p*j] = ((PetscReal) ((i+1)*(j+1))) / ((PetscReal) (p+q)); 40 T[i+p*j] = ((PetscReal) ((p-i)+j)) / ((PetscReal) (p*q)); 41 } 42 } 43 44 /* Test KAIJ when both S & T are not NULL */ 45 46 /* Create KAIJ matrix TA */ 47 CHKERRQ(MatCreateKAIJ(A,p,q,S,T,&TA)); 48 CHKERRQ(MatGetLocalSize(A,&m,&n)); 49 CHKERRQ(MatGetSize(A,&M,&N)); 50 51 CHKERRQ(MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B)); 52 53 /* Test MatKAIJGetScaledIdentity() */ 54 CHKERRQ(MatKAIJGetScaledIdentity(TA,&flg)); 55 PetscCheck(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 1: MatKAIJGetScaledIdentity()"); 56 /* Test MatMult() */ 57 CHKERRQ(MatMultEqual(TA,B,10,&flg)); 58 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 1: MatMult() for KAIJ matrix"); 59 /* Test MatMultAdd() */ 60 CHKERRQ(MatMultAddEqual(TA,B,10,&flg)); 61 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 1: MatMultAdd() for KAIJ matrix"); 62 63 CHKERRQ(MatDestroy(&TA)); 64 CHKERRQ(MatDestroy(&B)); 65 66 /* Test KAIJ when S is NULL */ 67 68 /* Create KAIJ matrix TA */ 69 CHKERRQ(MatCreateKAIJ(A,p,q,NULL,T,&TA)); 70 CHKERRQ(MatGetLocalSize(A,&m,&n)); 71 CHKERRQ(MatGetSize(A,&M,&N)); 72 73 CHKERRQ(MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B)); 74 75 /* Test MatKAIJGetScaledIdentity() */ 76 CHKERRQ(MatKAIJGetScaledIdentity(TA,&flg)); 77 PetscCheck(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 2: MatKAIJGetScaledIdentity()"); 78 /* Test MatMult() */ 79 CHKERRQ(MatMultEqual(TA,B,10,&flg)); 80 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 2: MatMult() for KAIJ matrix"); 81 /* Test MatMultAdd() */ 82 CHKERRQ(MatMultAddEqual(TA,B,10,&flg)); 83 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 2: MatMultAdd() for KAIJ matrix"); 84 85 CHKERRQ(MatDestroy(&TA)); 86 CHKERRQ(MatDestroy(&B)); 87 88 /* Test KAIJ when T is NULL */ 89 90 /* Create KAIJ matrix TA */ 91 CHKERRQ(MatCreateKAIJ(A,p,q,S,NULL,&TA)); 92 CHKERRQ(MatGetLocalSize(A,&m,&n)); 93 CHKERRQ(MatGetSize(A,&M,&N)); 94 95 CHKERRQ(MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B)); 96 97 /* Test MatKAIJGetScaledIdentity() */ 98 CHKERRQ(MatKAIJGetScaledIdentity(TA,&flg)); 99 PetscCheck(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 3: MatKAIJGetScaledIdentity()"); 100 /* Test MatMult() */ 101 CHKERRQ(MatMultEqual(TA,B,10,&flg)); 102 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 3: MatMult() for KAIJ matrix"); 103 /* Test MatMultAdd() */ 104 CHKERRQ(MatMultAddEqual(TA,B,10,&flg)); 105 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 3: MatMultAdd() for KAIJ matrix"); 106 107 CHKERRQ(MatDestroy(&TA)); 108 CHKERRQ(MatDestroy(&B)); 109 110 /* Test KAIJ when T is is an identity matrix */ 111 112 if (p == q) { 113 for (i=0; i<p; i++) { 114 for (j=0; j<q; j++) { 115 if (i==j) T[i+j*p] = 1.0; 116 else T[i+j*p] = 0.0; 117 } 118 } 119 120 /* Create KAIJ matrix TA */ 121 CHKERRQ(MatCreateKAIJ(A,p,q,S,T,&TA)); 122 CHKERRQ(MatGetLocalSize(A,&m,&n)); 123 CHKERRQ(MatGetSize(A,&M,&N)); 124 125 CHKERRQ(MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B)); 126 127 /* Test MatKAIJGetScaledIdentity() */ 128 CHKERRQ(MatKAIJGetScaledIdentity(TA,&flg)); 129 PetscCheck(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 4: MatKAIJGetScaledIdentity()"); 130 /* Test MatMult() */ 131 CHKERRQ(MatMultEqual(TA,B,10,&flg)); 132 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 4: MatMult() for KAIJ matrix"); 133 /* Test MatMultAdd() */ 134 CHKERRQ(MatMultAddEqual(TA,B,10,&flg)); 135 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 4: MatMultAdd() for KAIJ matrix"); 136 137 CHKERRQ(MatDestroy(&TA)); 138 CHKERRQ(MatDestroy(&B)); 139 140 CHKERRQ(MatCreateKAIJ(A,p,q,NULL,T,&TA)); 141 /* Test MatKAIJGetScaledIdentity() */ 142 CHKERRQ(MatKAIJGetScaledIdentity(TA,&flg)); 143 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 5: MatKAIJGetScaledIdentity()"); 144 CHKERRQ(MatDestroy(&TA)); 145 146 for (i=0; i<p; i++) { 147 for (j=0; j<q; j++) { 148 if (i==j) S[i+j*p] = T[i+j*p] = 2.0; 149 else S[i+j*p] = T[i+j*p] = 0.0; 150 } 151 } 152 CHKERRQ(MatCreateKAIJ(A,p,q,S,T,&TA)); 153 /* Test MatKAIJGetScaledIdentity() */ 154 CHKERRQ(MatKAIJGetScaledIdentity(TA,&flg)); 155 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 6: MatKAIJGetScaledIdentity()"); 156 CHKERRQ(MatDestroy(&TA)); 157 } 158 159 /* Done with all tests */ 160 161 CHKERRQ(PetscFree2(S,T)); 162 CHKERRQ(MatDestroy(&A)); 163 ierr = PetscFinalize(); 164 return ierr; 165 } 166 167 /*TEST 168 169 build: 170 requires: !complex 171 172 test: 173 requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 174 output_file: output/ex176.out 175 nsize: {{1 2 3 4}} 176 args: -f ${DATAFILESPATH}/matrices/small -p {{2 3 7}} -q {{3 7}} -viewer_binary_skip_info 177 178 TEST*/ 179