Actual source code: ex17.c
1: /*$Id: ex17.c,v 1.21 2001/08/07 21:30:08 bsmith Exp $*/
3: static char help[] = "Tests the use of MatSolveTranspose().\n\n";
5: #include petscmat.h
9: int main(int argc,char **args)
10: {
11: Mat C,A;
12: int i,j,m = 5,n = 5,I,J,ierr;
13: PetscScalar v,five = 5.0,one = 1.0,mone = -1.0;
14: IS isrow,row,col;
15: Vec x,u,b;
16: PetscReal norm;
18: PetscInitialize(&argc,&args,(char *)0,help);
19: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
20: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
22: MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,PETSC_NULL,&C);
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; I = j + n*i;
28: if (i>0) {J = I - n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
29: if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
30: if (j>0) {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
31: if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
32: v = 4.0; MatSetValues(C,1,&I,1,&I,&v,INSERT_VALUES);
33: }
34: }
35: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
36: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
38: ISCreateStride(PETSC_COMM_SELF,(m*n)/2,0,2,&isrow);
39: MatZeroRows(C,isrow,&five);
41: VecCreateSeq(PETSC_COMM_SELF,m*n,&u);
42: VecDuplicate(u,&x);
43: VecDuplicate(u,&b);
44: VecSet(&one,u);
46: MatMultTranspose(C,u,b);
48: /* Set default ordering to be Quotient Minimum Degree; also read
49: orderings from the options database */
50: MatGetOrdering(C,MATORDERING_QMD,&row,&col);
52: MatLUFactorSymbolic(C,row,col,PETSC_NULL,&A);
53: MatLUFactorNumeric(C,&A);
54: MatSolveTranspose(A,b,x);
56: ISView(row,PETSC_VIEWER_STDOUT_SELF);
57: VecAXPY(&mone,u,x);
58: VecNorm(x,NORM_2,&norm);
59: PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",norm);
61: ISDestroy(row);
62: ISDestroy(col);
63: ISDestroy(isrow);
64: VecDestroy(u);
65: VecDestroy(x);
66: VecDestroy(b);
67: MatDestroy(C);
68: MatDestroy(A);
69: PetscFinalize();
70: return 0;
71: }