xref: /petsc/src/mat/tests/ex17.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscScalar    v,five = 5.0,one = 1.0;
11c4762a1bSJed Brown   IS             isrow,row,col;
12c4762a1bSJed Brown   Vec            x,u,b;
13c4762a1bSJed Brown   PetscReal      norm;
14c4762a1bSJed Brown   MatFactorInfo  info;
15c4762a1bSJed Brown 
16*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
19c4762a1bSJed Brown 
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
24c4762a1bSJed Brown   for (i=0; i<m; i++) {
25c4762a1bSJed Brown     for (j=0; j<n; j++) {
26c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
275f80ce2aSJacob Faibussowitsch       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
285f80ce2aSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
295f80ce2aSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
305f80ce2aSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
315f80ce2aSJacob Faibussowitsch       v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
32c4762a1bSJed Brown     }
33c4762a1bSJed Brown   }
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
36c4762a1bSJed Brown 
375f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_SELF,(m*n)/2,0,2,&isrow));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroRowsIS(C,isrow,five,0,0));
39c4762a1bSJed Brown 
405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m*n,&u));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(u,&x));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(u,&b));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(u,one));
44c4762a1bSJed Brown 
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(C,u,b));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Set default ordering to be Quotient Minimum Degree; also read
48c4762a1bSJed Brown      orderings from the options database */
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(C,MATORDERINGQMD,&row,&col));
50c4762a1bSJed Brown 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatFactorInfoInitialize(&info));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&A));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactorSymbolic(A,C,row,col,&info));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactorNumeric(A,C,&info));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSolveTranspose(A,b,x));
56c4762a1bSJed Brown 
575f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(row,PETSC_VIEWER_STDOUT_SELF));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(x,-1.0,u));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(x,NORM_2,&norm));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",(double)norm));
61c4762a1bSJed Brown 
625f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&row));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&col));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrow));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
70*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
71*b122ec5aSJacob Faibussowitsch   return 0;
72c4762a1bSJed Brown }
73c4762a1bSJed Brown 
74c4762a1bSJed Brown /*TEST
75c4762a1bSJed Brown 
76c4762a1bSJed Brown    test:
77c4762a1bSJed Brown 
78c4762a1bSJed Brown TEST*/
79