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