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 23*06594da7SStefano Zampini static PetscErrorCode SNESQNGetMatrix_Private(SNES snes, Mat *B) 24*06594da7SStefano Zampini { 25*06594da7SStefano Zampini SNES_QN *qn = (SNES_QN *)snes->data; 26*06594da7SStefano Zampini const char *optionsprefix; 27*06594da7SStefano Zampini 28*06594da7SStefano Zampini PetscFunctionBegin; 29*06594da7SStefano Zampini if (!qn->B) { 30*06594da7SStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &qn->B)); 31*06594da7SStefano Zampini PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 32*06594da7SStefano Zampini PetscCall(MatSetOptionsPrefix(qn->B, optionsprefix)); 33*06594da7SStefano Zampini PetscCall(MatAppendOptionsPrefix(qn->B, "qn_")); 34*06594da7SStefano Zampini switch (qn->type) { 35*06594da7SStefano Zampini case SNES_QN_BROYDEN: 36*06594da7SStefano Zampini PetscCall(MatSetType(qn->B, MATLMVMBROYDEN)); 37*06594da7SStefano Zampini qn->scale_type = SNES_QN_SCALE_NONE; 38*06594da7SStefano Zampini break; 39*06594da7SStefano Zampini case SNES_QN_BADBROYDEN: 40*06594da7SStefano Zampini PetscCall(MatSetType(qn->B, MATLMVMBADBROYDEN)); 41*06594da7SStefano Zampini qn->scale_type = SNES_QN_SCALE_NONE; 42*06594da7SStefano Zampini break; 43*06594da7SStefano Zampini default: 44*06594da7SStefano Zampini PetscCall(MatSetType(qn->B, MATLMVMBFGS)); 45*06594da7SStefano Zampini switch (qn->scale_type) { 46*06594da7SStefano Zampini case SNES_QN_SCALE_NONE: 47*06594da7SStefano Zampini PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_NONE)); 48*06594da7SStefano Zampini break; 49*06594da7SStefano Zampini case SNES_QN_SCALE_SCALAR: 50*06594da7SStefano Zampini case SNES_QN_SCALE_DEFAULT: 51*06594da7SStefano Zampini PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_SCALAR)); 52*06594da7SStefano Zampini break; 53*06594da7SStefano Zampini case SNES_QN_SCALE_JACOBIAN: 54*06594da7SStefano Zampini PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_USER)); 55*06594da7SStefano Zampini break; 56*06594da7SStefano Zampini case SNES_QN_SCALE_DIAGONAL: 57*06594da7SStefano Zampini default: 58*06594da7SStefano Zampini break; 59*06594da7SStefano Zampini } 60*06594da7SStefano Zampini break; 61*06594da7SStefano Zampini } 62*06594da7SStefano Zampini } 63*06594da7SStefano Zampini *B = qn->B; 64*06594da7SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 65*06594da7SStefano Zampini } 66*06594da7SStefano Zampini 67d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_QN(SNES snes) 68d71ae5a4SJacob Faibussowitsch { 699f83bee8SJed Brown SNES_QN *qn = (SNES_QN *)snes->data; 70*06594da7SStefano Zampini Vec X, F, W; 7147f26062SPeter Brune Vec Y, D, Dold; 72b8840d6bSPeter Brune PetscInt i, i_r; 734279555eSSatish Balay PetscReal fnorm, xnorm, ynorm; 74422a814eSBarry Smith SNESLineSearchReason lssucceed; 75*06594da7SStefano Zampini PetscBool badstep, powell, periodic, restart; 761c6523dcSPeter Brune PetscScalar DolddotD, DolddotDold; 77b7281c8aSPeter Brune SNESConvergedReason reason; 784279555eSSatish Balay #if defined(PETSC_USE_INFO) 794279555eSSatish Balay PetscReal gnorm; 804279555eSSatish Balay #endif 810ecc9e1dSPeter Brune 826e111a19SKarl Rupp PetscFunctionBegin; 830b121fc5SBarry 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); 84c579b300SPatrick Farrell 859566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 869f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 873af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 88b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 893af51624SPeter Brune 90*06594da7SStefano Zampini /* directions */ 91*06594da7SStefano Zampini W = snes->work[0]; 92335fb975SPeter Brune D = snes->work[1]; 93335fb975SPeter Brune Dold = snes->work[2]; 944b11644fSPeter Brune 954b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 964b11644fSPeter Brune 979566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 984b11644fSPeter Brune snes->iter = 0; 994b11644fSPeter Brune snes->norm = 0.; 1009566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 10147f26062SPeter Brune 102efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 103*06594da7SStefano Zampini /* we need to sample the preconditioned function only once here, 104*06594da7SStefano Zampini since linesearch will always use the preconditioned function */ 1059566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, F)); 1069566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 107*06594da7SStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 108b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 1093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 110b7281c8aSPeter Brune } 11147f26062SPeter Brune } else { 112*06594da7SStefano Zampini if (!snes->vec_func_init_set) PetscCall(SNESComputeFunction(snes, X, F)); 113*06594da7SStefano Zampini else snes->vec_func_init_set = PETSC_FALSE; 114*06594da7SStefano Zampini } 1159566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 116422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 117b28a06ddSPeter Brune 1189566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 1194b11644fSPeter Brune snes->norm = fnorm; 1209566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1219566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 1229566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 1234b11644fSPeter Brune 1244b11644fSPeter Brune /* test convergence */ 1252d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 1263ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 12770d3b23bSPeter Brune 128*06594da7SStefano Zampini restart = PETSC_TRUE; 129*06594da7SStefano Zampini for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 130*06594da7SStefano Zampini PetscScalar gdx; 131*06594da7SStefano Zampini 132*06594da7SStefano Zampini /* general purpose update */ 133*06594da7SStefano Zampini PetscTryTypeMethod(snes, update, snes->iter); 134*06594da7SStefano Zampini 135*06594da7SStefano Zampini if (snes->npc) { 136*06594da7SStefano Zampini if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 137*06594da7SStefano Zampini PetscCall(SNESApplyNPC(snes, X, F, D)); 138*06594da7SStefano Zampini PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 139*06594da7SStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 140*06594da7SStefano Zampini snes->reason = SNES_DIVERGED_INNER; 141*06594da7SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 142*06594da7SStefano Zampini } 143*06594da7SStefano Zampini } else if (snes->npcside == PC_RIGHT) { 1449566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0)); 1459566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 1469566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0)); 1479566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 148*06594da7SStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 149ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 151ddd40ce5SPeter Brune } 1529566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 1539566063dSJacob Faibussowitsch PetscCall(VecCopy(F, D)); 154*06594da7SStefano Zampini } else { 155*06594da7SStefano Zampini PetscCall(VecCopy(F, D)); 156*06594da7SStefano Zampini } 157*06594da7SStefano Zampini } else { 158*06594da7SStefano Zampini PetscCall(VecCopy(F, D)); 1593cf07b75SPeter Brune } 1603cf07b75SPeter Brune 161f8e15203SPeter Brune /* scale the initial update */ 162*06594da7SStefano Zampini if (qn->scale_type == SNES_QN_SCALE_JACOBIAN && restart) { 1639566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 16407b62357SFande Kong SNESCheckJacobianDomainerror(snes); 1659566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 1669566063dSJacob Faibussowitsch PetscCall(MatLMVMSetJ0KSP(qn->B, snes->ksp)); 1670ecc9e1dSPeter Brune } 1680ecc9e1dSPeter Brune 16992f76d53SAlp Dener /* update QN approx and calculate step */ 1709566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(qn->B, X, D)); 1719566063dSJacob Faibussowitsch PetscCall(MatSolve(qn->B, D, Y)); 172*06594da7SStefano Zampini PetscCall(VecDot(F, Y, &gdx)); 173*06594da7SStefano Zampini if (PetscAbsScalar(gdx) <= 0) { 174*06594da7SStefano Zampini /* Step is not descent or solve was not successful 175*06594da7SStefano Zampini Use steepest descent direction */ 176*06594da7SStefano Zampini PetscCall(VecCopy(F, Y)); 177*06594da7SStefano Zampini PetscCall(MatLMVMReset(qn->B, PETSC_FALSE)); 178*06594da7SStefano Zampini } 17992f76d53SAlp Dener 18070d3b23bSPeter Brune /* line search for lambda */ 1819371c9d4SSatish Balay ynorm = 1; 1824279555eSSatish Balay #if defined(PETSC_USE_INFO) 1839371c9d4SSatish Balay gnorm = fnorm; 1844279555eSSatish Balay #endif 1859566063dSJacob Faibussowitsch PetscCall(VecCopy(D, Dold)); 1869566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y)); 1879f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 1889566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed)); 1899566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm)); 19092f76d53SAlp Dener badstep = PETSC_FALSE; 191422a814eSBarry Smith if (lssucceed) { 1929f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1939f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 1949f3a0142SPeter Brune break; 1959f3a0142SPeter Brune } 19692f76d53SAlp Dener badstep = PETSC_TRUE; 1970ecc9e1dSPeter Brune } 1983af51624SPeter Brune 19988d374b2SPeter Brune /* convergence monitoring */ 2009566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)fnorm, (double)gnorm, (double)ynorm, (int)lssucceed)); 201b21d5a53SPeter Brune 202*06594da7SStefano Zampini PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 203*06594da7SStefano Zampini snes->iter = i + 1; 20471dbe336SPeter Brune snes->norm = fnorm; 205c1e67a49SFande Kong snes->xnorm = xnorm; 206c1e67a49SFande Kong snes->ynorm = ynorm; 207*06594da7SStefano Zampini PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 208360c497dSPeter Brune 2099566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter)); 21092f76d53SAlp Dener 2114b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 2122d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 2132d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 2143ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 215*06594da7SStefano Zampini 216efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 2179566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, D)); 2189566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 219*06594da7SStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 220b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 2213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 222b7281c8aSPeter Brune } 22388d374b2SPeter Brune } else { 2249566063dSJacob Faibussowitsch PetscCall(VecCopy(F, D)); 22588d374b2SPeter Brune } 22601fe78eaSStefano Zampini 22792f76d53SAlp Dener /* restart conditions */ 2280c777b0cSPeter Brune powell = PETSC_FALSE; 2296bdcc836SBarry Smith if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) { 2306bdcc836SBarry Smith /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */ 23123c918e7SPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 232574a8f54SStefano Zampini PetscCall(MatMult(snes->jacobian, Dold, W)); 23323c918e7SPeter Brune } else { 2349566063dSJacob Faibussowitsch PetscCall(VecCopy(Dold, W)); 23523c918e7SPeter Brune } 2369566063dSJacob Faibussowitsch PetscCall(VecDotBegin(W, Dold, &DolddotDold)); 2379566063dSJacob Faibussowitsch PetscCall(VecDotBegin(W, D, &DolddotD)); 2389566063dSJacob Faibussowitsch PetscCall(VecDotEnd(W, Dold, &DolddotDold)); 2399566063dSJacob Faibussowitsch PetscCall(VecDotEnd(W, D, &DolddotD)); 2400c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma * PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 2410c777b0cSPeter Brune } 2420c777b0cSPeter Brune periodic = PETSC_FALSE; 243b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 244b8840d6bSPeter Brune if (i_r > qn->m - 1) periodic = PETSC_TRUE; 2450c777b0cSPeter Brune } 246*06594da7SStefano Zampini 2470c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 248*06594da7SStefano Zampini restart = PETSC_FALSE; 24992f76d53SAlp Dener if (badstep || powell || periodic) { 250*06594da7SStefano Zampini restart = PETSC_TRUE; 25192f76d53SAlp Dener if (qn->monflg) { 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(qn->monitor, ((PetscObject)snes)->tablevel + 2)); 2536bdcc836SBarry Smith if (powell) { 25463a3b9bcSJacob 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)); 2556bdcc836SBarry Smith } else { 25663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %" PetscInt_FMT "\n", i_r)); 2576bdcc836SBarry Smith } 2589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(qn->monitor, ((PetscObject)snes)->tablevel + 2)); 2595ba6227bSPeter Brune } 2605ba6227bSPeter Brune i_r = -1; 2619566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(qn->B, PETSC_FALSE)); 2620ecc9e1dSPeter Brune } 2635ba6227bSPeter Brune } 2643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2654b11644fSPeter Brune } 2664b11644fSPeter Brune 267d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_QN(SNES snes) 268d71ae5a4SJacob Faibussowitsch { 2699f83bee8SJed Brown SNES_QN *qn = (SNES_QN *)snes->data; 270fced5a79SAsbjørn Nilsen Riseth DM dm; 27192f76d53SAlp Dener PetscInt n, N; 272335fb975SPeter Brune 2734b11644fSPeter Brune PetscFunctionBegin; 274fced5a79SAsbjørn Nilsen Riseth if (!snes->vec_sol) { 2759566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 2769566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol)); 277fced5a79SAsbjørn Nilsen Riseth } 278*06594da7SStefano Zampini PetscCall(SNESSetWorkVecs(snes, 3)); 279*06594da7SStefano Zampini PetscCall(SNESQNGetMatrix_Private(snes, &qn->B)); 28092f76d53SAlp Dener 28148a46eb9SPierre Jolivet if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) PetscCall(SNESSetUpMatrices(snes)); 282*06594da7SStefano Zampini if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 28392f76d53SAlp Dener 28460c08b40SPeter Brune /* set method defaults */ 2851efc8c45SPeter Brune if (qn->scale_type == SNES_QN_SCALE_DEFAULT) { 28660c08b40SPeter Brune if (qn->type == SNES_QN_BADBROYDEN) { 28760c08b40SPeter Brune qn->scale_type = SNES_QN_SCALE_NONE; 28860c08b40SPeter Brune } else { 28992f76d53SAlp Dener qn->scale_type = SNES_QN_SCALE_SCALAR; 29060c08b40SPeter Brune } 29160c08b40SPeter Brune } 2921efc8c45SPeter Brune if (qn->restart_type == SNES_QN_RESTART_DEFAULT) { 29360c08b40SPeter Brune if (qn->type == SNES_QN_LBFGS) { 29460c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 29560c08b40SPeter Brune } else { 29660c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_PERIODIC; 29760c08b40SPeter Brune } 29860c08b40SPeter Brune } 29992f76d53SAlp Dener /* Set up the LMVM matrix */ 300*06594da7SStefano Zampini PetscCall(SNESQNGetMatrix_Private(snes, &qn->B)); 3019566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 3029566063dSJacob Faibussowitsch PetscCall(VecGetSize(snes->vec_sol, &N)); 3039566063dSJacob Faibussowitsch PetscCall(MatSetSizes(qn->B, n, n, N, N)); 3049566063dSJacob Faibussowitsch PetscCall(MatSetUp(qn->B)); 3059566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(qn->B, PETSC_TRUE)); 3069566063dSJacob Faibussowitsch PetscCall(MatLMVMSetHistorySize(qn->B, qn->m)); 3079566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(qn->B, snes->vec_sol, snes->vec_func)); 3083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3094b11644fSPeter Brune } 3104b11644fSPeter Brune 311d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_QN(SNES snes) 312d71ae5a4SJacob Faibussowitsch { 313*06594da7SStefano Zampini SNES_QN *qn = (SNES_QN *)snes->data; 3140adebc6cSBarry Smith 3154b11644fSPeter Brune PetscFunctionBegin; 3169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&qn->B)); 3173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3184b11644fSPeter Brune } 3194b11644fSPeter Brune 320d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_QN(SNES snes) 321d71ae5a4SJacob Faibussowitsch { 3224b11644fSPeter Brune PetscFunctionBegin; 3239566063dSJacob Faibussowitsch PetscCall(SNESReset_QN(snes)); 3249566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 3252e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", NULL)); 3262e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", NULL)); 3272e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", NULL)); 3283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3294b11644fSPeter Brune } 3304b11644fSPeter Brune 331d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_QN(SNES snes, PetscOptionItems *PetscOptionsObject) 332d71ae5a4SJacob Faibussowitsch { 3332150357eSBarry Smith SNES_QN *qn = (SNES_QN *)snes->data; 33492f76d53SAlp Dener PetscBool flg; 335f1c6b773SPeter Brune SNESLineSearch linesearch; 3362150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 3372150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 3381efc8c45SPeter Brune SNESQNType qtype = qn->type; 3392150357eSBarry Smith 3404b11644fSPeter Brune PetscFunctionBegin; 341d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES QN options"); 3429566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, NULL)); 3439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL)); 3449566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", qn->monflg, &qn->monflg, NULL)); 3459566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_qn_scale_type", "Scaling type", "SNESQNSetScaleType", SNESQNScaleTypes, (PetscEnum)stype, (PetscEnum *)&stype, &flg)); 3469566063dSJacob Faibussowitsch if (flg) PetscCall(SNESQNSetScaleType(snes, stype)); 34788f769c5SPeter Brune 3489566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_qn_restart_type", "Restart type", "SNESQNSetRestartType", SNESQNRestartTypes, (PetscEnum)rtype, (PetscEnum *)&rtype, &flg)); 3499566063dSJacob Faibussowitsch if (flg) PetscCall(SNESQNSetRestartType(snes, rtype)); 35088f769c5SPeter Brune 3519566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_qn_type", "Quasi-Newton update type", "", SNESQNTypes, (PetscEnum)qtype, (PetscEnum *)&qtype, &flg)); 3529566063dSJacob Faibussowitsch if (flg) PetscCall(SNESQNSetType(snes, qtype)); 353d0609cedSBarry Smith PetscOptionsHeadEnd(); 354*06594da7SStefano Zampini PetscCall(SNESQNGetMatrix_Private(snes, &qn->B)); 355*06594da7SStefano Zampini PetscCall(MatSetFromOptions(qn->B)); 3569e764e56SPeter Brune if (!snes->linesearch) { 3579566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 358ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 35960c08b40SPeter Brune if (qn->type == SNES_QN_LBFGS) { 3609566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP)); 36160c08b40SPeter Brune } else if (qn->type == SNES_QN_BROYDEN) { 3629566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC)); 36360c08b40SPeter Brune } else { 3649566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2)); 36560c08b40SPeter Brune } 3669e764e56SPeter Brune } 367ec786807SJed Brown } 36848a46eb9SPierre Jolivet if (qn->monflg) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes), &qn->monitor)); 3693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3704b11644fSPeter Brune } 3714b11644fSPeter Brune 372d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer) 373d71ae5a4SJacob Faibussowitsch { 3745cd83419SPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 3755cd83419SPeter Brune PetscBool iascii; 3765cd83419SPeter Brune 3775cd83419SPeter Brune PetscFunctionBegin; 3789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 3795cd83419SPeter Brune if (iascii) { 3809566063dSJacob 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])); 38163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Stored subspace size: %" PetscInt_FMT "\n", qn->m)); 382*06594da7SStefano Zampini PetscCall(MatView(qn->B, viewer)); 3835cd83419SPeter Brune } 3843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3855cd83419SPeter Brune } 38692c02d66SPeter Brune 3870c777b0cSPeter Brune /*@ 388f6dfbefdSBarry Smith SNESQNSetRestartType - Sets the restart type for `SNESQN`. 3890c777b0cSPeter Brune 390c3339decSBarry Smith Logically Collective 3910c777b0cSPeter Brune 3920c777b0cSPeter Brune Input Parameters: 3930c777b0cSPeter Brune + snes - the iterative context 394ceaaa498SBarry Smith - rtype - restart type, see `SNESQNRestartType` 3950c777b0cSPeter Brune 396f6dfbefdSBarry Smith Options Database Keys: 3970c777b0cSPeter Brune + -snes_qn_restart_type <powell,periodic,none> - set the restart type 39884c577daSBarry Smith - -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic 3990c777b0cSPeter Brune 4000c777b0cSPeter Brune Level: intermediate 4010c777b0cSPeter Brune 402420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESQN`, `SNESQNRestartType`, `SNES_QN_RESTART_NONE`, `SNES_QN_RESTART_POWELL`, `SNES_QN_RESTART_PERIODIC`, 403ceaaa498SBarry Smith `SNESQNType`, `SNESQNScaleType` 4040c777b0cSPeter Brune @*/ 405d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 406d71ae5a4SJacob Faibussowitsch { 4070c777b0cSPeter Brune PetscFunctionBegin; 4080c777b0cSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 409cac4c232SBarry Smith PetscTryMethod(snes, "SNESQNSetRestartType_C", (SNES, SNESQNRestartType), (snes, rtype)); 4103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4110c777b0cSPeter Brune } 4120c777b0cSPeter Brune 4130c777b0cSPeter Brune /*@ 414f6dfbefdSBarry Smith SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in `SNESQN`. 4150c777b0cSPeter Brune 416c3339decSBarry Smith Logically Collective 4170c777b0cSPeter Brune 4180c777b0cSPeter Brune Input Parameters: 419f6dfbefdSBarry Smith + snes - the nonlinear solver context 420ceaaa498SBarry Smith - stype - scale type, see `SNESQNScaleType` 4210c777b0cSPeter Brune 4223c7db156SBarry Smith Options Database Key: 42367b8a455SSatish Balay . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type 4240c777b0cSPeter Brune 4250c777b0cSPeter Brune Level: intermediate 4260c777b0cSPeter Brune 427420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESQN`, `SNESLineSearch`, `SNESQNScaleType`, `SNESSetJacobian()`, `SNESQNType`, `SNESQNRestartType` 4280c777b0cSPeter Brune @*/ 429d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 430d71ae5a4SJacob Faibussowitsch { 4310c777b0cSPeter Brune PetscFunctionBegin; 4320c777b0cSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 433cac4c232SBarry Smith PetscTryMethod(snes, "SNESQNSetScaleType_C", (SNES, SNESQNScaleType), (snes, stype)); 4343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4350c777b0cSPeter Brune } 4360c777b0cSPeter Brune 43766976f2fSJacob Faibussowitsch static PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 438d71ae5a4SJacob Faibussowitsch { 4390c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 4406e111a19SKarl Rupp 4410c777b0cSPeter Brune PetscFunctionBegin; 4420c777b0cSPeter Brune qn->scale_type = stype; 4437044b745SBarry Smith if (stype == SNES_QN_SCALE_JACOBIAN) snes->usesksp = PETSC_TRUE; 4443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4450c777b0cSPeter Brune } 4460c777b0cSPeter Brune 44766976f2fSJacob Faibussowitsch static PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 448d71ae5a4SJacob Faibussowitsch { 4490c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 4506e111a19SKarl Rupp 4510c777b0cSPeter Brune PetscFunctionBegin; 4520c777b0cSPeter Brune qn->restart_type = rtype; 4533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4540c777b0cSPeter Brune } 4550c777b0cSPeter Brune 4561efc8c45SPeter Brune /*@ 457f6dfbefdSBarry Smith SNESQNSetType - Sets the quasi-Newton variant to be used in `SNESQN`. 4581efc8c45SPeter Brune 459c3339decSBarry Smith Logically Collective 4601efc8c45SPeter Brune 4611efc8c45SPeter Brune Input Parameters: 4621efc8c45SPeter Brune + snes - the iterative context 463ceaaa498SBarry Smith - qtype - variant type, see `SNESQNType` 4641efc8c45SPeter Brune 4653c7db156SBarry Smith Options Database Key: 46667b8a455SSatish Balay . -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type 4671efc8c45SPeter Brune 468ceaaa498SBarry Smith Level: intermediate 4691efc8c45SPeter Brune 470420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESQN`, `SNES_QN_LBFGS`, `SNES_QN_BROYDEN`, `SNES_QN_BADBROYDEN`, `SNESQNType`, `SNESQNScaleType`, `TAOLMVM`, `TAOBLMVM` 4711efc8c45SPeter Brune @*/ 472d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype) 473d71ae5a4SJacob Faibussowitsch { 4741efc8c45SPeter Brune PetscFunctionBegin; 4751efc8c45SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 476cac4c232SBarry Smith PetscTryMethod(snes, "SNESQNSetType_C", (SNES, SNESQNType), (snes, qtype)); 4773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4781efc8c45SPeter Brune } 4791efc8c45SPeter Brune 48066976f2fSJacob Faibussowitsch static PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype) 481d71ae5a4SJacob Faibussowitsch { 4821efc8c45SPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 4831efc8c45SPeter Brune 4841efc8c45SPeter Brune PetscFunctionBegin; 4851efc8c45SPeter Brune qn->type = qtype; 4863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4871efc8c45SPeter Brune } 4881efc8c45SPeter Brune 4894b11644fSPeter Brune /*MC 4904b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 4914b11644fSPeter Brune 492f6dfbefdSBarry Smith Options Database Keys: 49384c577daSBarry Smith + -snes_qn_m <m> - Number of past states saved for the L-Broyden methods. 494f6dfbefdSBarry Smith . -snes_qn_restart_type <powell,periodic,none> - set the restart type 4956bdcc836SBarry Smith . -snes_qn_powell_gamma - Angle condition for restart. 4961867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 49784c577daSBarry Smith . -snes_qn_type <lbfgs,broyden,badbroyden> - QN type 49892f76d53SAlp Dener . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - scaling performed on inner Jacobian 49972835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 50084c577daSBarry Smith - -snes_qn_monitor - Monitors the quasi-newton Jacobian. 5016cc8130cSPeter Brune 5024b11644fSPeter Brune Level: beginner 5034b11644fSPeter Brune 504f6dfbefdSBarry Smith Notes: 505f6dfbefdSBarry Smith This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 506f6dfbefdSBarry Smith previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 507f6dfbefdSBarry Smith updates. 5086cc8130cSPeter Brune 509f6dfbefdSBarry Smith When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 510f6dfbefdSBarry Smith these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 511f6dfbefdSBarry Smith iteration as the current iteration's values when constructing the approximate Jacobian. The second, composed, 512f6dfbefdSBarry Smith perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 513f6dfbefdSBarry Smith 514f6dfbefdSBarry Smith Uses left nonlinear preconditioning by default. 515f6dfbefdSBarry Smith 516420bcc1bSBarry Smith References: 517420bcc1bSBarry Smith + * - Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995. 518420bcc1bSBarry Smith . * - R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods, 519420bcc1bSBarry Smith Technical Report, Northwestern University, June 1992. 520420bcc1bSBarry Smith . * - Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE 521420bcc1bSBarry Smith Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 522420bcc1bSBarry Smith . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 523420bcc1bSBarry Smith SIAM Review, 57(4), 2015 524420bcc1bSBarry Smith . * - Griewank, Andreas. "Broyden updating, the good and the bad!." Doc. Math (2012): 301-315. 525420bcc1bSBarry Smith . * - Gilbert, Jean Charles, and Claude Lemar{\'e}chal. "Some numerical experiments with variable-storage quasi-Newton algorithms." 526420bcc1bSBarry Smith Mathematical programming 45.1-3 (1989): 407-435. 527420bcc1bSBarry Smith - * - Dener A., Munson T. "Accelerating Limited-Memory Quasi-Newton Convergence for Large-Scale Optimization" 528420bcc1bSBarry Smith Computational Science - ICCS 2019. ICCS 2019. Lecture Notes in Computer Science, vol 11538. Springer, Cham 529420bcc1bSBarry Smith 530420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESQNRestartType`, `SNESQNSetRestartType()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, 531420bcc1bSBarry Smith `SNESQNScaleType`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()` 5324b11644fSPeter Brune M*/ 533d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) 534d71ae5a4SJacob Faibussowitsch { 5359f83bee8SJed Brown SNES_QN *qn; 5364b11644fSPeter Brune 5374b11644fSPeter Brune PetscFunctionBegin; 5384b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 5394b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 5404b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 5414b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 5425cd83419SPeter Brune snes->ops->view = SNESView_QN; 5434b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 5444b11644fSPeter Brune 545efd4aadfSBarry Smith snes->npcside = PC_LEFT; 54647f26062SPeter Brune 547efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 54842f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 54942f4f86dSBarry Smith 5504fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 5514fc747eaSLawrence Mitchell 55288976e71SPeter Brune if (!snes->tolerancesset) { 5530e444f03SPeter Brune snes->max_funcs = 30000; 5540e444f03SPeter Brune snes->max_its = 10000; 55588976e71SPeter Brune } 5560e444f03SPeter Brune 5574dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&qn)); 5584b11644fSPeter Brune snes->data = (void *)qn; 5590ecc9e1dSPeter Brune qn->m = 10; 560b21d5a53SPeter Brune qn->scaling = 1.0; 5610298fd71SBarry Smith qn->monitor = NULL; 56292f76d53SAlp Dener qn->monflg = PETSC_FALSE; 563b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 56460c08b40SPeter Brune qn->scale_type = SNES_QN_SCALE_DEFAULT; 56560c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_DEFAULT; 566b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 567ea630c6eSPeter Brune 5689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", SNESQNSetScaleType_QN)); 5699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", SNESQNSetRestartType_QN)); 5709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", SNESQNSetType_QN)); 5713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5724b11644fSPeter Brune } 573