Actual source code: ex103.c
1:
2: static char help[] = "Tests PLAPACK interface.\n\n";
4: #include petscmat.h
8: int main(int argc,char **args)
9: {
10: Mat C,C1,F;
11: Vec u,x,b;
13: PetscMPIInt rank,nproc;
14: PetscInt i,M = 10,m,n,nfact,nsolve;
15: PetscScalar *array,rval;
16: PetscReal norm,tol=1.e-12;
17: IS perm,iperm;
18: MatFactorInfo info;
19: PetscRandom rand;
20: PetscTruth flg;
22: PetscInitialize(&argc,&args,(char *)0,help);
23: MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
24: MPI_Comm_size(PETSC_COMM_WORLD, &nproc);
26: /* Create matrix and vectors */
27: PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
28: MatCreate(PETSC_COMM_WORLD,&C);
29: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,M,M);
30: MatSetType(C,MATDENSE);
31: MatSetFromOptions(C);
32:
33: MatGetLocalSize(C,&m,&n);
34: if (m != n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Matrix local size m %d must equal n %d",m,n);
36: VecCreate(PETSC_COMM_WORLD,&x);
37: VecSetSizes(x,n,PETSC_DECIDE);
38: VecSetFromOptions(x);
39: VecDuplicate(x,&b);
40: VecDuplicate(x,&u); /* save the true solution */
42: /* Assembly */
43: PetscRandomCreate(PETSC_COMM_WORLD,&rand);
44: PetscRandomSetFromOptions(rand);
45: MatGetArray(C,&array);
46: for (i=0; i<m*M; i++){
47: PetscRandomGetValue(rand,&rval);
48: array[i] = rval;
49: }
50: MatRestoreArray(C,&array);
51: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
52: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
53: /*if (!rank) {printf("main, C: \n");}
54: MatView(C,PETSC_VIEWER_STDOUT_WORLD); */
56: /* Test MatDuplicate() */
57: MatDuplicate(C,MAT_COPY_VALUES,&C1);
58: MatEqual(C,C1,&flg);
59: if (!flg){
60: SETERRQ(PETSC_ERR_ARG_WRONG,"Duplicate C1 != C");
61: }
63: /* Test LU Factorization */
64: MatGetOrdering(C1,MATORDERING_NATURAL,&perm,&iperm);
65: MatGetFactor(C1,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&F);
67: MatLUFactorSymbolic(F,C1,perm,iperm,&info);
69: for (nfact = 0; nfact < 2; nfact++){
70: if (!rank) printf(" LU nfact %d\n",nfact);
71: MatLUFactorNumeric(F,C1,&info);
73: /* Test MatSolve() */
74: for (nsolve = 0; nsolve < 5; nsolve++){
75: VecGetArray(x,&array);
76: for (i=0; i<m; i++){
77: PetscRandomGetValue(rand,&rval);
78: array[i] = rval;
79: }
80: VecRestoreArray(x,&array);
81: VecCopy(x,u);
82: MatMult(C,x,b);
84: MatSolve(F,b,x);
86: /* Check the error */
87: VecAXPY(u,-1.0,x); /* u <- (-1.0)x + u */
88: VecNorm(u,NORM_2,&norm);
89: if (norm > tol){
90: if (!rank){
91: PetscPrintf(PETSC_COMM_SELF,"Error: Norm of error %g, LU nfact %d\n",norm,nfact);
92: }
93: }
94: }
95: }
96: MatDestroy(C1);
97: MatDestroy(F);
99: /* Test Cholesky Factorization */
100: MatTranspose(C,MAT_INITIAL_MATRIX,&C1); /* C1 = C^T */
101: MatAXPY(C,1.0,C1,SAME_NONZERO_PATTERN); /* make C symmetric: C <- C + C^T */
102: MatShift(C,M); /* make C positive definite */
103: MatDestroy(C1);
104:
105: MatSetOption(C,MAT_SYMMETRIC,PETSC_TRUE);
106: MatSetOption(C,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);
108: MatGetFactor(C,MAT_SOLVER_PETSC,MAT_FACTOR_CHOLESKY,&F);
109: MatCholeskyFactorSymbolic(F,C,perm,&info);
110: for (nfact = 0; nfact < 2; nfact++){
111: if (!rank) printf(" Cholesky nfact %d\n",nfact);
112: MatCholeskyFactorNumeric(F,C,&info);
114: /* Test MatSolve() */
115: for (nsolve = 0; nsolve < 5; nsolve++){
116: VecGetArray(x,&array);
117: for (i=0; i<m; i++){
118: PetscRandomGetValue(rand,&rval);
119: array[i] = rval;
120: }
121: VecRestoreArray(x,&array);
122: VecCopy(x,u);
123: MatMult(C,x,b);
125: MatSolve(F,b,x);
127: /* Check the error */
128: VecAXPY(u,-1.0,x); /* u <- (-1.0)x + u */
129: VecNorm(u,NORM_2,&norm);
130: if (norm > tol){
131: if (!rank){
132: PetscPrintf(PETSC_COMM_SELF,"Error: Norm of error %g, Cholesky nfact %d\n",norm,nfact);
133: }
134: }
135: }
136: }
137: MatDestroy(F);
139: /* Free data structures */
140: PetscRandomDestroy(rand);
141: ISDestroy(perm);
142: ISDestroy(iperm);
143: VecDestroy(x);
144: VecDestroy(b);
145: VecDestroy(u);
146: MatDestroy(C);
148: PetscFinalize();
149: return 0;
150: }