1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPIBAIJ format.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown #define IMAX 15 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat A,B,C,At,Bt; 9c4762a1bSJed Brown PetscViewer fd; 10c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 11c4762a1bSJed Brown PetscRandom rand; 12c4762a1bSJed Brown Vec xx,yy,s1,s2; 13c4762a1bSJed Brown PetscReal s1norm,s2norm,rnorm,tol=1.e-10; 14c4762a1bSJed Brown PetscInt rstart,rend,rows[2],cols[2],m,n,i,j,M,N,ct,row,ncols1,ncols2,bs; 15c4762a1bSJed Brown PetscMPIInt rank,size; 16c4762a1bSJed Brown const PetscInt *cols1,*cols2; 17c4762a1bSJed Brown PetscScalar vals1[4],vals2[4],v; 18c4762a1bSJed Brown const PetscScalar *v1,*v2; 19c4762a1bSJed Brown PetscBool flg; 20c4762a1bSJed Brown 21*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 22*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 23*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* Check out if MatLoad() works */ 26*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg)); 2728b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Input file not specified"); 28*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 29*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 30*9566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATBAIJ)); 31*9566063dSJacob Faibussowitsch PetscCall(MatLoad(A,fd)); 32c4762a1bSJed Brown 33*9566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B)); 34*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 35c4762a1bSJed Brown 36*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 37*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 38*9566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,&n)); 39*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&xx)); 40*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(xx,m,PETSC_DECIDE)); 41*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(xx)); 42*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s1)); 43*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s2)); 44*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&yy)); 45*9566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A,&bs)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* Test MatNorm() */ 48*9566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_FROBENIUS,&s1norm)); 49*9566063dSJacob Faibussowitsch PetscCall(MatNorm(B,NORM_FROBENIUS,&s2norm)); 50c4762a1bSJed Brown rnorm = PetscAbsScalar(s2norm-s1norm)/s2norm; 51c4762a1bSJed Brown if (rnorm>tol) { 52*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatNorm_FROBENIUS()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n",rank,(double)s1norm,(double)s2norm,bs)); 53c4762a1bSJed Brown } 54*9566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_INFINITY,&s1norm)); 55*9566063dSJacob Faibussowitsch PetscCall(MatNorm(B,NORM_INFINITY,&s2norm)); 56c4762a1bSJed Brown rnorm = PetscAbsScalar(s2norm-s1norm)/s2norm; 57c4762a1bSJed Brown if (rnorm>tol) { 58*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatNorm_INFINITY()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n",rank,(double)s1norm,(double)s2norm,bs)); 59c4762a1bSJed Brown } 60*9566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_1,&s1norm)); 61*9566063dSJacob Faibussowitsch PetscCall(MatNorm(B,NORM_1,&s2norm)); 62c4762a1bSJed Brown rnorm = PetscAbsScalar(s2norm-s1norm)/s2norm; 63c4762a1bSJed Brown if (rnorm>tol) { 64*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatNorm_NORM_1()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n",rank,(double)s1norm,(double)s2norm,bs)); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* Test MatMult() */ 68c4762a1bSJed Brown for (i=0; i<IMAX; i++) { 69*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rand)); 70*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,xx,s1)); 71*9566063dSJacob Faibussowitsch PetscCall(MatMult(B,xx,s2)); 72*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2,-1.0,s1)); 73*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&rnorm)); 74c4762a1bSJed Brown if (rnorm>tol) { 75*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatMult - Norm2=%16.14e bs = %" PetscInt_FMT "\n",rank,(double)rnorm,bs)); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* test MatMultAdd() */ 80c4762a1bSJed Brown for (i=0; i<IMAX; i++) { 81*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rand)); 82*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(yy,rand)); 83*9566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A,xx,yy,s1)); 84*9566063dSJacob Faibussowitsch PetscCall(MatMultAdd(B,xx,yy,s2)); 85*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2,-1.0,s1)); 86*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&rnorm)); 87c4762a1bSJed Brown if (rnorm>tol) { 88*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatMultAdd - Norm2=%16.14e bs = %" PetscInt_FMT "\n",rank,(double)rnorm,bs)); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Test MatMultTranspose() */ 93c4762a1bSJed Brown for (i=0; i<IMAX; i++) { 94*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rand)); 95*9566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,xx,s1)); 96*9566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B,xx,s2)); 97*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 98*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 99c4762a1bSJed Brown rnorm = s2norm-s1norm; 100c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 101*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatMultTranspose - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",rank,(double)s1norm,(double)s2norm,bs)); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown } 104c4762a1bSJed Brown /* Test MatMultTransposeAdd() */ 105c4762a1bSJed Brown for (i=0; i<IMAX; i++) { 106*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rand)); 107*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(yy,rand)); 108*9566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A,xx,yy,s1)); 109*9566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(B,xx,yy,s2)); 110*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 111*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 112c4762a1bSJed Brown rnorm = s2norm-s1norm; 113c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 114*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatMultTransposeAdd - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",rank,(double)s1norm,(double)s2norm,bs)); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* Check MatGetValues() */ 119*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 120*9566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&M,&N)); 121c4762a1bSJed Brown 122c4762a1bSJed Brown for (i=0; i<IMAX; i++) { 123c4762a1bSJed Brown /* Create random row numbers ad col numbers */ 124*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&v)); 125c4762a1bSJed Brown cols[0] = (int)(PetscRealPart(v)*N); 126*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&v)); 127c4762a1bSJed Brown cols[1] = (int)(PetscRealPart(v)*N); 128*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&v)); 129c4762a1bSJed Brown rows[0] = rstart + (int)(PetscRealPart(v)*m); 130*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&v)); 131c4762a1bSJed Brown rows[1] = rstart + (int)(PetscRealPart(v)*m); 132c4762a1bSJed Brown 133*9566063dSJacob Faibussowitsch PetscCall(MatGetValues(A,2,rows,2,cols,vals1)); 134*9566063dSJacob Faibussowitsch PetscCall(MatGetValues(B,2,rows,2,cols,vals2)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown for (j=0; j<4; j++) { 137c4762a1bSJed Brown if (vals1[j] != vals2[j]) { 138*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d]: Error: MatGetValues rstart = %2" PetscInt_FMT " row = %2" PetscInt_FMT " col = %2" PetscInt_FMT " val1 = %e val2 = %e bs = %" PetscInt_FMT "\n",rank,rstart,rows[j/2],cols[j%2],(double)PetscRealPart(vals1[j]),(double)PetscRealPart(vals2[j]),bs)); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown } 141c4762a1bSJed Brown } 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* Test MatGetRow()/ MatRestoreRow() */ 144c4762a1bSJed Brown for (ct=0; ct<100; ct++) { 145*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&v)); 1468fffc762SJacob Faibussowitsch row = rstart + (PetscInt)(PetscRealPart(v)*m); 147*9566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,row,&ncols1,&cols1,&v1)); 148*9566063dSJacob Faibussowitsch PetscCall(MatGetRow(B,row,&ncols2,&cols2,&v2)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown for (i=0,j=0; i<ncols1 && j<ncols2; j++) { 151c4762a1bSJed Brown while (cols2[j] != cols1[i]) i++; 1522c71b3e2SJacob Faibussowitsch PetscCheckFalse(v1[i] != v2[j],PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - vals incorrect."); 153c4762a1bSJed Brown } 1542c71b3e2SJacob Faibussowitsch PetscCheckFalse(j<ncols2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - cols incorrect"); 155c4762a1bSJed Brown 156*9566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,row,&ncols1,&cols1,&v1)); 157*9566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(B,row,&ncols2,&cols2,&v2)); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* Test MatConvert() */ 161*9566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&C)); 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* See if MatMult Says both are same */ 164c4762a1bSJed Brown for (i=0; i<IMAX; i++) { 165*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rand)); 166*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,xx,s1)); 167*9566063dSJacob Faibussowitsch PetscCall(MatMult(C,xx,s2)); 168*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 169*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 170c4762a1bSJed Brown rnorm = s2norm-s1norm; 171c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 172*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] Error in MatConvert: MatMult - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",rank,(double)s1norm,(double)s2norm,bs)); 173c4762a1bSJed Brown } 174c4762a1bSJed Brown } 175*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* Test MatTranspose() */ 178*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&At)); 179*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(B,MAT_INITIAL_MATRIX,&Bt)); 180c4762a1bSJed Brown for (i=0; i<IMAX; i++) { 181*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rand)); 182*9566063dSJacob Faibussowitsch PetscCall(MatMult(At,xx,s1)); 183*9566063dSJacob Faibussowitsch PetscCall(MatMult(Bt,xx,s2)); 184*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 185*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 186c4762a1bSJed Brown rnorm = s2norm-s1norm; 187c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 188*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] Error in MatConvert:MatMult - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",rank,(double)s1norm,(double)s2norm,bs)); 189c4762a1bSJed Brown } 190c4762a1bSJed Brown } 191*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&At)); 192*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bt)); 193c4762a1bSJed Brown 194*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 195*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 196*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 197*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&yy)); 198*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 199*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 200*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 201*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 202b122ec5aSJacob Faibussowitsch return 0; 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown /*TEST 206c4762a1bSJed Brown 207c4762a1bSJed Brown build: 208c4762a1bSJed Brown requires: !complex 209c4762a1bSJed Brown 210c4762a1bSJed Brown test: 211dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 212c4762a1bSJed Brown nsize: 3 213c4762a1bSJed Brown args: -matload_block_size 1 -f ${DATAFILESPATH}/matrices/small 214c4762a1bSJed Brown 215c4762a1bSJed Brown test: 216c4762a1bSJed Brown suffix: 2 217dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 218c4762a1bSJed Brown nsize: 3 219c4762a1bSJed Brown args: -matload_block_size 2 -f ${DATAFILESPATH}/matrices/small 220c4762a1bSJed Brown output_file: output/ex53_1.out 221c4762a1bSJed Brown 222c4762a1bSJed Brown test: 223c4762a1bSJed Brown suffix: 3 224dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 225c4762a1bSJed Brown nsize: 3 226c4762a1bSJed Brown args: -matload_block_size 4 -f ${DATAFILESPATH}/matrices/small 227c4762a1bSJed Brown output_file: output/ex53_1.out 228c4762a1bSJed Brown 229c4762a1bSJed Brown test: 230c4762a1bSJed Brown TODO: Matrix row/column sizes are not compatible with block size 231c4762a1bSJed Brown suffix: 4 232dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 233c4762a1bSJed Brown nsize: 3 234c4762a1bSJed Brown args: -matload_block_size 5 -f ${DATAFILESPATH}/matrices/small 235c4762a1bSJed Brown output_file: output/ex53_1.out 236c4762a1bSJed Brown 237c4762a1bSJed Brown test: 238c4762a1bSJed Brown suffix: 5 239dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 240c4762a1bSJed Brown nsize: 3 241c4762a1bSJed Brown args: -matload_block_size 6 -f ${DATAFILESPATH}/matrices/small 242c4762a1bSJed Brown output_file: output/ex53_1.out 243c4762a1bSJed Brown 244c4762a1bSJed Brown test: 245c4762a1bSJed Brown TODO: Matrix row/column sizes are not compatible with block size 246c4762a1bSJed Brown suffix: 6 247dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 248c4762a1bSJed Brown nsize: 3 249c4762a1bSJed Brown args: -matload_block_size 7 -f ${DATAFILESPATH}/matrices/small 250c4762a1bSJed Brown output_file: output/ex53_1.out 251c4762a1bSJed Brown 252c4762a1bSJed Brown test: 253c4762a1bSJed Brown TODO: Matrix row/column sizes are not compatible with block size 254c4762a1bSJed Brown suffix: 7 255dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 256c4762a1bSJed Brown nsize: 3 257c4762a1bSJed Brown args: -matload_block_size 8 -f ${DATAFILESPATH}/matrices/small 258c4762a1bSJed Brown output_file: output/ex53_1.out 259c4762a1bSJed Brown 260c4762a1bSJed Brown test: 261c4762a1bSJed Brown suffix: 8 262dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 263c4762a1bSJed Brown nsize: 4 264c4762a1bSJed Brown args: -matload_block_size 3 -f ${DATAFILESPATH}/matrices/small 265c4762a1bSJed Brown output_file: output/ex53_1.out 266c4762a1bSJed Brown 267c4762a1bSJed Brown TEST*/ 268