Actual source code: ex145.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: static char help[] = "Tests LU, Cholesky factorization and MatMatSolve() for an Elemental dense matrix.\n\n";

  4: #include <petscmat.h>

  8: int main(int argc,char **argv)
  9: {
 10:   Mat            A,F,B,X,C,Aher,G;
 11:   Vec            b,x,c,d,e;
 13:   PetscInt       m = 5,n,p,i,j,nrows,ncols;
 14:   PetscScalar    *v,*barray,rval;
 15:   PetscReal      norm,tol=1.e-12;
 16:   PetscMPIInt    size,rank;
 17:   PetscRandom    rand;
 18:   const PetscInt *rows,*cols;
 19:   IS             isrows,iscols;
 20:   PetscBool      mats_view=PETSC_FALSE;
 21:   MatFactorInfo  finfo;

 23:   PetscInitialize(&argc,&argv,(char*) 0,help);
 24:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 25:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 27:   PetscRandomCreate(PETSC_COMM_WORLD,&rand);
 28:   PetscRandomSetFromOptions(rand);

 30:   /* Get local dimensions of matrices */
 31:   PetscOptionsGetInt(NULL,"-m",&m,NULL);
 32:   n    = m;
 33:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 34:   p    = m/2;
 35:   PetscOptionsGetInt(NULL,"-p",&p,NULL);
 36:   PetscOptionsHasName(NULL,"-mats_view",&mats_view);

 38:   /* Create matrix A */
 39:   PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix A\n");
 40:   MatCreate(PETSC_COMM_WORLD,&A);
 41:   MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);
 42:   MatSetType(A,MATELEMENTAL);
 43:   MatSetFromOptions(A);
 44:   MatSetUp(A);
 45:   /* Set local matrix entries */
 46:   MatGetOwnershipIS(A,&isrows,&iscols);
 47:   ISGetLocalSize(isrows,&nrows);
 48:   ISGetIndices(isrows,&rows);
 49:   ISGetLocalSize(iscols,&ncols);
 50:   ISGetIndices(iscols,&cols);
 51:   PetscMalloc1(nrows*ncols,&v);
 52:   for (i=0; i<nrows; i++) {
 53:     for (j=0; j<ncols; j++) {
 54:       PetscRandomGetValue(rand,&rval);
 55:       v[i*ncols+j] = rval;
 56:     }
 57:   }
 58:   MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES);
 59:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 60:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 61:   ISRestoreIndices(isrows,&rows);
 62:   ISRestoreIndices(iscols,&cols);
 63:   ISDestroy(&isrows);
 64:   ISDestroy(&iscols);
 65:   PetscFree(v);
 66:   if (mats_view) {
 67:     PetscPrintf(PETSC_COMM_WORLD, "A: nrows %d, m %d; ncols %d, n %d\n",nrows,m,ncols,n);
 68:     MatView(A,PETSC_VIEWER_STDOUT_WORLD);
 69:   }

 71:   /* Create rhs matrix B */
 72:   PetscPrintf(PETSC_COMM_WORLD," Create rhs matrix B\n");
 73:   MatCreate(PETSC_COMM_WORLD,&B);
 74:   MatSetSizes(B,m,p,PETSC_DECIDE,PETSC_DECIDE);
 75:   MatSetType(B,MATELEMENTAL);
 76:   MatSetFromOptions(B);
 77:   MatSetUp(B);
 78:   MatGetOwnershipIS(B,&isrows,&iscols);
 79:   ISGetLocalSize(isrows,&nrows);
 80:   ISGetIndices(isrows,&rows);
 81:   ISGetLocalSize(iscols,&ncols);
 82:   ISGetIndices(iscols,&cols);
 83:   PetscMalloc1(nrows*ncols,&v);
 84:   for (i=0; i<nrows; i++) {
 85:     for (j=0; j<ncols; j++) {
 86:       PetscRandomGetValue(rand,&rval);
 87:       v[i*ncols+j] = rval;
 88:     }
 89:   }
 90:   MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES);
 91:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 92:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 93:   ISRestoreIndices(isrows,&rows);
 94:   ISRestoreIndices(iscols,&cols);
 95:   ISDestroy(&isrows);
 96:   ISDestroy(&iscols);
 97:   PetscFree(v);
 98:   if (mats_view) {
 99:     PetscPrintf(PETSC_COMM_WORLD, "B: nrows %d, m %d; ncols %d, p %d\n",nrows,m,ncols,p);
100:     MatView(B,PETSC_VIEWER_STDOUT_WORLD);
101:   }

103:   /* Create rhs vector b and solution x (same size as b) */
104:   VecCreate(PETSC_COMM_WORLD,&b);
105:   VecSetSizes(b,m,PETSC_DECIDE);
106:   VecSetFromOptions(b);
107:   VecGetArray(b,&barray);
108:   for (j=0; j<m; j++) {
109:     PetscRandomGetValue(rand,&rval);
110:     barray[j] = rval;
111:   }
112:   VecRestoreArray(b,&barray);
113:   VecAssemblyBegin(b);
114:   VecAssemblyEnd(b);
115:   if (mats_view) {
116:     PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] b: m %d\n",rank,m);
117:     PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
118:     VecView(b,PETSC_VIEWER_STDOUT_WORLD);
119:   }
120:   VecDuplicate(b,&x);

122:   /* Create matrix X - same size as B */
123:   PetscPrintf(PETSC_COMM_WORLD," Create solution matrix X\n");
124:   MatCreate(PETSC_COMM_WORLD,&X);
125:   MatSetSizes(X,m,p,PETSC_DECIDE,PETSC_DECIDE);
126:   MatSetType(X,MATELEMENTAL);
127:   MatSetFromOptions(X);
128:   MatSetUp(X);
129:   MatAssemblyBegin(X,MAT_FINAL_ASSEMBLY);
130:   MatAssemblyEnd(X,MAT_FINAL_ASSEMBLY);

132:   /* Cholesky factorization */
133:   /*------------------------*/
134:   PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix Aher\n");
135:   MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&Aher);
136:   MatAXPY(Aher,1.0,A,SAME_NONZERO_PATTERN); /* Aher = A + A^T */
137:   if (!rank) { /* add 100.0 to diagonals of Aher to make it spd */
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:   MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);
248:   MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN);
249:   MatNorm(C,NORM_1,&norm);
250:   if (norm > tol) {
251:     PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*X - B| for LU %g\n",(double)norm);
252:   }

254:   /* Out-place LU */
255:   MatGetFactor(A,MATSOLVERELEMENTAL,MAT_FACTOR_LU,&F);
256:   MatLUFactorSymbolic(F,A,0,0,&finfo);
257:   MatLUFactorNumeric(F,A,&finfo);
258:   if (mats_view) {
259:     MatView(F,PETSC_VIEWER_STDOUT_WORLD);
260:   }
261:   MatSolve(F,b,x);
262:   MatMatSolve(F,B,X);
263:   MatDestroy(&F);

265:   /* Free space */
266:   MatDestroy(&A);
267:   MatDestroy(&Aher);
268:   MatDestroy(&B);
269:   MatDestroy(&C);
270:   MatDestroy(&X);
271:   VecDestroy(&b);
272:   VecDestroy(&c);
273:   VecDestroy(&d);
274:   VecDestroy(&e);
275:   VecDestroy(&x);
276:   PetscRandomDestroy(&rand);
277:   PetscFinalize();
278:   return 0;
279: }