1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests the various sequential routines in MATSEQSBAIJ format.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown PetscMPIInt size; 9c4762a1bSJed Brown PetscErrorCode ierr; 10c4762a1bSJed Brown Vec x,y,b,s1,s2; 11c4762a1bSJed Brown Mat A; /* linear system matrix */ 12c4762a1bSJed Brown Mat sA,sB,sFactor,B,C; /* symmetric matrices */ 13c4762a1bSJed Brown PetscInt n,mbs=16,bs=1,nz=3,prob=1,i,j,k1,k2,col[3],lf,block, row,Ii,J,n1,inc; 14c4762a1bSJed Brown PetscReal norm1,norm2,rnorm,tol=10*PETSC_SMALL; 15c4762a1bSJed Brown PetscScalar neg_one=-1.0,four=4.0,value[3]; 16c4762a1bSJed Brown IS perm, iscol; 17c4762a1bSJed Brown PetscRandom rdm; 18c4762a1bSJed Brown PetscBool doIcc=PETSC_TRUE,equal; 19c4762a1bSJed Brown MatInfo minfo1,minfo2; 20c4762a1bSJed Brown MatFactorInfo factinfo; 21c4762a1bSJed Brown MatType type; 22c4762a1bSJed Brown 23c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 24ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 25c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 26c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr); 27c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL);CHKERRQ(ierr); 28c4762a1bSJed Brown 29c4762a1bSJed Brown n = mbs*bs; 30c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&A);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = MatSetType(A,MATSEQBAIJ);CHKERRQ(ierr); 33c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = MatSeqBAIJSetPreallocation(A,bs,nz,NULL);CHKERRQ(ierr); 35c4762a1bSJed Brown 36c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&sA);CHKERRQ(ierr); 37c4762a1bSJed Brown ierr = MatSetSizes(sA,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = MatSetType(sA,MATSEQSBAIJ);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = MatSetFromOptions(sA);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = MatGetType(sA,&type);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)sA,MATSEQSBAIJ,&doIcc);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = MatSeqSBAIJSetPreallocation(sA,bs,nz,NULL);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = MatSetOption(sA,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 46c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&Ii,&J);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = MatGetOwnershipRange(sA,&i,&j);CHKERRQ(ierr); 48c4762a1bSJed Brown if (i-Ii || j-J) { 49c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ format\n");CHKERRQ(ierr); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Assemble matrix */ 53c4762a1bSJed Brown if (bs == 1) { 54c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL);CHKERRQ(ierr); 55c4762a1bSJed Brown if (prob == 1) { /* tridiagonal matrix */ 56c4762a1bSJed Brown value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; 57c4762a1bSJed Brown for (i=1; i<n-1; i++) { 58c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 59c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1; 63c4762a1bSJed Brown 64c4762a1bSJed Brown value[0]= 0.1; value[1]=-1; value[2]=2; 65c4762a1bSJed Brown 66c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 68c4762a1bSJed Brown 69c4762a1bSJed Brown i = 0; 70c4762a1bSJed Brown col[0] = n-1; col[1] = 1; col[2] = 0; 71c4762a1bSJed Brown value[0] = 0.1; value[1] = -1.0; value[2] = 2; 72c4762a1bSJed Brown 73c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 75c4762a1bSJed Brown 76c4762a1bSJed Brown } else if (prob == 2) { /* matrix for the five point stencil */ 77c4762a1bSJed Brown n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001); 78c4762a1bSJed Brown if (n1*n1 - n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive integer!"); 79c4762a1bSJed Brown for (i=0; i<n1; i++) { 80c4762a1bSJed Brown for (j=0; j<n1; j++) { 81c4762a1bSJed Brown Ii = j + n1*i; 82c4762a1bSJed Brown if (i>0) { 83c4762a1bSJed Brown J = Ii - n1; 84c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown if (i<n1-1) { 88c4762a1bSJed Brown J = Ii + n1; 89c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 91c4762a1bSJed Brown } 92c4762a1bSJed Brown if (j>0) { 93c4762a1bSJed Brown J = Ii - 1; 94c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown if (j<n1-1) { 98c4762a1bSJed Brown J = Ii + 1; 99c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES);CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = MatSetValues(sA,1,&Ii,1,&Ii,&four,INSERT_VALUES);CHKERRQ(ierr); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown } 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown } else { /* bs > 1 */ 109c4762a1bSJed Brown for (block=0; block<n/bs; block++) { 110c4762a1bSJed Brown /* diagonal blocks */ 111c4762a1bSJed Brown value[0] = -1.0; value[1] = 4.0; value[2] = -1.0; 112c4762a1bSJed Brown for (i=1+block*bs; i<bs-1+block*bs; i++) { 113c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 114c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 116c4762a1bSJed Brown } 117c4762a1bSJed Brown i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs; 118c4762a1bSJed Brown 119c4762a1bSJed Brown value[0]=-1.0; value[1]=4.0; 120c4762a1bSJed Brown 121c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 123c4762a1bSJed Brown 124c4762a1bSJed Brown i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs; 125c4762a1bSJed Brown 126c4762a1bSJed Brown value[0]=4.0; value[1] = -1.0; 127c4762a1bSJed Brown 128c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 129c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown /* off-diagonal blocks */ 132c4762a1bSJed Brown value[0]=-1.0; 133c4762a1bSJed Brown for (i=0; i<(n/bs-1)*bs; i++) { 134c4762a1bSJed Brown col[0]=i+bs; 135c4762a1bSJed Brown 136c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 138c4762a1bSJed Brown 139c4762a1bSJed Brown col[0]=i; row=i+bs; 140c4762a1bSJed Brown 141c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);CHKERRQ(ierr); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown } 145c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 147c4762a1bSJed Brown 148c4762a1bSJed Brown ierr = MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* Test MatGetInfo() of A and sA */ 152c4762a1bSJed Brown ierr = MatGetInfo(A,MAT_LOCAL,&minfo1);CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = MatGetInfo(sA,MAT_LOCAL,&minfo2);CHKERRQ(ierr); 154c4762a1bSJed Brown i = (int) (minfo1.nz_used - minfo2.nz_used); 155c4762a1bSJed Brown j = (int) (minfo1.nz_allocated - minfo2.nz_allocated); 156c4762a1bSJed Brown k1 = (int) (minfo1.nz_allocated - minfo1.nz_used); 157c4762a1bSJed Brown k2 = (int) (minfo2.nz_allocated - minfo2.nz_used); 158c4762a1bSJed Brown if (i < 0 || j < 0 || k1 < 0 || k2 < 0) { 159c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error (compare A and sA): MatGetInfo()\n");CHKERRQ(ierr); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* Test MatDuplicate() */ 163c4762a1bSJed Brown ierr = MatNorm(A,NORM_FROBENIUS,&norm1);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = MatDuplicate(sA,MAT_COPY_VALUES,&sB);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = MatEqual(sA,sB,&equal);CHKERRQ(ierr); 166c4762a1bSJed Brown if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDuplicate()"); 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* Test MatNorm() */ 169c4762a1bSJed Brown ierr = MatNorm(A,NORM_FROBENIUS,&norm1);CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = MatNorm(sB,NORM_FROBENIUS,&norm2);CHKERRQ(ierr); 171c4762a1bSJed Brown rnorm = PetscAbsReal(norm1-norm2)/norm2; 172c4762a1bSJed Brown if (rnorm > tol) { 173*8fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS, NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2);CHKERRQ(ierr); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown ierr = MatNorm(A,NORM_INFINITY,&norm1);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = MatNorm(sB,NORM_INFINITY,&norm2);CHKERRQ(ierr); 177c4762a1bSJed Brown rnorm = PetscAbsReal(norm1-norm2)/norm2; 178c4762a1bSJed Brown if (rnorm > tol) { 179*8fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2);CHKERRQ(ierr); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown ierr = MatNorm(A,NORM_1,&norm1);CHKERRQ(ierr); 182c4762a1bSJed Brown ierr = MatNorm(sB,NORM_1,&norm2);CHKERRQ(ierr); 183c4762a1bSJed Brown rnorm = PetscAbsReal(norm1-norm2)/norm2; 184c4762a1bSJed Brown if (rnorm > tol) { 185*8fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2);CHKERRQ(ierr); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */ 189c4762a1bSJed Brown ierr = MatGetInfo(A,MAT_LOCAL,&minfo1);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = MatGetInfo(sB,MAT_LOCAL,&minfo2);CHKERRQ(ierr); 191c4762a1bSJed Brown i = (int) (minfo1.nz_used - minfo2.nz_used); 192c4762a1bSJed Brown j = (int) (minfo1.nz_allocated - minfo2.nz_allocated); 193c4762a1bSJed Brown k1 = (int) (minfo1.nz_allocated - minfo1.nz_used); 194c4762a1bSJed Brown k2 = (int) (minfo2.nz_allocated - minfo2.nz_used); 195c4762a1bSJed Brown if (i < 0 || j < 0 || k1 < 0 || k2 < 0) { 196c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error(compare A and sB): MatGetInfo()\n");CHKERRQ(ierr); 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199c4762a1bSJed Brown ierr = MatGetSize(A,&Ii,&J);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = MatGetSize(sB,&i,&j);CHKERRQ(ierr); 201c4762a1bSJed Brown if (i-Ii || j-J) { 2022f613bf5SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()\n");CHKERRQ(ierr); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown ierr = MatGetBlockSize(A, &Ii);CHKERRQ(ierr); 206c4762a1bSJed Brown ierr = MatGetBlockSize(sB, &i);CHKERRQ(ierr); 207c4762a1bSJed Brown if (i-Ii) { 208c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()\n");CHKERRQ(ierr); 209c4762a1bSJed Brown } 210c4762a1bSJed Brown 211c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF,&rdm);CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 213c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,n,&x);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = VecDuplicate(x,&s1);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = VecDuplicate(x,&s2);CHKERRQ(ierr); 216c4762a1bSJed Brown ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = VecDuplicate(x,&b);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */ 221c4762a1bSJed Brown #if !defined(PETSC_USE_COMPLEX) 222c4762a1bSJed Brown /* Scaling matrix with complex numbers results non-spd matrix, 223c4762a1bSJed Brown causing crash of MatForwardSolve() and MatBackwardSolve() */ 224c4762a1bSJed Brown ierr = MatDiagonalScale(A,x,x);CHKERRQ(ierr); 225c4762a1bSJed Brown ierr = MatDiagonalScale(sB,x,x);CHKERRQ(ierr); 226c4762a1bSJed Brown ierr = MatMultEqual(A,sB,10,&equal);CHKERRQ(ierr); 227c4762a1bSJed Brown if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDiagonalScale"); 228c4762a1bSJed Brown 229c4762a1bSJed Brown ierr = MatGetDiagonal(A,s1);CHKERRQ(ierr); 230c4762a1bSJed Brown ierr = MatGetDiagonal(sB,s2);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = VecAXPY(s2,neg_one,s1);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = VecNorm(s2,NORM_1,&norm1);CHKERRQ(ierr); 233c4762a1bSJed Brown if (norm1>tol) { 234c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal(), ||s1-s2||=%g\n",(double)norm1);CHKERRQ(ierr); 235c4762a1bSJed Brown } 236c4762a1bSJed Brown 237c4762a1bSJed Brown { 238c4762a1bSJed Brown PetscScalar alpha=0.1; 239c4762a1bSJed Brown ierr = MatScale(A,alpha);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = MatScale(sB,alpha);CHKERRQ(ierr); 241c4762a1bSJed Brown } 242c4762a1bSJed Brown #endif 243c4762a1bSJed Brown 244c4762a1bSJed Brown /* Test MatGetRowMaxAbs() */ 245c4762a1bSJed Brown ierr = MatGetRowMaxAbs(A,s1,NULL);CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = MatGetRowMaxAbs(sB,s2,NULL);CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = VecNorm(s1,NORM_1,&norm1);CHKERRQ(ierr); 248c4762a1bSJed Brown ierr = VecNorm(s2,NORM_1,&norm2);CHKERRQ(ierr); 249c4762a1bSJed Brown norm1 -= norm2; 250c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 251c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error:MatGetRowMaxAbs() \n");CHKERRQ(ierr); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* Test MatMult() */ 255c4762a1bSJed Brown for (i=0; i<40; i++) { 256c4762a1bSJed Brown ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = MatMult(A,x,s1);CHKERRQ(ierr); 258c4762a1bSJed Brown ierr = MatMult(sB,x,s2);CHKERRQ(ierr); 259c4762a1bSJed Brown ierr = VecNorm(s1,NORM_1,&norm1);CHKERRQ(ierr); 260c4762a1bSJed Brown ierr = VecNorm(s2,NORM_1,&norm2);CHKERRQ(ierr); 261c4762a1bSJed Brown norm1 -= norm2; 262c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 263c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatMult(), norm1-norm2: %g\n",(double)norm1);CHKERRQ(ierr); 264c4762a1bSJed Brown } 265c4762a1bSJed Brown } 266c4762a1bSJed Brown 267c4762a1bSJed Brown /* MatMultAdd() */ 268c4762a1bSJed Brown for (i=0; i<40; i++) { 269c4762a1bSJed Brown ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = VecSetRandom(y,rdm);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr); 272c4762a1bSJed Brown ierr = MatMultAdd(sB,x,y,s2);CHKERRQ(ierr); 273c4762a1bSJed Brown ierr = VecNorm(s1,NORM_1,&norm1);CHKERRQ(ierr); 274c4762a1bSJed Brown ierr = VecNorm(s2,NORM_1,&norm2);CHKERRQ(ierr); 275c4762a1bSJed Brown norm1 -= norm2; 276c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 277c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd(), norm1-norm2: %g\n",(double)norm1);CHKERRQ(ierr); 278c4762a1bSJed Brown } 279c4762a1bSJed Brown } 280c4762a1bSJed Brown 281c20d7725SJed Brown /* Test MatMatMult() for sbaij and dense matrices */ 282c4762a1bSJed Brown ierr = MatCreateSeqDense(PETSC_COMM_SELF,n,5*n,NULL,&B);CHKERRQ(ierr); 283c4762a1bSJed Brown ierr = MatSetRandom(B,rdm);CHKERRQ(ierr); 284c4762a1bSJed Brown ierr = MatMatMult(sA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr); 285c4762a1bSJed Brown ierr = MatMatMultEqual(sA,B,C,5*n,&equal);CHKERRQ(ierr); 286c4762a1bSJed Brown if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult()"); 287c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 288c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 289c4762a1bSJed Brown 290c4762a1bSJed Brown /* Test MatCholeskyFactor(), MatICCFactor() with natural ordering */ 291c4762a1bSJed Brown ierr = MatGetOrdering(A,MATORDERINGNATURAL,&perm,&iscol);CHKERRQ(ierr); 292c4762a1bSJed Brown ierr = ISDestroy(&iscol);CHKERRQ(ierr); 293c4762a1bSJed Brown norm1 = tol; 294c4762a1bSJed Brown inc = bs; 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* initialize factinfo */ 297c4762a1bSJed Brown ierr = PetscMemzero(&factinfo,sizeof(MatFactorInfo));CHKERRQ(ierr); 298c4762a1bSJed Brown 299c4762a1bSJed Brown for (lf=-1; lf<10; lf += inc) { 300c4762a1bSJed Brown if (lf==-1) { /* Cholesky factor of sB (duplicate sA) */ 301c4762a1bSJed Brown factinfo.fill = 5.0; 302c4762a1bSJed Brown 303c4762a1bSJed Brown ierr = MatGetFactor(sB,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sFactor);CHKERRQ(ierr); 304c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(sFactor,sB,perm,&factinfo);CHKERRQ(ierr); 305c4762a1bSJed Brown } else if (!doIcc) break; 306c4762a1bSJed Brown else { /* incomplete Cholesky factor */ 307c4762a1bSJed Brown factinfo.fill = 5.0; 308c4762a1bSJed Brown factinfo.levels = lf; 309c4762a1bSJed Brown 310c4762a1bSJed Brown ierr = MatGetFactor(sB,MATSOLVERPETSC,MAT_FACTOR_ICC,&sFactor);CHKERRQ(ierr); 311c4762a1bSJed Brown ierr = MatICCFactorSymbolic(sFactor,sB,perm,&factinfo);CHKERRQ(ierr); 312c4762a1bSJed Brown } 313c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(sFactor,sB,&factinfo);CHKERRQ(ierr); 314c4762a1bSJed Brown /* MatView(sFactor, PETSC_VIEWER_DRAW_WORLD); */ 315c4762a1bSJed Brown 316c4762a1bSJed Brown /* test MatGetDiagonal on numeric factor */ 317c4762a1bSJed Brown /* 318c4762a1bSJed Brown if (lf == -1) { 319c4762a1bSJed Brown ierr = MatGetDiagonal(sFactor,s1);CHKERRQ(ierr); 320c4762a1bSJed Brown printf(" in ex74.c, diag: \n"); 321c4762a1bSJed Brown ierr = VecView(s1,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 322c4762a1bSJed Brown } 323c4762a1bSJed Brown */ 324c4762a1bSJed Brown 325c4762a1bSJed Brown ierr = MatMult(sB,x,b);CHKERRQ(ierr); 326c4762a1bSJed Brown 327c4762a1bSJed Brown /* test MatForwardSolve() and MatBackwardSolve() */ 328c4762a1bSJed Brown if (lf == -1) { 329c4762a1bSJed Brown ierr = MatForwardSolve(sFactor,b,s1);CHKERRQ(ierr); 330c4762a1bSJed Brown ierr = MatBackwardSolve(sFactor,s1,s2);CHKERRQ(ierr); 331c4762a1bSJed Brown ierr = VecAXPY(s2,neg_one,x);CHKERRQ(ierr); 332c4762a1bSJed Brown ierr = VecNorm(s2,NORM_2,&norm2);CHKERRQ(ierr); 333c4762a1bSJed Brown if (10*norm1 < norm2) { 334*8fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g, bs=%" PetscInt_FMT "\n",(double)norm2,bs);CHKERRQ(ierr); 335c4762a1bSJed Brown } 336c4762a1bSJed Brown } 337c4762a1bSJed Brown 338c4762a1bSJed Brown /* test MatSolve() */ 339c4762a1bSJed Brown ierr = MatSolve(sFactor,b,y);CHKERRQ(ierr); 340c4762a1bSJed Brown ierr = MatDestroy(&sFactor);CHKERRQ(ierr); 341c4762a1bSJed Brown /* Check the error */ 342c4762a1bSJed Brown ierr = VecAXPY(y,neg_one,x);CHKERRQ(ierr); 343c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); 344c4762a1bSJed Brown if (10*norm1 < norm2 && lf-inc != -1) { 345*8fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF,"lf=%" PetscInt_FMT ", %" PetscInt_FMT ", Norm of error=%g, %g\n",lf-inc,lf,(double)norm1,(double)norm2);CHKERRQ(ierr); 346c4762a1bSJed Brown } 347c4762a1bSJed Brown norm1 = norm2; 348c4762a1bSJed Brown if (norm2 < tol && lf != -1) break; 349c4762a1bSJed Brown } 350c4762a1bSJed Brown 351c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 352c4762a1bSJed Brown ierr = MatGetFactor(sA,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&sFactor);CHKERRQ(ierr); 353c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(sFactor,sA,NULL,NULL);CHKERRQ(ierr); 354c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(sFactor,sA,NULL);CHKERRQ(ierr); 355c4762a1bSJed Brown for (i=0; i<10; i++) { 356c4762a1bSJed Brown ierr = VecSetRandom(b,rdm);CHKERRQ(ierr); 357c4762a1bSJed Brown ierr = MatSolve(sFactor,b,y);CHKERRQ(ierr); 358c4762a1bSJed Brown /* Check the error */ 359c4762a1bSJed Brown ierr = MatMult(sA,y,x);CHKERRQ(ierr); 360c4762a1bSJed Brown ierr = VecAXPY(x,neg_one,b);CHKERRQ(ierr); 361c4762a1bSJed Brown ierr = VecNorm(x,NORM_2,&norm2);CHKERRQ(ierr); 362c4762a1bSJed Brown if (norm2>tol) { 363c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error:MatSolve(), norm2: %g\n",(double)norm2);CHKERRQ(ierr); 364c4762a1bSJed Brown } 365c4762a1bSJed Brown } 366c4762a1bSJed Brown ierr = MatDestroy(&sFactor);CHKERRQ(ierr); 367c4762a1bSJed Brown #endif 368c4762a1bSJed Brown 369c4762a1bSJed Brown ierr = ISDestroy(&perm);CHKERRQ(ierr); 370c4762a1bSJed Brown 371c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 372c4762a1bSJed Brown ierr = MatDestroy(&sB);CHKERRQ(ierr); 373c4762a1bSJed Brown ierr = MatDestroy(&sA);CHKERRQ(ierr); 374c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 375c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 376c4762a1bSJed Brown ierr = VecDestroy(&s1);CHKERRQ(ierr); 377c4762a1bSJed Brown ierr = VecDestroy(&s2);CHKERRQ(ierr); 378c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 379c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 380c4762a1bSJed Brown 381c4762a1bSJed Brown ierr = PetscFinalize(); 382c4762a1bSJed Brown return ierr; 383c4762a1bSJed Brown } 384c4762a1bSJed Brown 385c4762a1bSJed Brown /*TEST 386c4762a1bSJed Brown 387c4762a1bSJed Brown test: 388c4762a1bSJed Brown args: -bs {{1 2 3 4 5 6 7 8}} 389c4762a1bSJed Brown 390c4762a1bSJed Brown TEST*/ 391