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