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