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;
17: PetscTruth flg;
18: IS perm,iperm;
19: MatFactorInfo info;
20: PetscRandom rand;
22: PetscInitialize(&argc,&args,(char *)0,help);
23: MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
24: MPI_Comm_size(PETSC_COMM_WORLD, &nproc);
26: #ifdef PETSC_HAVE_PLAPACK
27: /* Create matrix and vectors */
28: PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
29: MatCreate(PETSC_COMM_WORLD,&C);
30: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,M,M);
31: MatSetType(C,MATPLAPACK);
32: MatSetFromOptions(C);
33:
34: MatGetLocalSize(C,&m,&n);
35: if (m != n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Matrix local size m %d must equal n %d",m,n);
37: VecCreate(PETSC_COMM_WORLD,&x);
38: VecSetSizes(x,n,PETSC_DECIDE);
39: VecSetFromOptions(x);
40: VecDuplicate(x,&b);
41: VecDuplicate(x,&u); /* save the true solution */
43: /* Assembly */
44: PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&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);
59: /* Test LU Factorization */
60: MatGetOrdering(C1,MATORDERING_NATURAL,&perm,&iperm);
61: MatLUFactorSymbolic(C1,perm,iperm,&info,&F);
62: for (nfact = 0; nfact < 2; nfact++){
63: if (!rank) printf(" LU nfact %d\n",nfact);
64: MatLUFactorNumeric(C1,&info,&F);
66: /* Test MatSolve() */
67: for (nsolve = 0; nsolve < 5; nsolve++){
68: VecGetArray(x,&array);
69: for (i=0; i<m; i++){
70: PetscRandomGetValue(rand,&rval);
71: array[i] = rval;
72: }
73: VecRestoreArray(x,&array);
74: VecCopy(x,u);
75: MatMult(C,x,b);
77: MatSolve(F,b,x);
79: /* Check the error */
80: VecAXPY(u,-1.0,x); /* u <- (-1.0)x + u */
81: VecNorm(u,NORM_2,&norm);
82: if (!rank){
83: PetscPrintf(PETSC_COMM_SELF,"Norm of error %A\n",norm);
84: }
85: }
86: }
87: MatDestroy(C1);
88: MatDestroy(F);
90: /* Test Cholesky Factorization */
91: MatTranspose(C,&C1); /* C1 = C^T */
92: MatAXPY(C,1.0,C1,SAME_NONZERO_PATTERN); /* make C symmetric: C <- C + C^T */
93: MatShift(C,M); /* make C positive definite */
94: MatDestroy(C1);
95:
96: MatSetOption(C,MAT_SYMMETRIC);
97: MatSetOption(C,MAT_SYMMETRY_ETERNAL);
99: MatDuplicate(C,MAT_COPY_VALUES,&C1);
100: MatCholeskyFactorSymbolic(C,perm,&info,&F);
101: for (nfact = 0; nfact < 2; nfact++){
102: if (!rank) printf(" Cholesky nfact %d\n",nfact);
103: MatCholeskyFactorNumeric(C1,&info,&F);
105: /* Test MatSolve() */
106: for (nsolve = 0; nsolve < 5; nsolve++){
107: VecGetArray(x,&array);
108: for (i=0; i<m; i++){
109: PetscRandomGetValue(rand,&rval);
110: array[i] = rval;
111: }
112: VecRestoreArray(x,&array);
113: VecCopy(x,u);
114: MatMult(C,x,b);
116: MatSolve(F,b,x);
118: /* Check the error */
119: VecAXPY(u,-1.0,x); /* u <- (-1.0)x + u */
120: VecNorm(u,NORM_2,&norm);
121: if (!rank){
122: PetscPrintf(PETSC_COMM_SELF,"Norm of error %A\n",norm);
123: }
124: }
125: }
126: MatDestroy(C1);
127: MatDestroy(F);
129: /* Free data structures */
130: PetscRandomDestroy(rand);
131: ISDestroy(perm);
132: ISDestroy(iperm);
133: VecDestroy(x);
134: VecDestroy(b);
135: VecDestroy(u);
136: MatDestroy(C);
138: #else
139: if (!rank) printf("This example needs PLAPLAPACK\n");
140: #endif
141: PetscFinalize();
142: return 0;
143: }