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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 235f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 242c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown n = mbs*bs; 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSEQBAIJ)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetPreallocation(A,bs,nz,NULL)); 34c4762a1bSJed Brown 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&sA)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(sA,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(sA,MATSEQSBAIJ)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(sA)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetType(sA,&type)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)sA,MATSEQSBAIJ,&doIcc)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqSBAIJSetPreallocation(sA,bs,nz,NULL)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(sA,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&Ii,&J)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(sA,&i,&j)); 47c4762a1bSJed Brown if (i-Ii || j-J) { 485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ format\n")); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Assemble matrix */ 52c4762a1bSJed Brown if (bs == 1) { 535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL)); 54c4762a1bSJed Brown if (prob == 1) { /* tridiagonal matrix */ 55c4762a1bSJed Brown value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; 56c4762a1bSJed Brown for (i=1; i<n-1; i++) { 57c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1; 62c4762a1bSJed Brown 63c4762a1bSJed Brown value[0]= 0.1; value[1]=-1; value[2]=2; 64c4762a1bSJed Brown 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown i = 0; 69c4762a1bSJed Brown col[0] = n-1; col[1] = 1; col[2] = 0; 70c4762a1bSJed Brown value[0] = 0.1; value[1] = -1.0; value[2] = 2; 71c4762a1bSJed Brown 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 74c4762a1bSJed Brown 75c4762a1bSJed Brown } else if (prob == 2) { /* matrix for the five point stencil */ 76c4762a1bSJed Brown n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001); 772c71b3e2SJacob Faibussowitsch PetscCheckFalse(n1*n1 - n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive integer!"); 78c4762a1bSJed Brown for (i=0; i<n1; i++) { 79c4762a1bSJed Brown for (j=0; j<n1; j++) { 80c4762a1bSJed Brown Ii = j + n1*i; 81c4762a1bSJed Brown if (i>0) { 82c4762a1bSJed Brown J = Ii - n1; 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 85c4762a1bSJed Brown } 86c4762a1bSJed Brown if (i<n1-1) { 87c4762a1bSJed Brown J = Ii + n1; 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown if (j>0) { 92c4762a1bSJed Brown J = Ii - 1; 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown if (j<n1-1) { 97c4762a1bSJed Brown J = Ii + 1; 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 100c4762a1bSJed Brown } 1015f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&Ii,1,&Ii,&four,INSERT_VALUES)); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown } 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107c4762a1bSJed Brown } else { /* bs > 1 */ 108c4762a1bSJed Brown for (block=0; block<n/bs; block++) { 109c4762a1bSJed Brown /* diagonal blocks */ 110c4762a1bSJed Brown value[0] = -1.0; value[1] = 4.0; value[2] = -1.0; 111c4762a1bSJed Brown for (i=1+block*bs; i<bs-1+block*bs; i++) { 112c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES)); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs; 117c4762a1bSJed Brown 118c4762a1bSJed Brown value[0]=-1.0; value[1]=4.0; 119c4762a1bSJed Brown 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES)); 122c4762a1bSJed Brown 123c4762a1bSJed Brown i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs; 124c4762a1bSJed Brown 125c4762a1bSJed Brown value[0]=4.0; value[1] = -1.0; 126c4762a1bSJed Brown 1275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES)); 129c4762a1bSJed Brown } 130c4762a1bSJed Brown /* off-diagonal blocks */ 131c4762a1bSJed Brown value[0]=-1.0; 132c4762a1bSJed Brown for (i=0; i<(n/bs-1)*bs; i++) { 133c4762a1bSJed Brown col[0]=i+bs; 134c4762a1bSJed Brown 1355f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,col,value,INSERT_VALUES)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown col[0]=i; row=i+bs; 139c4762a1bSJed Brown 1405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES)); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown } 1445f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 146c4762a1bSJed Brown 1475f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* Test MatGetInfo() of A and sA */ 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(A,MAT_LOCAL,&minfo1)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(sA,MAT_LOCAL,&minfo2)); 153c4762a1bSJed Brown i = (int) (minfo1.nz_used - minfo2.nz_used); 154c4762a1bSJed Brown j = (int) (minfo1.nz_allocated - minfo2.nz_allocated); 155c4762a1bSJed Brown k1 = (int) (minfo1.nz_allocated - minfo1.nz_used); 156c4762a1bSJed Brown k2 = (int) (minfo2.nz_allocated - minfo2.nz_used); 157c4762a1bSJed Brown if (i < 0 || j < 0 || k1 < 0 || k2 < 0) { 1585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error (compare A and sA): MatGetInfo()\n")); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* Test MatDuplicate() */ 1625f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_FROBENIUS,&norm1)); 1635f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(sA,MAT_COPY_VALUES,&sB)); 1645f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(sA,sB,&equal)); 16528b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDuplicate()"); 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* Test MatNorm() */ 1685f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_FROBENIUS,&norm1)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(sB,NORM_FROBENIUS,&norm2)); 170c4762a1bSJed Brown rnorm = PetscAbsReal(norm1-norm2)/norm2; 171c4762a1bSJed Brown if (rnorm > tol) { 1725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS, NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2)); 173c4762a1bSJed Brown } 1745f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_INFINITY,&norm1)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(sB,NORM_INFINITY,&norm2)); 176c4762a1bSJed Brown rnorm = PetscAbsReal(norm1-norm2)/norm2; 177c4762a1bSJed Brown if (rnorm > tol) { 1785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2)); 179c4762a1bSJed Brown } 1805f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_1,&norm1)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(sB,NORM_1,&norm2)); 182c4762a1bSJed Brown rnorm = PetscAbsReal(norm1-norm2)/norm2; 183c4762a1bSJed Brown if (rnorm > tol) { 1845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2)); 185c4762a1bSJed Brown } 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */ 1885f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(A,MAT_LOCAL,&minfo1)); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(sB,MAT_LOCAL,&minfo2)); 190c4762a1bSJed Brown i = (int) (minfo1.nz_used - minfo2.nz_used); 191c4762a1bSJed Brown j = (int) (minfo1.nz_allocated - minfo2.nz_allocated); 192c4762a1bSJed Brown k1 = (int) (minfo1.nz_allocated - minfo1.nz_used); 193c4762a1bSJed Brown k2 = (int) (minfo2.nz_allocated - minfo2.nz_used); 194c4762a1bSJed Brown if (i < 0 || j < 0 || k1 < 0 || k2 < 0) { 1955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error(compare A and sB): MatGetInfo()\n")); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 1985f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&Ii,&J)); 1995f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(sB,&i,&j)); 200c4762a1bSJed Brown if (i-Ii || j-J) { 2015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()\n")); 202c4762a1bSJed Brown } 203c4762a1bSJed Brown 2045f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(A, &Ii)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(sB, &i)); 206c4762a1bSJed Brown if (i-Ii) { 2075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()\n")); 208c4762a1bSJed Brown } 209c4762a1bSJed Brown 2105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); 2115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rdm)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&x)); 2135f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&s1)); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&s2)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&y)); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&b)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rdm)); 218c4762a1bSJed Brown 219c4762a1bSJed Brown /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */ 220c4762a1bSJed Brown #if !defined(PETSC_USE_COMPLEX) 221c4762a1bSJed Brown /* Scaling matrix with complex numbers results non-spd matrix, 222c4762a1bSJed Brown causing crash of MatForwardSolve() and MatBackwardSolve() */ 2235f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(A,x,x)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(sB,x,x)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,sB,10,&equal)); 22628b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDiagonalScale"); 227c4762a1bSJed Brown 2285f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(A,s1)); 2295f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(sB,s2)); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(s2,neg_one,s1)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_1,&norm1)); 232c4762a1bSJed Brown if (norm1>tol) { 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal(), ||s1-s2||=%g\n",(double)norm1)); 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 236c4762a1bSJed Brown { 237c4762a1bSJed Brown PetscScalar alpha=0.1; 2385f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,alpha)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(sB,alpha)); 240c4762a1bSJed Brown } 241c4762a1bSJed Brown #endif 242c4762a1bSJed Brown 243c4762a1bSJed Brown /* Test MatGetRowMaxAbs() */ 2445f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMaxAbs(A,s1,NULL)); 2455f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMaxAbs(sB,s2,NULL)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s1,NORM_1,&norm1)); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_1,&norm2)); 248c4762a1bSJed Brown norm1 -= norm2; 249c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 2505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatGetRowMaxAbs() \n")); 251c4762a1bSJed Brown } 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* Test MatMult() */ 254c4762a1bSJed Brown for (i=0; i<40; i++) { 2555f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rdm)); 2565f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,s1)); 2575f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(sB,x,s2)); 2585f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s1,NORM_1,&norm1)); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_1,&norm2)); 260c4762a1bSJed Brown norm1 -= norm2; 261c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 2625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMult(), norm1-norm2: %g\n",(double)norm1)); 263c4762a1bSJed Brown } 264c4762a1bSJed Brown } 265c4762a1bSJed Brown 266c4762a1bSJed Brown /* MatMultAdd() */ 267c4762a1bSJed Brown for (i=0; i<40; i++) { 2685f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rdm)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(y,rdm)); 2705f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(A,x,y,s1)); 2715f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(sB,x,y,s2)); 2725f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s1,NORM_1,&norm1)); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_1,&norm2)); 274c4762a1bSJed Brown norm1 -= norm2; 275c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 2765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd(), norm1-norm2: %g\n",(double)norm1)); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown } 279c4762a1bSJed Brown 280c20d7725SJed Brown /* Test MatMatMult() for sbaij and dense matrices */ 2815f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,n,5*n,NULL,&B)); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(B,rdm)); 2835f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(sA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 2845f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(sA,B,C,5*n,&equal)); 28528b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult()"); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 288c4762a1bSJed Brown 289c4762a1bSJed Brown /* Test MatCholeskyFactor(), MatICCFactor() with natural ordering */ 2905f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(A,MATORDERINGNATURAL,&perm,&iscol)); 2915f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscol)); 292c4762a1bSJed Brown norm1 = tol; 293c4762a1bSJed Brown inc = bs; 294c4762a1bSJed Brown 295c4762a1bSJed Brown /* initialize factinfo */ 2965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemzero(&factinfo,sizeof(MatFactorInfo))); 297c4762a1bSJed Brown 298c4762a1bSJed Brown for (lf=-1; lf<10; lf += inc) { 299c4762a1bSJed Brown if (lf==-1) { /* Cholesky factor of sB (duplicate sA) */ 300c4762a1bSJed Brown factinfo.fill = 5.0; 301c4762a1bSJed Brown 3025f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(sB,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sFactor)); 3035f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorSymbolic(sFactor,sB,perm,&factinfo)); 304c4762a1bSJed Brown } else if (!doIcc) break; 305c4762a1bSJed Brown else { /* incomplete Cholesky factor */ 306c4762a1bSJed Brown factinfo.fill = 5.0; 307c4762a1bSJed Brown factinfo.levels = lf; 308c4762a1bSJed Brown 3095f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(sB,MATSOLVERPETSC,MAT_FACTOR_ICC,&sFactor)); 3105f80ce2aSJacob Faibussowitsch CHKERRQ(MatICCFactorSymbolic(sFactor,sB,perm,&factinfo)); 311c4762a1bSJed Brown } 3125f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorNumeric(sFactor,sB,&factinfo)); 313c4762a1bSJed Brown /* MatView(sFactor, PETSC_VIEWER_DRAW_WORLD); */ 314c4762a1bSJed Brown 315c4762a1bSJed Brown /* test MatGetDiagonal on numeric factor */ 316c4762a1bSJed Brown /* 317c4762a1bSJed Brown if (lf == -1) { 3185f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(sFactor,s1)); 319c4762a1bSJed Brown printf(" in ex74.c, diag: \n"); 3205f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(s1,PETSC_VIEWER_STDOUT_SELF)); 321c4762a1bSJed Brown } 322c4762a1bSJed Brown */ 323c4762a1bSJed Brown 3245f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(sB,x,b)); 325c4762a1bSJed Brown 326c4762a1bSJed Brown /* test MatForwardSolve() and MatBackwardSolve() */ 327c4762a1bSJed Brown if (lf == -1) { 3285f80ce2aSJacob Faibussowitsch CHKERRQ(MatForwardSolve(sFactor,b,s1)); 3295f80ce2aSJacob Faibussowitsch CHKERRQ(MatBackwardSolve(sFactor,s1,s2)); 3305f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(s2,neg_one,x)); 3315f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_2,&norm2)); 332c4762a1bSJed Brown if (10*norm1 < norm2) { 3335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g, bs=%" PetscInt_FMT "\n",(double)norm2,bs)); 334c4762a1bSJed Brown } 335c4762a1bSJed Brown } 336c4762a1bSJed Brown 337c4762a1bSJed Brown /* test MatSolve() */ 3385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(sFactor,b,y)); 3395f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sFactor)); 340c4762a1bSJed Brown /* Check the error */ 3415f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,neg_one,x)); 3425f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm2)); 343c4762a1bSJed Brown if (10*norm1 < norm2 && lf-inc != -1) { 3445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"lf=%" PetscInt_FMT ", %" PetscInt_FMT ", Norm of error=%g, %g\n",lf-inc,lf,(double)norm1,(double)norm2)); 345c4762a1bSJed Brown } 346c4762a1bSJed Brown norm1 = norm2; 347c4762a1bSJed Brown if (norm2 < tol && lf != -1) break; 348c4762a1bSJed Brown } 349c4762a1bSJed Brown 350c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 3515f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(sA,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&sFactor)); 3525f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorSymbolic(sFactor,sA,NULL,NULL)); 3535f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorNumeric(sFactor,sA,NULL)); 354c4762a1bSJed Brown for (i=0; i<10; i++) { 3555f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(b,rdm)); 3565f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(sFactor,b,y)); 357c4762a1bSJed Brown /* Check the error */ 3585f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(sA,y,x)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(x,neg_one,b)); 3605f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(x,NORM_2,&norm2)); 361c4762a1bSJed Brown if (norm2>tol) { 3625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolve(), norm2: %g\n",(double)norm2)); 363c4762a1bSJed Brown } 364c4762a1bSJed Brown } 3655f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sFactor)); 366c4762a1bSJed Brown #endif 367c4762a1bSJed Brown 3685f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&perm)); 369c4762a1bSJed Brown 3705f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 3715f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sB)); 3725f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sA)); 3735f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 3745f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 3755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&s1)); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&s2)); 3775f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 3785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rdm)); 379c4762a1bSJed Brown 380*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 381*b122ec5aSJacob Faibussowitsch return 0; 382c4762a1bSJed Brown } 383c4762a1bSJed Brown 384c4762a1bSJed Brown /*TEST 385c4762a1bSJed Brown 386c4762a1bSJed Brown test: 387c4762a1bSJed Brown args: -bs {{1 2 3 4 5 6 7 8}} 388c4762a1bSJed Brown 389c4762a1bSJed Brown TEST*/ 390