xref: /petsc/src/snes/impls/qn/qn.c (revision 6a6fc655e6b1b6c3ff84be1691b1c5064ef5f4c2)
10c777b0cSPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/
24b11644fSPeter Brune 
38e231d97SPeter Brune #define H(i,j)  qn->dXdFmat[i*qn->m + j]
48e231d97SPeter Brune 
5*6a6fc655SJed Brown const char *const SNESQNCompositionTypes[] =  {"SEQUENTIAL","COMPOSED","SNESQNCompositionType","SNES_QN_",0};
6*6a6fc655SJed Brown const char *const SNESQNScaleTypes[] =      {"NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0};
7*6a6fc655SJed Brown const char *const SNESQNRestartTypes[] =      {"NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0};
86bf1b2e5SPeter Brune 
94b11644fSPeter Brune typedef struct {
104b11644fSPeter Brune   Vec          *dX;              /* The change in X */
116bf1b2e5SPeter Brune   Vec          *dF;              /* The change in F */
128e231d97SPeter Brune   PetscInt     m;                /* The number of kept previous steps */
138e231d97SPeter Brune   PetscScalar  *alpha, *beta;
148e231d97SPeter Brune   PetscScalar  *dXtdF, *dFtdX, *YtdX;
1518aa0c0cSPeter Brune   PetscBool    singlereduction;  /* Aggregated reduction implementation */
168e231d97SPeter Brune   PetscScalar  *dXdFmat;         /* A matrix of values for dX_i dot dF_j */
1744f7e39eSPeter Brune   PetscViewer  monitor;
186bf1b2e5SPeter Brune   PetscReal    powell_gamma;     /* Powell angle restart condition */
196bf1b2e5SPeter Brune   PetscReal    powell_downhill;  /* Powell descent restart condition */
20b21d5a53SPeter Brune   PetscReal    scaling;          /* scaling of H0 */
210c777b0cSPeter Brune   PetscInt     restart_periodic; /* the maximum iterations between restart */
2288d374b2SPeter Brune 
230c777b0cSPeter Brune   SNESQNCompositionType composition_type; /* determine if the composition is done sequentially or as a composition */
240c777b0cSPeter Brune   SNESQNScaleType       scale_type;       /* determine if the composition is done sequentially or as a composition */
250c777b0cSPeter Brune   SNESQNRestartType     restart_type;     /* determine the frequency and type of restart conditions */
269f83bee8SJed Brown } SNES_QN;
274b11644fSPeter Brune 
284b11644fSPeter Brune #undef __FUNCT__
298c40d5fbSBarry Smith #define __FUNCT__ "SNESQNApplyJinv_Private"
308c40d5fbSBarry Smith PetscErrorCode SNESQNApplyJinv_Private(SNES snes, PetscInt it, Vec D, Vec Y) {
314b11644fSPeter Brune 
324b11644fSPeter Brune   PetscErrorCode ierr;
334b11644fSPeter Brune 
349f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*)snes->data;
354b11644fSPeter Brune 
36335fb975SPeter Brune   Vec Yin = snes->work[3];
370ecc9e1dSPeter Brune 
384b11644fSPeter Brune   Vec *dX = qn->dX;
396bf1b2e5SPeter Brune   Vec *dF = qn->dF;
404b11644fSPeter Brune 
41bd052dfeSPeter Brune   PetscScalar *alpha    = qn->alpha;
42bd052dfeSPeter Brune   PetscScalar *beta     = qn->beta;
438e231d97SPeter Brune   PetscScalar *dXtdF    = qn->dXtdF;
448e231d97SPeter Brune   PetscScalar *YtdX     = qn->YtdX;
45bd052dfeSPeter Brune 
460ecc9e1dSPeter Brune   /* ksp thing for jacobian scaling */
470ecc9e1dSPeter Brune   KSPConvergedReason kspreason;
480ecc9e1dSPeter Brune   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
490ecc9e1dSPeter Brune 
508e231d97SPeter Brune   PetscInt k, i, j, g, lits;
514b11644fSPeter Brune   PetscInt m = qn->m;
524b11644fSPeter Brune   PetscScalar t;
534b11644fSPeter Brune   PetscInt l = m;
544b11644fSPeter Brune 
558e231d97SPeter Brune   Mat jac, jac_pre;
568e231d97SPeter Brune 
574b11644fSPeter Brune   PetscFunctionBegin;
584b11644fSPeter Brune 
593af51624SPeter Brune   ierr = VecCopy(D, Y);CHKERRQ(ierr);
604b11644fSPeter Brune 
615ba6227bSPeter Brune   if (it < m) l = it;
624b11644fSPeter Brune 
6318aa0c0cSPeter Brune   if (qn->singlereduction) {
648e231d97SPeter Brune     ierr = VecMDot(Y, l, qn->dX, YtdX);CHKERRQ(ierr);
658e231d97SPeter Brune   }
664b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
674b11644fSPeter Brune   for (i = 0; i < l; i++) {
68b21d5a53SPeter Brune     k = (it - i - 1) % l;
6918aa0c0cSPeter Brune     if (qn->singlereduction) {
708e231d97SPeter Brune       /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */
718e231d97SPeter Brune       t = YtdX[k];
728e231d97SPeter Brune       for (j = 0; j < i; j++) {
738e231d97SPeter Brune         g = (it - j - 1) % l;
748e231d97SPeter Brune         t += -alpha[g]*H(g, k);
758e231d97SPeter Brune       }
768e231d97SPeter Brune       alpha[k] = t / H(k, k);
778e231d97SPeter Brune     } else {
783af51624SPeter Brune       ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr);
798e231d97SPeter Brune       alpha[k] = t / dXtdF[k];
808e231d97SPeter Brune     }
8144f7e39eSPeter Brune     if (qn->monitor) {
823af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
835ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha:        %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr);
843af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
8544f7e39eSPeter Brune     }
866bf1b2e5SPeter Brune     ierr = VecAXPY(Y, -alpha[k], dF[k]);CHKERRQ(ierr);
874b11644fSPeter Brune   }
884b11644fSPeter Brune 
890c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
908e231d97SPeter Brune     ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
918e231d97SPeter Brune     ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr);
920ecc9e1dSPeter Brune     ierr = SNES_KSPSolve(snes,snes->ksp,Y,Yin);CHKERRQ(ierr);
930ecc9e1dSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
940ecc9e1dSPeter Brune     if (kspreason < 0) {
950ecc9e1dSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
960ecc9e1dSPeter Brune         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
970ecc9e1dSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
980ecc9e1dSPeter Brune         PetscFunctionReturn(0);
990ecc9e1dSPeter Brune       }
1000ecc9e1dSPeter Brune     }
1010ecc9e1dSPeter Brune     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
1020ecc9e1dSPeter Brune     snes->linear_its += lits;
1030ecc9e1dSPeter Brune     ierr = VecCopy(Yin, Y);CHKERRQ(ierr);
1040ecc9e1dSPeter Brune   } else {
105b21d5a53SPeter Brune     ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
1060ecc9e1dSPeter Brune   }
10718aa0c0cSPeter Brune   if (qn->singlereduction) {
1088e231d97SPeter Brune     ierr = VecMDot(Y, l, qn->dF, YtdX);CHKERRQ(ierr);
1098e231d97SPeter Brune   }
1104b11644fSPeter Brune   /* inward recursion starting at the first update and working forward */
111bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
112b21d5a53SPeter Brune     k = (it + i - l) % l;
11318aa0c0cSPeter Brune     if (qn->singlereduction) {
1148e231d97SPeter Brune       t = YtdX[k];
1158e231d97SPeter Brune       for (j = 0; j < i; j++) {
1168e231d97SPeter Brune         g = (it + j - l) % l;
1178e231d97SPeter Brune         t += (alpha[g] - beta[g])*H(k, g);
1188e231d97SPeter Brune       }
1198e231d97SPeter Brune       beta[k] = t / H(k, k);
1208e231d97SPeter Brune     } else {
1216bf1b2e5SPeter Brune       ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
1228e231d97SPeter Brune       beta[k] = t / dXtdF[k];
1238e231d97SPeter Brune     }
1243af51624SPeter Brune     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);
12544f7e39eSPeter Brune     if (qn->monitor) {
1263af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1275ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
1283af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
12944f7e39eSPeter Brune     }
1304b11644fSPeter Brune   }
1314b11644fSPeter Brune   PetscFunctionReturn(0);
1324b11644fSPeter Brune }
1334b11644fSPeter Brune 
1344b11644fSPeter Brune #undef __FUNCT__
1354b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN"
1364b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
1374b11644fSPeter Brune {
1384b11644fSPeter Brune 
1394b11644fSPeter Brune   PetscErrorCode ierr;
1409f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*) snes->data;
1414b11644fSPeter Brune 
14215f5eeeaSPeter Brune   Vec X, Xold;
143335fb975SPeter Brune   Vec F, B;
144335fb975SPeter Brune   Vec Y, FPC, D, Dold;
1453af51624SPeter Brune   SNESConvergedReason reason;
1468e231d97SPeter Brune   PetscInt i, i_r, k, l, j;
1474b11644fSPeter Brune 
148335fb975SPeter Brune   PetscReal fnorm, xnorm, ynorm, gnorm;
1494b11644fSPeter Brune   PetscInt m = qn->m;
1500c777b0cSPeter Brune   PetscBool lssucceed,powell,periodic;
1514b11644fSPeter Brune 
152bd052dfeSPeter Brune   Vec *dX = qn->dX;
1536bf1b2e5SPeter Brune   Vec *dF = qn->dF;
1548e231d97SPeter Brune   PetscScalar *dXtdF = qn->dXtdF;
1558e231d97SPeter Brune   PetscScalar *dFtdX = qn->dFtdX;
15688d374b2SPeter Brune   PetscScalar DolddotD, DolddotDold, DdotD, YdotD, a;
1574b11644fSPeter Brune 
1580ecc9e1dSPeter Brune   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
1590ecc9e1dSPeter Brune 
1604b11644fSPeter Brune   /* basically just a regular newton's method except for the application of the jacobian */
1614b11644fSPeter Brune   PetscFunctionBegin;
1624b11644fSPeter Brune 
1639f3a0142SPeter Brune   X             = snes->vec_sol;        /* solution vector */
1649f3a0142SPeter Brune   F             = snes->vec_func;       /* residual vector */
1653af51624SPeter Brune   Y             = snes->vec_sol_update; /* search direction generated by J^-1D*/
1663af51624SPeter Brune   B             = snes->vec_rhs;
167335fb975SPeter Brune   Xold          = snes->work[0];
1683af51624SPeter Brune 
1693af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
170335fb975SPeter Brune   D             = snes->work[1];
171335fb975SPeter Brune   Dold          = snes->work[2];
1724b11644fSPeter Brune 
1734b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
1744b11644fSPeter Brune 
1754b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1764b11644fSPeter Brune   snes->iter = 0;
1774b11644fSPeter Brune   snes->norm = 0.;
1784b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
179e4ed7901SPeter Brune   if (!snes->vec_func_init_set){
18015f5eeeaSPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1814b11644fSPeter Brune     if (snes->domainerror) {
1824b11644fSPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1834b11644fSPeter Brune       PetscFunctionReturn(0);
1844b11644fSPeter Brune     }
185e4ed7901SPeter Brune   } else {
186e4ed7901SPeter Brune     snes->vec_func_init_set = PETSC_FALSE;
187e4ed7901SPeter Brune   }
188e4ed7901SPeter Brune 
189e4ed7901SPeter Brune   if (!snes->norm_init_set) {
19015f5eeeaSPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1914b11644fSPeter Brune     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
192e4ed7901SPeter Brune   } else {
193e4ed7901SPeter Brune     fnorm = snes->norm_init;
194e4ed7901SPeter Brune     snes->norm_init_set = PETSC_FALSE;
195e4ed7901SPeter Brune   }
196e4ed7901SPeter Brune 
1974b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1984b11644fSPeter Brune   snes->norm = fnorm;
1994b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2004b11644fSPeter Brune   SNESLogConvHistory(snes,fnorm,0);
2014b11644fSPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
2024b11644fSPeter Brune 
2034b11644fSPeter Brune   /* set parameter for default relative tolerance convergence test */
2044b11644fSPeter Brune    snes->ttol = fnorm*snes->rtol;
2054b11644fSPeter Brune   /* test convergence */
2064b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2074b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
20870d3b23bSPeter Brune 
20988d374b2SPeter Brune   /* composed solve -- either sequential or composed */
21088d374b2SPeter Brune   if (snes->pc) {
2110c777b0cSPeter Brune     if (qn->composition_type == SNES_QN_SEQUENTIAL) {
212217b9c2eSPeter Brune       ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
213217b9c2eSPeter Brune       ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
21488d374b2SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
21588d374b2SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
21688d374b2SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
21788d374b2SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
21888d374b2SPeter Brune         PetscFunctionReturn(0);
21988d374b2SPeter Brune       }
22088d374b2SPeter Brune       ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
22188d374b2SPeter Brune       ierr = VecCopy(FPC, F);CHKERRQ(ierr);
22288d374b2SPeter Brune       ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
2236bf1b2e5SPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
22488d374b2SPeter Brune     } else {
22588d374b2SPeter Brune       ierr = VecCopy(X, Y);CHKERRQ(ierr);
226217b9c2eSPeter Brune       ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
227217b9c2eSPeter Brune       ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
22888d374b2SPeter Brune       ierr = SNESSolve(snes->pc, B, Y);CHKERRQ(ierr);
22988d374b2SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
23088d374b2SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
23188d374b2SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
23288d374b2SPeter Brune         PetscFunctionReturn(0);
23388d374b2SPeter Brune       }
23488d374b2SPeter Brune       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);
23588d374b2SPeter Brune     }
23688d374b2SPeter Brune   } else {
23788d374b2SPeter Brune     ierr = VecCopy(F, Y);CHKERRQ(ierr);
23888d374b2SPeter Brune   }
23988d374b2SPeter Brune   ierr = VecCopy(Y, D);CHKERRQ(ierr);
2403af51624SPeter Brune 
241f8e15203SPeter Brune   /* scale the initial update */
2420c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
2430ecc9e1dSPeter Brune     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
2440ecc9e1dSPeter Brune   }
2450ecc9e1dSPeter Brune 
2465ba6227bSPeter Brune   for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
2478c40d5fbSBarry Smith     ierr = SNESQNApplyJinv_Private(snes, i_r, D, Y);CHKERRQ(ierr);
24870d3b23bSPeter Brune     /* line search for lambda */
24970d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
25088d374b2SPeter Brune     ierr = VecCopy(D, Dold);CHKERRQ(ierr);
25115f5eeeaSPeter Brune     ierr = VecCopy(X, Xold);CHKERRQ(ierr);
252f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
2539f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
2549f3a0142SPeter Brune     if (snes->domainerror) {
2559f3a0142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2569f3a0142SPeter Brune       PetscFunctionReturn(0);
2579f3a0142SPeter Brune       }
258f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr);
2599f3a0142SPeter Brune     if (!lssucceed) {
2609f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
2619f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2629f3a0142SPeter Brune         break;
2639f3a0142SPeter Brune       }
2649f3a0142SPeter Brune     }
265f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
2660c777b0cSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
267f1c6b773SPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr);
2680ecc9e1dSPeter Brune     }
2693af51624SPeter Brune 
27088d374b2SPeter Brune     /* convergence monitoring */
271b21d5a53SPeter Brune     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr);
272b21d5a53SPeter Brune 
273360c497dSPeter Brune     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
274360c497dSPeter Brune     ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr);
275360c497dSPeter Brune 
2768409ca45SMatthew G Knepley     SNESLogConvHistory(snes,snes->norm,snes->iter);
2778409ca45SMatthew G Knepley     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
2784b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
279d484d688SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2804b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
2814b11644fSPeter Brune 
28288d374b2SPeter Brune 
28388d374b2SPeter Brune     if (snes->pc) {
2840c777b0cSPeter Brune       if (qn->composition_type == SNES_QN_SEQUENTIAL) {
285217b9c2eSPeter Brune         ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
286217b9c2eSPeter Brune         ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
28788d374b2SPeter Brune         ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
28888d374b2SPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
28988d374b2SPeter Brune         if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
29088d374b2SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
29188d374b2SPeter Brune           PetscFunctionReturn(0);
29288d374b2SPeter Brune         }
29388d374b2SPeter Brune         ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
29488d374b2SPeter Brune         ierr = VecCopy(FPC, F);CHKERRQ(ierr);
29588d374b2SPeter Brune         ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
29688d374b2SPeter Brune         ierr = VecCopy(F, D);CHKERRQ(ierr);
29788d374b2SPeter Brune       } else {
29888d374b2SPeter Brune         ierr = VecCopy(X, D);CHKERRQ(ierr);
299217b9c2eSPeter Brune         ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
300217b9c2eSPeter Brune         ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
30188d374b2SPeter Brune         ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr);
30288d374b2SPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
30388d374b2SPeter Brune         if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
30488d374b2SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
30588d374b2SPeter Brune           PetscFunctionReturn(0);
30688d374b2SPeter Brune         }
30788d374b2SPeter Brune         ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr);
30888d374b2SPeter Brune       }
30988d374b2SPeter Brune     } else {
31088d374b2SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
31188d374b2SPeter Brune     }
31288d374b2SPeter Brune 
3130c777b0cSPeter Brune     powell = PETSC_FALSE;
3140c777b0cSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_POWELL) {
3156bf1b2e5SPeter Brune       /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
3168e231d97SPeter Brune       ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
3178e231d97SPeter Brune       ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr);
3188e231d97SPeter Brune       ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr);
3198e231d97SPeter Brune       ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr);
3208e231d97SPeter Brune       ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
3218e231d97SPeter Brune       ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr);
3228e231d97SPeter Brune       ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr);
3238e231d97SPeter Brune       ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr);
3240c777b0cSPeter Brune       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
3250c777b0cSPeter Brune     }
3260c777b0cSPeter Brune     periodic = PETSC_FALSE;
3270c777b0cSPeter Brune     if (qn->restart_type != SNES_QN_RESTART_NONE) {
3280c777b0cSPeter Brune       if ((i_r > qn->restart_periodic - 1 && qn->restart_periodic > 0)) periodic = PETSC_TRUE;
3290c777b0cSPeter Brune     }
3300c777b0cSPeter Brune 
3310c777b0cSPeter Brune     /* restart if either powell or periodic restart is satisfied. */
3320c777b0cSPeter Brune     if (powell || periodic) {
3335ba6227bSPeter Brune       if (qn->monitor) {
3345ba6227bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
335516fe3c3SPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e| or i_r = %d\n", PetscRealPart(DolddotD), qn->powell_gamma, PetscRealPart(DolddotDold), i_r);CHKERRQ(ierr);
3365ba6227bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
3375ba6227bSPeter Brune       }
3385ba6227bSPeter Brune       i_r = -1;
3395ba6227bSPeter Brune       /* general purpose update */
3405ba6227bSPeter Brune       if (snes->ops->update) {
3415ba6227bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
3425ba6227bSPeter Brune       }
3430c777b0cSPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
3440ecc9e1dSPeter Brune         ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
3450ecc9e1dSPeter Brune       }
3465ba6227bSPeter Brune     } else {
347bd052dfeSPeter Brune       /* set the differences */
3485ba6227bSPeter Brune       k = i_r % m;
3498e231d97SPeter Brune       l = m;
3508e231d97SPeter Brune       if (i_r + 1 < m) l = i_r + 1;
35188d374b2SPeter Brune       ierr = VecCopy(D, dF[k]);CHKERRQ(ierr);
35288d374b2SPeter Brune       ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr);
35315f5eeeaSPeter Brune       ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
35415f5eeeaSPeter Brune       ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
35518aa0c0cSPeter Brune       if (qn->singlereduction) {
3568e231d97SPeter Brune         ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr);
3578e231d97SPeter Brune         ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr);
3588e231d97SPeter Brune         for (j = 0; j < l; j++) {
3598e231d97SPeter Brune           H(k, j) = dFtdX[j];
3608e231d97SPeter Brune           H(j, k) = dXtdF[j];
3618e231d97SPeter Brune         }
3628e231d97SPeter Brune         /* copy back over to make the computation of alpha and beta easier */
3638e231d97SPeter Brune         for (j = 0; j < l; j++) {
3648e231d97SPeter Brune           dXtdF[j] = H(j, j);
3658e231d97SPeter Brune         }
3668e231d97SPeter Brune       } else {
3678e231d97SPeter Brune         ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr);
3688e231d97SPeter Brune       }
3696bf1b2e5SPeter Brune       /* set scaling to be shanno scaling */
3700c777b0cSPeter Brune       if (qn->scale_type == SNES_QN_SCALE_SHANNO) {
3716bf1b2e5SPeter Brune         ierr = VecDot(dF[k], dF[k], &a);CHKERRQ(ierr);
3728e231d97SPeter Brune         qn->scaling = PetscRealPart(dXtdF[k]) / PetscRealPart(a);
3730ecc9e1dSPeter Brune       }
37470d3b23bSPeter Brune       /* general purpose update */
37570d3b23bSPeter Brune       if (snes->ops->update) {
37670d3b23bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
37770d3b23bSPeter Brune       }
3785ba6227bSPeter Brune     }
3794b11644fSPeter Brune   }
3804b11644fSPeter Brune   if (i == snes->max_its) {
3814b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
3824b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
3834b11644fSPeter Brune   }
3844b11644fSPeter Brune   PetscFunctionReturn(0);
3854b11644fSPeter Brune }
3864b11644fSPeter Brune 
3874b11644fSPeter Brune 
3884b11644fSPeter Brune #undef __FUNCT__
3894b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
3904b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
3914b11644fSPeter Brune {
3929f83bee8SJed Brown   SNES_QN        *qn = (SNES_QN*)snes->data;
3934b11644fSPeter Brune   PetscErrorCode ierr;
394335fb975SPeter Brune 
3954b11644fSPeter Brune   PetscFunctionBegin;
3964b11644fSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr);
3976bf1b2e5SPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr);
3988e231d97SPeter Brune   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha,
3998e231d97SPeter Brune                       qn->m, PetscScalar, &qn->beta,
4008e231d97SPeter Brune                       qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr);
4018e231d97SPeter Brune 
40218aa0c0cSPeter Brune   if (qn->singlereduction) {
4038e231d97SPeter Brune     ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat,
4048e231d97SPeter Brune                         qn->m, PetscScalar, &qn->dFtdX,
4058e231d97SPeter Brune                         qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr);
4068e231d97SPeter Brune   }
407335fb975SPeter Brune   ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr);
408335fb975SPeter Brune 
409335fb975SPeter Brune   /* set up the line search */
4100c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
4118e231d97SPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
4128e231d97SPeter Brune   }
4134b11644fSPeter Brune   PetscFunctionReturn(0);
4144b11644fSPeter Brune }
4154b11644fSPeter Brune 
4164b11644fSPeter Brune #undef __FUNCT__
4174b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
4184b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
4194b11644fSPeter Brune {
4204b11644fSPeter Brune   PetscErrorCode ierr;
4219f83bee8SJed Brown   SNES_QN *qn;
4224b11644fSPeter Brune   PetscFunctionBegin;
4234b11644fSPeter Brune   if (snes->data) {
4249f83bee8SJed Brown     qn = (SNES_QN*)snes->data;
4254b11644fSPeter Brune     if (qn->dX) {
4264b11644fSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr);
4274b11644fSPeter Brune     }
4286bf1b2e5SPeter Brune     if (qn->dF) {
4296bf1b2e5SPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr);
4304b11644fSPeter Brune     }
43118aa0c0cSPeter Brune     if (qn->singlereduction) {
4328e231d97SPeter Brune       ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr);
4338e231d97SPeter Brune     }
4348e231d97SPeter Brune     ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr);
4354b11644fSPeter Brune   }
4364b11644fSPeter Brune   PetscFunctionReturn(0);
4374b11644fSPeter Brune }
4384b11644fSPeter Brune 
4394b11644fSPeter Brune #undef __FUNCT__
4404b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
4414b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
4424b11644fSPeter Brune {
4434b11644fSPeter Brune   PetscErrorCode ierr;
4444b11644fSPeter Brune   PetscFunctionBegin;
4454b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
4464b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
4479c05af5eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr);
4484b11644fSPeter Brune   PetscFunctionReturn(0);
4494b11644fSPeter Brune }
4504b11644fSPeter Brune 
4514b11644fSPeter Brune #undef __FUNCT__
4524b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
4534b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
4544b11644fSPeter Brune {
4554b11644fSPeter Brune 
4564b11644fSPeter Brune   PetscErrorCode ierr;
4579f83bee8SJed Brown   SNES_QN    *qn;
45844f7e39eSPeter Brune   PetscBool  monflg = PETSC_FALSE;
459f1c6b773SPeter Brune   SNESLineSearch linesearch;
4604b11644fSPeter Brune   PetscFunctionBegin;
4614b11644fSPeter Brune 
4629f83bee8SJed Brown   qn = (SNES_QN*)snes->data;
4634b11644fSPeter Brune 
4644b11644fSPeter Brune   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
4655b4627e8SPeter Brune   ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,PETSC_NULL);CHKERRQ(ierr);
4660c777b0cSPeter Brune   ierr = PetscOptionsInt("-snes_qn_restart","Maximum number of iterations between restarts","SNESQN",qn->restart_periodic,&qn->restart_periodic, PETSC_NULL);CHKERRQ(ierr);
4675b4627e8SPeter Brune   ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr);
4685b4627e8SPeter Brune   ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance",        "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr);
4695b4627e8SPeter Brune   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr);
4704d116c7dSPeter Brune   ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions",           "SNESQN", qn->singlereduction, &qn->singlereduction, PETSC_NULL);CHKERRQ(ierr);
4710c777b0cSPeter Brune   ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)qn->scale_type,(PetscEnum*)&qn->scale_type,PETSC_NULL);CHKERRQ(ierr);
4720c777b0cSPeter Brune   ierr = PetscOptionsEnum("-snes_qn_composition_type","Composition type","SNESQNSetCompositionType",SNESQNCompositionTypes,
4730c777b0cSPeter Brune                           (PetscEnum)qn->composition_type,(PetscEnum*)&qn->composition_type,PETSC_NULL);CHKERRQ(ierr);
4740c777b0cSPeter Brune   ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,
4750c777b0cSPeter Brune                           (PetscEnum)qn->restart_type,(PetscEnum*)&qn->restart_type,PETSC_NULL);CHKERRQ(ierr);
4764b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
4779e764e56SPeter Brune   if (!snes->linesearch) {
478f1c6b773SPeter Brune     ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
4790c777b0cSPeter Brune     if (!snes->pc || qn->composition_type == SNES_QN_SEQUENTIAL) {
4801a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
4819e764e56SPeter Brune     } else {
4821a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
4839e764e56SPeter Brune     }
4849e764e56SPeter Brune   }
48544f7e39eSPeter Brune   if (monflg) {
48644f7e39eSPeter Brune     qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
48744f7e39eSPeter Brune   }
4884b11644fSPeter Brune   PetscFunctionReturn(0);
4894b11644fSPeter Brune }
4904b11644fSPeter Brune 
49192c02d66SPeter Brune 
4920c777b0cSPeter Brune #undef __FUNCT__
4930c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType"
4940c777b0cSPeter Brune /*@
4950c777b0cSPeter Brune     SNESQNSetRestartType - Sets the restart type for SNESQN.
4960c777b0cSPeter Brune 
4970c777b0cSPeter Brune     Logically Collective on SNES
4980c777b0cSPeter Brune 
4990c777b0cSPeter Brune     Input Parameters:
5000c777b0cSPeter Brune +   snes - the iterative context
5010c777b0cSPeter Brune -   rtype - restart type
5020c777b0cSPeter Brune 
5030c777b0cSPeter Brune     Options Database:
5040c777b0cSPeter Brune +   -snes_qn_restart_type<powell,periodic,none> - set the restart type
5050c777b0cSPeter Brune -   -snes_qn_restart[30] - sets the number of iterations before restart for periodic
5060c777b0cSPeter Brune 
5070c777b0cSPeter Brune     Level: intermediate
5080c777b0cSPeter Brune 
5090c777b0cSPeter Brune     SNESQNRestartTypes:
5100c777b0cSPeter Brune +   SNES_QN_RESTART_NONE - never restart
5110c777b0cSPeter Brune .   SNES_QN_RESTART_POWELL - restart based upon descent criteria
5120c777b0cSPeter Brune -   SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations
5130c777b0cSPeter Brune 
5140c777b0cSPeter Brune     Notes:
5150c777b0cSPeter Brune     The default line search used is the L2 line search and it requires two additional function evaluations.
5160c777b0cSPeter Brune 
5170c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch
5180c777b0cSPeter Brune @*/
5190c777b0cSPeter Brune PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) {
5200c777b0cSPeter Brune   PetscErrorCode ierr;
5210c777b0cSPeter Brune   PetscFunctionBegin;
5220c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
5230c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr);
5240c777b0cSPeter Brune   PetscFunctionReturn(0);
5250c777b0cSPeter Brune }
5260c777b0cSPeter Brune 
5270c777b0cSPeter Brune #undef __FUNCT__
5280c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType"
5290c777b0cSPeter Brune /*@
5300c777b0cSPeter Brune     SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN.
5310c777b0cSPeter Brune 
5320c777b0cSPeter Brune     Logically Collective on SNES
5330c777b0cSPeter Brune 
5340c777b0cSPeter Brune     Input Parameters:
5350c777b0cSPeter Brune +   snes - the iterative context
5360c777b0cSPeter Brune -   stype - scale type
5370c777b0cSPeter Brune 
5380c777b0cSPeter Brune     Options Database:
5390c777b0cSPeter Brune .   -snes_qn_scale_type<shanno,none,linesearch,jacobian>
5400c777b0cSPeter Brune 
5410c777b0cSPeter Brune     Level: intermediate
5420c777b0cSPeter Brune 
5430c777b0cSPeter Brune     SNESQNSelectTypes:
5440c777b0cSPeter Brune +   SNES_QN_SCALE_NONE - don't scale the problem
5450c777b0cSPeter Brune .   SNES_QN_SCALE_SHANNO - use shanno scaling
5460c777b0cSPeter Brune .   SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda
5470c777b0cSPeter Brune -   SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian.
5480c777b0cSPeter Brune 
5490c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch
5500c777b0cSPeter Brune @*/
5510c777b0cSPeter Brune 
5520c777b0cSPeter Brune PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) {
5530c777b0cSPeter Brune   PetscErrorCode ierr;
5540c777b0cSPeter Brune   PetscFunctionBegin;
5550c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
5560c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr);
5570c777b0cSPeter Brune   PetscFunctionReturn(0);
5580c777b0cSPeter Brune }
5590c777b0cSPeter Brune 
5600c777b0cSPeter Brune #undef __FUNCT__
5610c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetCompositionType"
5620c777b0cSPeter Brune /*@
5630c777b0cSPeter Brune     SNESQNSetCompositionType - Sets the composition type
5640c777b0cSPeter Brune 
5650c777b0cSPeter Brune     Logically Collective on SNES
5660c777b0cSPeter Brune 
5670c777b0cSPeter Brune     Input Parameters:
5680c777b0cSPeter Brune +   snes - the iterative context
5690c777b0cSPeter Brune -   stype - composition type
5700c777b0cSPeter Brune 
5710c777b0cSPeter Brune     Options Database:
5720c777b0cSPeter Brune .   -snes_qn_composition_type<sequential, composed>
5730c777b0cSPeter Brune 
5740c777b0cSPeter Brune     Level: intermediate
5750c777b0cSPeter Brune 
5760c777b0cSPeter Brune     SNESQNSelectTypes:
5770c777b0cSPeter Brune +   SNES_QN_COMPOSITION_SEQUENTIAL - Solve the system with X = PC(X) and D = F(PC(X))
5780c777b0cSPeter Brune -   SNES_QN_COMPOSITION_COMPOSED - solve the system with X = X and D = PC(X) - X
5790c777b0cSPeter Brune 
5800c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch
5810c777b0cSPeter Brune @*/
5820c777b0cSPeter Brune 
5830c777b0cSPeter Brune PetscErrorCode SNESQNSetCompositionType(SNES snes, SNESQNCompositionType ctype) {
5840c777b0cSPeter Brune   PetscErrorCode ierr;
5850c777b0cSPeter Brune   PetscFunctionBegin;
5860c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
5870c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetCompositionType_C",(SNES,SNESQNCompositionType),(snes,ctype));CHKERRQ(ierr);
5880c777b0cSPeter Brune   PetscFunctionReturn(0);
5890c777b0cSPeter Brune }
5900c777b0cSPeter Brune 
5910c777b0cSPeter Brune EXTERN_C_BEGIN
5920c777b0cSPeter Brune #undef __FUNCT__
5930c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN"
5940c777b0cSPeter Brune PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) {
5950c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
5960c777b0cSPeter Brune   PetscFunctionBegin;
5970c777b0cSPeter Brune   qn->scale_type = stype;
5980c777b0cSPeter Brune   PetscFunctionReturn(0);
5990c777b0cSPeter Brune }
6000c777b0cSPeter Brune 
6010c777b0cSPeter Brune #undef __FUNCT__
6020c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN"
6030c777b0cSPeter Brune PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) {
6040c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
6050c777b0cSPeter Brune   PetscFunctionBegin;
6060c777b0cSPeter Brune   qn->restart_type = rtype;
6070c777b0cSPeter Brune   PetscFunctionReturn(0);
6080c777b0cSPeter Brune }
6090c777b0cSPeter Brune 
6100c777b0cSPeter Brune #undef __FUNCT__
6110c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetCompositionType_QN"
6120c777b0cSPeter Brune 
6130c777b0cSPeter Brune PetscErrorCode SNESQNSetCompositionType_QN(SNES snes, SNESQNCompositionType ctype) {
6140c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
6150c777b0cSPeter Brune   PetscFunctionBegin;
6160c777b0cSPeter Brune   qn->composition_type = ctype;
6170c777b0cSPeter Brune   PetscFunctionReturn(0);
6180c777b0cSPeter Brune }
6190c777b0cSPeter Brune EXTERN_C_END
6200c777b0cSPeter Brune 
6214b11644fSPeter Brune /* -------------------------------------------------------------------------- */
6224b11644fSPeter Brune /*MC
6234b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
6244b11644fSPeter Brune 
6256cc8130cSPeter Brune       Options Database:
6266cc8130cSPeter Brune 
6276cc8130cSPeter Brune +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
6281867fe5bSPeter Brune .     -snes_qn_powell_angle - Angle condition for restart.
6291867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
6301867fe5bSPeter Brune .     -snes_qn_composition <sequential, composed>- Type of composition.
63172835e02SPeter Brune .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
63244f7e39eSPeter Brune -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
6336cc8130cSPeter Brune 
63444f7e39eSPeter Brune       Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to
6356cc8130cSPeter Brune       form the approximate inverse Jacobian using a series of multiplicative rank-one updates.  This will eventually be
6366cc8130cSPeter Brune       generalized to implement several limited-memory Broyden methods.
6376cc8130cSPeter Brune 
6381867fe5bSPeter Brune       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
6391867fe5bSPeter Brune       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
6401867fe5bSPeter Brune       iteration as the current iteration's values when constructing the approximate jacobian.  The second, composed,
6411867fe5bSPeter Brune       perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
6421867fe5bSPeter Brune 
6436cc8130cSPeter Brune       References:
6446cc8130cSPeter Brune 
6456cc8130cSPeter Brune       L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed,
6466cc8130cSPeter Brune       International Journal of Computer Mathematics, vol. 86, 2009.
6476cc8130cSPeter Brune 
6484b11644fSPeter Brune 
6494b11644fSPeter Brune       Level: beginner
6504b11644fSPeter Brune 
6514b11644fSPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
6526cc8130cSPeter Brune 
6534b11644fSPeter Brune M*/
6544b11644fSPeter Brune EXTERN_C_BEGIN
6554b11644fSPeter Brune #undef __FUNCT__
6564b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
6574b11644fSPeter Brune PetscErrorCode  SNESCreate_QN(SNES snes)
6584b11644fSPeter Brune {
6594b11644fSPeter Brune 
6604b11644fSPeter Brune   PetscErrorCode ierr;
6619f83bee8SJed Brown   SNES_QN *qn;
6624b11644fSPeter Brune 
6634b11644fSPeter Brune   PetscFunctionBegin;
6644b11644fSPeter Brune   snes->ops->setup           = SNESSetUp_QN;
6654b11644fSPeter Brune   snes->ops->solve           = SNESSolve_QN;
6664b11644fSPeter Brune   snes->ops->destroy         = SNESDestroy_QN;
6674b11644fSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_QN;
6684b11644fSPeter Brune   snes->ops->view            = 0;
6694b11644fSPeter Brune   snes->ops->reset           = SNESReset_QN;
6704b11644fSPeter Brune 
67142f4f86dSBarry Smith   snes->usespc          = PETSC_TRUE;
67242f4f86dSBarry Smith   snes->usesksp         = PETSC_FALSE;
67342f4f86dSBarry Smith 
67488976e71SPeter Brune   if (!snes->tolerancesset) {
6750e444f03SPeter Brune     snes->max_funcs = 30000;
6760e444f03SPeter Brune     snes->max_its   = 10000;
67788976e71SPeter Brune   }
6780e444f03SPeter Brune 
6799f83bee8SJed Brown   ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr);
6804b11644fSPeter Brune   snes->data = (void *) qn;
6810ecc9e1dSPeter Brune   qn->m               = 10;
682b21d5a53SPeter Brune   qn->scaling         = 1.0;
6834b11644fSPeter Brune   qn->dX              = PETSC_NULL;
6846bf1b2e5SPeter Brune   qn->dF              = PETSC_NULL;
6858e231d97SPeter Brune   qn->dXtdF           = PETSC_NULL;
6868e231d97SPeter Brune   qn->dFtdX           = PETSC_NULL;
6878e231d97SPeter Brune   qn->dXdFmat         = PETSC_NULL;
68844f7e39eSPeter Brune   qn->monitor         = PETSC_NULL;
68918aa0c0cSPeter Brune   qn->singlereduction = PETSC_FALSE;
69088d374b2SPeter Brune   qn->powell_gamma    = 0.9;
69188d374b2SPeter Brune   qn->powell_downhill = 0.2;
6920c777b0cSPeter Brune   qn->composition_type= SNES_QN_SEQUENTIAL;
6930c777b0cSPeter Brune   qn->scale_type      = SNES_QN_SCALE_SHANNO;
6940c777b0cSPeter Brune   qn->restart_type    = SNES_QN_RESTART_POWELL;
6950c777b0cSPeter Brune   qn->restart_periodic= -1;
696ea630c6eSPeter Brune 
6974b11644fSPeter Brune   PetscFunctionReturn(0);
6984b11644fSPeter Brune }
6998e231d97SPeter Brune 
7004b11644fSPeter Brune EXTERN_C_END
701