Actual source code: mlCheck.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: mlCheck.c,v 1.2 2000/01/10 03:20:33 knepley Exp $";
3: #endif
4: /*
5: Defines error checking routines for the multilevel preconditioner
6: */
7: #include src/sles/pc/pcimpl.h
8: #include ml.h
10: #undef __FUNCT__
12: /*
13: PCValidQ_Multilevel - Validates the ML data structure
15: Collective on PC
17: Input Parameter:
18: . pc - The PC
20: Level: advanced
22: .keywords PC, Valid
23: .seealso PCDebug_Multilevel()
24: */
25: int PCValidQ_Multilevel(PC pc)
26: {
27: PC_Multilevel *ml = (PC_Multilevel *) pc->data;
28: int numRows, numCols;
29: int level, part;
30: int ierr;
33: if (pc->setupcalled < 2)
34: PetscFunctionReturn(1);
37: /* Check dimensions */
38: if (ml->locB != PETSC_NULL) {
39: MatGetSize(ml->locB, &numRows, &numCols);
40: if ((ml->numRows != numRows) || (ml->numCols != numCols)) {
41: PetscLogInfo(pc, "PCValidQ_Multilevel: Invalid dimensions (%d,%d) for factorization", numRows, numCols);
42: PetscFunctionReturn(1);
43: }
44: }
45: if (ml->numLevels < 0) {
46: PetscLogInfo(pc, "PCValidQ_Multilevel: Invalid number of levels %d for factorization", ml->numLevels);
47: PetscFunctionReturn(1);
48: }
49: /* Check thresholds */
50: if (ml->QRthresh < 0) {
51: PetscLogInfo(pc, "PCValidQ_Multilevel: Invalid threshold %d for final QR factorization", ml->QRthresh);
52: PetscFunctionReturn(1);
53: }
54: if (ml->zeroTol < PETSC_MACHINE_EPSILON) {
55: PetscLogInfo(pc, "PCValidQ_Multilevel: Numeric tolerance %g less than machine epsilon", ml->zeroTol);
56: PetscFunctionReturn(1);
57: }
59: /* Check meshes */
60: if (ml->numMeshes != ml->numLevels + 1)
61: PetscLogInfo(pc, "PCValidQ_Multilevel: Invalid number %d of hierarchical meshes", ml->numMeshes);
65: for(level = 0; level <= ml->numLevels; level++) {
70: }
72: /* Bail out if no factorization was done */
73: if (ml->numLevels == 0)
74: return(0);
76: /* Check numbering */
82: for(level = 0; level < ml->numLevels; level++) {
83: if (ml->numPartitions[level] == 0)
84: continue;
92: if (ml->numPartitionRows[level][ml->numPartitions[level]] > 0)
94: if (ml->numPartitionRows[level][ml->numPartitions[level]+1] > 0)
96: for(part = 0; part < ml->numPartitions[level]; part++) {
97: if (ml->numPartitionCols[level][part] > 0)
99: if (ml->numPartitionRows[level][part] > 0)
101: }
102: }
104: /* Check factorization */
106: for(level = 0; level < ml->numLevels; level++) {
107: if (ml->numPartitions[level] == 0)
108: continue;
110: for(part = 0; part < ml->numPartitions[level]; part++) {
111: numRows = ml->numPartitionRows[level][part];
112: numCols = ml->numPartitionCols[level][part];
114: if (numRows > 0) {
116: if (PCMultiLevelDoQR_Private(pc, numRows, numCols) == PETSC_TRUE) {
119: }
120: }
121: if (numCols > 0) {
124: }
125: }
126: }
128: /* Check boundary gradients */
135: for(level = 0; level < ml->numLevels; level++) {
136: if (ml->numPartitions[level] == 0)
137: continue;
142: }
144: return(0);
145: }
147: #undef __FUNCT__
149: int PCDebugPrintMat_Multilevel_Private(PC pc, VarOrdering rowOrder, VarOrdering colOrder, int (*applyM)(PC, GVec, GVec))
150: {
151: PC_Multilevel *ml = (PC_Multilevel *) pc->data;
152: Grid grid = ml->grid;
153: PetscScalar zero = 0.0;
154: GVec rowVec, colVec;
155: PetscScalar *rowArray, *colArray;
156: int rank;
157: int tempLevels, colVars, locRowVars;
158: int *firstCol;
159: int row, col;
160: PetscTruth opt;
161: int ierr;
164: /* Setup storage */
165: MPI_Comm_rank(pc->comm, &rank);
166: GVecCreateRectangular(grid, rowOrder, &rowVec);
167: GVecCreateRectangular(grid, colOrder, &colVec);
168: VecGetArray(rowVec, &rowArray);
169: VecGetArray(colVec, &colArray);
171: /* Allow only the first k levels to be used */
172: tempLevels = ml->numLevels;
173: PetscOptionsGetInt(pc->prefix, "-pc_ml_debug_level", &tempLevels, &opt);
175: /* Print out the matrix M */
176: VarOrderingGetSize(colOrder, &colVars);
177: VarOrderingGetLocalSize(rowOrder, &locRowVars);
178: VarOrderingGetFirst(colOrder, &firstCol);
179: for(col = 0; col < colVars; col++)
180: {
181: /* Start with e_col */
182: VecSet(&zero, colVec);
183: if ((col >= firstCol[rank]) && (col < firstCol[rank+1]))
184: colArray[col-firstCol[rank]] = 1.0;
185: /* M_col = M e_col */
186: (*applyM)(pc, colVec, rowVec);
187: /* Print out as a row so we get the transpose */
188: for(row = 0; row < locRowVars; row++)
189: PetscSynchronizedPrintf(pc->comm, "%g ", rowArray[row]);
190: PetscSynchronizedFlush(pc->comm);
191: PetscPrintf(pc->comm, "n");
192: }
194: /* Cleanup stroage */
195: VecRestoreArray(rowVec, &rowArray);
196: VecRestoreArray(colVec, &colArray);
197: GVecDestroy(rowVec);
198: GVecDestroy(colVec);
199: return(0);
200: }
202: #undef __FUNCT__
204: /*
205: PCDebugPrintMat_Multilevel - Prints out the dense form of the level matrices.
207: Collective on PC
209: Input Parameter:
210: . pc - The PC
212: Level: advanced
214: .keywords PC, debug
215: .seealso PCDebug_Multilevel()
216: */
217: int PCDebugPrintMat_Multilevel(PC pc, const char mat[])
218: {
219: PC_Multilevel *ml = (PC_Multilevel *) pc->data;
220: PetscTruth isP, isB, isV, isD, isall;
221: int ierr;
224: PetscStrcasecmp(mat, "P", &isP);
225: PetscStrcasecmp(mat, "B", &isB);
226: PetscStrcasecmp(mat, "V", &isV);
227: PetscStrcasecmp(mat, "D", &isD);
228: PetscStrcasecmp(mat, "all", &isall);
229: if (isP == PETSC_TRUE) {
230: PetscPrintf(pc->comm, "P^Tn");
231: PCDebugPrintMat_Multilevel_Private(pc, ml->tOrder, ml->tOrder, PCMultiLevelApplyP);
232: } else if (isB == PETSC_TRUE) {
233: PetscPrintf(pc->comm, "B^Tn");
234: PCDebugPrintMat_Multilevel_Private(pc, ml->tOrder, ml->sOrder, PCMultiLevelApplyGradient);
235: } else if (isV == PETSC_TRUE) {
236: PetscPrintf(PETSC_COMM_WORLD, "V^Tn");
237: PCDebugPrintMat_Multilevel_Private(pc, ml->sOrder, ml->sOrder, PCMultiLevelApplyV);
238: } else if (isD == PETSC_TRUE) {
239: PetscPrintf(PETSC_COMM_WORLD, "D^{-T}n");
240: PCDebugPrintMat_Multilevel_Private(pc, ml->sOrder, ml->sOrder, PCMultiLevelApplyDInv);
241: } else if (isall == PETSC_TRUE) {
242: PCDebugPrintMat_Multilevel_Private(pc, ml->tOrder, ml->tOrder, PCMultiLevelApplyP);
243: PCDebugPrintMat_Multilevel_Private(pc, ml->tOrder, ml->sOrder, PCMultiLevelApplyGradient);
244: PCDebugPrintMat_Multilevel_Private(pc, ml->sOrder, ml->sOrder, PCMultiLevelApplyV);
245: PCDebugPrintMat_Multilevel_Private(pc, ml->sOrder, ml->sOrder, PCMultiLevelApplyDInv);
246: } else {
247: SETERRQ1(PETSC_ERR_ARG_WRONG, "Unknown ML level matrix %s.", mat);
248: }
249: return(0);
250: }
252: #undef __FUNCT__
254: int PCDebug_Multilevel(PC pc)
255: {
256: PC_Multilevel *ml = (PC_Multilevel *) pc->data;
257: Grid grid = ml->grid;
258: PetscScalar zero = 0.0;
259: GVec testVec, testVec2, testVec3;
260: PetscScalar *testArray, *testArray2, *testArray3;
261: GVec singVec;
262: PetscScalar *singArray;
263: PetscScalar norm;
264: char printMat[256];
265: int rank;
266: int range;
267: int rowVars, colVars, locColVars;
268: int *firstRow, *firstCol;
269: int singCount, singRow;
270: int locSingCount, locSingRow;
271: int row, col, locCol;
272: PetscTruth opt;
273: int ierr;
276: /* Initialize testing */
277: MPI_Comm_rank(pc->comm, &rank);
278: GVecCreateRectangular(grid, ml->sOrder, &testVec);
279: GVecDuplicate(testVec, &singVec);
280: GVecCreateConstrained(grid, &testVec2);
281: GVecCreateConstrained(grid, &testVec3);
282: VecGetArray(testVec, &testArray);
283: VecGetArray(testVec2, &testArray2);
284: VecGetArray(testVec3, &testArray3);
285: VecGetArray(singVec, &singArray);
287: /* Print matrices */
288: PetscOptionsGetString(pc->prefix, "-pc_ml_print_mat", printMat, 255, &opt);
289: if (opt == PETSC_TRUE) {
290: PCDebugPrintMat_Multilevel(pc, printMat);
291: }
293: VarOrderingGetSize(ml->tOrder, &rowVars);
294: VarOrderingGetSize(ml->sOrder, &colVars);
295: VarOrderingGetLocalSize(ml->sOrder, &locColVars);
296: VarOrderingGetFirst(ml->tOrder, &firstRow);
297: VarOrderingGetFirst(ml->sOrder, &firstCol);
299: /* Check that B^T P = 0 for numNullSpace rows */
300: for(row = 0, range = 0; row < rowVars; row++)
301: {
302: /* Start with e_row */
303: VecSet(&zero, testVec2);
304: if ((row >= firstRow[rank]) && (row < firstRow[rank+1]))
305: testArray2[row-firstRow[rank]] = 1.0;
306: /* P e_row */
307: PCMultiLevelApplyP(pc, testVec2, testVec2);
308: /* B^T P e_row */
309: PCMultiLevelApplyGradientTrans(pc, testVec2, testVec);
310: /* Check for nonzeros only in the range */
311: VecNorm(testVec, NORM_2, &norm);
312: if (norm > ml->zeroTol)
313: range++;
314: }
315: if (range != ml->globalRank) {
316: PetscLogInfo(pc, "PCDebug_Multilevel: P_2 is not a null space for B^T, %d != %d range vectorsn", range, ml->globalRank);
317: PetscFunctionReturn(1);
318: }
320: /* Test to see that D^{-1} P^T_1 B Z = I */
321: for(col = 0; col < colVars; col++)
322: {
323: /* Start with e_col */
324: locCol = col - firstCol[rank];
325: VecSet(&zero, testVec);
326: if ((col >= firstCol[rank]) && (col < firstCol[rank+1]))
327: testArray[locCol] = 1.0;
329: /* Z e_col */
330: PCMultiLevelApplyV(pc, testVec, testVec);
331: /* B Z e_col */
332: PCMultiLevelApplyGradient(pc, testVec, testVec2);
333: /* P^T B Z e_col */
334: PCMultiLevelApplyPTrans(pc, testVec2, testVec2);
335: /* Scatter to a column vector */
336: VecScatterBegin(testVec2, testVec, INSERT_VALUES, SCATTER_FORWARD, ml->rangeScatter);
337: VecScatterEnd(testVec2, testVec, INSERT_VALUES, SCATTER_FORWARD, ml->rangeScatter);
338: /* D^{-1} P^T B Z e_col */
339: PCMultiLevelApplyDInv(pc, testVec, testVec);
340: /* Check the row */
341: for(row = 0, locSingCount = 0, locSingRow = -1; row < locColVars; row++)
342: if (testArray[row] > ml->zeroTol) {
343: locSingCount++;
344: locSingRow = row;
345: }
346: MPI_Allreduce(&locSingCount, &singCount, 1, MPI_INT, MPI_SUM, pc->comm);
347: MPI_Allreduce(&locSingRow, &singRow, 1, MPI_INT, MPI_MAX, pc->comm);
348: if (singCount > 1) {
349: PetscLogInfo(pc, "PCDebug_Multilevel: Invalid column %d in P^T B Z, %d nonzerosn", col, singCount);
350: PetscFunctionReturn(1);
351: } else if (singCount == 0) {
352: PetscLogInfo(pc, "PCDebug_Multilevel: B is rank deficient in column %dn", col);
353: }
354: /* Check the singular value */
355: if (locSingRow > -1)
356: {
357: singRow += firstCol[rank];
358: if (singRow != col) {
359: PetscLogInfo(pc, "PCDebug_Multilevel: Invalid ordering in P^T B Z, value in column %d and row %dn", col, singRow);
360: PetscFunctionReturn(1);
361: }
362: if (PetscAbsScalar(testArray[locSingRow] - 1.0) > ml->zeroTol) {
363: PetscLogInfo(pc, "PCDebug_Multilevel: Invalid singular value in column %dn", col);
364: PetscFunctionReturn(1);
365: }
366: }
367: }
369: /* Cleanup testing */
370: VecRestoreArray(testVec, &testArray);
371: VecRestoreArray(testVec2, &testArray2);
372: VecRestoreArray(testVec3, &testArray3);
373: VecRestoreArray(singVec, &singArray);
374: GVecDestroy(testVec);
375: GVecDestroy(testVec2);
376: GVecDestroy(testVec3);
377: GVecDestroy(singVec);
379: return(0);
380: }