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