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 2306594da7SStefano Zampini static PetscErrorCode SNESQNGetMatrix_Private(SNES snes, Mat *B) 2406594da7SStefano Zampini { 2506594da7SStefano Zampini SNES_QN *qn = (SNES_QN *)snes->data; 2606594da7SStefano Zampini const char *optionsprefix; 2706594da7SStefano Zampini 2806594da7SStefano Zampini PetscFunctionBegin; 2906594da7SStefano Zampini if (!qn->B) { 3006594da7SStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &qn->B)); 3106594da7SStefano Zampini PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 3206594da7SStefano Zampini PetscCall(MatSetOptionsPrefix(qn->B, optionsprefix)); 3306594da7SStefano Zampini PetscCall(MatAppendOptionsPrefix(qn->B, "qn_")); 3406594da7SStefano Zampini switch (qn->type) { 3506594da7SStefano Zampini case SNES_QN_BROYDEN: 3606594da7SStefano Zampini PetscCall(MatSetType(qn->B, MATLMVMBROYDEN)); 3706594da7SStefano Zampini qn->scale_type = SNES_QN_SCALE_NONE; 3806594da7SStefano Zampini break; 3906594da7SStefano Zampini case SNES_QN_BADBROYDEN: 4006594da7SStefano Zampini PetscCall(MatSetType(qn->B, MATLMVMBADBROYDEN)); 4106594da7SStefano Zampini qn->scale_type = SNES_QN_SCALE_NONE; 4206594da7SStefano Zampini break; 4306594da7SStefano Zampini default: 4406594da7SStefano Zampini PetscCall(MatSetType(qn->B, MATLMVMBFGS)); 4506594da7SStefano Zampini switch (qn->scale_type) { 4606594da7SStefano Zampini case SNES_QN_SCALE_NONE: 475839f191SToby Isaac case SNES_QN_SCALE_JACOBIAN: 4806594da7SStefano Zampini PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_NONE)); 4906594da7SStefano Zampini break; 5006594da7SStefano Zampini case SNES_QN_SCALE_SCALAR: 5106594da7SStefano Zampini case SNES_QN_SCALE_DEFAULT: 5206594da7SStefano Zampini PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_SCALAR)); 5306594da7SStefano Zampini break; 5406594da7SStefano Zampini default: 5506594da7SStefano Zampini break; 5606594da7SStefano Zampini } 5706594da7SStefano Zampini break; 5806594da7SStefano Zampini } 5906594da7SStefano Zampini } 6006594da7SStefano Zampini *B = qn->B; 6106594da7SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 6206594da7SStefano Zampini } 6306594da7SStefano Zampini 64d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_QN(SNES snes) 65d71ae5a4SJacob Faibussowitsch { 669f83bee8SJed Brown SNES_QN *qn = (SNES_QN *)snes->data; 6706594da7SStefano Zampini Vec X, F, W; 6847f26062SPeter Brune Vec Y, D, Dold; 69b8840d6bSPeter Brune PetscInt i, i_r; 704279555eSSatish Balay PetscReal fnorm, xnorm, ynorm; 71422a814eSBarry Smith SNESLineSearchReason lssucceed; 7206594da7SStefano Zampini PetscBool badstep, powell, periodic, restart; 731c6523dcSPeter Brune PetscScalar DolddotD, DolddotDold; 74b7281c8aSPeter Brune SNESConvergedReason reason; 754279555eSSatish Balay #if defined(PETSC_USE_INFO) 764279555eSSatish Balay PetscReal gnorm; 774279555eSSatish Balay #endif 780ecc9e1dSPeter Brune 796e111a19SKarl Rupp PetscFunctionBegin; 800b121fc5SBarry 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); 81c579b300SPatrick Farrell 829566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 839f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 843af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 85b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 863af51624SPeter Brune 8706594da7SStefano Zampini /* directions */ 8806594da7SStefano Zampini W = snes->work[0]; 89335fb975SPeter Brune D = snes->work[1]; 90335fb975SPeter Brune Dold = snes->work[2]; 914b11644fSPeter Brune 924b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 934b11644fSPeter Brune 949566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 954b11644fSPeter Brune snes->iter = 0; 964b11644fSPeter Brune snes->norm = 0.; 979566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 9847f26062SPeter Brune 99efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 10006594da7SStefano Zampini /* we need to sample the preconditioned function only once here, 10106594da7SStefano Zampini since linesearch will always use the preconditioned function */ 1029566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, F)); 1039566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 10406594da7SStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 105b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 1063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 107b7281c8aSPeter Brune } 10847f26062SPeter Brune } else { 10906594da7SStefano Zampini if (!snes->vec_func_init_set) PetscCall(SNESComputeFunction(snes, X, F)); 11006594da7SStefano Zampini else snes->vec_func_init_set = PETSC_FALSE; 11106594da7SStefano Zampini } 1129566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 113422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 114b28a06ddSPeter Brune 1159566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 1164b11644fSPeter Brune snes->norm = fnorm; 1179566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1189566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 1199566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 1204b11644fSPeter Brune 1214b11644fSPeter Brune /* test convergence */ 1222d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 1233ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 12470d3b23bSPeter Brune 12506594da7SStefano Zampini restart = PETSC_TRUE; 12606594da7SStefano Zampini for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 12706594da7SStefano Zampini PetscScalar gdx; 12806594da7SStefano Zampini 12906594da7SStefano Zampini /* general purpose update */ 13006594da7SStefano Zampini PetscTryTypeMethod(snes, update, snes->iter); 13106594da7SStefano Zampini 13206594da7SStefano Zampini if (snes->npc) { 13306594da7SStefano Zampini if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 13406594da7SStefano Zampini PetscCall(SNESApplyNPC(snes, X, F, D)); 13506594da7SStefano Zampini PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 13606594da7SStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 13706594da7SStefano Zampini snes->reason = SNES_DIVERGED_INNER; 13806594da7SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 13906594da7SStefano Zampini } 14006594da7SStefano Zampini } else if (snes->npcside == PC_RIGHT) { 1419566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0)); 1429566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 1439566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0)); 1449566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 14506594da7SStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 146ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 148ddd40ce5SPeter Brune } 1499566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 1509566063dSJacob Faibussowitsch PetscCall(VecCopy(F, D)); 15106594da7SStefano Zampini } else { 15206594da7SStefano Zampini PetscCall(VecCopy(F, D)); 15306594da7SStefano Zampini } 15406594da7SStefano Zampini } else { 15506594da7SStefano Zampini PetscCall(VecCopy(F, D)); 1563cf07b75SPeter Brune } 1573cf07b75SPeter Brune 158f8e15203SPeter Brune /* scale the initial update */ 15906594da7SStefano Zampini if (qn->scale_type == SNES_QN_SCALE_JACOBIAN && restart) { 1609566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 16107b62357SFande Kong SNESCheckJacobianDomainerror(snes); 1629566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 1639566063dSJacob Faibussowitsch PetscCall(MatLMVMSetJ0KSP(qn->B, snes->ksp)); 1640ecc9e1dSPeter Brune } 1650ecc9e1dSPeter Brune 16692f76d53SAlp Dener /* update QN approx and calculate step */ 1679566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(qn->B, X, D)); 1689566063dSJacob Faibussowitsch PetscCall(MatSolve(qn->B, D, Y)); 16906594da7SStefano Zampini PetscCall(VecDot(F, Y, &gdx)); 17006594da7SStefano Zampini if (PetscAbsScalar(gdx) <= 0) { 17106594da7SStefano Zampini /* Step is not descent or solve was not successful 17206594da7SStefano Zampini Use steepest descent direction */ 17306594da7SStefano Zampini PetscCall(VecCopy(F, Y)); 17406594da7SStefano Zampini PetscCall(MatLMVMReset(qn->B, PETSC_FALSE)); 17506594da7SStefano Zampini } 17692f76d53SAlp Dener 17770d3b23bSPeter Brune /* line search for lambda */ 1789371c9d4SSatish Balay ynorm = 1; 1794279555eSSatish Balay #if defined(PETSC_USE_INFO) 1809371c9d4SSatish Balay gnorm = fnorm; 1814279555eSSatish Balay #endif 1829566063dSJacob Faibussowitsch PetscCall(VecCopy(D, Dold)); 1839566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y)); 1849f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 1859566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed)); 1869566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm)); 18792f76d53SAlp Dener badstep = PETSC_FALSE; 188422a814eSBarry Smith if (lssucceed) { 1899f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1909f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 1919f3a0142SPeter Brune break; 1929f3a0142SPeter Brune } 19392f76d53SAlp Dener badstep = PETSC_TRUE; 1940ecc9e1dSPeter Brune } 1953af51624SPeter Brune 19688d374b2SPeter Brune /* convergence monitoring */ 1979566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)fnorm, (double)gnorm, (double)ynorm, (int)lssucceed)); 198b21d5a53SPeter Brune 19906594da7SStefano Zampini PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 20006594da7SStefano Zampini snes->iter = i + 1; 20171dbe336SPeter Brune snes->norm = fnorm; 202c1e67a49SFande Kong snes->xnorm = xnorm; 203c1e67a49SFande Kong snes->ynorm = ynorm; 20406594da7SStefano Zampini PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 205360c497dSPeter Brune 2069566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter)); 20792f76d53SAlp Dener 2084b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 2092d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 2102d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 2113ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 21206594da7SStefano Zampini 213efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 2149566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, D)); 2159566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 21606594da7SStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 217b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 219b7281c8aSPeter Brune } 22088d374b2SPeter Brune } else { 2219566063dSJacob Faibussowitsch PetscCall(VecCopy(F, D)); 22288d374b2SPeter Brune } 22301fe78eaSStefano Zampini 22492f76d53SAlp Dener /* restart conditions */ 2250c777b0cSPeter Brune powell = PETSC_FALSE; 2266bdcc836SBarry Smith if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) { 2276bdcc836SBarry Smith /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */ 22823c918e7SPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 229574a8f54SStefano Zampini PetscCall(MatMult(snes->jacobian, Dold, W)); 23023c918e7SPeter Brune } else { 2319566063dSJacob Faibussowitsch PetscCall(VecCopy(Dold, W)); 23223c918e7SPeter Brune } 2339566063dSJacob Faibussowitsch PetscCall(VecDotBegin(W, Dold, &DolddotDold)); 2349566063dSJacob Faibussowitsch PetscCall(VecDotBegin(W, D, &DolddotD)); 2359566063dSJacob Faibussowitsch PetscCall(VecDotEnd(W, Dold, &DolddotDold)); 2369566063dSJacob Faibussowitsch PetscCall(VecDotEnd(W, D, &DolddotD)); 2370c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma * PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 2380c777b0cSPeter Brune } 2390c777b0cSPeter Brune periodic = PETSC_FALSE; 240b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 241b8840d6bSPeter Brune if (i_r > qn->m - 1) periodic = PETSC_TRUE; 2420c777b0cSPeter Brune } 24306594da7SStefano Zampini 2440c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 24506594da7SStefano Zampini restart = PETSC_FALSE; 24692f76d53SAlp Dener if (badstep || powell || periodic) { 24706594da7SStefano Zampini restart = PETSC_TRUE; 24892f76d53SAlp Dener if (qn->monflg) { 2499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(qn->monitor, ((PetscObject)snes)->tablevel + 2)); 2506bdcc836SBarry Smith if (powell) { 25163a3b9bcSJacob 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)); 2526bdcc836SBarry Smith } else { 25363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %" PetscInt_FMT "\n", i_r)); 2546bdcc836SBarry Smith } 2559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(qn->monitor, ((PetscObject)snes)->tablevel + 2)); 2565ba6227bSPeter Brune } 2575ba6227bSPeter Brune i_r = -1; 2589566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(qn->B, PETSC_FALSE)); 2590ecc9e1dSPeter Brune } 2605ba6227bSPeter Brune } 2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2624b11644fSPeter Brune } 2634b11644fSPeter Brune 264d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_QN(SNES snes) 265d71ae5a4SJacob Faibussowitsch { 2669f83bee8SJed Brown SNES_QN *qn = (SNES_QN *)snes->data; 267fced5a79SAsbjørn Nilsen Riseth DM dm; 26892f76d53SAlp Dener PetscInt n, N; 269335fb975SPeter Brune 2704b11644fSPeter Brune PetscFunctionBegin; 271fced5a79SAsbjørn Nilsen Riseth if (!snes->vec_sol) { 2729566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 2739566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol)); 274fced5a79SAsbjørn Nilsen Riseth } 27506594da7SStefano Zampini PetscCall(SNESSetWorkVecs(snes, 3)); 27606594da7SStefano Zampini PetscCall(SNESQNGetMatrix_Private(snes, &qn->B)); 27792f76d53SAlp Dener 27848a46eb9SPierre Jolivet if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) PetscCall(SNESSetUpMatrices(snes)); 27906594da7SStefano Zampini if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 28092f76d53SAlp Dener 28160c08b40SPeter Brune /* set method defaults */ 2821efc8c45SPeter Brune if (qn->scale_type == SNES_QN_SCALE_DEFAULT) { 28360c08b40SPeter Brune if (qn->type == SNES_QN_BADBROYDEN) { 28460c08b40SPeter Brune qn->scale_type = SNES_QN_SCALE_NONE; 28560c08b40SPeter Brune } else { 28692f76d53SAlp Dener qn->scale_type = SNES_QN_SCALE_SCALAR; 28760c08b40SPeter Brune } 28860c08b40SPeter Brune } 2891efc8c45SPeter Brune if (qn->restart_type == SNES_QN_RESTART_DEFAULT) { 29060c08b40SPeter Brune if (qn->type == SNES_QN_LBFGS) { 29160c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 29260c08b40SPeter Brune } else { 29360c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_PERIODIC; 29460c08b40SPeter Brune } 29560c08b40SPeter Brune } 29692f76d53SAlp Dener /* Set up the LMVM matrix */ 29706594da7SStefano Zampini PetscCall(SNESQNGetMatrix_Private(snes, &qn->B)); 2989566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 2999566063dSJacob Faibussowitsch PetscCall(VecGetSize(snes->vec_sol, &N)); 3009566063dSJacob Faibussowitsch PetscCall(MatSetSizes(qn->B, n, n, N, N)); 3019566063dSJacob Faibussowitsch PetscCall(MatSetUp(qn->B)); 3029566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(qn->B, PETSC_TRUE)); 3039566063dSJacob Faibussowitsch PetscCall(MatLMVMSetHistorySize(qn->B, qn->m)); 3049566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(qn->B, snes->vec_sol, snes->vec_func)); 3053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3064b11644fSPeter Brune } 3074b11644fSPeter Brune 308d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_QN(SNES snes) 309d71ae5a4SJacob Faibussowitsch { 31006594da7SStefano Zampini SNES_QN *qn = (SNES_QN *)snes->data; 3110adebc6cSBarry Smith 3124b11644fSPeter Brune PetscFunctionBegin; 3139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&qn->B)); 3143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3154b11644fSPeter Brune } 3164b11644fSPeter Brune 317d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_QN(SNES snes) 318d71ae5a4SJacob Faibussowitsch { 3194b11644fSPeter Brune PetscFunctionBegin; 3209566063dSJacob Faibussowitsch PetscCall(SNESReset_QN(snes)); 3219566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 3222e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", NULL)); 3232e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", NULL)); 3242e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", NULL)); 3253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3264b11644fSPeter Brune } 3274b11644fSPeter Brune 328ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_QN(SNES snes, PetscOptionItems PetscOptionsObject) 329d71ae5a4SJacob Faibussowitsch { 3302150357eSBarry Smith SNES_QN *qn = (SNES_QN *)snes->data; 33192f76d53SAlp Dener PetscBool flg; 332f1c6b773SPeter Brune SNESLineSearch linesearch; 3332150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 3342150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 3351efc8c45SPeter Brune SNESQNType qtype = qn->type; 3362150357eSBarry Smith 3374b11644fSPeter Brune PetscFunctionBegin; 338d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES QN options"); 3399566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, NULL)); 3409566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL)); 3419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", qn->monflg, &qn->monflg, NULL)); 3429566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_qn_scale_type", "Scaling type", "SNESQNSetScaleType", SNESQNScaleTypes, (PetscEnum)stype, (PetscEnum *)&stype, &flg)); 3439566063dSJacob Faibussowitsch if (flg) PetscCall(SNESQNSetScaleType(snes, stype)); 34488f769c5SPeter Brune 3459566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_qn_restart_type", "Restart type", "SNESQNSetRestartType", SNESQNRestartTypes, (PetscEnum)rtype, (PetscEnum *)&rtype, &flg)); 3469566063dSJacob Faibussowitsch if (flg) PetscCall(SNESQNSetRestartType(snes, rtype)); 34788f769c5SPeter Brune 3489566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_qn_type", "Quasi-Newton update type", "", SNESQNTypes, (PetscEnum)qtype, (PetscEnum *)&qtype, &flg)); 3499566063dSJacob Faibussowitsch if (flg) PetscCall(SNESQNSetType(snes, qtype)); 350d0609cedSBarry Smith PetscOptionsHeadEnd(); 35106594da7SStefano Zampini PetscCall(SNESQNGetMatrix_Private(snes, &qn->B)); 35206594da7SStefano Zampini PetscCall(MatSetFromOptions(qn->B)); 3539e764e56SPeter Brune if (!snes->linesearch) { 3549566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 355ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 35660c08b40SPeter Brune if (qn->type == SNES_QN_LBFGS) { 3579566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP)); 35860c08b40SPeter Brune } else if (qn->type == SNES_QN_BROYDEN) { 3599566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC)); 36060c08b40SPeter Brune } else { 361a99ef635SJonas Heinzmann PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHSECANT)); 36260c08b40SPeter Brune } 3639e764e56SPeter Brune } 364ec786807SJed Brown } 36548a46eb9SPierre Jolivet if (qn->monflg) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes), &qn->monitor)); 3663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3674b11644fSPeter Brune } 3684b11644fSPeter Brune 369d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer) 370d71ae5a4SJacob Faibussowitsch { 3715cd83419SPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 372*9f196a02SMartin Diehl PetscBool isascii; 3735cd83419SPeter Brune 3745cd83419SPeter Brune PetscFunctionBegin; 375*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 376*9f196a02SMartin Diehl if (isascii) { 3779566063dSJacob 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])); 37863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Stored subspace size: %" PetscInt_FMT "\n", qn->m)); 37906594da7SStefano Zampini PetscCall(MatView(qn->B, viewer)); 3805cd83419SPeter Brune } 3813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3825cd83419SPeter Brune } 38392c02d66SPeter Brune 3840c777b0cSPeter Brune /*@ 385f6dfbefdSBarry Smith SNESQNSetRestartType - Sets the restart type for `SNESQN`. 3860c777b0cSPeter Brune 387c3339decSBarry Smith Logically Collective 3880c777b0cSPeter Brune 3890c777b0cSPeter Brune Input Parameters: 3900c777b0cSPeter Brune + snes - the iterative context 391ceaaa498SBarry Smith - rtype - restart type, see `SNESQNRestartType` 3920c777b0cSPeter Brune 393f6dfbefdSBarry Smith Options Database Keys: 3940c777b0cSPeter Brune + -snes_qn_restart_type <powell,periodic,none> - set the restart type 39584c577daSBarry Smith - -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic 3960c777b0cSPeter Brune 3970c777b0cSPeter Brune Level: intermediate 3980c777b0cSPeter Brune 399420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESQN`, `SNESQNRestartType`, `SNES_QN_RESTART_NONE`, `SNES_QN_RESTART_POWELL`, `SNES_QN_RESTART_PERIODIC`, 400ceaaa498SBarry Smith `SNESQNType`, `SNESQNScaleType` 4010c777b0cSPeter Brune @*/ 402d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 403d71ae5a4SJacob Faibussowitsch { 4040c777b0cSPeter Brune PetscFunctionBegin; 4050c777b0cSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 406cac4c232SBarry Smith PetscTryMethod(snes, "SNESQNSetRestartType_C", (SNES, SNESQNRestartType), (snes, rtype)); 4073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4080c777b0cSPeter Brune } 4090c777b0cSPeter Brune 4100c777b0cSPeter Brune /*@ 411f6dfbefdSBarry Smith SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in `SNESQN`. 4120c777b0cSPeter Brune 413c3339decSBarry Smith Logically Collective 4140c777b0cSPeter Brune 4150c777b0cSPeter Brune Input Parameters: 416f6dfbefdSBarry Smith + snes - the nonlinear solver context 417ceaaa498SBarry Smith - stype - scale type, see `SNESQNScaleType` 4180c777b0cSPeter Brune 4193c7db156SBarry Smith Options Database Key: 42067b8a455SSatish Balay . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type 4210c777b0cSPeter Brune 4220c777b0cSPeter Brune Level: intermediate 4230c777b0cSPeter Brune 424420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESQN`, `SNESLineSearch`, `SNESQNScaleType`, `SNESSetJacobian()`, `SNESQNType`, `SNESQNRestartType` 4250c777b0cSPeter Brune @*/ 426d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 427d71ae5a4SJacob Faibussowitsch { 4280c777b0cSPeter Brune PetscFunctionBegin; 4290c777b0cSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 430cac4c232SBarry Smith PetscTryMethod(snes, "SNESQNSetScaleType_C", (SNES, SNESQNScaleType), (snes, stype)); 4313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4320c777b0cSPeter Brune } 4330c777b0cSPeter Brune 43466976f2fSJacob Faibussowitsch static PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 435d71ae5a4SJacob Faibussowitsch { 4360c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 4376e111a19SKarl Rupp 4380c777b0cSPeter Brune PetscFunctionBegin; 4390c777b0cSPeter Brune qn->scale_type = stype; 4407044b745SBarry Smith if (stype == SNES_QN_SCALE_JACOBIAN) snes->usesksp = PETSC_TRUE; 4413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4420c777b0cSPeter Brune } 4430c777b0cSPeter Brune 44466976f2fSJacob Faibussowitsch static PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 445d71ae5a4SJacob Faibussowitsch { 4460c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 4476e111a19SKarl Rupp 4480c777b0cSPeter Brune PetscFunctionBegin; 4490c777b0cSPeter Brune qn->restart_type = rtype; 4503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4510c777b0cSPeter Brune } 4520c777b0cSPeter Brune 4531efc8c45SPeter Brune /*@ 454f6dfbefdSBarry Smith SNESQNSetType - Sets the quasi-Newton variant to be used in `SNESQN`. 4551efc8c45SPeter Brune 456c3339decSBarry Smith Logically Collective 4571efc8c45SPeter Brune 4581efc8c45SPeter Brune Input Parameters: 4591efc8c45SPeter Brune + snes - the iterative context 460ceaaa498SBarry Smith - qtype - variant type, see `SNESQNType` 4611efc8c45SPeter Brune 4623c7db156SBarry Smith Options Database Key: 46367b8a455SSatish Balay . -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type 4641efc8c45SPeter Brune 465ceaaa498SBarry Smith Level: intermediate 4661efc8c45SPeter Brune 467420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESQN`, `SNES_QN_LBFGS`, `SNES_QN_BROYDEN`, `SNES_QN_BADBROYDEN`, `SNESQNType`, `SNESQNScaleType`, `TAOLMVM`, `TAOBLMVM` 4681efc8c45SPeter Brune @*/ 469d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype) 470d71ae5a4SJacob Faibussowitsch { 4711efc8c45SPeter Brune PetscFunctionBegin; 4721efc8c45SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 473cac4c232SBarry Smith PetscTryMethod(snes, "SNESQNSetType_C", (SNES, SNESQNType), (snes, qtype)); 4743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4751efc8c45SPeter Brune } 4761efc8c45SPeter Brune 47766976f2fSJacob Faibussowitsch static PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype) 478d71ae5a4SJacob Faibussowitsch { 4791efc8c45SPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 4801efc8c45SPeter Brune 4811efc8c45SPeter Brune PetscFunctionBegin; 4821efc8c45SPeter Brune qn->type = qtype; 4833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4841efc8c45SPeter Brune } 4851efc8c45SPeter Brune 4864b11644fSPeter Brune /*MC 4874b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 4884b11644fSPeter Brune 489f6dfbefdSBarry Smith Options Database Keys: 49084c577daSBarry Smith + -snes_qn_m <m> - Number of past states saved for the L-Broyden methods. 491f6dfbefdSBarry Smith . -snes_qn_restart_type <powell,periodic,none> - set the restart type 4926bdcc836SBarry Smith . -snes_qn_powell_gamma - Angle condition for restart. 4931867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 49484c577daSBarry Smith . -snes_qn_type <lbfgs,broyden,badbroyden> - QN type 49592f76d53SAlp Dener . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - scaling performed on inner Jacobian 49672835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 49784c577daSBarry Smith - -snes_qn_monitor - Monitors the quasi-newton Jacobian. 4986cc8130cSPeter Brune 4994b11644fSPeter Brune Level: beginner 5004b11644fSPeter Brune 501f6dfbefdSBarry Smith Notes: 5021d27aa22SBarry Smith This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of $F(x) = b$ using 5031d27aa22SBarry Smith previous change in $F(x)$ and $x$ to form the approximate inverse Jacobian using a series of multiplicative rank-one 504f6dfbefdSBarry Smith updates. 5056cc8130cSPeter Brune 506f6dfbefdSBarry Smith When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 507f6dfbefdSBarry Smith these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 508f6dfbefdSBarry Smith iteration as the current iteration's values when constructing the approximate Jacobian. The second, composed, 5091d27aa22SBarry Smith perturbs the problem the Jacobian represents to be $P(x, b) - x = 0 $, where $P(x, b)$ is the preconditioner. 510f6dfbefdSBarry Smith 511f6dfbefdSBarry Smith Uses left nonlinear preconditioning by default. 512f6dfbefdSBarry Smith 5131d27aa22SBarry Smith See {cite}`byrd1994representations`, {cite}`kelley95`, {cite}`brown1985experiments`, {cite}`griewank2012broyden`, {cite}`gilbert1989some`, 5141d27aa22SBarry Smith {cite}`dener2019accelerating`, and {cite}`bruneknepleysmithtu15` 515420bcc1bSBarry Smith 516420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESQNRestartType`, `SNESQNSetRestartType()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, 517420bcc1bSBarry Smith `SNESQNScaleType`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()` 5184b11644fSPeter Brune M*/ 519d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) 520d71ae5a4SJacob Faibussowitsch { 5219f83bee8SJed Brown SNES_QN *qn; 5224b11644fSPeter Brune 5234b11644fSPeter Brune PetscFunctionBegin; 5244b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 5254b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 5264b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 5274b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 5285cd83419SPeter Brune snes->ops->view = SNESView_QN; 5294b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 5304b11644fSPeter Brune 531efd4aadfSBarry Smith snes->npcside = PC_LEFT; 53247f26062SPeter Brune 533efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 53442f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 53542f4f86dSBarry Smith 5364fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 5374fc747eaSLawrence Mitchell 53877e5a1f9SBarry Smith PetscCall(SNESParametersInitialize(snes)); 53977e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, max_funcs, 30000); 54077e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, max_its, 10000); 5410e444f03SPeter Brune 5424dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&qn)); 5434b11644fSPeter Brune snes->data = (void *)qn; 5440ecc9e1dSPeter Brune qn->m = 10; 545b21d5a53SPeter Brune qn->scaling = 1.0; 5460298fd71SBarry Smith qn->monitor = NULL; 54792f76d53SAlp Dener qn->monflg = PETSC_FALSE; 548b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 54960c08b40SPeter Brune qn->scale_type = SNES_QN_SCALE_DEFAULT; 55060c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_DEFAULT; 551b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 552ea630c6eSPeter Brune 5539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", SNESQNSetScaleType_QN)); 5549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", SNESQNSetRestartType_QN)); 5559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", SNESQNSetType_QN)); 5563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5574b11644fSPeter Brune } 558