Actual source code: gsnes.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: gsnes.c,v 1.11 2000/07/16 05:40:27 knepley Exp $";
3: #endif
5: #include src/snes/snesimpl.h
6: #include src/gvec/gvecimpl.h
7: #include gsolver.h
9: /*
10: This file provides routines for grid SNES objects. They support solution of
11: systems of nonlinear equations generated from a mesh and discretization
12: specified by a Grid object.
13: */
15: int GSNESCreateStructures_Private(GSNES snes, void *funcCtx, void *jacCtx)
16: {
17: Grid grid;
18: GMat J;
19: GVec f;
20: PetscTruth isMatrixFree;
21: int ierr;
24: GSNESGetGrid(snes, &grid);
25: GVecCreateConstrained(grid, &f);
26: GridIsMatrixFree(grid, &isMatrixFree);
27: if (isMatrixFree == PETSC_TRUE) {
28: if (grid->explicitConstraints == PETSC_TRUE) {
29: GMatCreateMF(grid, grid->constraintOrder, grid->constraintOrder, &J);
30: } else {
31: GMatCreateMF(grid, grid->order, grid->order, &J);
32: }
33: } else {
34: if (grid->explicitConstraints == PETSC_TRUE) {
35: GMatCreateRectangular(grid, grid->constraintOrder, grid->constraintOrder, &J);
36: } else {
37: GMatCreate(grid, &J);
38: }
39: }
40: SNESSetFunction(snes, f, (int (*)(SNES, Vec, Vec, void *)) GSNESEvaluateRhs, funcCtx);
41: SNESSetJacobian(snes, J, J, (int (*)(SNES, Vec, Mat *, Mat *, MatStructure *, void *)) GSNESEvaluateJacobian, jacCtx);
42: PetscLogObjectParent(grid, f);
43: PetscLogObjectParent(grid, J);
44: return(0);
45: }
47: int GSNESDestroyStructures_Private(GSNES snes)
48: {
49: GMat J;
50: GVec f;
51: int ierr;
54: SNESGetFunction(snes, &f, PETSC_NULL, PETSC_NULL);
55: SNESGetJacobian(snes, &J, PETSC_NULL, PETSC_NULL, PETSC_NULL);
56: GMatDestroy(J);
57: GVecDestroy(f);
58: return(0);
59: }
61: /*@C
62: GSNESDestroy - Destroys a grid SNES.
64: Collective on GSNES
66: Input Parameter:
67: . snes - The nonlinear solver context
69: Level: beginner
71: .keywords: grid SNES, destroy
72: .seealso: SNESDestroy(), GSNESDuplicate()
73: @*/
74: int GSNESDestroy(GSNES snes)
75: {
80: if (--snes->refct > 0)
81: return(0);
82: GSNESDestroyStructures_Private(snes);
83: SNESDestroy(snes);
84: return(0);
85: }
87: /*@
88: GSNESView - Views a grid SNES.
90: Collective on GSNES
92: Input Parameters:
93: + snes - The grid SNES
94: - viewer - [Optional] A visualization context
96: Options Database Key:
97: $ -snes_view : calls SNESView() at end of SNESSolve()
99: Notes:
100: The available visualization contexts include
101: $ VIEWER_STDOUT_SELF - standard output (default)
102: $ VIEWER_STDOUT_WORLD - synchronized standard
103: $ output where only the first processor opens
104: $ the file. All other processors send their
105: $ data to the first processor to print.
107: The user can open alternative vistualization contexts with
108: $ PetscViewerFileOpenASCII() - output vector to a specified file
110: Level: beginner
112: .keywords: grid SNES, view, visualize
113: .seealso: SNESView()
114: @*/
115: int GSNESView(GSNES snes, PetscViewer viewer)
116: {
117: Grid grid;
118: FILE *fd;
119: int numActiveFields;
120: int field, func, op;
121: PetscTruth isascii;
122: int ierr;
126: if (!viewer) {
127: viewer = PETSC_VIEWER_STDOUT_SELF;
128: } else {
130: }
132: GSNESGetGrid(snes, &grid);
133: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &isascii);
134: if (isascii == PETSC_TRUE) {
135: GridView(grid, viewer);
136: PetscViewerFlush(viewer);
137: PetscViewerASCIIGetPointer(viewer, &fd);
138: PetscFPrintf(grid->comm, fd, "GSNES Object:n");
139: GridGetNumActiveFields(grid, &numActiveFields);
140: if (numActiveFields == 1) {
141: PetscFPrintf(snes->comm, fd, " %d active field:n", numActiveFields);
142: } else {
143: PetscFPrintf(snes->comm, fd, " %d active fields:n", numActiveFields);
144: }
145: for(field = 0; field < grid->numFields; field++) {
146: if (grid->fields[field].isActive == PETSC_FALSE) continue;
147: if (grid->fields[field].name) {
148: PetscFPrintf(grid->comm, fd, " %s field with %d components:n", grid->fields[field].name, grid->fields[field].numComp);
149: } else {
150: PetscFPrintf(grid->comm, fd, " field %d with %d components:n", field, grid->fields[field].numComp);
151: }
152: for(func = 0; func < grid->numRhsFuncs; func++) {
153: if (grid->rhsFuncs[func].field == field) {
154: PetscFPrintf(grid->comm, fd, " Rhs function with scalar multiplier %gn", PetscRealPart(grid->rhsFuncs[func].alpha));
155: }
156: }
157: for(op = 0; op < grid->numRhsOps; op++) {
158: if (grid->rhsOps[op].field == field) {
159: if (grid->rhsOps[op].nonlinearOp != PETSC_NULL) {
160: PetscFPrintf(grid->comm, fd, " Rhs nonlinear operator with scalar multiplier %gn",
161: PetscRealPart(grid->rhsOps[field].alpha));
162: } else {
163: PetscFPrintf(grid->comm, fd, " Rhs operator %d with scalar multiplier %gn",
164: grid->rhsOps[op].op, PetscRealPart(grid->rhsOps[op].alpha));
165: }
166: }
167: }
168: for(op = 0; op < grid->numMatOps; op++) {
169: if (grid->matOps[op].field == field) {
170: PetscFPrintf(grid->comm, fd, " Jacobian operator %d with scalar multiplier %gn",
171: grid->matOps[op].op, PetscRealPart(grid->matOps[op].alpha));
172: }
173: }
174: }
175: SNESView(snes, viewer);
176: }
178: return(0);
179: }
181: /*@C
182: GSNESDuplicate - Duplicates a grid SNES.
184: Collective on GSNES
186: Input Parameter:
187: . snes - The GSNES
189: Output Parameter:
190: . newsnes - The new GSNES
192: Level: beginner
194: .keywords: grid SNES, duplicate
195: .seealso: SNESDuplicate(), GSNESDestroy()
196: @*/
197: int GSNESDuplicate(GSNES snes, GSNES *newsnes)
198: {
202: SETERRQ(PETSC_ERR_SUP, " ");
203: #if 0
204: return(0);
205: #endif
206: }
208: /*@C
209: GSNESDuplicateMonitors - Duplicates the monitors of a grid SNES in another.
211: Collective on GSNES
213: Input Parameter:
214: . snes - The GSNES to copy from
216: Output Parameter:
217: . newsnes - The GSNES to copy to
219: Level: intermediate
221: .keywords: grid snes, duplicate, monitor
222: .seealso: GSNESDuplicate(), GSNESDestroy()
223: @*/
224: int GSNESDuplicateMonitors(GSNES snes, GSNES newsnes)
225: {
226: int mon, mon2;
233: for(mon = 0; mon < snes->numbermonitors; mon++)
234: {
235: for(mon2 = 0; mon2 < newsnes->numbermonitors; mon2++)
236: if (newsnes->monitor[mon] == snes->monitor[mon])
237: break;
238: if (mon2 == newsnes->numbermonitors) {
239: SNESSetMonitor(newsnes, snes->monitor[mon], snes->monitorcontext[mon], PETSC_NULL);
240: }
241: }
242: return(0);
243: }
245: /*@C
246: GSNESReallocate - Reallocates the storage for a grid SNES.
248: Collective on GSNES
250: Input Parameter:
251: . snes - The GSNES
253: Level: advanced
255: .keywords: grid snes, reallocate
256: .seealso: GridReform(), GSNESDuplicate(), GSNESDestroy()
257: @*/
258: int GSNESReallocate(GSNES snes)
259: {
260: void *funcCtx;
261: void *jacCtx;
262: int ierr;
266: /* Preserve contexts */
267: SNESGetFunction(snes, PETSC_NULL, &funcCtx, PETSC_NULL);
268: SNESGetJacobian(snes, PETSC_NULL, PETSC_NULL, &jacCtx, PETSC_NULL);
269: GSNESDestroyStructures_Private(snes);
270: GSNESCreateStructures_Private(snes, funcCtx, jacCtx);
271: return(0);
272: }
274: /*@
275: GSNESGetGrid - Gets the grid from a grid SNES.
277: Not collective
279: Input Parameter:
280: . snes - The grid SNES
282: Output Parameter:
283: . grid - The grid context
285: Level: intermediate
287: .keywords: grid SNES, grid, get
288: .seealso: GVecGetGrid(), GMatGetGrid()
289: @*/
290: int GSNESGetGrid(GSNES snes, Grid *grid)
291: {
298: PetscObjectQuery((PetscObject) snes, "Grid", (PetscObject *) grid);
300: return(0);
301: }
303: /*@
304: GSNESEvaluateRhs - Constructs the Rhs vector for Newton's equation.
306: Collective on GSNES
308: Input Parameters:
309: + snes - The grid SNES
310: . x - The current iterate
311: - ctx - The optional user context
313: Output Parameter:
314: . f - The function value
316: Note:
317: This function initializes the vector f to zero before calculation.
319: Level: advanced
321: .keywords: grid SNES, rhs
322: .seealso: GridEvaluateRhs()
323: @*/
324: int GSNESEvaluateRhs(GSNES snes, GVec x, GVec f, PetscObject ctx)
325: {
326: Grid grid;
327: PetscScalar zero = 0.0;
328: int ierr;
332: GSNESGetGrid(snes, &grid);
333: /* Initialize vector */
334: VecSet(&zero, f);
335: /* Form function */
336: GridEvaluateRhs(grid, x, f, ctx);
337: return(0);
338: }
340: /*@
341: GSNESEvaluateRhsFunction - Constructs only the constant portion of the nonlinear equation,
342: that is the portion arising from functions independent of the solution variable x.
344: Collective on GSNES
346: Input Parameters:
347: + snes - The grid SNES
348: . x - The current iterate
349: - ctx - The optional user context
351: Output Parameter:
352: . f - The function value
354: Notes:
355: This function actually constructs -f since f is usually written on the right, and the user
356: cannot just multiply by -1 after boundary conditions are imposed.
358: Level: advanced
360: .keywords: grid SNES, rhs, function
361: .seealso: GSNESEvaluateRhsOperator(), GSNESEvaluateRhs(), GSNESEvaluateJacobian()
362: @*/
363: int GSNESEvaluateRhsFunction(GSNES snes, GVec x, GVec f, void *ctx)
364: {
365: Grid grid;
366: PetscScalar zero = 0.0;
367: int ierr;
371: GSNESGetGrid(snes, &grid);
372: /* Initialize vector */
373: VecSet(&zero, f);
374: /* Form function */
375: GridScaleRhs(grid, -1.0);
376: GridEvaluateRhsFunction(grid, x, f, ctx);
377: GridScaleRhs(grid, -1.0);
378: return(0);
379: }
381: /*@
382: GSNESEvaluateRhsOperator - Constructs only the portion of the nonlinear equation dependent of the solution variable x.
384: Collective on GSNES
386: Input Parameters:
387: + snes - The grid SNES
388: . x - The current iterate
389: - ctx - The optional user context
391: Output Parameter:
392: . f - The function value
394: Level: advanced
396: .keywords: grid SNES, rhs, operator
397: .seealso: GSNESEvaluateRhsLinearOperator(), GSNESEvaluateRhsNonlinearOperator(), GSNESEvaluateRhsFunction(),
398: GSNESEvaluateRhs(), GSNESEvaluateJacobian()
399: @*/
400: int GSNESEvaluateRhsOperator(GSNES snes, GVec x, GVec f, void *ctx)
401: {
402: Grid grid;
403: PetscTruth reduceElement;
404: PetscScalar zero = 0.0;
405: int ierr;
409: GSNESGetGrid(snes, &grid);
410: /* Initialize vector */
411: VecSet(&zero, f);
412: /* Turn off boundary values */
413: GridGetReduceElement(grid, &reduceElement);
414: GridSetReduceElement(grid, PETSC_FALSE);
415: /* Apply operator */
416: GridEvaluateRhsOperator(grid, x, f, PETSC_TRUE, PETSC_TRUE, ctx);
417: /* Restore boundary values */
418: GridSetReduceElement(grid, reduceElement);
419: return(0);
420: }
422: /*@
423: GSNESEvaluateRhsLinearOperator - Constructs only the portion of the nonlinear equation
424: linearly dependent on the solution variable x.
426: Collective on GSNES
428: Input Parameters:
429: + snes - The grid SNES
430: . x - The current iterate
431: - ctx - The optional user context
433: Output Parameter:
434: . f - The function value
436: Level: advanced
438: .keywords: grid SNES, rhs, operator, linear
439: .seealso: GSNESEvaluateRhsOperator(), GSNESEvaluateRhsFunction(), GSNESEvaluateRhs(), GSNESEvaluateJacobian()
440: @*/
441: int GSNESEvaluateRhsLinearOperator(GSNES snes, GVec x, GVec f, void *ctx)
442: {
443: Grid grid;
444: PetscTruth reduceElement;
445: PetscScalar zero = 0.0;
446: int ierr;
450: GSNESGetGrid(snes, &grid);
451: /* Initialize vector */
452: VecSet(&zero, f);
453: /* Turn off boundary values */
454: GridGetReduceElement(grid, &reduceElement);
455: GridSetReduceElement(grid, PETSC_FALSE);
456: /* Apply operator */
457: GridEvaluateRhsOperator(grid, x, f, PETSC_TRUE, PETSC_FALSE, ctx);
458: /* Restore boundary values */
459: GridSetReduceElement(grid, reduceElement);
460: return(0);
461: }
463: /*@
464: GSNESEvaluateRhsNonlinearOperator - Constructs only the portion of the nonlinear equation
465: nonlinearly dependent on the solution variable x.
467: Collective on GSNES
469: Input Parameters:
470: + snes - The grid SNES
471: . n - The number of input vectors
472: . vecs - The input vectors
473: - ctx - The optional user context
475: Output Parameter:
476: . f - The function value
478: Level: advanced
480: Note:
481: If there is only one input argument, it will be used as every argument for
482: multilinear functions.
484: .keywords: grid SNES, rhs, operator, nonlinear
485: .seealso: GSNESEvaluateRhsOperator(), GSNESEvaluateRhsFunction(), GSNESEvaluateRhs(), GSNESEvaluateJacobian()
486: @*/
487: int GSNESEvaluateRhsNonlinearOperator(GSNES snes, int n, GVec vecs[], GVec f, void *ctx)
488: {
489: Grid grid;
490: NonlinearOperator op;
491: int field;
492: PetscScalar alpha;
493: PetscTruth isALE;
494: PetscTruth reduceElement;
495: PetscScalar zero = 0.0;
496: int ierr;
500: GSNESGetGrid(snes, &grid);
501: /* Initialize vector */
502: VecSet(&zero, f);
503: /* Turn off boundary values */
504: GridGetReduceElement(grid, &reduceElement);
505: GridSetReduceElement(grid, PETSC_FALSE);
506: /* Apply operator */
507: if (n == 1) {
508: GridEvaluateRhsOperator(grid, vecs[0], f, PETSC_FALSE, PETSC_TRUE, ctx);
509: } else {
510: GridGetNonlinearOperatorStart(grid, &op, &field, &alpha, &isALE);
511: while(op != PETSC_NULL) {
512: GVecEvaluateNonlinearOperatorGalerkin(f, n, vecs, 1, &field, op, alpha, isALE, ctx);
513: GridGetNonlinearOperatorNext(grid, &op, &field, &alpha, &isALE);
514: }
515: }
516: /* Restore boundary values */
517: GridSetReduceElement(grid, reduceElement);
518: return(0);
519: }
521: /*@
522: GSNESEvaluateJacobian - Constructs the system matrix for Newton's equation
524: Collective on GSNES
526: Input Parameters:
527: + snes - The grid SNES
528: . x - The current iterate
529: . flag - The matrix structure flag
530: - ctx - The optional user context
532: Output Parameters:
533: + J - The Jacobian matrix
534: - M - The preconditioner for the Jacobian matrix, usually the same as J
536: Level: advanced
538: .keywords: grid SNES, jacobian
539: .seealso: GSNESEvaluateJacobianMF(), GSNESEvaluateRhs(), GridEvaluateSystemMatrix()
540: @*/
541: int GSNESEvaluateJacobian(GSNES snes, GVec x, GMat *J, GMat *M, MatStructure *flag, PetscObject ctx)
542: {
543: Grid grid;
544: int ierr;
548: GSNESGetGrid(snes, &grid);
549: /* Initialize matrix */
550: MatZeroEntries(*J);
551: /* Form function */
552: GridEvaluateSystemMatrix(grid, x, J, M, flag, ctx);
553: /* Apply boundary conditions */
554: GMatSetBoundary(*J, 1.0, ctx);
555: return(0);
556: }
558: /*@
559: GSNESEvaluateJacobianMF - Constructs the application of the system matrix for Newton's equation
561: Collective on GSNES
563: Input Parameters:
564: + snes - The grid SNES
565: . x - The current iterate
566: . y - The vector to which to apply the Jacobian
567: - ctx - The optional user context
569: Output Parameter:
570: . f - The vector J(x) y
572: Level: advanced
574: .keywords: grid SNES< jacobian, mf
575: .seealso: GSNESEvaluateJacobian(), GSNESEvaluateRhs(), GridEvaluateSystemMatrix(), GVecEvaluateJacobian()
576: @*/
577: int GSNESEvaluateJacobianMF(GSNES snes, GVec x, GVec y, GVec f, void *ctx)
578: {
579: Grid grid;
580: PetscScalar zero = 0.0;
581: int ierr;
585: GSNESGetGrid(snes, &grid);
586: /* Initialize vector */
587: VecSet(&zero, f);
588: /* Form function */
589: GVecEvaluateJacobian(f, x, y, ctx);
590: return(0);
591: }
593: /*@
594: GSNESRhsBC - This sets zero boundary conditions for the Rhs.
596: Collective on GSNES
598: Input Parameters:
599: + snes - The SNES
600: - ctx - The user-supplied context
602: Output Parameter:
603: . rhs - The modified Rhs
605: Level: advanced
607: .keywords: grid SNES, set, boundary conditions
608: .seealso GSNESSetSolutionBC()
609: @*/
610: int GSNESRhsBC(GSNES snes, GVec rhs, void *ctx)
611: {
617: GVecSetBoundaryZero(rhs, ctx);
618: return(0);
619: }
621: /*@
622: GSNESSolutionBC - This sets the boundary conditions registered with the grid
623: on the solution vector.
625: Collective on GSNES
627: Input Parameters:
628: + snes - The SNES
629: - ctx - The user-supplied context
631: Output Parameter:
632: . sol - The modified solution
634: Level: advanced
636: .keywords: GSNES, set, boundary conditions
637: .seealso: GSNESSetRhsBC()
638: @*/
639: int GSNESSolutionBC(GSNES snes, GVec sol, void *ctx)
640: {
646: GVecSetBoundary(sol, ctx);
647: return(0);
648: }
650: /*@
651: GSNESUpdate - This is a general purpose function which is called at the beginning of every iterate.
653: Collective on GSNES
655: Input Parameters:
656: + snes - The nonlinear solver context
657: - step - The current step
659: Level: advanced
661: .keywords: GSNES, update
662: .seealso: SNESSolve()
663: @*/
664: int GSNESUpdate(GSNES snes, int step)
665: {
666: Grid grid;
667: Mesh mesh;
668: MeshMover mover;
669: int ierr;
672: /* Calculate mesh velocity */
673: GSNESGetGrid(snes, &grid);
674: GridGetMesh(grid, &mesh);
675: MeshGetMover(mesh, &mover);
676: MeshMoverCalcNodeVelocities(mover, PETSC_TRUE);
677: return(0);
678: }
680: /*@
681: GSNESCreate - Creates a grid SNES associated with a particular discretization context.
683: Collective on Grid
685: Input Parameter:
686: + grid - The discretization context
687: - ctx - The optional user context
689: Output Parameter:
690: . snes - The nonlinear solver context
692: Level: beginner
694: .keywords: GSNES, create
695: .seealso: GSNESDestroy()
696: @*/
697: int GSNESCreate(Grid grid, void *ctx, GSNES *snes)
698: {
704: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
705: GSolverInitializePackage(PETSC_NULL);
706: #endif
708: /* Create the SNES */
709: SNESCreate(grid->comm,snes);
710: PetscObjectCompose((PetscObject) *snes, "Grid", (PetscObject) grid);
711: (*snes)->isGSNES = PETSC_TRUE;
713: /* Set boundary conditions */
714: SNESSetRhsBC(*snes, GSNESRhsBC);
715: SNESSetSolutionBC(*snes, GSNESSolutionBC);
716: if (grid->reduceSystem == PETSC_TRUE) {
717: /* Save space and time by reducing the system right in construction rather than using an explicit matvec and subtraction */
718: grid->reduceElement = PETSC_TRUE;
719: /* Since we solve J x = -F, we need a BC multiplier of -1 */
720: GridSetBCMultiplier(grid, -1.0);
721: }
723: /* Set update function */
724: SNESSetUpdate(*snes, GSNESUpdate);
726: /* Setup SNES data structures */
727: GSNESCreateStructures_Private(*snes, ctx, ctx);
728: return(0);
729: }
731: int GSNESDefaultComputeJacobianWithColoring(SNES snes, GVec x, GMat *J, GMat *M, MatStructure *flag, void *ctx)
732: {
733: SETERRQ(PETSC_ERR_SUP, " ");
734: #if 0
735: Grid grid;
736: int ierr;
743: GMatGetGrid(*M, &grid);
745: if (!grid->fdcoloring) {
746: GMatFDColoringCreate(*M);
747: }
748: SNESDefaultComputeJacobianWithColoring(snes, x, J, M, flag, grid->fdcoloring);
749: return(0);
750: #endif
751: }
753: /*@
754: GSNESSolutionMonitor - Monitors solution at each GSNES iteration.
756: Collective on GSNES
758: Input Parameters:
759: + snes - The nonlinear solver context
760: . it - The number of iterations so far
761: . rnorm - The current (approximate) residual norm
762: - ctx - The viewer
764: Level: advanced
766: .keywords: grid SNES, monitor, solution
767: .seealso: GSNESResidualMonitor(), GSNESErrorMonitor(), SNESDefaultMonitor(), SNESSetMonitor()
768: @*/
769: int GSNESSolutionMonitor(GSNES snes, int it, PetscReal rnorm, void *ctx)
770: {
771: PetscViewer viewer = (PetscViewer) ctx;
772: Vec x;
773: int ierr;
777: SNESGetSolution(snes, &x);
778: VecView(x, viewer);
779: return(0);
780: }
782: /*@
783: GSNESResidualMonitor - Monitors residual at each SNES iteration.
785: Collective on GSNES
787: Input Parameters:
788: + snes - The nonlinear solver context
789: . it - The number of iterations so far
790: . rnorm - The current (approximate) residual norm
791: - ctx - The viewer
793: Level: advanced
795: .keywords: grid SNES, monitor, residual
796: .seealso: GSNESSolutionMonitor(), GSNESErrorMonitor(), SNESDefaultMonitor(), SNESSetMonitor()
797: @*/
798: int GSNESResidualMonitor(GSNES snes, int it, PetscReal rnorm, void *ctx)
799: {
800: PetscViewer viewer = (PetscViewer) ctx;
801: Vec x;
802: int ierr;
806: SNESGetFunction(snes, &x, PETSC_NULL, PETSC_NULL);
807: VecView(x, viewer);
808: return(0);
809: }
811: /*@
812: GSNESErrorMonitor - Displays the error at each iteration.
814: Collective on GSNES
816: Input Parameters:
817: + snes - The nonlinear context
818: . it - The number of iterations so far
819: . rnorm - The current (approximate) residual norm
820: - monCtx - A GSNESErrorMonitorCtx
822: Notes:
823: The final argument to SNESSetMonitor() with this routine must be a pointer to a GSNESErrorMonitorCtx.
825: Level: advanced
827: .keywords: grid SNES, monitor, residual
828: .seealso: GSNESSolutionMonitor(), GSNESResidualMonitor(), SNESDefaultMonitor(), SNESSetMonitor()
829: @*/
830: int GSNESErrorMonitor(GSNES snes, int it, PetscReal rnorm, void *monCtx)
831: {
832: GSNESErrorMonitorCtx *errorCtx = (GSNESErrorMonitorCtx *) monCtx;
833: void *ctx = errorCtx->ctx;
834: PetscScalar minusOne = -1.0;
835: GVec x, e;
836: PetscReal norm_2, norm_max, sol_norm_2;
837: int ierr;
841: if (it == 0) {
842: /* Build solution at the beginning */
843: (*errorCtx->solutionFunc)(errorCtx->solution, ctx);
844: }
845: SNESGetSolution(snes, &x);
846: GVecGetWorkGVec(errorCtx->solution, &e);
847: VecWAXPY(&minusOne, x, errorCtx->solution, e);
848: VecView(e, errorCtx->error_viewer);
850: /* Compute 2-norm and max-norm of error */
851: if (errorCtx->norm_error_viewer) {
852: VecNorm(e, NORM_2, &norm_2);
853: VecNorm(e, NORM_MAX, &norm_max);
854: VecNorm(errorCtx->solution, NORM_2, &sol_norm_2);
855: PetscViewerASCIIPrintf(errorCtx->norm_error_viewer,
856: "Iteration %d residual norm %g error 2-norm %g error max-norm %g exact sol norm-2 %gn",
857: it, rnorm, norm_2, norm_max, sol_norm_2);
858: }
859: GVecRestoreWorkGVec(errorCtx->solution, &e);
860: return(0);
861: }
863: EXTERN_C_BEGIN
864: int GSNESOptionsChecker_Private(GSNES snes)
865: {
866: char *prefix;
867: int loc[4], nmax;
868: MPI_Comm comm;
869: PetscViewer viewer;
870: PetscDraw draw;
871: PetscTruth opt;
872: int ierr;
875: SNESGetOptionsPrefix(snes, &prefix);
876: PetscObjectGetComm((PetscObject) snes, &comm);
878: nmax = 4;
879: loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
880: PetscOptionsGetIntArray(prefix, "-gvec_snes_solutionmonitor", loc, &nmax, &opt);
881: if (opt) {
882: PetscViewerDrawOpen(comm, 0, 0, loc[0], loc[1], loc[2], loc[3], &viewer);
883: PetscViewerDrawGetDraw(viewer, 0, &draw);
884: PetscDrawSetTitle(draw, "Approx. Solution");
885: PetscLogObjectParent(snes, (PetscObject) viewer);
886: SNESSetMonitor(snes, GSNESSolutionMonitor, (void *) viewer, PETSC_NULL);
887: }
888: nmax = 4;
889: loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
890: PetscOptionsGetIntArray(prefix, "-gvec_snes_residualmonitor", loc, &nmax, &opt);
891: if (opt) {
892: PetscViewerDrawOpen(comm, 0, 0, loc[0], loc[1], loc[2], loc[3], &viewer);
893: PetscViewerDrawGetDraw(viewer, 0, &draw);
894: PetscDrawSetTitle(draw, "Residual");
895: PetscLogObjectParent(snes, (PetscObject) viewer);
896: SNESSetMonitor(snes, GSNESResidualMonitor, (void *) viewer, PETSC_NULL);
897: }
898: #if 0
899: nmax = 4;
900: loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
901: PetscOptionsGetIntArray(prefix, "-gvec_snes_errormonitor", loc, &nmax, &opt);
902: if (opt) {
903: PetscViewerDrawOpen(comm, 0, 0, loc[0], loc[1], loc[2], loc[3], &viewer);
904: PetscViewerDrawGetDraw(viewer, 0, &draw);
905: PetscDrawSetTitle(draw, "Error");
906: PetscLogObjectParent(snes, (PetscObject) viewer);
907: SNESSetMonitor(snes, GSNESErrorMonitor, (void *) viewer, PETSC_NULL);
908: }
909: #endif
910: /*--------------------------------------------------------------------- */
911: PetscOptionsHasName(prefix, "-gvec_snes_fd", &opt);
912: if (opt) {
913: GMat A,B;
914: SNESGetJacobian(snes, &A, &B, PETSC_NULL, PETSC_NULL);
915: SNESSetJacobian(snes, A, B, GSNESDefaultComputeJacobianWithColoring, 0);
916: PetscLogInfo(snes, "GSNESSetFromOptions: Setting gvec finite difference Jacobian matrixn");
917: }
918: /*--------------------------------------------------------------------- */
919: PetscOptionsHasName(PETSC_NULL, "-help", &opt);
920: if (opt) {
921: char *pprefix;
922: int len = 2;
923: if (prefix != PETSC_NULL) {
924: PetscStrlen(prefix, &len);
925: }
926: PetscMalloc((len+2) * sizeof(char), &pprefix);
927: PetscStrcpy(pprefix, "-");
928: if (prefix != PETSC_NULL)
929: PetscStrcat(pprefix, prefix);
930: PetscPrintf(comm," Additional SNES Monitor options for grid vectorsn");
931: PetscPrintf(comm," %sgvec_snes_solutionmonitorn", pprefix);
932: PetscPrintf(comm," %sgvec_snes_residualmonitorn", pprefix);
933: PetscPrintf(comm," %sgvec_snes_errormonitor -- not working yetn", pprefix);
934: PetscFree(pprefix);
935: }
936: return(0);
937: }
938: EXTERN_C_END