1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests the use of MatSolveTranspose().\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat C,A; 9c4762a1bSJed Brown PetscInt i,j,m = 5,n = 5,Ii,J; 10c4762a1bSJed Brown PetscScalar v,five = 5.0,one = 1.0; 11c4762a1bSJed Brown IS isrow,row,col; 12c4762a1bSJed Brown Vec x,u,b; 13c4762a1bSJed Brown PetscReal norm; 14c4762a1bSJed Brown MatFactorInfo info; 15c4762a1bSJed Brown 16*327415f7SBarry Smith PetscFunctionBeginUser; 179566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 20c4762a1bSJed Brown 219566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C)); 229566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* create the matrix for the five point stencil, YET AGAIN*/ 25c4762a1bSJed Brown for (i=0; i<m; i++) { 26c4762a1bSJed Brown for (j=0; j<n; j++) { 27c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 289566063dSJacob Faibussowitsch if (i>0) {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 299566063dSJacob Faibussowitsch if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 309566063dSJacob Faibussowitsch if (j>0) {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 319566063dSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 329566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 33c4762a1bSJed Brown } 34c4762a1bSJed Brown } 359566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 369566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,(m*n)/2,0,2,&isrow)); 399566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(C,isrow,five,0,0)); 40c4762a1bSJed Brown 419566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,m*n,&u)); 429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&x)); 439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&b)); 449566063dSJacob Faibussowitsch PetscCall(VecSet(u,one)); 45c4762a1bSJed Brown 469566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(C,u,b)); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Set default ordering to be Quotient Minimum Degree; also read 49c4762a1bSJed Brown orderings from the options database */ 509566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C,MATORDERINGQMD,&row,&col)); 51c4762a1bSJed Brown 529566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 539566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&A)); 549566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(A,C,row,col,&info)); 559566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(A,C,&info)); 569566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(A,b,x)); 57c4762a1bSJed Brown 589566063dSJacob Faibussowitsch PetscCall(ISView(row,PETSC_VIEWER_STDOUT_SELF)); 599566063dSJacob Faibussowitsch PetscCall(VecAXPY(x,-1.0,u)); 609566063dSJacob Faibussowitsch PetscCall(VecNorm(x,NORM_2,&norm)); 619566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",(double)norm)); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 649566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col)); 659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 719566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 72b122ec5aSJacob Faibussowitsch return 0; 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown /*TEST 76c4762a1bSJed Brown 77c4762a1bSJed Brown test: 78c4762a1bSJed Brown 79c4762a1bSJed Brown TEST*/ 80