1 2 static char help[] = "Tests the use of MatSolveTranspose().\n\n"; 3 4 #include <petscmat.h> 5 6 int main(int argc,char **args) 7 { 8 Mat C,A; 9 PetscInt i,j,m = 5,n = 5,Ii,J; 10 PetscErrorCode ierr; 11 PetscScalar v,five = 5.0,one = 1.0; 12 IS isrow,row,col; 13 Vec x,u,b; 14 PetscReal norm; 15 MatFactorInfo info; 16 17 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 18 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 19 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 20 21 CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C)); 22 CHKERRQ(MatSetUp(C)); 23 24 /* create the matrix for the five point stencil, YET AGAIN*/ 25 for (i=0; i<m; i++) { 26 for (j=0; j<n; j++) { 27 v = -1.0; Ii = j + n*i; 28 if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 29 if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 30 if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 31 if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 32 v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 33 } 34 } 35 CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 36 CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 37 38 CHKERRQ(ISCreateStride(PETSC_COMM_SELF,(m*n)/2,0,2,&isrow)); 39 CHKERRQ(MatZeroRowsIS(C,isrow,five,0,0)); 40 41 CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m*n,&u)); 42 CHKERRQ(VecDuplicate(u,&x)); 43 CHKERRQ(VecDuplicate(u,&b)); 44 CHKERRQ(VecSet(u,one)); 45 46 CHKERRQ(MatMultTranspose(C,u,b)); 47 48 /* Set default ordering to be Quotient Minimum Degree; also read 49 orderings from the options database */ 50 CHKERRQ(MatGetOrdering(C,MATORDERINGQMD,&row,&col)); 51 52 CHKERRQ(MatFactorInfoInitialize(&info)); 53 CHKERRQ(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&A)); 54 CHKERRQ(MatLUFactorSymbolic(A,C,row,col,&info)); 55 CHKERRQ(MatLUFactorNumeric(A,C,&info)); 56 CHKERRQ(MatSolveTranspose(A,b,x)); 57 58 CHKERRQ(ISView(row,PETSC_VIEWER_STDOUT_SELF)); 59 CHKERRQ(VecAXPY(x,-1.0,u)); 60 CHKERRQ(VecNorm(x,NORM_2,&norm)); 61 CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",(double)norm)); 62 63 CHKERRQ(ISDestroy(&row)); 64 CHKERRQ(ISDestroy(&col)); 65 CHKERRQ(ISDestroy(&isrow)); 66 CHKERRQ(VecDestroy(&u)); 67 CHKERRQ(VecDestroy(&x)); 68 CHKERRQ(VecDestroy(&b)); 69 CHKERRQ(MatDestroy(&C)); 70 CHKERRQ(MatDestroy(&A)); 71 ierr = PetscFinalize(); 72 return ierr; 73 } 74 75 /*TEST 76 77 test: 78 79 TEST*/ 80