xref: /petsc/src/mat/tests/ex17.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
11c4762a1bSJed Brown   PetscScalar    v,five = 5.0,one = 1.0;
12c4762a1bSJed Brown   IS             isrow,row,col;
13c4762a1bSJed Brown   Vec            x,u,b;
14c4762a1bSJed Brown   PetscReal      norm;
15c4762a1bSJed Brown   MatFactorInfo  info;
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
20c4762a1bSJed Brown 
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
28*5f80ce2aSJacob Faibussowitsch       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
29*5f80ce2aSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
30*5f80ce2aSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
31*5f80ce2aSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
32*5f80ce2aSJacob Faibussowitsch       v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
33c4762a1bSJed Brown     }
34c4762a1bSJed Brown   }
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
37c4762a1bSJed Brown 
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_SELF,(m*n)/2,0,2,&isrow));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroRowsIS(C,isrow,five,0,0));
40c4762a1bSJed Brown 
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m*n,&u));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(u,&x));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(u,&b));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(u,one));
45c4762a1bSJed Brown 
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 */
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(C,MATORDERINGQMD,&row,&col));
51c4762a1bSJed Brown 
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatFactorInfoInitialize(&info));
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&A));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactorSymbolic(A,C,row,col,&info));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactorNumeric(A,C,&info));
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSolveTranspose(A,b,x));
57c4762a1bSJed Brown 
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(row,PETSC_VIEWER_STDOUT_SELF));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(x,-1.0,u));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(x,NORM_2,&norm));
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",(double)norm));
62c4762a1bSJed Brown 
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&row));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&col));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrow));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
71c4762a1bSJed Brown   ierr = PetscFinalize();
72c4762a1bSJed Brown   return ierr;
73c4762a1bSJed Brown }
74c4762a1bSJed Brown 
75c4762a1bSJed Brown /*TEST
76c4762a1bSJed Brown 
77c4762a1bSJed Brown    test:
78c4762a1bSJed Brown 
79c4762a1bSJed Brown TEST*/
80