Actual source code: ex145.c
2: static char help[] = "Tests LU, Cholesky factorization and MatMatSolve() for an Elemental dense matrix.\n\n";
4: #include <petscmat.h>
6: int main(int argc,char **argv)
7: {
8: Mat A,F,B,X,C,Aher,G;
9: Vec b,x,c,d,e;
10: PetscInt m = 5,n,p,i,j,nrows,ncols;
11: PetscScalar *v,*barray,rval;
12: PetscReal norm,tol=1.e-11;
13: PetscMPIInt size,rank;
14: PetscRandom rand;
15: const PetscInt *rows,*cols;
16: IS isrows,iscols;
17: PetscBool mats_view=PETSC_FALSE;
18: MatFactorInfo finfo;
20: PetscInitialize(&argc,&argv,(char*) 0,help);
21: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
22: MPI_Comm_size(PETSC_COMM_WORLD,&size);
24: PetscRandomCreate(PETSC_COMM_WORLD,&rand);
25: PetscRandomSetFromOptions(rand);
27: /* Get local dimensions of matrices */
28: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
29: n = m;
30: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
31: p = m/2;
32: PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL);
33: PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view);
35: /* Create matrix A */
36: PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix A\n");
37: MatCreate(PETSC_COMM_WORLD,&A);
38: MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);
39: MatSetType(A,MATELEMENTAL);
40: MatSetFromOptions(A);
41: MatSetUp(A);
42: /* Set local matrix entries */
43: MatGetOwnershipIS(A,&isrows,&iscols);
44: ISGetLocalSize(isrows,&nrows);
45: ISGetIndices(isrows,&rows);
46: ISGetLocalSize(iscols,&ncols);
47: ISGetIndices(iscols,&cols);
48: PetscMalloc1(nrows*ncols,&v);
49: for (i=0; i<nrows; i++) {
50: for (j=0; j<ncols; j++) {
51: PetscRandomGetValue(rand,&rval);
52: v[i*ncols+j] = rval;
53: }
54: }
55: MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES);
56: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
57: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
58: ISRestoreIndices(isrows,&rows);
59: ISRestoreIndices(iscols,&cols);
60: ISDestroy(&isrows);
61: ISDestroy(&iscols);
62: PetscFree(v);
63: if (mats_view) {
64: PetscPrintf(PETSC_COMM_WORLD, "A: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", n %" PetscInt_FMT "\n",nrows,m,ncols,n);
65: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
66: }
68: /* Create rhs matrix B */
69: PetscPrintf(PETSC_COMM_WORLD," Create rhs matrix B\n");
70: MatCreate(PETSC_COMM_WORLD,&B);
71: MatSetSizes(B,m,p,PETSC_DECIDE,PETSC_DECIDE);
72: MatSetType(B,MATELEMENTAL);
73: MatSetFromOptions(B);
74: MatSetUp(B);
75: MatGetOwnershipIS(B,&isrows,&iscols);
76: ISGetLocalSize(isrows,&nrows);
77: ISGetIndices(isrows,&rows);
78: ISGetLocalSize(iscols,&ncols);
79: ISGetIndices(iscols,&cols);
80: PetscMalloc1(nrows*ncols,&v);
81: for (i=0; i<nrows; i++) {
82: for (j=0; j<ncols; j++) {
83: PetscRandomGetValue(rand,&rval);
84: v[i*ncols+j] = rval;
85: }
86: }
87: MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES);
88: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
89: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
90: ISRestoreIndices(isrows,&rows);
91: ISRestoreIndices(iscols,&cols);
92: ISDestroy(&isrows);
93: ISDestroy(&iscols);
94: PetscFree(v);
95: if (mats_view) {
96: PetscPrintf(PETSC_COMM_WORLD, "B: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", p %" PetscInt_FMT "\n",nrows,m,ncols,p);
97: MatView(B,PETSC_VIEWER_STDOUT_WORLD);
98: }
100: /* Create rhs vector b and solution x (same size as b) */
101: VecCreate(PETSC_COMM_WORLD,&b);
102: VecSetSizes(b,m,PETSC_DECIDE);
103: VecSetFromOptions(b);
104: VecGetArray(b,&barray);
105: for (j=0; j<m; j++) {
106: PetscRandomGetValue(rand,&rval);
107: barray[j] = rval;
108: }
109: VecRestoreArray(b,&barray);
110: VecAssemblyBegin(b);
111: VecAssemblyEnd(b);
112: if (mats_view) {
113: PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] b: m %" PetscInt_FMT "\n",rank,m);
114: PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
115: VecView(b,PETSC_VIEWER_STDOUT_WORLD);
116: }
117: VecDuplicate(b,&x);
119: /* Create matrix X - same size as B */
120: PetscPrintf(PETSC_COMM_WORLD," Create solution matrix X\n");
121: MatCreate(PETSC_COMM_WORLD,&X);
122: MatSetSizes(X,m,p,PETSC_DECIDE,PETSC_DECIDE);
123: MatSetType(X,MATELEMENTAL);
124: MatSetFromOptions(X);
125: MatSetUp(X);
126: MatAssemblyBegin(X,MAT_FINAL_ASSEMBLY);
127: MatAssemblyEnd(X,MAT_FINAL_ASSEMBLY);
129: /* Cholesky factorization */
130: /*------------------------*/
131: PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix Aher\n");
132: MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&Aher);
133: MatAXPY(Aher,1.0,A,SAME_NONZERO_PATTERN); /* Aher = A + A^T */
134: if (rank == 0) { /* add 100.0 to diagonals of Aher to make it spd */
136: /* TODO: Replace this with a call to El::ShiftDiagonal( A, 100.),
137: or at least pre-allocate the right amount of space */
138: PetscInt M,N;
139: MatGetSize(Aher,&M,&N);
140: for (i=0; i<M; i++) {
141: rval = 100.0;
142: MatSetValues(Aher,1,&i,1,&i,&rval,ADD_VALUES);
143: }
144: }
145: MatAssemblyBegin(Aher,MAT_FINAL_ASSEMBLY);
146: MatAssemblyEnd(Aher,MAT_FINAL_ASSEMBLY);
147: if (mats_view) {
148: PetscPrintf(PETSC_COMM_WORLD, "Aher:\n");
149: MatView(Aher,PETSC_VIEWER_STDOUT_WORLD);
150: }
152: /* Cholesky factorization */
153: /*------------------------*/
154: PetscPrintf(PETSC_COMM_WORLD," Test Cholesky Solver \n");
155: /* In-place Cholesky */
156: /* Create matrix factor G, then copy Aher to G */
157: MatCreate(PETSC_COMM_WORLD,&G);
158: MatSetSizes(G,m,n,PETSC_DECIDE,PETSC_DECIDE);
159: MatSetType(G,MATELEMENTAL);
160: MatSetFromOptions(G);
161: MatSetUp(G);
162: MatAssemblyBegin(G,MAT_FINAL_ASSEMBLY);
163: MatAssemblyEnd(G,MAT_FINAL_ASSEMBLY);
164: MatCopy(Aher,G,SAME_NONZERO_PATTERN);
166: /* Only G = U^T * U is implemented for now */
167: MatCholeskyFactor(G,0,0);
168: if (mats_view) {
169: PetscPrintf(PETSC_COMM_WORLD, "Cholesky Factor G:\n");
170: MatView(G,PETSC_VIEWER_STDOUT_WORLD);
171: }
173: /* Solve U^T * U x = b and U^T * U X = B */
174: MatSolve(G,b,x);
175: MatMatSolve(G,B,X);
176: MatDestroy(&G);
178: /* Out-place Cholesky */
179: MatGetFactor(Aher,MATSOLVERELEMENTAL,MAT_FACTOR_CHOLESKY,&G);
180: MatCholeskyFactorSymbolic(G,Aher,0,&finfo);
181: MatCholeskyFactorNumeric(G,Aher,&finfo);
182: if (mats_view) {
183: MatView(G,PETSC_VIEWER_STDOUT_WORLD);
184: }
185: MatSolve(G,b,x);
186: MatMatSolve(G,B,X);
187: MatDestroy(&G);
189: /* Check norm(Aher*x - b) */
190: VecCreate(PETSC_COMM_WORLD,&c);
191: VecSetSizes(c,m,PETSC_DECIDE);
192: VecSetFromOptions(c);
193: MatMult(Aher,x,c);
194: VecAXPY(c,-1.0,b);
195: VecNorm(c,NORM_1,&norm);
196: if (norm > tol) {
197: PetscPrintf(PETSC_COMM_WORLD,"Warning: |Aher*x - b| for Cholesky %g\n",(double)norm);
198: }
200: /* Check norm(Aher*X - B) */
201: MatMatMult(Aher,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);
202: MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN);
203: MatNorm(C,NORM_1,&norm);
204: if (norm > tol) {
205: PetscPrintf(PETSC_COMM_WORLD,"Warning: |Aher*X - B| for Cholesky %g\n",(double)norm);
206: }
208: /* LU factorization */
209: /*------------------*/
210: PetscPrintf(PETSC_COMM_WORLD," Test LU Solver \n");
211: /* In-place LU */
212: /* Create matrix factor F, then copy A to F */
213: MatCreate(PETSC_COMM_WORLD,&F);
214: MatSetSizes(F,m,n,PETSC_DECIDE,PETSC_DECIDE);
215: MatSetType(F,MATELEMENTAL);
216: MatSetFromOptions(F);
217: MatSetUp(F);
218: MatAssemblyBegin(F,MAT_FINAL_ASSEMBLY);
219: MatAssemblyEnd(F,MAT_FINAL_ASSEMBLY);
220: MatCopy(A,F,SAME_NONZERO_PATTERN);
221: /* Create vector d to test MatSolveAdd() */
222: VecDuplicate(x,&d);
223: VecCopy(x,d);
225: /* PF=LU or F=LU factorization - perms is ignored by Elemental;
226: set finfo.dtcol !0 or 0 to enable/disable partial pivoting */
227: finfo.dtcol = 0.1;
228: MatLUFactor(F,0,0,&finfo);
230: /* Solve LUX = PB or LUX = B */
231: MatSolveAdd(F,b,d,x);
232: MatMatSolve(F,B,X);
233: MatDestroy(&F);
235: /* Check norm(A*X - B) */
236: VecCreate(PETSC_COMM_WORLD,&e);
237: VecSetSizes(e,m,PETSC_DECIDE);
238: VecSetFromOptions(e);
239: MatMult(A,x,c);
240: MatMult(A,d,e);
241: VecAXPY(c,-1.0,e);
242: VecAXPY(c,-1.0,b);
243: VecNorm(c,NORM_1,&norm);
244: if (norm > tol) {
245: PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*x - b| for LU %g\n",(double)norm);
246: }
247: /* Reuse product C; replace Aher with A */
248: MatProductReplaceMats(A,NULL,NULL,C);
249: MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);
250: MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN);
251: MatNorm(C,NORM_1,&norm);
252: if (norm > tol) {
253: PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*X - B| for LU %g\n",(double)norm);
254: }
256: /* Out-place LU */
257: MatGetFactor(A,MATSOLVERELEMENTAL,MAT_FACTOR_LU,&F);
258: MatLUFactorSymbolic(F,A,0,0,&finfo);
259: MatLUFactorNumeric(F,A,&finfo);
260: if (mats_view) {
261: MatView(F,PETSC_VIEWER_STDOUT_WORLD);
262: }
263: MatSolve(F,b,x);
264: MatMatSolve(F,B,X);
265: MatDestroy(&F);
267: /* Free space */
268: MatDestroy(&A);
269: MatDestroy(&Aher);
270: MatDestroy(&B);
271: MatDestroy(&C);
272: MatDestroy(&X);
273: VecDestroy(&b);
274: VecDestroy(&c);
275: VecDestroy(&d);
276: VecDestroy(&e);
277: VecDestroy(&x);
278: PetscRandomDestroy(&rand);
279: PetscFinalize();
280: return 0;
281: }
283: /*TEST
285: build:
286: requires: elemental
288: test:
289: nsize: 2
290: output_file: output/ex145.out
292: test:
293: suffix: 2
294: nsize: 6
295: output_file: output/ex145.out
297: TEST*/