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