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 Vec x,y,b,s1,s2; 10c4762a1bSJed Brown Mat A; /* linear system matrix */ 11c4762a1bSJed Brown Mat sA,sB,sFactor,B,C; /* symmetric matrices */ 12c4762a1bSJed Brown PetscInt n,mbs=16,bs=1,nz=3,prob=1,i,j,k1,k2,col[3],lf,block, row,Ii,J,n1,inc; 13c4762a1bSJed Brown PetscReal norm1,norm2,rnorm,tol=10*PETSC_SMALL; 14c4762a1bSJed Brown PetscScalar neg_one=-1.0,four=4.0,value[3]; 15c4762a1bSJed Brown IS perm, iscol; 16c4762a1bSJed Brown PetscRandom rdm; 17c4762a1bSJed Brown PetscBool doIcc=PETSC_TRUE,equal; 18c4762a1bSJed Brown MatInfo minfo1,minfo2; 19c4762a1bSJed Brown MatFactorInfo factinfo; 20c4762a1bSJed Brown MatType type; 21c4762a1bSJed Brown 22*327415f7SBarry Smith PetscFunctionBeginUser; 239566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 25be096a46SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL)); 28c4762a1bSJed Brown 29c4762a1bSJed Brown n = mbs*bs; 309566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&A)); 319566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 329566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSEQBAIJ)); 339566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 349566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(A,bs,nz,NULL)); 35c4762a1bSJed Brown 369566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&sA)); 379566063dSJacob Faibussowitsch PetscCall(MatSetSizes(sA,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 389566063dSJacob Faibussowitsch PetscCall(MatSetType(sA,MATSEQSBAIJ)); 399566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(sA)); 409566063dSJacob Faibussowitsch PetscCall(MatGetType(sA,&type)); 419566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)sA,MATSEQSBAIJ,&doIcc)); 429566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(sA,bs,nz,NULL)); 439566063dSJacob Faibussowitsch PetscCall(MatSetOption(sA,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 469566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Ii,&J)); 479566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(sA,&i,&j)); 48c4762a1bSJed Brown if (i-Ii || j-J) { 499566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ format\n")); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Assemble matrix */ 53c4762a1bSJed Brown if (bs == 1) { 549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL)); 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; 599566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 609566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 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 669566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 679566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 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 739566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 749566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown } else if (prob == 2) { /* matrix for the five point stencil */ 77c4762a1bSJed Brown n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001); 78bc30f867SBarry Smith PetscCheck(n1*n1 == n,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; 849566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 859566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown if (i<n1-1) { 88c4762a1bSJed Brown J = Ii + n1; 899566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 909566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 91c4762a1bSJed Brown } 92c4762a1bSJed Brown if (j>0) { 93c4762a1bSJed Brown J = Ii - 1; 949566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 959566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown if (j<n1-1) { 98c4762a1bSJed Brown J = Ii + 1; 999566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 1009566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 101c4762a1bSJed Brown } 1029566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES)); 1039566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&Ii,1,&Ii,&four,INSERT_VALUES)); 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; 1149566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 1159566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 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 1219566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 1229566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES)); 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 1289566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 1299566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES)); 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 1369566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,col,value,INSERT_VALUES)); 1379566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES)); 138c4762a1bSJed Brown 139c4762a1bSJed Brown col[0]=i; row=i+bs; 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES)); 1429566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES)); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown } 1459566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1469566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 147c4762a1bSJed Brown 1489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY)); 1499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY)); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* Test MatGetInfo() of A and sA */ 1529566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A,MAT_LOCAL,&minfo1)); 1539566063dSJacob Faibussowitsch PetscCall(MatGetInfo(sA,MAT_LOCAL,&minfo2)); 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) { 1599566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error (compare A and sA): MatGetInfo()\n")); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* Test MatDuplicate() */ 1639566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_FROBENIUS,&norm1)); 1649566063dSJacob Faibussowitsch PetscCall(MatDuplicate(sA,MAT_COPY_VALUES,&sB)); 1659566063dSJacob Faibussowitsch PetscCall(MatEqual(sA,sB,&equal)); 16628b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDuplicate()"); 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* Test MatNorm() */ 1699566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_FROBENIUS,&norm1)); 1709566063dSJacob Faibussowitsch PetscCall(MatNorm(sB,NORM_FROBENIUS,&norm2)); 171c4762a1bSJed Brown rnorm = PetscAbsReal(norm1-norm2)/norm2; 172c4762a1bSJed Brown if (rnorm > tol) { 1739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS, NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2)); 174c4762a1bSJed Brown } 1759566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_INFINITY,&norm1)); 1769566063dSJacob Faibussowitsch PetscCall(MatNorm(sB,NORM_INFINITY,&norm2)); 177c4762a1bSJed Brown rnorm = PetscAbsReal(norm1-norm2)/norm2; 178c4762a1bSJed Brown if (rnorm > tol) { 1799566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2)); 180c4762a1bSJed Brown } 1819566063dSJacob Faibussowitsch PetscCall(MatNorm(A,NORM_1,&norm1)); 1829566063dSJacob Faibussowitsch PetscCall(MatNorm(sB,NORM_1,&norm2)); 183c4762a1bSJed Brown rnorm = PetscAbsReal(norm1-norm2)/norm2; 184c4762a1bSJed Brown if (rnorm > tol) { 1859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2)); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */ 1899566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A,MAT_LOCAL,&minfo1)); 1909566063dSJacob Faibussowitsch PetscCall(MatGetInfo(sB,MAT_LOCAL,&minfo2)); 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) { 1969566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error(compare A and sB): MatGetInfo()\n")); 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 1999566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&Ii,&J)); 2009566063dSJacob Faibussowitsch PetscCall(MatGetSize(sB,&i,&j)); 201c4762a1bSJed Brown if (i-Ii || j-J) { 2029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()\n")); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 2059566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &Ii)); 2069566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(sB, &i)); 207c4762a1bSJed Brown if (i-Ii) { 2089566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()\n")); 209c4762a1bSJed Brown } 210c4762a1bSJed Brown 2119566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); 2129566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 2139566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&x)); 2149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&s1)); 2159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&s2)); 2169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&y)); 2179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&b)); 2189566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rdm)); 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() */ 2249566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A,x,x)); 2259566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(sB,x,x)); 2269566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,sB,10,&equal)); 22728b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDiagonalScale"); 228c4762a1bSJed Brown 2299566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A,s1)); 2309566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(sB,s2)); 2319566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2,neg_one,s1)); 2329566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_1,&norm1)); 233c4762a1bSJed Brown if (norm1>tol) { 2349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal(), ||s1-s2||=%g\n",(double)norm1)); 235c4762a1bSJed Brown } 236c4762a1bSJed Brown 237c4762a1bSJed Brown { 238c4762a1bSJed Brown PetscScalar alpha=0.1; 2399566063dSJacob Faibussowitsch PetscCall(MatScale(A,alpha)); 2409566063dSJacob Faibussowitsch PetscCall(MatScale(sB,alpha)); 241c4762a1bSJed Brown } 242c4762a1bSJed Brown #endif 243c4762a1bSJed Brown 244c4762a1bSJed Brown /* Test MatGetRowMaxAbs() */ 2459566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(A,s1,NULL)); 2469566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(sB,s2,NULL)); 2479566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_1,&norm1)); 2489566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_1,&norm2)); 249c4762a1bSJed Brown norm1 -= norm2; 250c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 2519566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatGetRowMaxAbs() \n")); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* Test MatMult() */ 255c4762a1bSJed Brown for (i=0; i<40; i++) { 2569566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rdm)); 2579566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,s1)); 2589566063dSJacob Faibussowitsch PetscCall(MatMult(sB,x,s2)); 2599566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_1,&norm1)); 2609566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_1,&norm2)); 261c4762a1bSJed Brown norm1 -= norm2; 262c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 2639566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatMult(), norm1-norm2: %g\n",(double)norm1)); 264c4762a1bSJed Brown } 265c4762a1bSJed Brown } 266c4762a1bSJed Brown 267c4762a1bSJed Brown /* MatMultAdd() */ 268c4762a1bSJed Brown for (i=0; i<40; i++) { 2699566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rdm)); 2709566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y,rdm)); 2719566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A,x,y,s1)); 2729566063dSJacob Faibussowitsch PetscCall(MatMultAdd(sB,x,y,s2)); 2739566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_1,&norm1)); 2749566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_1,&norm2)); 275c4762a1bSJed Brown norm1 -= norm2; 276c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 2779566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd(), norm1-norm2: %g\n",(double)norm1)); 278c4762a1bSJed Brown } 279c4762a1bSJed Brown } 280c4762a1bSJed Brown 281c20d7725SJed Brown /* Test MatMatMult() for sbaij and dense matrices */ 2829566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,5*n,NULL,&B)); 2839566063dSJacob Faibussowitsch PetscCall(MatSetRandom(B,rdm)); 2849566063dSJacob Faibussowitsch PetscCall(MatMatMult(sA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 2859566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(sA,B,C,5*n,&equal)); 28628b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult()"); 2879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 289c4762a1bSJed Brown 290c4762a1bSJed Brown /* Test MatCholeskyFactor(), MatICCFactor() with natural ordering */ 2919566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A,MATORDERINGNATURAL,&perm,&iscol)); 2929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 293c4762a1bSJed Brown norm1 = tol; 294c4762a1bSJed Brown inc = bs; 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* initialize factinfo */ 2979566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&factinfo,sizeof(MatFactorInfo))); 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 3039566063dSJacob Faibussowitsch PetscCall(MatGetFactor(sB,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sFactor)); 3049566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(sFactor,sB,perm,&factinfo)); 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 3109566063dSJacob Faibussowitsch PetscCall(MatGetFactor(sB,MATSOLVERPETSC,MAT_FACTOR_ICC,&sFactor)); 3119566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(sFactor,sB,perm,&factinfo)); 312c4762a1bSJed Brown } 3139566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(sFactor,sB,&factinfo)); 314c4762a1bSJed Brown /* MatView(sFactor, PETSC_VIEWER_DRAW_WORLD); */ 315c4762a1bSJed Brown 316c4762a1bSJed Brown /* test MatGetDiagonal on numeric factor */ 317c4762a1bSJed Brown /* 318c4762a1bSJed Brown if (lf == -1) { 3199566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(sFactor,s1)); 320c4762a1bSJed Brown printf(" in ex74.c, diag: \n"); 3219566063dSJacob Faibussowitsch PetscCall(VecView(s1,PETSC_VIEWER_STDOUT_SELF)); 322c4762a1bSJed Brown } 323c4762a1bSJed Brown */ 324c4762a1bSJed Brown 3259566063dSJacob Faibussowitsch PetscCall(MatMult(sB,x,b)); 326c4762a1bSJed Brown 327c4762a1bSJed Brown /* test MatForwardSolve() and MatBackwardSolve() */ 328c4762a1bSJed Brown if (lf == -1) { 3299566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(sFactor,b,s1)); 3309566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(sFactor,s1,s2)); 3319566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2,neg_one,x)); 3329566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&norm2)); 333c4762a1bSJed Brown if (10*norm1 < norm2) { 3349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g, bs=%" PetscInt_FMT "\n",(double)norm2,bs)); 335c4762a1bSJed Brown } 336c4762a1bSJed Brown } 337c4762a1bSJed Brown 338c4762a1bSJed Brown /* test MatSolve() */ 3399566063dSJacob Faibussowitsch PetscCall(MatSolve(sFactor,b,y)); 3409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sFactor)); 341c4762a1bSJed Brown /* Check the error */ 3429566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,neg_one,x)); 3439566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&norm2)); 344c4762a1bSJed Brown if (10*norm1 < norm2 && lf-inc != -1) { 3459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"lf=%" PetscInt_FMT ", %" PetscInt_FMT ", Norm of error=%g, %g\n",lf-inc,lf,(double)norm1,(double)norm2)); 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) 3529566063dSJacob Faibussowitsch PetscCall(MatGetFactor(sA,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&sFactor)); 3539566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(sFactor,sA,NULL,NULL)); 3549566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(sFactor,sA,NULL)); 355c4762a1bSJed Brown for (i=0; i<10; i++) { 3569566063dSJacob Faibussowitsch PetscCall(VecSetRandom(b,rdm)); 3579566063dSJacob Faibussowitsch PetscCall(MatSolve(sFactor,b,y)); 358c4762a1bSJed Brown /* Check the error */ 3599566063dSJacob Faibussowitsch PetscCall(MatMult(sA,y,x)); 3609566063dSJacob Faibussowitsch PetscCall(VecAXPY(x,neg_one,b)); 3619566063dSJacob Faibussowitsch PetscCall(VecNorm(x,NORM_2,&norm2)); 362c4762a1bSJed Brown if (norm2>tol) { 3639566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolve(), norm2: %g\n",(double)norm2)); 364c4762a1bSJed Brown } 365c4762a1bSJed Brown } 3669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sFactor)); 367c4762a1bSJed Brown #endif 368c4762a1bSJed Brown 3699566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 370c4762a1bSJed Brown 3719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sB)); 3739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 3749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 3769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 3779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 3789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 3799566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 380c4762a1bSJed Brown 3819566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 382b122ec5aSJacob Faibussowitsch return 0; 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