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; 245f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 252c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL)); 28c4762a1bSJed Brown 29c4762a1bSJed Brown n = mbs*bs; 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSEQBAIJ)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetPreallocation(A,bs,nz,NULL)); 35c4762a1bSJed Brown 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&sA)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(sA,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(sA,MATSEQSBAIJ)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(sA)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetType(sA,&type)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)sA,MATSEQSBAIJ,&doIcc)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqSBAIJSetPreallocation(sA,bs,nz,NULL)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(sA,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&Ii,&J)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(sA,&i,&j)); 48c4762a1bSJed Brown if (i-Ii || j-J) { 495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ format\n")); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Assemble matrix */ 53c4762a1bSJed Brown if (bs == 1) { 545f80ce2aSJacob Faibussowitsch CHKERRQ(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; 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(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 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(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 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(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); 782c71b3e2SJacob Faibussowitsch PetscCheckFalse(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; 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown if (i<n1-1) { 88c4762a1bSJed Brown J = Ii + n1; 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 91c4762a1bSJed Brown } 92c4762a1bSJed Brown if (j>0) { 93c4762a1bSJed Brown J = Ii - 1; 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown if (j<n1-1) { 98c4762a1bSJed Brown J = Ii + 1; 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES)); 101c4762a1bSJed Brown } 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(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; 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(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 1215f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(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 1285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(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 1365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,col,value,INSERT_VALUES)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES)); 138c4762a1bSJed Brown 139c4762a1bSJed Brown col[0]=i; row=i+bs; 140c4762a1bSJed Brown 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES)); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown } 1455f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 147c4762a1bSJed Brown 1485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY)); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* Test MatGetInfo() of A and sA */ 1525f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(A,MAT_LOCAL,&minfo1)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 1595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error (compare A and sA): MatGetInfo()\n")); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* Test MatDuplicate() */ 1635f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_FROBENIUS,&norm1)); 1645f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(sA,MAT_COPY_VALUES,&sB)); 1655f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(sA,sB,&equal)); 166*28b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDuplicate()"); 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* Test MatNorm() */ 1695f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_FROBENIUS,&norm1)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(sB,NORM_FROBENIUS,&norm2)); 171c4762a1bSJed Brown rnorm = PetscAbsReal(norm1-norm2)/norm2; 172c4762a1bSJed Brown if (rnorm > tol) { 1735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS, NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2)); 174c4762a1bSJed Brown } 1755f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_INFINITY,&norm1)); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(sB,NORM_INFINITY,&norm2)); 177c4762a1bSJed Brown rnorm = PetscAbsReal(norm1-norm2)/norm2; 178c4762a1bSJed Brown if (rnorm > tol) { 1795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n",(double)norm1,(double)norm2)); 180c4762a1bSJed Brown } 1815f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(A,NORM_1,&norm1)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(sB,NORM_1,&norm2)); 183c4762a1bSJed Brown rnorm = PetscAbsReal(norm1-norm2)/norm2; 184c4762a1bSJed Brown if (rnorm > tol) { 1855f80ce2aSJacob Faibussowitsch CHKERRQ(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() */ 1895f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(A,MAT_LOCAL,&minfo1)); 1905f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 1965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error(compare A and sB): MatGetInfo()\n")); 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 1995f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&Ii,&J)); 2005f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(sB,&i,&j)); 201c4762a1bSJed Brown if (i-Ii || j-J) { 2025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()\n")); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 2055f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(A, &Ii)); 2065f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(sB, &i)); 207c4762a1bSJed Brown if (i-Ii) { 2085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()\n")); 209c4762a1bSJed Brown } 210c4762a1bSJed Brown 2115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rdm)); 2135f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&x)); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&s1)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&s2)); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&y)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&b)); 2185f80ce2aSJacob Faibussowitsch CHKERRQ(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() */ 2245f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(A,x,x)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(sB,x,x)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,sB,10,&equal)); 227*28b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDiagonalScale"); 228c4762a1bSJed Brown 2295f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(A,s1)); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(sB,s2)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(s2,neg_one,s1)); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_1,&norm1)); 233c4762a1bSJed Brown if (norm1>tol) { 2345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal(), ||s1-s2||=%g\n",(double)norm1)); 235c4762a1bSJed Brown } 236c4762a1bSJed Brown 237c4762a1bSJed Brown { 238c4762a1bSJed Brown PetscScalar alpha=0.1; 2395f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,alpha)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(sB,alpha)); 241c4762a1bSJed Brown } 242c4762a1bSJed Brown #endif 243c4762a1bSJed Brown 244c4762a1bSJed Brown /* Test MatGetRowMaxAbs() */ 2455f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMaxAbs(A,s1,NULL)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowMaxAbs(sB,s2,NULL)); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s1,NORM_1,&norm1)); 2485f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_1,&norm2)); 249c4762a1bSJed Brown norm1 -= norm2; 250c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 2515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatGetRowMaxAbs() \n")); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* Test MatMult() */ 255c4762a1bSJed Brown for (i=0; i<40; i++) { 2565f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rdm)); 2575f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,s1)); 2585f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(sB,x,s2)); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s1,NORM_1,&norm1)); 2605f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_1,&norm2)); 261c4762a1bSJed Brown norm1 -= norm2; 262c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 2635f80ce2aSJacob Faibussowitsch CHKERRQ(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++) { 2695f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rdm)); 2705f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(y,rdm)); 2715f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(A,x,y,s1)); 2725f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(sB,x,y,s2)); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s1,NORM_1,&norm1)); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_1,&norm2)); 275c4762a1bSJed Brown norm1 -= norm2; 276c4762a1bSJed Brown if (norm1<-tol || norm1>tol) { 2775f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 2825f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,n,5*n,NULL,&B)); 2835f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(B,rdm)); 2845f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(sA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 2855f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(sA,B,C,5*n,&equal)); 286*28b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult()"); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 2885f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 289c4762a1bSJed Brown 290c4762a1bSJed Brown /* Test MatCholeskyFactor(), MatICCFactor() with natural ordering */ 2915f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(A,MATORDERINGNATURAL,&perm,&iscol)); 2925f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscol)); 293c4762a1bSJed Brown norm1 = tol; 294c4762a1bSJed Brown inc = bs; 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* initialize factinfo */ 2975f80ce2aSJacob Faibussowitsch CHKERRQ(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 3035f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(sB,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sFactor)); 3045f80ce2aSJacob Faibussowitsch CHKERRQ(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 3105f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(sB,MATSOLVERPETSC,MAT_FACTOR_ICC,&sFactor)); 3115f80ce2aSJacob Faibussowitsch CHKERRQ(MatICCFactorSymbolic(sFactor,sB,perm,&factinfo)); 312c4762a1bSJed Brown } 3135f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 3195f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(sFactor,s1)); 320c4762a1bSJed Brown printf(" in ex74.c, diag: \n"); 3215f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(s1,PETSC_VIEWER_STDOUT_SELF)); 322c4762a1bSJed Brown } 323c4762a1bSJed Brown */ 324c4762a1bSJed Brown 3255f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(sB,x,b)); 326c4762a1bSJed Brown 327c4762a1bSJed Brown /* test MatForwardSolve() and MatBackwardSolve() */ 328c4762a1bSJed Brown if (lf == -1) { 3295f80ce2aSJacob Faibussowitsch CHKERRQ(MatForwardSolve(sFactor,b,s1)); 3305f80ce2aSJacob Faibussowitsch CHKERRQ(MatBackwardSolve(sFactor,s1,s2)); 3315f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(s2,neg_one,x)); 3325f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_2,&norm2)); 333c4762a1bSJed Brown if (10*norm1 < norm2) { 3345f80ce2aSJacob Faibussowitsch CHKERRQ(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() */ 3395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(sFactor,b,y)); 3405f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sFactor)); 341c4762a1bSJed Brown /* Check the error */ 3425f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,neg_one,x)); 3435f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm2)); 344c4762a1bSJed Brown if (10*norm1 < norm2 && lf-inc != -1) { 3455f80ce2aSJacob Faibussowitsch CHKERRQ(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) 3525f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(sA,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&sFactor)); 3535f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorSymbolic(sFactor,sA,NULL,NULL)); 3545f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorNumeric(sFactor,sA,NULL)); 355c4762a1bSJed Brown for (i=0; i<10; i++) { 3565f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(b,rdm)); 3575f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(sFactor,b,y)); 358c4762a1bSJed Brown /* Check the error */ 3595f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(sA,y,x)); 3605f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(x,neg_one,b)); 3615f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(x,NORM_2,&norm2)); 362c4762a1bSJed Brown if (norm2>tol) { 3635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolve(), norm2: %g\n",(double)norm2)); 364c4762a1bSJed Brown } 365c4762a1bSJed Brown } 3665f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sFactor)); 367c4762a1bSJed Brown #endif 368c4762a1bSJed Brown 3695f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&perm)); 370c4762a1bSJed Brown 3715f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 3725f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sB)); 3735f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sA)); 3745f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 3755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&s1)); 3775f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&s2)); 3785f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 3795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rdm)); 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