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 234b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 244b11644fSPeter Brune { 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; 30335fb975SPeter Brune PetscReal fnorm,xnorm,ynorm,gnorm; 31422a814eSBarry Smith SNESLineSearchReason lssucceed; 3292f76d53SAlp Dener PetscBool badstep,powell,periodic; 331c6523dcSPeter Brune PetscScalar DolddotD,DolddotDold; 34b7281c8aSPeter Brune SNESConvergedReason reason; 350ecc9e1dSPeter Brune 3684c577daSBarry Smith /* basically just a regular newton's method except for the application of the Jacobian */ 374b11644fSPeter Brune 386e111a19SKarl Rupp PetscFunctionBegin; 390b121fc5SBarry 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); 40c579b300SPatrick Farrell 419566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite)); 429f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 433af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 440a3905e1SPeter Brune W = snes->work[3]; 45b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 46335fb975SPeter Brune Xold = snes->work[0]; 473af51624SPeter Brune 483af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 49335fb975SPeter Brune D = snes->work[1]; 50335fb975SPeter Brune Dold = snes->work[2]; 514b11644fSPeter Brune 524b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 534b11644fSPeter Brune 549566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 554b11644fSPeter Brune snes->iter = 0; 564b11644fSPeter Brune snes->norm = 0.; 579566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 5847f26062SPeter Brune 59efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 609566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes,X,NULL,F)); 619566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc,&reason)); 62b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 63b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 64b7281c8aSPeter Brune PetscFunctionReturn(0); 65b7281c8aSPeter Brune } 669566063dSJacob Faibussowitsch PetscCall(VecNorm(F,NORM_2,&fnorm)); 6747f26062SPeter Brune } else { 68e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 699566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 701aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 71c1c75074SPeter Brune 729566063dSJacob Faibussowitsch PetscCall(VecNorm(F,NORM_2,&fnorm)); 73422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 7447f26062SPeter Brune } 75efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 769566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes,X,F,D)); 779566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc,&reason)); 78b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 79b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 80b7281c8aSPeter Brune PetscFunctionReturn(0); 81b7281c8aSPeter Brune } 8247f26062SPeter Brune } else { 839566063dSJacob Faibussowitsch PetscCall(VecCopy(F,D)); 8447f26062SPeter Brune } 85b28a06ddSPeter Brune 869566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 874b11644fSPeter Brune snes->norm = fnorm; 889566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 899566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,fnorm,0)); 909566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,0,fnorm)); 914b11644fSPeter Brune 924b11644fSPeter Brune /* test convergence */ 939566063dSJacob Faibussowitsch PetscCall((*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 944b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 9570d3b23bSPeter Brune 96efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_RIGHT) { 979566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve,snes->npc,X,0,0)); 989566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc,snes->vec_rhs,X)); 999566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve,snes->npc,X,0,0)); 1009566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc,&reason)); 101ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 102ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 103ddd40ce5SPeter Brune PetscFunctionReturn(0); 104ddd40ce5SPeter Brune } 1059566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes,F,&fnorm)); 1069566063dSJacob Faibussowitsch PetscCall(VecCopy(F,D)); 1073cf07b75SPeter Brune } 1083cf07b75SPeter Brune 10901fe78eaSStefano Zampini /* general purpose update */ 11001fe78eaSStefano Zampini if (snes->ops->update) { 1119566063dSJacob Faibussowitsch PetscCall((*snes->ops->update)(snes, snes->iter)); 11201fe78eaSStefano Zampini } 11301fe78eaSStefano Zampini 114f8e15203SPeter Brune /* scale the initial update */ 1150c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 1169566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre)); 11707b62357SFande Kong SNESCheckJacobianDomainerror(snes); 1189566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre)); 1199566063dSJacob Faibussowitsch PetscCall(MatLMVMSetJ0KSP(qn->B, snes->ksp)); 1200ecc9e1dSPeter Brune } 1210ecc9e1dSPeter Brune 1225ba6227bSPeter Brune for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 12392f76d53SAlp Dener /* update QN approx and calculate step */ 1249566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(qn->B, X, D)); 1259566063dSJacob Faibussowitsch PetscCall(MatSolve(qn->B, D, Y)); 12692f76d53SAlp Dener 12770d3b23bSPeter Brune /* line search for lambda */ 12870d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 1299566063dSJacob Faibussowitsch PetscCall(VecCopy(D, Dold)); 1309566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xold)); 1319566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y)); 1329f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 1339566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed)); 1349566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm)); 13592f76d53SAlp Dener badstep = PETSC_FALSE; 136422a814eSBarry Smith if (lssucceed) { 1379f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1389f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 1399f3a0142SPeter Brune break; 1409f3a0142SPeter Brune } 14192f76d53SAlp Dener badstep = PETSC_TRUE; 1420ecc9e1dSPeter Brune } 1433af51624SPeter Brune 14488d374b2SPeter Brune /* convergence monitoring */ 1459566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed)); 146b21d5a53SPeter Brune 147efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_RIGHT) { 1489566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve,snes->npc,X,0,0)); 1499566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc,snes->vec_rhs,X)); 1509566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve,snes->npc,X,0,0)); 1519566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc,&reason)); 152ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 153ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 154ddd40ce5SPeter Brune PetscFunctionReturn(0); 155ddd40ce5SPeter Brune } 1569566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes,F,&fnorm)); 157b28a06ddSPeter Brune } 158b28a06ddSPeter Brune 1599566063dSJacob Faibussowitsch PetscCall(SNESSetIterationNumber(snes, i+1)); 16071dbe336SPeter Brune snes->norm = fnorm; 161c1e67a49SFande Kong snes->xnorm = xnorm; 162c1e67a49SFande Kong snes->ynorm = ynorm; 163360c497dSPeter Brune 1649566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,snes->iter)); 1659566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,snes->iter,snes->norm)); 16692f76d53SAlp Dener 1674b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1689566063dSJacob Faibussowitsch PetscCall((*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP)); 1694b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 170efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 1719566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes,X,F,D)); 1729566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc,&reason)); 173b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 174b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 175b7281c8aSPeter Brune PetscFunctionReturn(0); 176b7281c8aSPeter Brune } 17788d374b2SPeter Brune } else { 1789566063dSJacob Faibussowitsch PetscCall(VecCopy(F, D)); 17988d374b2SPeter Brune } 18001fe78eaSStefano Zampini 18101fe78eaSStefano Zampini /* general purpose update */ 18201fe78eaSStefano Zampini if (snes->ops->update) { 1839566063dSJacob Faibussowitsch PetscCall((*snes->ops->update)(snes, snes->iter)); 18401fe78eaSStefano Zampini } 18501fe78eaSStefano Zampini 18692f76d53SAlp Dener /* restart conditions */ 1870c777b0cSPeter Brune powell = PETSC_FALSE; 1886bdcc836SBarry Smith if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) { 1896bdcc836SBarry Smith /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */ 19023c918e7SPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 1919566063dSJacob Faibussowitsch PetscCall(MatMult(snes->jacobian_pre,Dold,W)); 19223c918e7SPeter Brune } else { 1939566063dSJacob Faibussowitsch PetscCall(VecCopy(Dold,W)); 19423c918e7SPeter Brune } 1959566063dSJacob Faibussowitsch PetscCall(VecDotBegin(W, Dold, &DolddotDold)); 1969566063dSJacob Faibussowitsch PetscCall(VecDotBegin(W, D, &DolddotD)); 1979566063dSJacob Faibussowitsch PetscCall(VecDotEnd(W, Dold, &DolddotDold)); 1989566063dSJacob Faibussowitsch PetscCall(VecDotEnd(W, D, &DolddotD)); 1990c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 2000c777b0cSPeter Brune } 2010c777b0cSPeter Brune periodic = PETSC_FALSE; 202b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 203b8840d6bSPeter Brune if (i_r>qn->m-1) periodic = PETSC_TRUE; 2040c777b0cSPeter Brune } 2050c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 20692f76d53SAlp Dener if (badstep || powell || periodic) { 20792f76d53SAlp Dener if (qn->monflg) { 2089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2)); 2096bdcc836SBarry Smith if (powell) { 21063a3b9bcSJacob 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)); 2116bdcc836SBarry Smith } else { 21263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %" PetscInt_FMT "\n", i_r)); 2136bdcc836SBarry Smith } 2149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2)); 2155ba6227bSPeter Brune } 2165ba6227bSPeter Brune i_r = -1; 2170c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 2189566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre)); 21907b62357SFande Kong SNESCheckJacobianDomainerror(snes); 2200ecc9e1dSPeter Brune } 2219566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(qn->B, PETSC_FALSE)); 2220ecc9e1dSPeter Brune } 2235ba6227bSPeter Brune } 2244b11644fSPeter Brune if (i == snes->max_its) { 22563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", snes->max_its)); 2264b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 2274b11644fSPeter Brune } 2284b11644fSPeter Brune PetscFunctionReturn(0); 2294b11644fSPeter Brune } 2304b11644fSPeter Brune 2314b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 2324b11644fSPeter Brune { 2339f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 234fced5a79SAsbjørn Nilsen Riseth DM dm; 23592f76d53SAlp Dener PetscInt n, N; 236335fb975SPeter Brune 2374b11644fSPeter Brune PetscFunctionBegin; 238fced5a79SAsbjørn Nilsen Riseth 239fced5a79SAsbjørn Nilsen Riseth if (!snes->vec_sol) { 2409566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 2419566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm,&snes->vec_sol)); 242fced5a79SAsbjørn Nilsen Riseth } 2439566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes,4)); 24492f76d53SAlp Dener 24592f76d53SAlp Dener if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 2469566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 24792f76d53SAlp Dener } 24892f76d53SAlp Dener if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;} 24992f76d53SAlp Dener 25060c08b40SPeter Brune /* set method defaults */ 2511efc8c45SPeter Brune if (qn->scale_type == SNES_QN_SCALE_DEFAULT) { 25260c08b40SPeter Brune if (qn->type == SNES_QN_BADBROYDEN) { 25360c08b40SPeter Brune qn->scale_type = SNES_QN_SCALE_NONE; 25460c08b40SPeter Brune } else { 25592f76d53SAlp Dener qn->scale_type = SNES_QN_SCALE_SCALAR; 25660c08b40SPeter Brune } 25760c08b40SPeter Brune } 2581efc8c45SPeter Brune if (qn->restart_type == SNES_QN_RESTART_DEFAULT) { 25960c08b40SPeter Brune if (qn->type == SNES_QN_LBFGS) { 26060c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 26160c08b40SPeter Brune } else { 26260c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_PERIODIC; 26360c08b40SPeter Brune } 26460c08b40SPeter Brune } 26560c08b40SPeter Brune 26692f76d53SAlp Dener /* Set up the LMVM matrix */ 26792f76d53SAlp Dener switch (qn->type) { 26892f76d53SAlp Dener case SNES_QN_BROYDEN: 2699566063dSJacob Faibussowitsch PetscCall(MatSetType(qn->B, MATLMVMBROYDEN)); 27092f76d53SAlp Dener qn->scale_type = SNES_QN_SCALE_NONE; 27192f76d53SAlp Dener break; 27292f76d53SAlp Dener case SNES_QN_BADBROYDEN: 2739566063dSJacob Faibussowitsch PetscCall(MatSetType(qn->B, MATLMVMBADBROYDEN)); 27492f76d53SAlp Dener qn->scale_type = SNES_QN_SCALE_NONE; 27592f76d53SAlp Dener break; 27692f76d53SAlp Dener default: 2779566063dSJacob Faibussowitsch PetscCall(MatSetType(qn->B, MATLMVMBFGS)); 27892f76d53SAlp Dener switch (qn->scale_type) { 27992f76d53SAlp Dener case SNES_QN_SCALE_NONE: 2809566063dSJacob Faibussowitsch PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_NONE)); 28192f76d53SAlp Dener break; 28292f76d53SAlp Dener case SNES_QN_SCALE_SCALAR: 2839566063dSJacob Faibussowitsch PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_SCALAR)); 28492f76d53SAlp Dener break; 28592f76d53SAlp Dener case SNES_QN_SCALE_JACOBIAN: 2869566063dSJacob Faibussowitsch PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_USER)); 28792f76d53SAlp Dener break; 28892f76d53SAlp Dener case SNES_QN_SCALE_DIAGONAL: 28992f76d53SAlp Dener case SNES_QN_SCALE_DEFAULT: 29092f76d53SAlp Dener default: 29192f76d53SAlp Dener break; 2928e231d97SPeter Brune } 29392f76d53SAlp Dener break; 29492f76d53SAlp Dener } 2959566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 2969566063dSJacob Faibussowitsch PetscCall(VecGetSize(snes->vec_sol, &N)); 2979566063dSJacob Faibussowitsch PetscCall(MatSetSizes(qn->B, n, n, N, N)); 2989566063dSJacob Faibussowitsch PetscCall(MatSetUp(qn->B)); 2999566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(qn->B, PETSC_TRUE)); 3009566063dSJacob Faibussowitsch PetscCall(MatLMVMSetHistorySize(qn->B, qn->m)); 3019566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(qn->B, snes->vec_sol, snes->vec_func)); 3024b11644fSPeter Brune PetscFunctionReturn(0); 3034b11644fSPeter Brune } 3044b11644fSPeter Brune 3054b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 3064b11644fSPeter Brune { 3079f83bee8SJed Brown SNES_QN *qn; 3080adebc6cSBarry Smith 3094b11644fSPeter Brune PetscFunctionBegin; 3104b11644fSPeter Brune if (snes->data) { 3119f83bee8SJed Brown qn = (SNES_QN*)snes->data; 3129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&qn->B)); 3134b11644fSPeter Brune } 3144b11644fSPeter Brune PetscFunctionReturn(0); 3154b11644fSPeter Brune } 3164b11644fSPeter Brune 3174b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 3184b11644fSPeter Brune { 3194b11644fSPeter Brune PetscFunctionBegin; 3209566063dSJacob Faibussowitsch PetscCall(SNESReset_QN(snes)); 3219566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 322*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",NULL)); 323*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",NULL)); 324*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",NULL)); 3254b11644fSPeter Brune PetscFunctionReturn(0); 3264b11644fSPeter Brune } 3274b11644fSPeter Brune 3284416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_QN(PetscOptionItems *PetscOptionsObject,SNES snes) 3294b11644fSPeter Brune { 3304b11644fSPeter Brune 3312150357eSBarry Smith SNES_QN *qn = (SNES_QN*)snes->data; 33292f76d53SAlp Dener PetscBool flg; 333f1c6b773SPeter Brune SNESLineSearch linesearch; 3342150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 3352150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 3361efc8c45SPeter Brune SNESQNType qtype = qn->type; 3372150357eSBarry Smith 3384b11644fSPeter Brune PetscFunctionBegin; 339d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"SNES QN options"); 3409566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL)); 3419566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL)); 3429566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", qn->monflg, &qn->monflg, NULL)); 3439566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg)); 3449566063dSJacob Faibussowitsch if (flg) PetscCall(SNESQNSetScaleType(snes,stype)); 34588f769c5SPeter Brune 3469566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg)); 3479566063dSJacob Faibussowitsch if (flg) PetscCall(SNESQNSetRestartType(snes,rtype)); 34888f769c5SPeter Brune 3499566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,(PetscEnum)qtype,(PetscEnum*)&qtype,&flg)); 3509566063dSJacob Faibussowitsch if (flg) PetscCall(SNESQNSetType(snes,qtype)); 3519566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(qn->B)); 352d0609cedSBarry Smith PetscOptionsHeadEnd(); 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 { 3619566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2)); 36260c08b40SPeter Brune } 3639e764e56SPeter Brune } 364ec786807SJed Brown } 36592f76d53SAlp Dener if (qn->monflg) { 3669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes), &qn->monitor)); 36744f7e39eSPeter Brune } 3684b11644fSPeter Brune PetscFunctionReturn(0); 3694b11644fSPeter Brune } 3704b11644fSPeter Brune 3715cd83419SPeter Brune static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer) 3725cd83419SPeter Brune { 3735cd83419SPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 3745cd83419SPeter Brune PetscBool iascii; 3755cd83419SPeter Brune 3765cd83419SPeter Brune PetscFunctionBegin; 3779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 3785cd83419SPeter Brune if (iascii) { 3799566063dSJacob 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])); 38063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Stored subspace size: %" PetscInt_FMT "\n", qn->m)); 3815cd83419SPeter Brune } 3825cd83419SPeter Brune PetscFunctionReturn(0); 3835cd83419SPeter Brune } 38492c02d66SPeter Brune 3850c777b0cSPeter Brune /*@ 3860c777b0cSPeter Brune SNESQNSetRestartType - Sets the restart type for SNESQN. 3870c777b0cSPeter Brune 3880c777b0cSPeter Brune Logically Collective on SNES 3890c777b0cSPeter Brune 3900c777b0cSPeter Brune Input Parameters: 3910c777b0cSPeter Brune + snes - the iterative context 3920c777b0cSPeter Brune - rtype - restart type 3930c777b0cSPeter Brune 3940c777b0cSPeter Brune Options Database: 3950c777b0cSPeter Brune + -snes_qn_restart_type <powell,periodic,none> - set the restart type 39684c577daSBarry Smith - -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic 3970c777b0cSPeter Brune 3980c777b0cSPeter Brune Level: intermediate 3990c777b0cSPeter Brune 4000c777b0cSPeter Brune SNESQNRestartTypes: 4010c777b0cSPeter Brune + SNES_QN_RESTART_NONE - never restart 4020c777b0cSPeter Brune . SNES_QN_RESTART_POWELL - restart based upon descent criteria 4030c777b0cSPeter Brune - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 4040c777b0cSPeter Brune 4050c777b0cSPeter Brune @*/ 4062150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 4072150357eSBarry Smith { 4080c777b0cSPeter Brune PetscFunctionBegin; 4090c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 410cac4c232SBarry Smith PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype)); 4110c777b0cSPeter Brune PetscFunctionReturn(0); 4120c777b0cSPeter Brune } 4130c777b0cSPeter Brune 4140c777b0cSPeter Brune /*@ 41584c577daSBarry Smith SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in SNESQN. 4160c777b0cSPeter Brune 4170c777b0cSPeter Brune Logically Collective on SNES 4180c777b0cSPeter Brune 4190c777b0cSPeter Brune Input Parameters: 4200c777b0cSPeter Brune + snes - the iterative context 4210c777b0cSPeter Brune - stype - scale type 4220c777b0cSPeter Brune 4230c777b0cSPeter Brune Options Database: 42467b8a455SSatish Balay . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type 4250c777b0cSPeter Brune 4260c777b0cSPeter Brune Level: intermediate 4270c777b0cSPeter Brune 42884c577daSBarry Smith SNESQNScaleTypes: 4290c777b0cSPeter Brune + SNES_QN_SCALE_NONE - don't scale the problem 4307044b745SBarry Smith . SNES_QN_SCALE_SCALAR - use Shanno scaling 43192f76d53SAlp Dener . SNES_QN_SCALE_DIAGONAL - scale with a diagonalized BFGS formula (see Gilbert and Lemarechal 1989), available 432a01a0525SBarry Smith - SNES_QN_SCALE_JACOBIAN - scale by solving a linear system coming from the Jacobian you provided with SNESSetJacobian() computed at the first iteration 433a01a0525SBarry Smith of QN and at ever restart. 4340c777b0cSPeter Brune 435db781477SPatrick Sanan .seealso: `SNES`, `SNESQN`, `SNESLineSearch`, `SNESQNScaleType`, `SNESSetJacobian()` 4360c777b0cSPeter Brune @*/ 4370c777b0cSPeter Brune 4382150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 4392150357eSBarry Smith { 4400c777b0cSPeter Brune PetscFunctionBegin; 4410c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 442cac4c232SBarry Smith PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype)); 4430c777b0cSPeter Brune PetscFunctionReturn(0); 4440c777b0cSPeter Brune } 4450c777b0cSPeter Brune 4460adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 4470adebc6cSBarry Smith { 4480c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 4496e111a19SKarl Rupp 4500c777b0cSPeter Brune PetscFunctionBegin; 4510c777b0cSPeter Brune qn->scale_type = stype; 4527044b745SBarry Smith if (stype == SNES_QN_SCALE_JACOBIAN) snes->usesksp = PETSC_TRUE; 4530c777b0cSPeter Brune PetscFunctionReturn(0); 4540c777b0cSPeter Brune } 4550c777b0cSPeter Brune 4560adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 4570adebc6cSBarry Smith { 4580c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 4596e111a19SKarl Rupp 4600c777b0cSPeter Brune PetscFunctionBegin; 4610c777b0cSPeter Brune qn->restart_type = rtype; 4620c777b0cSPeter Brune PetscFunctionReturn(0); 4630c777b0cSPeter Brune } 4640c777b0cSPeter Brune 4651efc8c45SPeter Brune /*@ 4661efc8c45SPeter Brune SNESQNSetType - Sets the quasi-Newton variant to be used in SNESQN. 4671efc8c45SPeter Brune 4681efc8c45SPeter Brune Logically Collective on SNES 4691efc8c45SPeter Brune 4701efc8c45SPeter Brune Input Parameters: 4711efc8c45SPeter Brune + snes - the iterative context 4721efc8c45SPeter Brune - qtype - variant type 4731efc8c45SPeter Brune 4741efc8c45SPeter Brune Options Database: 47567b8a455SSatish Balay . -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type 4761efc8c45SPeter Brune 4771efc8c45SPeter Brune Level: beginner 4781efc8c45SPeter Brune 4791efc8c45SPeter Brune SNESQNTypes: 4801efc8c45SPeter Brune + SNES_QN_LBFGS - LBFGS variant 4811efc8c45SPeter Brune . SNES_QN_BROYDEN - Broyden variant 4821efc8c45SPeter Brune - SNES_QN_BADBROYDEN - Bad Broyden variant 4831efc8c45SPeter Brune 4841efc8c45SPeter Brune @*/ 4851efc8c45SPeter Brune 4861efc8c45SPeter Brune PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype) 4871efc8c45SPeter Brune { 4881efc8c45SPeter Brune PetscFunctionBegin; 4891efc8c45SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 490cac4c232SBarry Smith PetscTryMethod(snes,"SNESQNSetType_C",(SNES,SNESQNType),(snes,qtype)); 4911efc8c45SPeter Brune PetscFunctionReturn(0); 4921efc8c45SPeter Brune } 4931efc8c45SPeter Brune 4941efc8c45SPeter Brune PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype) 4951efc8c45SPeter Brune { 4961efc8c45SPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 4971efc8c45SPeter Brune 4981efc8c45SPeter Brune PetscFunctionBegin; 4991efc8c45SPeter Brune qn->type = qtype; 5001efc8c45SPeter Brune PetscFunctionReturn(0); 5011efc8c45SPeter Brune } 5021efc8c45SPeter Brune 5034b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 5044b11644fSPeter Brune /*MC 5054b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 5064b11644fSPeter Brune 5076cc8130cSPeter Brune Options Database: 5086cc8130cSPeter Brune 50984c577daSBarry Smith + -snes_qn_m <m> - Number of past states saved for the L-Broyden methods. 51084c577daSBarry Smith + -snes_qn_restart_type <powell,periodic,none> - set the restart type 5116bdcc836SBarry Smith . -snes_qn_powell_gamma - Angle condition for restart. 5121867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 51384c577daSBarry Smith . -snes_qn_type <lbfgs,broyden,badbroyden> - QN type 51492f76d53SAlp Dener . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - scaling performed on inner Jacobian 51572835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 51684c577daSBarry Smith - -snes_qn_monitor - Monitors the quasi-newton Jacobian. 5176cc8130cSPeter Brune 51895452b02SPatrick Sanan Notes: 51995452b02SPatrick Sanan This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 520b8840d6bSPeter Brune previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 521b8840d6bSPeter Brune updates. 5226cc8130cSPeter Brune 5231867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 5241867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 52584c577daSBarry Smith iteration as the current iteration's values when constructing the approximate Jacobian. The second, composed, 52684c577daSBarry Smith perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 5271867fe5bSPeter Brune 5282d547940SBarry Smith Uses left nonlinear preconditioning by default. 5292d547940SBarry Smith 5306cc8130cSPeter Brune References: 531606c0280SSatish Balay + * - Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995. 532606c0280SSatish Balay . * - R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods, 5330a8e8854SPeter Brune Technical Report, Northwestern University, June 1992. 534606c0280SSatish Balay . * - Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE 5350a8e8854SPeter Brune Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 536606c0280SSatish Balay . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 5374f02bc6aSBarry Smith SIAM Review, 57(4), 2015 538606c0280SSatish Balay . * - Griewank, Andreas. "Broyden updating, the good and the bad!." Doc. Math (2012): 301-315. 539606c0280SSatish Balay . * - Gilbert, Jean Charles, and Claude Lemar{\'e}chal. "Some numerical experiments with variable-storage quasi-Newton algorithms." 54092f76d53SAlp Dener Mathematical programming 45.1-3 (1989): 407-435. 541606c0280SSatish Balay - * - Dener A., Munson T. "Accelerating Limited-Memory Quasi-Newton Convergence for Large-Scale Optimization" 542110fc3b0SBarry Smith Computational Science - ICCS 2019. ICCS 2019. Lecture Notes in Computer Science, vol 11538. Springer, Cham 5434b11644fSPeter Brune 5444b11644fSPeter Brune Level: beginner 5454b11644fSPeter Brune 546db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR` 5476cc8130cSPeter Brune 5484b11644fSPeter Brune M*/ 5498cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) 5504b11644fSPeter Brune { 5519f83bee8SJed Brown SNES_QN *qn; 55292f76d53SAlp Dener const char *optionsprefix; 5534b11644fSPeter Brune 5544b11644fSPeter Brune PetscFunctionBegin; 5554b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 5564b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 5574b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 5584b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 5595cd83419SPeter Brune snes->ops->view = SNESView_QN; 5604b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 5614b11644fSPeter Brune 562efd4aadfSBarry Smith snes->npcside= PC_LEFT; 56347f26062SPeter Brune 564efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 56542f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 56642f4f86dSBarry Smith 5674fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 5684fc747eaSLawrence Mitchell 56988976e71SPeter Brune if (!snes->tolerancesset) { 5700e444f03SPeter Brune snes->max_funcs = 30000; 5710e444f03SPeter Brune snes->max_its = 10000; 57288976e71SPeter Brune } 5730e444f03SPeter Brune 5749566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes,&qn)); 5754b11644fSPeter Brune snes->data = (void*) qn; 5760ecc9e1dSPeter Brune qn->m = 10; 577b21d5a53SPeter Brune qn->scaling = 1.0; 5780298fd71SBarry Smith qn->monitor = NULL; 57992f76d53SAlp Dener qn->monflg = PETSC_FALSE; 580b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 58160c08b40SPeter Brune qn->scale_type = SNES_QN_SCALE_DEFAULT; 58260c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_DEFAULT; 583b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 584ea630c6eSPeter Brune 5859566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &qn->B)); 5869566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 5879566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(qn->B, optionsprefix)); 58892f76d53SAlp Dener 5899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN)); 5909566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN)); 5919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",SNESQNSetType_QN)); 5924b11644fSPeter Brune PetscFunctionReturn(0); 5934b11644fSPeter Brune } 594