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