Actual source code: fasfunc.c
petsc-3.13.6 2020-09-29
1: #include <../src/snes/impls/fas/fasimpls.h>
4: /* -------------- functions called on the fine level -------------- */
6: /*@
7: SNESFASSetType - Sets the update and correction type used for FAS.
9: Logically Collective
11: Input Parameters:
12: + snes - FAS context
13: - fastype - SNES_FAS_ADDITIVE, SNES_FAS_MULTIPLICATIVE, SNES_FAS_FULL, or SNES_FAS_KASKADE
15: Level: intermediate
17: .seealso: PCMGSetType()
18: @*/
19: PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype)
20: {
21: SNES_FAS *fas;
27: fas = (SNES_FAS*)snes->data;
28: fas->fastype = fastype;
29: if (fas->next) {
30: SNESFASSetType(fas->next, fastype);
31: }
32: return(0);
33: }
36: /*@
37: SNESFASGetType - Sets the update and correction type used for FAS.
39: Logically Collective
41: Input Parameters:
42: . snes - FAS context
44: Output Parameters:
45: . fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE
47: Level: intermediate
49: .seealso: PCMGSetType()
50: @*/
51: PetscErrorCode SNESFASGetType(SNES snes,SNESFASType *fastype)
52: {
53: SNES_FAS *fas;
58: fas = (SNES_FAS*)snes->data;
59: *fastype = fas->fastype;
60: return(0);
61: }
63: /*@C
64: SNESFASSetLevels - Sets the number of levels to use with FAS.
65: Must be called before any other FAS routine.
67: Input Parameters:
68: + snes - the snes context
69: . levels - the number of levels
70: - comms - optional communicators for each level; this is to allow solving the coarser
71: problems on smaller sets of processors.
73: Level: intermediate
75: Notes:
76: If the number of levels is one then the multigrid uses the -fas_levels prefix
77: for setting the level options rather than the -fas_coarse prefix.
79: .seealso: SNESFASGetLevels()
80: @*/
81: PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm *comms)
82: {
84: PetscInt i;
85: const char *optionsprefix;
86: char tprefix[128];
87: SNES_FAS *fas;
88: SNES prevsnes;
89: MPI_Comm comm;
93: fas = (SNES_FAS*)snes->data;
94: PetscObjectGetComm((PetscObject)snes,&comm);
95: if (levels == fas->levels) {
96: if (!comms) return(0);
97: }
98: /* user has changed the number of levels; reset */
99: (*snes->ops->reset)(snes);
100: /* destroy any coarser levels if necessary */
101: SNESDestroy(&fas->next);
102: fas->next = NULL;
103: fas->previous = NULL;
104: prevsnes = snes;
105: /* setup the finest level */
106: SNESGetOptionsPrefix(snes, &optionsprefix);
107: PetscObjectComposedDataSetInt((PetscObject) snes, PetscMGLevelId, levels-1);
108: for (i = levels - 1; i >= 0; i--) {
109: if (comms) comm = comms[i];
110: fas->level = i;
111: fas->levels = levels;
112: fas->fine = snes;
113: fas->next = NULL;
114: if (i > 0) {
115: SNESCreate(comm, &fas->next);
116: SNESGetOptionsPrefix(fas->fine, &optionsprefix);
117: PetscSNPrintf(tprefix,sizeof(tprefix),"fas_levels_%d_cycle_",(int)fas->level);
118: SNESAppendOptionsPrefix(fas->next,optionsprefix);
119: SNESAppendOptionsPrefix(fas->next,tprefix);
120: SNESSetType(fas->next, SNESFAS);
121: SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);
122: PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);
123: PetscObjectComposedDataSetInt((PetscObject) fas->next, PetscMGLevelId, i-1);
125: ((SNES_FAS*)fas->next->data)->previous = prevsnes;
127: prevsnes = fas->next;
128: fas = (SNES_FAS*)prevsnes->data;
129: }
130: }
131: return(0);
132: }
135: /*@
136: SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids
138: Input Parameter:
139: . snes - the nonlinear solver context
141: Output parameter:
142: . levels - the number of levels
144: Level: advanced
146: .seealso: SNESFASSetLevels(), PCMGGetLevels()
147: @*/
148: PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels)
149: {
150: SNES_FAS *fas;
155: fas = (SNES_FAS*)snes->data;
156: *levels = fas->levels;
157: return(0);
158: }
161: /*@
162: SNESFASGetCycleSNES - Gets the SNES corresponding to a particular
163: level of the FAS hierarchy.
165: Input Parameters:
166: + snes - the multigrid context
167: level - the level to get
168: - lsnes - whether to use the nonlinear smoother or not
170: Level: advanced
172: .seealso: SNESFASSetLevels(), SNESFASGetLevels()
173: @*/
174: PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes)
175: {
176: SNES_FAS *fas;
177: PetscInt i;
182: fas = (SNES_FAS*)snes->data;
183: if (level > fas->levels-1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels);
184: if (fas->level != fas->levels - 1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESFASGetCycleSNES may only be called on the finest-level SNES.",level,fas->level);
186: *lsnes = snes;
187: for (i = fas->level; i > level; i--) {
188: *lsnes = fas->next;
189: fas = (SNES_FAS*)(*lsnes)->data;
190: }
191: if (fas->level != level) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
192: return(0);
193: }
195: /*@
196: SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
197: use on all levels.
199: Logically Collective on SNES
201: Input Parameters:
202: + snes - the multigrid context
203: - n - the number of smoothing steps
205: Options Database Key:
206: . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
208: Level: advanced
210: .seealso: SNESFASSetNumberSmoothDown()
211: @*/
212: PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n)
213: {
214: SNES_FAS *fas;
219: fas = (SNES_FAS*)snes->data;
220: fas->max_up_it = n;
221: if (!fas->smoothu && fas->level != 0) {
222: SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
223: }
224: if (fas->smoothu) {
225: SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);
226: }
227: if (fas->next) {
228: SNESFASSetNumberSmoothUp(fas->next, n);
229: }
230: return(0);
231: }
233: /*@
234: SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
235: use on all levels.
237: Logically Collective on SNES
239: Input Parameters:
240: + snes - the multigrid context
241: - n - the number of smoothing steps
243: Options Database Key:
244: . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
246: Level: advanced
248: .seealso: SNESFASSetNumberSmoothUp()
249: @*/
250: PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n)
251: {
252: SNES_FAS *fas;
257: fas = (SNES_FAS*)snes->data;
258: if (!fas->smoothd) {
259: SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);
260: }
261: SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);
263: fas->max_down_it = n;
264: if (fas->next) {
265: SNESFASSetNumberSmoothDown(fas->next, n);
266: }
267: return(0);
268: }
271: /*@
272: SNESFASSetContinuation - Sets the FAS cycle to default to exact Newton solves on the upsweep
274: Logically Collective on SNES
276: Input Parameters:
277: + snes - the multigrid context
278: - n - the number of smoothing steps
280: Options Database Key:
281: . -snes_fas_continuation - sets continuation to true
283: Level: advanced
285: Notes:
286: This sets the prefix on the upsweep smoothers to -fas_continuation
288: .seealso: SNESFAS
289: @*/
290: PetscErrorCode SNESFASSetContinuation(SNES snes,PetscBool continuation)
291: {
292: const char *optionsprefix;
293: char tprefix[128];
294: SNES_FAS *fas;
299: fas = (SNES_FAS*)snes->data;
300: SNESGetOptionsPrefix(fas->fine, &optionsprefix);
301: if (!fas->smoothu) {
302: SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
303: }
304: PetscStrncpy(tprefix,"fas_levels_continuation_",sizeof(tprefix));
305: SNESSetOptionsPrefix(fas->smoothu, optionsprefix);
306: SNESAppendOptionsPrefix(fas->smoothu, tprefix);
307: SNESSetType(fas->smoothu,SNESNEWTONLS);
308: SNESSetTolerances(fas->smoothu,fas->fine->abstol,fas->fine->rtol,fas->fine->stol,50,100);
309: fas->continuation = continuation;
310: if (fas->next) {
311: SNESFASSetContinuation(fas->next,continuation);
312: }
313: return(0);
314: }
317: /*@
318: SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use SNESFASSetCyclesOnLevel() for more
319: complicated cycling.
321: Logically Collective on SNES
323: Input Parameters:
324: + snes - the multigrid context
325: - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
327: Options Database Key:
328: . -snes_fas_cycles 1 or 2
330: Level: advanced
332: .seealso: SNESFASSetCyclesOnLevel()
333: @*/
334: PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles)
335: {
336: SNES_FAS *fas;
338: PetscBool isFine;
342: SNESFASCycleIsFine(snes, &isFine);
343: fas = (SNES_FAS*)snes->data;
344: fas->n_cycles = cycles;
345: if (!isFine) {
346: SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);
347: }
348: if (fas->next) {
349: SNESFASSetCycles(fas->next, cycles);
350: }
351: return(0);
352: }
355: /*@
356: SNESFASSetMonitor - Sets the method-specific cycle monitoring
358: Logically Collective on SNES
360: Input Parameters:
361: + snes - the FAS context
362: . vf - viewer and format structure (may be NULL if flg is FALSE)
363: - flg - monitor or not
365: Level: advanced
367: .seealso: SNESFASSetCyclesOnLevel()
368: @*/
369: PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg)
370: {
371: SNES_FAS *fas;
373: PetscBool isFine;
374: PetscInt i, levels;
375: SNES levelsnes;
379: SNESFASCycleIsFine(snes, &isFine);
380: fas = (SNES_FAS*)snes->data;
381: levels = fas->levels;
382: if (isFine) {
383: for (i = 0; i < levels; i++) {
384: SNESFASGetCycleSNES(snes, i, &levelsnes);
385: fas = (SNES_FAS*)levelsnes->data;
386: if (flg) {
387: /* set the monitors for the upsmoother and downsmoother */
388: SNESMonitorCancel(levelsnes);
389: /* Only register destroy on finest level */
390: SNESMonitorSet(levelsnes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorDefault,vf,(!i ? (PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy : NULL));
391: } else if (i != fas->levels - 1) {
392: /* unset the monitors on the coarse levels */
393: SNESMonitorCancel(levelsnes);
394: }
395: }
396: }
397: return(0);
398: }
400: /*@
401: SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels
403: Logically Collective on SNES
405: Input Parameters:
406: + snes - the FAS context
407: - flg - monitor or not
409: Level: advanced
411: .seealso: SNESFASSetMonitor()
412: @*/
413: PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
414: {
415: SNES_FAS *fas;
417: PetscBool isFine;
418: PetscInt i, levels;
419: SNES levelsnes;
420: char eventname[128];
424: SNESFASCycleIsFine(snes, &isFine);
425: fas = (SNES_FAS*)snes->data;
426: levels = fas->levels;
427: if (isFine) {
428: for (i = 0; i < levels; i++) {
429: SNESFASGetCycleSNES(snes, i, &levelsnes);
430: fas = (SNES_FAS*)levelsnes->data;
431: if (flg) {
432: PetscSNPrintf(eventname,sizeof(eventname),"FASSetup %d",(int)i);
433: PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);
434: PetscSNPrintf(eventname,sizeof(eventname),"FASSmooth %d",(int)i);
435: PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);
436: PetscSNPrintf(eventname,sizeof(eventname),"FASResid %d",(int)i);
437: PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);
438: PetscSNPrintf(eventname,sizeof(eventname),"FASInterp %d",(int)i);
439: PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);
440: } else {
441: fas->eventsmoothsetup = 0;
442: fas->eventsmoothsolve = 0;
443: fas->eventresidual = 0;
444: fas->eventinterprestrict = 0;
445: }
446: }
447: }
448: return(0);
449: }
451: /*
452: Creates the default smoother type.
454: This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level.
456: */
457: PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth)
458: {
459: SNES_FAS *fas;
460: const char *optionsprefix;
461: char tprefix[128];
463: SNES nsmooth;
468: fas = (SNES_FAS*)snes->data;
469: SNESGetOptionsPrefix(fas->fine, &optionsprefix);
470: /* create the default smoother */
471: SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth);
472: if (fas->level == 0) {
473: PetscStrncpy(tprefix,"fas_coarse_",sizeof(tprefix));
474: SNESAppendOptionsPrefix(nsmooth, optionsprefix);
475: SNESAppendOptionsPrefix(nsmooth, tprefix);
476: SNESSetType(nsmooth, SNESNEWTONLS);
477: SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);
478: } else {
479: PetscSNPrintf(tprefix,sizeof(tprefix),"fas_levels_%d_",(int)fas->level);
480: SNESAppendOptionsPrefix(nsmooth, optionsprefix);
481: SNESAppendOptionsPrefix(nsmooth, tprefix);
482: SNESSetType(nsmooth, SNESNRICHARDSON);
483: SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);
484: }
485: PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);
486: PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth);
487: PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);
488: PetscObjectComposedDataSetInt((PetscObject) nsmooth, PetscMGLevelId, fas->level);
489: *smooth = nsmooth;
490: return(0);
491: }
493: /* ------------- Functions called on a particular level ----------------- */
495: /*@
496: SNESFASCycleSetCycles - Sets the number of cycles on a particular level.
498: Logically Collective on SNES
500: Input Parameters:
501: + snes - the multigrid context
502: . level - the level to set the number of cycles on
503: - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
505: Level: advanced
507: .seealso: SNESFASSetCycles()
508: @*/
509: PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
510: {
511: SNES_FAS *fas;
516: fas = (SNES_FAS*)snes->data;
517: fas->n_cycles = cycles;
518: SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);
519: return(0);
520: }
523: /*@
524: SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.
526: Logically Collective on SNES
528: Input Parameters:
529: . snes - the multigrid context
531: Output Parameters:
532: . smooth - the smoother
534: Level: advanced
536: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown()
537: @*/
538: PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
539: {
540: SNES_FAS *fas;
545: fas = (SNES_FAS*)snes->data;
546: *smooth = fas->smoothd;
547: return(0);
548: }
549: /*@
550: SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.
552: Logically Collective on SNES
554: Input Parameters:
555: . snes - the multigrid context
557: Output Parameters:
558: . smoothu - the smoother
560: Notes:
561: Returns the downsmoother if no up smoother is available. This enables transparent
562: default behavior in the process of the solve.
564: Level: advanced
566: .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown()
567: @*/
568: PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
569: {
570: SNES_FAS *fas;
575: fas = (SNES_FAS*)snes->data;
576: if (!fas->smoothu) *smoothu = fas->smoothd;
577: else *smoothu = fas->smoothu;
578: return(0);
579: }
581: /*@
582: SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level.
584: Logically Collective on SNES
586: Input Parameters:
587: . snes - the multigrid context
589: Output Parameters:
590: . smoothd - the smoother
592: Level: advanced
594: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
595: @*/
596: PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
597: {
598: SNES_FAS *fas;
603: fas = (SNES_FAS*)snes->data;
604: *smoothd = fas->smoothd;
605: return(0);
606: }
609: /*@
610: SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level
612: Logically Collective on SNES
614: Input Parameters:
615: . snes - the multigrid context
617: Output Parameters:
618: . correction - the coarse correction on this level
620: Notes:
621: Returns NULL on the coarsest level.
623: Level: advanced
625: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
626: @*/
627: PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
628: {
629: SNES_FAS *fas;
634: fas = (SNES_FAS*)snes->data;
635: *correction = fas->next;
636: return(0);
637: }
639: /*@
640: SNESFASCycleGetInterpolation - Gets the interpolation on this level
642: Logically Collective on SNES
644: Input Parameters:
645: . snes - the multigrid context
647: Output Parameters:
648: . mat - the interpolation operator on this level
650: Level: developer
652: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
653: @*/
654: PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
655: {
656: SNES_FAS *fas;
661: fas = (SNES_FAS*)snes->data;
662: *mat = fas->interpolate;
663: return(0);
664: }
667: /*@
668: SNESFASCycleGetRestriction - Gets the restriction on this level
670: Logically Collective on SNES
672: Input Parameters:
673: . snes - the multigrid context
675: Output Parameters:
676: . mat - the restriction operator on this level
678: Level: developer
680: .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation()
681: @*/
682: PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
683: {
684: SNES_FAS *fas;
689: fas = (SNES_FAS*)snes->data;
690: *mat = fas->restrct;
691: return(0);
692: }
695: /*@
696: SNESFASCycleGetInjection - Gets the injection on this level
698: Logically Collective on SNES
700: Input Parameters:
701: . snes - the multigrid context
703: Output Parameters:
704: . mat - the restriction operator on this level
706: Level: developer
708: .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction()
709: @*/
710: PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
711: {
712: SNES_FAS *fas;
717: fas = (SNES_FAS*)snes->data;
718: *mat = fas->inject;
719: return(0);
720: }
722: /*@
723: SNESFASCycleGetRScale - Gets the injection on this level
725: Logically Collective on SNES
727: Input Parameters:
728: . snes - the multigrid context
730: Output Parameters:
731: . mat - the restriction operator on this level
733: Level: developer
735: .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale()
736: @*/
737: PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
738: {
739: SNES_FAS *fas;
744: fas = (SNES_FAS*)snes->data;
745: *vec = fas->rscale;
746: return(0);
747: }
749: /*@
750: SNESFASCycleIsFine - Determines if a given cycle is the fine level.
752: Logically Collective on SNES
754: Input Parameters:
755: . snes - the FAS context
757: Output Parameters:
758: . flg - indicates if this is the fine level or not
760: Level: advanced
762: .seealso: SNESFASSetLevels()
763: @*/
764: PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
765: {
766: SNES_FAS *fas;
771: fas = (SNES_FAS*)snes->data;
772: if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
773: else *flg = PETSC_FALSE;
774: return(0);
775: }
777: /* ---------- functions called on the finest level that return level-specific information ---------- */
779: /*@
780: SNESFASSetInterpolation - Sets the function to be used to calculate the
781: interpolation from l-1 to the lth level
783: Input Parameters:
784: + snes - the multigrid context
785: . mat - the interpolation operator
786: - level - the level (0 is coarsest) to supply [do not supply 0]
788: Level: advanced
790: Notes:
791: Usually this is the same matrix used also to set the restriction
792: for the same level.
794: One can pass in the interpolation matrix or its transpose; PETSc figures
795: out from the matrix size which one it is.
797: .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
798: @*/
799: PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
800: {
801: SNES_FAS *fas;
803: SNES levelsnes;
808: SNESFASGetCycleSNES(snes, level, &levelsnes);
809: fas = (SNES_FAS*)levelsnes->data;
810: PetscObjectReference((PetscObject)mat);
811: MatDestroy(&fas->interpolate);
812: fas->interpolate = mat;
813: return(0);
814: }
816: /*@
817: SNESFASGetInterpolation - Gets the matrix used to calculate the
818: interpolation from l-1 to the lth level
820: Input Parameters:
821: + snes - the multigrid context
822: - level - the level (0 is coarsest) to supply [do not supply 0]
824: Output Parameters:
825: . mat - the interpolation operator
827: Level: advanced
829: .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale()
830: @*/
831: PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
832: {
833: SNES_FAS *fas;
835: SNES levelsnes;
840: SNESFASGetCycleSNES(snes, level, &levelsnes);
841: fas = (SNES_FAS*)levelsnes->data;
842: *mat = fas->interpolate;
843: return(0);
844: }
846: /*@
847: SNESFASSetRestriction - Sets the function to be used to restrict the defect
848: from level l to l-1.
850: Input Parameters:
851: + snes - the multigrid context
852: . mat - the restriction matrix
853: - level - the level (0 is coarsest) to supply [Do not supply 0]
855: Level: advanced
857: Notes:
858: Usually this is the same matrix used also to set the interpolation
859: for the same level.
861: One can pass in the interpolation matrix or its transpose; PETSc figures
862: out from the matrix size which one it is.
864: If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
865: is used.
867: .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
868: @*/
869: PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
870: {
871: SNES_FAS *fas;
873: SNES levelsnes;
878: SNESFASGetCycleSNES(snes, level, &levelsnes);
879: fas = (SNES_FAS*)levelsnes->data;
880: PetscObjectReference((PetscObject)mat);
881: MatDestroy(&fas->restrct);
882: fas->restrct = mat;
883: return(0);
884: }
886: /*@
887: SNESFASGetRestriction - Gets the matrix used to calculate the
888: restriction from l to the l-1th level
890: Input Parameters:
891: + snes - the multigrid context
892: - level - the level (0 is coarsest) to supply [do not supply 0]
894: Output Parameters:
895: . mat - the interpolation operator
897: Level: advanced
899: .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
900: @*/
901: PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
902: {
903: SNES_FAS *fas;
905: SNES levelsnes;
910: SNESFASGetCycleSNES(snes, level, &levelsnes);
911: fas = (SNES_FAS*)levelsnes->data;
912: *mat = fas->restrct;
913: return(0);
914: }
917: /*@
918: SNESFASSetInjection - Sets the function to be used to inject the solution
919: from level l to l-1.
921: Input Parameters:
922: + snes - the multigrid context
923: . mat - the restriction matrix
924: - level - the level (0 is coarsest) to supply [Do not supply 0]
926: Level: advanced
928: Notes:
929: If you do not set this, the restriction and rscale is used to
930: project the solution instead.
932: .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
933: @*/
934: PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
935: {
936: SNES_FAS *fas;
938: SNES levelsnes;
943: SNESFASGetCycleSNES(snes, level, &levelsnes);
944: fas = (SNES_FAS*)levelsnes->data;
945: PetscObjectReference((PetscObject)mat);
946: MatDestroy(&fas->inject);
948: fas->inject = mat;
949: return(0);
950: }
953: /*@
954: SNESFASGetInjection - Gets the matrix used to calculate the
955: injection from l-1 to the lth level
957: Input Parameters:
958: + snes - the multigrid context
959: - level - the level (0 is coarsest) to supply [do not supply 0]
961: Output Parameters:
962: . mat - the injection operator
964: Level: advanced
966: .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale()
967: @*/
968: PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
969: {
970: SNES_FAS *fas;
972: SNES levelsnes;
977: SNESFASGetCycleSNES(snes, level, &levelsnes);
978: fas = (SNES_FAS*)levelsnes->data;
979: *mat = fas->inject;
980: return(0);
981: }
983: /*@
984: SNESFASSetRScale - Sets the scaling factor of the restriction
985: operator from level l to l-1.
987: Input Parameters:
988: + snes - the multigrid context
989: . rscale - the restriction scaling
990: - level - the level (0 is coarsest) to supply [Do not supply 0]
992: Level: advanced
994: Notes:
995: This is only used in the case that the injection is not set.
997: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
998: @*/
999: PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
1000: {
1001: SNES_FAS *fas;
1003: SNES levelsnes;
1008: SNESFASGetCycleSNES(snes, level, &levelsnes);
1009: fas = (SNES_FAS*)levelsnes->data;
1010: PetscObjectReference((PetscObject)rscale);
1011: VecDestroy(&fas->rscale);
1012: fas->rscale = rscale;
1013: return(0);
1014: }
1016: /*@
1017: SNESFASGetSmoother - Gets the default smoother on a level.
1019: Input Parameters:
1020: + snes - the multigrid context
1021: - level - the level (0 is coarsest) to supply
1023: Output Parameters:
1024: smooth - the smoother
1026: Level: advanced
1028: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1029: @*/
1030: PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
1031: {
1032: SNES_FAS *fas;
1034: SNES levelsnes;
1039: SNESFASGetCycleSNES(snes, level, &levelsnes);
1040: fas = (SNES_FAS*)levelsnes->data;
1041: if (!fas->smoothd) {
1042: SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1043: }
1044: *smooth = fas->smoothd;
1045: return(0);
1046: }
1048: /*@
1049: SNESFASGetSmootherDown - Gets the downsmoother on a level.
1051: Input Parameters:
1052: + snes - the multigrid context
1053: - level - the level (0 is coarsest) to supply
1055: Output Parameters:
1056: smooth - the smoother
1058: Level: advanced
1060: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1061: @*/
1062: PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
1063: {
1064: SNES_FAS *fas;
1066: SNES levelsnes;
1071: SNESFASGetCycleSNES(snes, level, &levelsnes);
1072: fas = (SNES_FAS*)levelsnes->data;
1073: /* if the user chooses to differentiate smoothers, create them both at this point */
1074: if (!fas->smoothd) {
1075: SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1076: }
1077: if (!fas->smoothu) {
1078: SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);
1079: }
1080: *smooth = fas->smoothd;
1081: return(0);
1082: }
1084: /*@
1085: SNESFASGetSmootherUp - Gets the upsmoother on a level.
1087: Input Parameters:
1088: + snes - the multigrid context
1089: - level - the level (0 is coarsest)
1091: Output Parameters:
1092: smooth - the smoother
1094: Level: advanced
1096: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1097: @*/
1098: PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
1099: {
1100: SNES_FAS *fas;
1102: SNES levelsnes;
1107: SNESFASGetCycleSNES(snes, level, &levelsnes);
1108: fas = (SNES_FAS*)levelsnes->data;
1109: /* if the user chooses to differentiate smoothers, create them both at this point */
1110: if (!fas->smoothd) {
1111: SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1112: }
1113: if (!fas->smoothu) {
1114: SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);
1115: }
1116: *smooth = fas->smoothu;
1117: return(0);
1118: }
1120: /*@
1121: SNESFASGetCoarseSolve - Gets the coarsest solver.
1123: Input Parameters:
1124: . snes - the multigrid context
1126: Output Parameters:
1127: . coarse - the coarse-level solver
1129: Level: advanced
1131: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1132: @*/
1133: PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse)
1134: {
1135: SNES_FAS *fas;
1137: SNES levelsnes;
1142: SNESFASGetCycleSNES(snes, 0, &levelsnes);
1143: fas = (SNES_FAS*)levelsnes->data;
1144: /* if the user chooses to differentiate smoothers, create them both at this point */
1145: if (!fas->smoothd) {
1146: SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1147: }
1148: *coarse = fas->smoothd;
1149: return(0);
1150: }
1152: /*@
1153: SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS
1155: Logically Collective on SNES
1157: Input Parameters:
1158: + snes - the multigrid context
1159: - swp - whether to downsweep or not
1161: Options Database Key:
1162: . -snes_fas_full_downsweep - Sets number of pre-smoothing steps
1164: Level: advanced
1166: .seealso: SNESFASSetNumberSmoothUp()
1167: @*/
1168: PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp)
1169: {
1170: SNES_FAS *fas;
1175: fas = (SNES_FAS*)snes->data;
1176: fas->full_downsweep = swp;
1177: if (fas->next) {
1178: SNESFASFullSetDownSweep(fas->next,swp);
1179: }
1180: return(0);
1181: }