xref: /petsc/src/snes/impls/qn/qn.c (revision ceaaa4989964adb3f5eb130cb04b8f49c83e49c9)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2fced5a79SAsbjørn Nilsen Riseth #include <petscdm.h>
34b11644fSPeter Brune 
48e231d97SPeter Brune #define H(i, j) qn->dXdFmat[i * qn->m + j]
58e231d97SPeter Brune 
69e5d0892SLisandro Dalcin const char *const SNESQNScaleTypes[]   = {"DEFAULT", "NONE", "SCALAR", "DIAGONAL", "JACOBIAN", "SNESQNScaleType", "SNES_QN_SCALING_", NULL};
79e5d0892SLisandro Dalcin const char *const SNESQNRestartTypes[] = {"DEFAULT", "NONE", "POWELL", "PERIODIC", "SNESQNRestartType", "SNES_QN_RESTART_", NULL};
89e5d0892SLisandro Dalcin const char *const SNESQNTypes[]        = {"LBFGS", "BROYDEN", "BADBROYDEN", "SNESQNType", "SNES_QN_", NULL};
9b8840d6bSPeter Brune 
104b11644fSPeter Brune typedef struct {
1192f76d53SAlp Dener   Mat               B;      /* Quasi-Newton approximation Matrix (MATLMVM) */
128e231d97SPeter Brune   PetscInt          m;      /* The number of kept previous steps */
135c7a0a03SPeter Brune   PetscReal        *lambda; /* The line search history of the method */
1492f76d53SAlp Dener   PetscBool         monflg;
1544f7e39eSPeter Brune   PetscViewer       monitor;
166bf1b2e5SPeter Brune   PetscReal         powell_gamma; /* Powell angle restart condition */
17b21d5a53SPeter Brune   PetscReal         scaling;      /* scaling of H0 */
18b8840d6bSPeter Brune   SNESQNType        type;         /* the type of quasi-newton method used */
1988f769c5SPeter Brune   SNESQNScaleType   scale_type;   /* the type of scaling used */
200c777b0cSPeter Brune   SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */
219f83bee8SJed Brown } SNES_QN;
224b11644fSPeter Brune 
23d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_QN(SNES snes)
24d71ae5a4SJacob Faibussowitsch {
259f83bee8SJed Brown   SNES_QN             *qn = (SNES_QN *)snes->data;
2615f5eeeaSPeter Brune   Vec                  X, Xold;
270a3905e1SPeter Brune   Vec                  F, W;
2847f26062SPeter Brune   Vec                  Y, D, Dold;
29b8840d6bSPeter Brune   PetscInt             i, i_r;
304279555eSSatish Balay   PetscReal            fnorm, xnorm, ynorm;
31422a814eSBarry Smith   SNESLineSearchReason lssucceed;
3292f76d53SAlp Dener   PetscBool            badstep, powell, periodic;
331c6523dcSPeter Brune   PetscScalar          DolddotD, DolddotDold;
34b7281c8aSPeter Brune   SNESConvergedReason  reason;
354279555eSSatish Balay #if defined(PETSC_USE_INFO)
364279555eSSatish Balay   PetscReal gnorm;
374279555eSSatish Balay #endif
380ecc9e1dSPeter Brune 
3984c577daSBarry Smith   /* basically just a regular newton's method except for the application of the Jacobian */
404b11644fSPeter Brune 
416e111a19SKarl Rupp   PetscFunctionBegin;
420b121fc5SBarry Smith   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
43c579b300SPatrick Farrell 
449566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
459f3a0142SPeter Brune   F    = snes->vec_func;       /* residual vector */
463af51624SPeter Brune   Y    = snes->vec_sol_update; /* search direction generated by J^-1D*/
470a3905e1SPeter Brune   W    = snes->work[3];
48b8840d6bSPeter Brune   X    = snes->vec_sol; /* solution vector */
49335fb975SPeter Brune   Xold = snes->work[0];
503af51624SPeter Brune 
513af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
52335fb975SPeter Brune   D    = snes->work[1];
53335fb975SPeter Brune   Dold = snes->work[2];
544b11644fSPeter Brune 
554b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
564b11644fSPeter Brune 
579566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
584b11644fSPeter Brune   snes->iter = 0;
594b11644fSPeter Brune   snes->norm = 0.;
609566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
6147f26062SPeter Brune 
62efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
639566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
649566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
65b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
66b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
673ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
68b7281c8aSPeter Brune     }
699566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
7047f26062SPeter Brune   } else {
71e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
729566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
731aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
74c1c75074SPeter Brune 
759566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
76422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
7747f26062SPeter Brune   }
78efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
799566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, F, D));
809566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
81b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
82b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
833ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
84b7281c8aSPeter Brune     }
8547f26062SPeter Brune   } else {
869566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, D));
8747f26062SPeter Brune   }
88b28a06ddSPeter Brune 
899566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
904b11644fSPeter Brune   snes->norm = fnorm;
919566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
929566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
939566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
944b11644fSPeter Brune 
954b11644fSPeter Brune   /* test convergence */
962d157150SStefano Zampini   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
973ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
9870d3b23bSPeter Brune 
99efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_RIGHT) {
1009566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0));
1019566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
1029566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0));
1039566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
104ddd40ce5SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
105ddd40ce5SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
1063ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
107ddd40ce5SPeter Brune     }
1089566063dSJacob Faibussowitsch     PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
1099566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, D));
1103cf07b75SPeter Brune   }
1113cf07b75SPeter Brune 
11201fe78eaSStefano Zampini   /* general purpose update */
113dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
11401fe78eaSStefano Zampini 
115f8e15203SPeter Brune   /* scale the initial update */
1160c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
1179566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
11807b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
1199566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
1209566063dSJacob Faibussowitsch     PetscCall(MatLMVMSetJ0KSP(qn->B, snes->ksp));
1210ecc9e1dSPeter Brune   }
1220ecc9e1dSPeter Brune 
1235ba6227bSPeter Brune   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
12492f76d53SAlp Dener     /* update QN approx and calculate step */
1259566063dSJacob Faibussowitsch     PetscCall(MatLMVMUpdate(qn->B, X, D));
1269566063dSJacob Faibussowitsch     PetscCall(MatSolve(qn->B, D, Y));
12792f76d53SAlp Dener 
12870d3b23bSPeter Brune     /* line search for lambda */
1299371c9d4SSatish Balay     ynorm = 1;
1304279555eSSatish Balay #if defined(PETSC_USE_INFO)
1319371c9d4SSatish Balay     gnorm = fnorm;
1324279555eSSatish Balay #endif
1339566063dSJacob Faibussowitsch     PetscCall(VecCopy(D, Dold));
1349566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Xold));
1359566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y));
1369f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
1379566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
1389566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm));
13992f76d53SAlp Dener     badstep = PETSC_FALSE;
140422a814eSBarry Smith     if (lssucceed) {
1419f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
1429f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
1439f3a0142SPeter Brune         break;
1449f3a0142SPeter Brune       }
14592f76d53SAlp Dener       badstep = PETSC_TRUE;
1460ecc9e1dSPeter Brune     }
1473af51624SPeter Brune 
14888d374b2SPeter Brune     /* convergence monitoring */
1499566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)fnorm, (double)gnorm, (double)ynorm, (int)lssucceed));
150b21d5a53SPeter Brune 
151efd4aadfSBarry Smith     if (snes->npc && snes->npcside == PC_RIGHT) {
1529566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0));
1539566063dSJacob Faibussowitsch       PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
1549566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0));
1559566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
156ddd40ce5SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
157ddd40ce5SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
1583ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
159ddd40ce5SPeter Brune       }
1609566063dSJacob Faibussowitsch       PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
161b28a06ddSPeter Brune     }
162b28a06ddSPeter Brune 
1639566063dSJacob Faibussowitsch     PetscCall(SNESSetIterationNumber(snes, i + 1));
16471dbe336SPeter Brune     snes->norm  = fnorm;
165c1e67a49SFande Kong     snes->xnorm = xnorm;
166c1e67a49SFande Kong     snes->ynorm = ynorm;
167360c497dSPeter Brune 
1689566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
16992f76d53SAlp Dener 
1704b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
1712d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
1722d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
1733ba16761SJacob Faibussowitsch     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
174efd4aadfSBarry Smith     if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
1759566063dSJacob Faibussowitsch       PetscCall(SNESApplyNPC(snes, X, F, D));
1769566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
177b7281c8aSPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
178b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
1793ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
180b7281c8aSPeter Brune       }
18188d374b2SPeter Brune     } else {
1829566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, D));
18388d374b2SPeter Brune     }
18401fe78eaSStefano Zampini 
18501fe78eaSStefano Zampini     /* general purpose update */
186dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
18701fe78eaSStefano Zampini 
18892f76d53SAlp Dener     /* restart conditions */
1890c777b0cSPeter Brune     powell = PETSC_FALSE;
1906bdcc836SBarry Smith     if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) {
1916bdcc836SBarry Smith       /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */
19223c918e7SPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
193574a8f54SStefano Zampini         PetscCall(MatMult(snes->jacobian, Dold, W));
19423c918e7SPeter Brune       } else {
1959566063dSJacob Faibussowitsch         PetscCall(VecCopy(Dold, W));
19623c918e7SPeter Brune       }
1979566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(W, Dold, &DolddotDold));
1989566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(W, D, &DolddotD));
1999566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(W, Dold, &DolddotDold));
2009566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(W, D, &DolddotD));
2010c777b0cSPeter Brune       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma * PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
2020c777b0cSPeter Brune     }
2030c777b0cSPeter Brune     periodic = PETSC_FALSE;
204b8840d6bSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
205b8840d6bSPeter Brune       if (i_r > qn->m - 1) periodic = PETSC_TRUE;
2060c777b0cSPeter Brune     }
2070c777b0cSPeter Brune     /* restart if either powell or periodic restart is satisfied. */
20892f76d53SAlp Dener     if (badstep || powell || periodic) {
20992f76d53SAlp Dener       if (qn->monflg) {
2109566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(qn->monitor, ((PetscObject)snes)->tablevel + 2));
2116bdcc836SBarry Smith         if (powell) {
21263a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(qn->monitor, "Powell restart! |%14.12e| > %6.4f*|%14.12e| i_r = %" PetscInt_FMT "\n", (double)PetscRealPart(DolddotD), (double)qn->powell_gamma, (double)PetscRealPart(DolddotDold), i_r));
2136bdcc836SBarry Smith         } else {
21463a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %" PetscInt_FMT "\n", i_r));
2156bdcc836SBarry Smith         }
2169566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(qn->monitor, ((PetscObject)snes)->tablevel + 2));
2175ba6227bSPeter Brune       }
2185ba6227bSPeter Brune       i_r = -1;
2190c777b0cSPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
2209566063dSJacob Faibussowitsch         PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
22107b62357SFande Kong         SNESCheckJacobianDomainerror(snes);
2220ecc9e1dSPeter Brune       }
2239566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(qn->B, PETSC_FALSE));
2240ecc9e1dSPeter Brune     }
2255ba6227bSPeter Brune   }
2263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2274b11644fSPeter Brune }
2284b11644fSPeter Brune 
229d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_QN(SNES snes)
230d71ae5a4SJacob Faibussowitsch {
2319f83bee8SJed Brown   SNES_QN *qn = (SNES_QN *)snes->data;
232fced5a79SAsbjørn Nilsen Riseth   DM       dm;
23392f76d53SAlp Dener   PetscInt n, N;
234335fb975SPeter Brune 
2354b11644fSPeter Brune   PetscFunctionBegin;
236fced5a79SAsbjørn Nilsen Riseth 
237fced5a79SAsbjørn Nilsen Riseth   if (!snes->vec_sol) {
2389566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes, &dm));
2399566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol));
240fced5a79SAsbjørn Nilsen Riseth   }
2419566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 4));
24292f76d53SAlp Dener 
24348a46eb9SPierre Jolivet   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) PetscCall(SNESSetUpMatrices(snes));
244ad540459SPierre Jolivet   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
24592f76d53SAlp Dener 
24660c08b40SPeter Brune   /* set method defaults */
2471efc8c45SPeter Brune   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
24860c08b40SPeter Brune     if (qn->type == SNES_QN_BADBROYDEN) {
24960c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_NONE;
25060c08b40SPeter Brune     } else {
25192f76d53SAlp Dener       qn->scale_type = SNES_QN_SCALE_SCALAR;
25260c08b40SPeter Brune     }
25360c08b40SPeter Brune   }
2541efc8c45SPeter Brune   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
25560c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
25660c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_POWELL;
25760c08b40SPeter Brune     } else {
25860c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_PERIODIC;
25960c08b40SPeter Brune     }
26060c08b40SPeter Brune   }
26160c08b40SPeter Brune 
26292f76d53SAlp Dener   /* Set up the LMVM matrix */
26392f76d53SAlp Dener   switch (qn->type) {
26492f76d53SAlp Dener   case SNES_QN_BROYDEN:
2659566063dSJacob Faibussowitsch     PetscCall(MatSetType(qn->B, MATLMVMBROYDEN));
26692f76d53SAlp Dener     qn->scale_type = SNES_QN_SCALE_NONE;
26792f76d53SAlp Dener     break;
26892f76d53SAlp Dener   case SNES_QN_BADBROYDEN:
2699566063dSJacob Faibussowitsch     PetscCall(MatSetType(qn->B, MATLMVMBADBROYDEN));
27092f76d53SAlp Dener     qn->scale_type = SNES_QN_SCALE_NONE;
27192f76d53SAlp Dener     break;
27292f76d53SAlp Dener   default:
2739566063dSJacob Faibussowitsch     PetscCall(MatSetType(qn->B, MATLMVMBFGS));
27492f76d53SAlp Dener     switch (qn->scale_type) {
275d71ae5a4SJacob Faibussowitsch     case SNES_QN_SCALE_NONE:
276d71ae5a4SJacob Faibussowitsch       PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_NONE));
277d71ae5a4SJacob Faibussowitsch       break;
278d71ae5a4SJacob Faibussowitsch     case SNES_QN_SCALE_SCALAR:
279d71ae5a4SJacob Faibussowitsch       PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_SCALAR));
280d71ae5a4SJacob Faibussowitsch       break;
281d71ae5a4SJacob Faibussowitsch     case SNES_QN_SCALE_JACOBIAN:
282d71ae5a4SJacob Faibussowitsch       PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_USER));
283d71ae5a4SJacob Faibussowitsch       break;
28492f76d53SAlp Dener     case SNES_QN_SCALE_DIAGONAL:
28592f76d53SAlp Dener     case SNES_QN_SCALE_DEFAULT:
286d71ae5a4SJacob Faibussowitsch     default:
287d71ae5a4SJacob Faibussowitsch       break;
2888e231d97SPeter Brune     }
28992f76d53SAlp Dener     break;
29092f76d53SAlp Dener   }
2919566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
2929566063dSJacob Faibussowitsch   PetscCall(VecGetSize(snes->vec_sol, &N));
2939566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(qn->B, n, n, N, N));
2949566063dSJacob Faibussowitsch   PetscCall(MatSetUp(qn->B));
2959566063dSJacob Faibussowitsch   PetscCall(MatLMVMReset(qn->B, PETSC_TRUE));
2969566063dSJacob Faibussowitsch   PetscCall(MatLMVMSetHistorySize(qn->B, qn->m));
2979566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(qn->B, snes->vec_sol, snes->vec_func));
2983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2994b11644fSPeter Brune }
3004b11644fSPeter Brune 
301d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_QN(SNES snes)
302d71ae5a4SJacob Faibussowitsch {
3039f83bee8SJed Brown   SNES_QN *qn;
3040adebc6cSBarry Smith 
3054b11644fSPeter Brune   PetscFunctionBegin;
3064b11644fSPeter Brune   if (snes->data) {
3079f83bee8SJed Brown     qn = (SNES_QN *)snes->data;
3089566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&qn->B));
3094b11644fSPeter Brune   }
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3114b11644fSPeter Brune }
3124b11644fSPeter Brune 
313d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_QN(SNES snes)
314d71ae5a4SJacob Faibussowitsch {
3154b11644fSPeter Brune   PetscFunctionBegin;
3169566063dSJacob Faibussowitsch   PetscCall(SNESReset_QN(snes));
3179566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
3182e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", NULL));
3192e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", NULL));
3202e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", NULL));
3213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3224b11644fSPeter Brune }
3234b11644fSPeter Brune 
324d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_QN(SNES snes, PetscOptionItems *PetscOptionsObject)
325d71ae5a4SJacob Faibussowitsch {
3262150357eSBarry Smith   SNES_QN          *qn = (SNES_QN *)snes->data;
32792f76d53SAlp Dener   PetscBool         flg;
328f1c6b773SPeter Brune   SNESLineSearch    linesearch;
3292150357eSBarry Smith   SNESQNRestartType rtype = qn->restart_type;
3302150357eSBarry Smith   SNESQNScaleType   stype = qn->scale_type;
3311efc8c45SPeter Brune   SNESQNType        qtype = qn->type;
3322150357eSBarry Smith 
3334b11644fSPeter Brune   PetscFunctionBegin;
334d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES QN options");
3359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, NULL));
3369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL));
3379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", qn->monflg, &qn->monflg, NULL));
3389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_scale_type", "Scaling type", "SNESQNSetScaleType", SNESQNScaleTypes, (PetscEnum)stype, (PetscEnum *)&stype, &flg));
3399566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetScaleType(snes, stype));
34088f769c5SPeter Brune 
3419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_restart_type", "Restart type", "SNESQNSetRestartType", SNESQNRestartTypes, (PetscEnum)rtype, (PetscEnum *)&rtype, &flg));
3429566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetRestartType(snes, rtype));
34388f769c5SPeter Brune 
3449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_type", "Quasi-Newton update type", "", SNESQNTypes, (PetscEnum)qtype, (PetscEnum *)&qtype, &flg));
3459566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetType(snes, qtype));
3469566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(qn->B));
347d0609cedSBarry Smith   PetscOptionsHeadEnd();
3489e764e56SPeter Brune   if (!snes->linesearch) {
3499566063dSJacob Faibussowitsch     PetscCall(SNESGetLineSearch(snes, &linesearch));
350ec786807SJed Brown     if (!((PetscObject)linesearch)->type_name) {
35160c08b40SPeter Brune       if (qn->type == SNES_QN_LBFGS) {
3529566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
35360c08b40SPeter Brune       } else if (qn->type == SNES_QN_BROYDEN) {
3549566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
35560c08b40SPeter Brune       } else {
3569566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
35760c08b40SPeter Brune       }
3589e764e56SPeter Brune     }
359ec786807SJed Brown   }
36048a46eb9SPierre Jolivet   if (qn->monflg) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes), &qn->monitor));
3613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3624b11644fSPeter Brune }
3634b11644fSPeter Brune 
364d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer)
365d71ae5a4SJacob Faibussowitsch {
3665cd83419SPeter Brune   SNES_QN  *qn = (SNES_QN *)snes->data;
3675cd83419SPeter Brune   PetscBool iascii;
3685cd83419SPeter Brune 
3695cd83419SPeter Brune   PetscFunctionBegin;
3709566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3715cd83419SPeter Brune   if (iascii) {
3729566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  type is %s, restart type is %s, scale type is %s\n", SNESQNTypes[qn->type], SNESQNRestartTypes[qn->restart_type], SNESQNScaleTypes[qn->scale_type]));
37363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Stored subspace size: %" PetscInt_FMT "\n", qn->m));
3745cd83419SPeter Brune   }
3753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3765cd83419SPeter Brune }
37792c02d66SPeter Brune 
3780c777b0cSPeter Brune /*@
379f6dfbefdSBarry Smith   SNESQNSetRestartType - Sets the restart type for `SNESQN`.
3800c777b0cSPeter Brune 
381c3339decSBarry Smith   Logically Collective
3820c777b0cSPeter Brune 
3830c777b0cSPeter Brune   Input Parameters:
3840c777b0cSPeter Brune + snes  - the iterative context
385*ceaaa498SBarry Smith - rtype - restart type, see `SNESQNRestartType`
3860c777b0cSPeter Brune 
387f6dfbefdSBarry Smith   Options Database Keys:
3880c777b0cSPeter Brune + -snes_qn_restart_type <powell,periodic,none> - set the restart type
38984c577daSBarry Smith - -snes_qn_m <m>                               - sets the number of stored updates and the restart period for periodic
3900c777b0cSPeter Brune 
3910c777b0cSPeter Brune   Level: intermediate
3920c777b0cSPeter Brune 
393*ceaaa498SBarry Smith .seealso: `SNESQN`, `SNESQNRestartType`, `SNES_QN_RESTART_NONE`, `SNES_QN_RESTART_POWELL`, `SNES_QN_RESTART_PERIODIC`,
394*ceaaa498SBarry Smith           `SNESQNType`, `SNESQNScaleType`
3950c777b0cSPeter Brune @*/
396d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
397d71ae5a4SJacob Faibussowitsch {
3980c777b0cSPeter Brune   PetscFunctionBegin;
3990c777b0cSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
400cac4c232SBarry Smith   PetscTryMethod(snes, "SNESQNSetRestartType_C", (SNES, SNESQNRestartType), (snes, rtype));
4013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4020c777b0cSPeter Brune }
4030c777b0cSPeter Brune 
4040c777b0cSPeter Brune /*@
405f6dfbefdSBarry Smith   SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in `SNESQN`.
4060c777b0cSPeter Brune 
407c3339decSBarry Smith   Logically Collective
4080c777b0cSPeter Brune 
4090c777b0cSPeter Brune   Input Parameters:
410f6dfbefdSBarry Smith + snes  - the nonlinear solver context
411*ceaaa498SBarry Smith - stype - scale type, see `SNESQNScaleType`
4120c777b0cSPeter Brune 
4133c7db156SBarry Smith   Options Database Key:
41467b8a455SSatish Balay . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type
4150c777b0cSPeter Brune 
4160c777b0cSPeter Brune   Level: intermediate
4170c777b0cSPeter Brune 
418*ceaaa498SBarry Smith .seealso: `SNES`, `SNESQN`, `SNESLineSearch`, `SNESQNScaleType`, `SNESSetJacobian()`, `SNESQNType`, `SNESQNRestartType`
4190c777b0cSPeter Brune @*/
420d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
421d71ae5a4SJacob Faibussowitsch {
4220c777b0cSPeter Brune   PetscFunctionBegin;
4230c777b0cSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
424cac4c232SBarry Smith   PetscTryMethod(snes, "SNESQNSetScaleType_C", (SNES, SNESQNScaleType), (snes, stype));
4253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4260c777b0cSPeter Brune }
4270c777b0cSPeter Brune 
428d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
429d71ae5a4SJacob Faibussowitsch {
4300c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
4316e111a19SKarl Rupp 
4320c777b0cSPeter Brune   PetscFunctionBegin;
4330c777b0cSPeter Brune   qn->scale_type = stype;
4347044b745SBarry Smith   if (stype == SNES_QN_SCALE_JACOBIAN) snes->usesksp = PETSC_TRUE;
4353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4360c777b0cSPeter Brune }
4370c777b0cSPeter Brune 
438d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
439d71ae5a4SJacob Faibussowitsch {
4400c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
4416e111a19SKarl Rupp 
4420c777b0cSPeter Brune   PetscFunctionBegin;
4430c777b0cSPeter Brune   qn->restart_type = rtype;
4443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4450c777b0cSPeter Brune }
4460c777b0cSPeter Brune 
4471efc8c45SPeter Brune /*@
448f6dfbefdSBarry Smith   SNESQNSetType - Sets the quasi-Newton variant to be used in `SNESQN`.
4491efc8c45SPeter Brune 
450c3339decSBarry Smith   Logically Collective
4511efc8c45SPeter Brune 
4521efc8c45SPeter Brune   Input Parameters:
4531efc8c45SPeter Brune + snes  - the iterative context
454*ceaaa498SBarry Smith - qtype - variant type, see `SNESQNType`
4551efc8c45SPeter Brune 
4563c7db156SBarry Smith   Options Database Key:
45767b8a455SSatish Balay . -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type
4581efc8c45SPeter Brune 
459*ceaaa498SBarry Smith   Level: intermediate
4601efc8c45SPeter Brune 
461*ceaaa498SBarry Smith .seealso: `SNESQN`, `SNES_QN_LBFGS`, `SNES_QN_BROYDEN`, `SNES_QN_BADBROYDEN`, `SNESQNType`,  `SNESQNScaleType`, `TAOLMVM`, `TAOBLMVM`
4621efc8c45SPeter Brune @*/
463d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype)
464d71ae5a4SJacob Faibussowitsch {
4651efc8c45SPeter Brune   PetscFunctionBegin;
4661efc8c45SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
467cac4c232SBarry Smith   PetscTryMethod(snes, "SNESQNSetType_C", (SNES, SNESQNType), (snes, qtype));
4683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4691efc8c45SPeter Brune }
4701efc8c45SPeter Brune 
471d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype)
472d71ae5a4SJacob Faibussowitsch {
4731efc8c45SPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
4741efc8c45SPeter Brune 
4751efc8c45SPeter Brune   PetscFunctionBegin;
4761efc8c45SPeter Brune   qn->type = qtype;
4773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4781efc8c45SPeter Brune }
4791efc8c45SPeter Brune 
4804b11644fSPeter Brune /*MC
4814b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
4824b11644fSPeter Brune 
483f6dfbefdSBarry Smith       Options Database Keys:
48484c577daSBarry Smith +     -snes_qn_m <m> - Number of past states saved for the L-Broyden methods.
485f6dfbefdSBarry Smith .     -snes_qn_restart_type <powell,periodic,none> - set the restart type
4866bdcc836SBarry Smith .     -snes_qn_powell_gamma - Angle condition for restart.
4871867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
48884c577daSBarry Smith .     -snes_qn_type <lbfgs,broyden,badbroyden> - QN type
48992f76d53SAlp Dener .     -snes_qn_scale_type <diagonal,none,scalar,jacobian> - scaling performed on inner Jacobian
49072835e02SPeter Brune .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
49184c577daSBarry Smith -     -snes_qn_monitor - Monitors the quasi-newton Jacobian.
4926cc8130cSPeter Brune 
4936cc8130cSPeter Brune       References:
494606c0280SSatish Balay +   * -   Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.
495606c0280SSatish Balay .   * -   R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods,
4960a8e8854SPeter Brune       Technical Report, Northwestern University, June 1992.
497606c0280SSatish Balay .   * -   Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
4980a8e8854SPeter Brune       Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
499606c0280SSatish Balay .   * -   Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
5004f02bc6aSBarry Smith        SIAM Review, 57(4), 2015
501606c0280SSatish Balay .   * -   Griewank, Andreas. "Broyden updating, the good and the bad!." Doc. Math (2012): 301-315.
502606c0280SSatish Balay .   * -   Gilbert, Jean Charles, and Claude Lemar{\'e}chal. "Some numerical experiments with variable-storage quasi-Newton algorithms."
50392f76d53SAlp Dener       Mathematical programming 45.1-3 (1989): 407-435.
504606c0280SSatish Balay -   * -   Dener A., Munson T. "Accelerating Limited-Memory Quasi-Newton Convergence for Large-Scale Optimization"
505110fc3b0SBarry Smith       Computational Science - ICCS 2019. ICCS 2019. Lecture Notes in Computer Science, vol 11538. Springer, Cham
5064b11644fSPeter Brune 
5074b11644fSPeter Brune       Level: beginner
5084b11644fSPeter Brune 
509f6dfbefdSBarry Smith       Notes:
510f6dfbefdSBarry Smith       This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
511f6dfbefdSBarry Smith       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
512f6dfbefdSBarry Smith       updates.
5136cc8130cSPeter Brune 
514f6dfbefdSBarry Smith       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
515f6dfbefdSBarry Smith       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
516f6dfbefdSBarry Smith       iteration as the current iteration's values when constructing the approximate Jacobian.  The second, composed,
517f6dfbefdSBarry Smith       perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
518f6dfbefdSBarry Smith 
519f6dfbefdSBarry Smith       Uses left nonlinear preconditioning by default.
520f6dfbefdSBarry Smith 
521f6dfbefdSBarry Smith .seealso: `SNESQNRestartType`, `SNESQNSetRestartType()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`,
522f6dfbefdSBarry Smith           `SNESQNScaleType`, `SNESQNSetScaleType()`, `SNESQNSetType`, `SNESQNSetType()`
5234b11644fSPeter Brune M*/
524d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
525d71ae5a4SJacob Faibussowitsch {
5269f83bee8SJed Brown   SNES_QN    *qn;
52792f76d53SAlp Dener   const char *optionsprefix;
5284b11644fSPeter Brune 
5294b11644fSPeter Brune   PetscFunctionBegin;
5304b11644fSPeter Brune   snes->ops->setup          = SNESSetUp_QN;
5314b11644fSPeter Brune   snes->ops->solve          = SNESSolve_QN;
5324b11644fSPeter Brune   snes->ops->destroy        = SNESDestroy_QN;
5334b11644fSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_QN;
5345cd83419SPeter Brune   snes->ops->view           = SNESView_QN;
5354b11644fSPeter Brune   snes->ops->reset          = SNESReset_QN;
5364b11644fSPeter Brune 
537efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
53847f26062SPeter Brune 
539efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
54042f4f86dSBarry Smith   snes->usesksp = PETSC_FALSE;
54142f4f86dSBarry Smith 
5424fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
5434fc747eaSLawrence Mitchell 
54488976e71SPeter Brune   if (!snes->tolerancesset) {
5450e444f03SPeter Brune     snes->max_funcs = 30000;
5460e444f03SPeter Brune     snes->max_its   = 10000;
54788976e71SPeter Brune   }
5480e444f03SPeter Brune 
5494dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&qn));
5504b11644fSPeter Brune   snes->data       = (void *)qn;
5510ecc9e1dSPeter Brune   qn->m            = 10;
552b21d5a53SPeter Brune   qn->scaling      = 1.0;
5530298fd71SBarry Smith   qn->monitor      = NULL;
55492f76d53SAlp Dener   qn->monflg       = PETSC_FALSE;
555b8840d6bSPeter Brune   qn->powell_gamma = 0.9999;
55660c08b40SPeter Brune   qn->scale_type   = SNES_QN_SCALE_DEFAULT;
55760c08b40SPeter Brune   qn->restart_type = SNES_QN_RESTART_DEFAULT;
558b8840d6bSPeter Brune   qn->type         = SNES_QN_LBFGS;
559ea630c6eSPeter Brune 
5609566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &qn->B));
5619566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
5629566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(qn->B, optionsprefix));
56392f76d53SAlp Dener 
5649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", SNESQNSetScaleType_QN));
5659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", SNESQNSetRestartType_QN));
5669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", SNESQNSetType_QN));
5673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5684b11644fSPeter Brune }
569