1 #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2 3 #define H(i,j) qn->dXdFmat[i*qn->m + j] 4 5 const char *const SNESQNCompositionTypes[] = {"SEQUENTIAL","COMPOSED","SNESQNCompositionType","SNES_QN_",0}; 6 const char *const SNESQNScaleTypes[] = {"NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0}; 7 const char *const SNESQNRestartTypes[] = {"NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0}; 8 9 typedef struct { 10 Vec *dX; /* The change in X */ 11 Vec *dF; /* The change in F */ 12 PetscInt m; /* The number of kept previous steps */ 13 PetscScalar *alpha, *beta; 14 PetscScalar *dXtdF, *dFtdX, *YtdX; 15 PetscBool singlereduction; /* Aggregated reduction implementation */ 16 PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */ 17 PetscViewer monitor; 18 PetscReal powell_gamma; /* Powell angle restart condition */ 19 PetscReal powell_downhill; /* Powell descent restart condition */ 20 PetscReal scaling; /* scaling of H0 */ 21 PetscInt restart_periodic; /* the maximum iterations between restart */ 22 23 SNESQNCompositionType composition_type; /* determine if the composition is done sequentially or as a composition */ 24 SNESQNScaleType scale_type; /* determine if the composition is done sequentially or as a composition */ 25 SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */ 26 } SNES_QN; 27 28 #undef __FUNCT__ 29 #define __FUNCT__ "SNESQNApplyJinv_Private" 30 PetscErrorCode SNESQNApplyJinv_Private(SNES snes, PetscInt it, Vec D, Vec Y) { 31 32 PetscErrorCode ierr; 33 34 SNES_QN *qn = (SNES_QN*)snes->data; 35 36 Vec Yin = snes->work[3]; 37 38 Vec *dX = qn->dX; 39 Vec *dF = qn->dF; 40 41 PetscScalar *alpha = qn->alpha; 42 PetscScalar *beta = qn->beta; 43 PetscScalar *dXtdF = qn->dXtdF; 44 PetscScalar *YtdX = qn->YtdX; 45 46 /* ksp thing for jacobian scaling */ 47 KSPConvergedReason kspreason; 48 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 49 50 PetscInt k, i, j, g, lits; 51 PetscInt m = qn->m; 52 PetscScalar t; 53 PetscInt l = m; 54 55 Mat jac, jac_pre; 56 57 PetscFunctionBegin; 58 59 ierr = VecCopy(D, Y);CHKERRQ(ierr); 60 61 if (it < m) l = it; 62 63 if (qn->singlereduction) { 64 ierr = VecMDot(Y, l, qn->dX, YtdX);CHKERRQ(ierr); 65 } 66 /* outward recursion starting at iteration k's update and working back */ 67 for (i = 0; i < l; i++) { 68 k = (it - i - 1) % l; 69 if (qn->singlereduction) { 70 /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 71 t = YtdX[k]; 72 for (j = 0; j < i; j++) { 73 g = (it - j - 1) % l; 74 t += -alpha[g]*H(g, k); 75 } 76 alpha[k] = t / H(k, k); 77 } else { 78 ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr); 79 alpha[k] = t / dXtdF[k]; 80 } 81 if (qn->monitor) { 82 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 83 ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 84 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 85 } 86 ierr = VecAXPY(Y, -alpha[k], dF[k]);CHKERRQ(ierr); 87 } 88 89 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 90 ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 91 ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 92 ierr = SNES_KSPSolve(snes,snes->ksp,Y,Yin);CHKERRQ(ierr); 93 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 94 if (kspreason < 0) { 95 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 96 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 97 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 98 PetscFunctionReturn(0); 99 } 100 } 101 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 102 snes->linear_its += lits; 103 ierr = VecCopy(Yin, Y);CHKERRQ(ierr); 104 } else { 105 ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 106 } 107 if (qn->singlereduction) { 108 ierr = VecMDot(Y, l, qn->dF, YtdX);CHKERRQ(ierr); 109 } 110 /* inward recursion starting at the first update and working forward */ 111 for (i = 0; i < l; i++) { 112 k = (it + i - l) % l; 113 if (qn->singlereduction) { 114 t = YtdX[k]; 115 for (j = 0; j < i; j++) { 116 g = (it + j - l) % l; 117 t += (alpha[g] - beta[g])*H(k, g); 118 } 119 beta[k] = t / H(k, k); 120 } else { 121 ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 122 beta[k] = t / dXtdF[k]; 123 } 124 ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]); 125 if (qn->monitor) { 126 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 127 ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 128 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 129 } 130 } 131 PetscFunctionReturn(0); 132 } 133 134 #undef __FUNCT__ 135 #define __FUNCT__ "SNESSolve_QN" 136 static PetscErrorCode SNESSolve_QN(SNES snes) 137 { 138 139 PetscErrorCode ierr; 140 SNES_QN *qn = (SNES_QN*) snes->data; 141 142 Vec X, Xold; 143 Vec F, B; 144 Vec Y, FPC, D, Dold; 145 SNESConvergedReason reason; 146 PetscInt i, i_r, k, l, j; 147 148 PetscReal fnorm, xnorm, ynorm, gnorm; 149 PetscInt m = qn->m; 150 PetscBool lssucceed,powell,periodic; 151 152 Vec *dX = qn->dX; 153 Vec *dF = qn->dF; 154 PetscScalar *dXtdF = qn->dXtdF; 155 PetscScalar *dFtdX = qn->dFtdX; 156 PetscScalar DolddotD, DolddotDold, DdotD, YdotD, a; 157 158 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 159 160 /* basically just a regular newton's method except for the application of the jacobian */ 161 PetscFunctionBegin; 162 163 X = snes->vec_sol; /* solution vector */ 164 F = snes->vec_func; /* residual vector */ 165 Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 166 B = snes->vec_rhs; 167 Xold = snes->work[0]; 168 169 /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 170 D = snes->work[1]; 171 Dold = snes->work[2]; 172 173 snes->reason = SNES_CONVERGED_ITERATING; 174 175 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 176 snes->iter = 0; 177 snes->norm = 0.; 178 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 179 if (!snes->vec_func_init_set){ 180 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 181 if (snes->domainerror) { 182 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 183 PetscFunctionReturn(0); 184 } 185 } else { 186 snes->vec_func_init_set = PETSC_FALSE; 187 } 188 189 if (!snes->norm_init_set) { 190 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 191 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 192 } else { 193 fnorm = snes->norm_init; 194 snes->norm_init_set = PETSC_FALSE; 195 } 196 197 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 198 snes->norm = fnorm; 199 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 200 SNESLogConvHistory(snes,fnorm,0); 201 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 202 203 /* set parameter for default relative tolerance convergence test */ 204 snes->ttol = fnorm*snes->rtol; 205 /* test convergence */ 206 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 207 if (snes->reason) PetscFunctionReturn(0); 208 209 /* composed solve -- either sequential or composed */ 210 if (snes->pc) { 211 if (qn->composition_type == SNES_QN_SEQUENTIAL) { 212 ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 213 ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 214 ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 215 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 216 if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 217 snes->reason = SNES_DIVERGED_INNER; 218 PetscFunctionReturn(0); 219 } 220 ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 221 ierr = VecCopy(FPC, F);CHKERRQ(ierr); 222 ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 223 ierr = VecCopy(F, Y);CHKERRQ(ierr); 224 } else { 225 ierr = VecCopy(X, Y);CHKERRQ(ierr); 226 ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 227 ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 228 ierr = SNESSolve(snes->pc, B, Y);CHKERRQ(ierr); 229 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 230 if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 231 snes->reason = SNES_DIVERGED_INNER; 232 PetscFunctionReturn(0); 233 } 234 ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); 235 } 236 } else { 237 ierr = VecCopy(F, Y);CHKERRQ(ierr); 238 } 239 ierr = VecCopy(Y, D);CHKERRQ(ierr); 240 241 /* scale the initial update */ 242 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 243 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 244 } 245 246 for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 247 ierr = SNESQNApplyJinv_Private(snes, i_r, D, Y);CHKERRQ(ierr); 248 /* line search for lambda */ 249 ynorm = 1; gnorm = fnorm; 250 ierr = VecCopy(D, Dold);CHKERRQ(ierr); 251 ierr = VecCopy(X, Xold);CHKERRQ(ierr); 252 ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 253 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 254 if (snes->domainerror) { 255 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 256 PetscFunctionReturn(0); 257 } 258 ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 259 if (!lssucceed) { 260 if (++snes->numFailures >= snes->maxFailures) { 261 snes->reason = SNES_DIVERGED_LINE_SEARCH; 262 break; 263 } 264 } 265 ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 266 if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 267 ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 268 } 269 270 /* convergence monitoring */ 271 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); 272 273 ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 274 ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 275 276 SNESLogConvHistory(snes,snes->norm,snes->iter); 277 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 278 /* set parameter for default relative tolerance convergence test */ 279 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 280 if (snes->reason) PetscFunctionReturn(0); 281 282 283 if (snes->pc) { 284 if (qn->composition_type == SNES_QN_SEQUENTIAL) { 285 ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 286 ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 287 ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 288 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 289 if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 290 snes->reason = SNES_DIVERGED_INNER; 291 PetscFunctionReturn(0); 292 } 293 ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 294 ierr = VecCopy(FPC, F);CHKERRQ(ierr); 295 ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 296 ierr = VecCopy(F, D);CHKERRQ(ierr); 297 } else { 298 ierr = VecCopy(X, D);CHKERRQ(ierr); 299 ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 300 ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 301 ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr); 302 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 303 if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 304 snes->reason = SNES_DIVERGED_INNER; 305 PetscFunctionReturn(0); 306 } 307 ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr); 308 } 309 } else { 310 ierr = VecCopy(F, D);CHKERRQ(ierr); 311 } 312 313 powell = PETSC_FALSE; 314 if (qn->restart_type == SNES_QN_RESTART_POWELL) { 315 /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 316 ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 317 ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr); 318 ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr); 319 ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr); 320 ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 321 ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr); 322 ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr); 323 ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr); 324 if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 325 } 326 periodic = PETSC_FALSE; 327 if (qn->restart_type != SNES_QN_RESTART_NONE) { 328 if ((i_r > qn->restart_periodic - 1 && qn->restart_periodic > 0)) periodic = PETSC_TRUE; 329 } 330 331 /* restart if either powell or periodic restart is satisfied. */ 332 if (powell || periodic) { 333 if (qn->monitor) { 334 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 335 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); 336 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 337 } 338 i_r = -1; 339 /* general purpose update */ 340 if (snes->ops->update) { 341 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 342 } 343 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 344 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 345 } 346 } else { 347 /* set the differences */ 348 k = i_r % m; 349 l = m; 350 if (i_r + 1 < m) l = i_r + 1; 351 ierr = VecCopy(D, dF[k]);CHKERRQ(ierr); 352 ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 353 ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 354 ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 355 if (qn->singlereduction) { 356 ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr); 357 ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr); 358 for (j = 0; j < l; j++) { 359 H(k, j) = dFtdX[j]; 360 H(j, k) = dXtdF[j]; 361 } 362 /* copy back over to make the computation of alpha and beta easier */ 363 for (j = 0; j < l; j++) { 364 dXtdF[j] = H(j, j); 365 } 366 } else { 367 ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 368 } 369 /* set scaling to be shanno scaling */ 370 if (qn->scale_type == SNES_QN_SCALE_SHANNO) { 371 ierr = VecDot(dF[k], dF[k], &a);CHKERRQ(ierr); 372 qn->scaling = PetscRealPart(dXtdF[k]) / PetscRealPart(a); 373 } 374 /* general purpose update */ 375 if (snes->ops->update) { 376 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 377 } 378 } 379 } 380 if (i == snes->max_its) { 381 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 382 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 383 } 384 PetscFunctionReturn(0); 385 } 386 387 388 #undef __FUNCT__ 389 #define __FUNCT__ "SNESSetUp_QN" 390 static PetscErrorCode SNESSetUp_QN(SNES snes) 391 { 392 SNES_QN *qn = (SNES_QN*)snes->data; 393 PetscErrorCode ierr; 394 395 PetscFunctionBegin; 396 ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 397 ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr); 398 ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, 399 qn->m, PetscScalar, &qn->beta, 400 qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr); 401 402 if (qn->singlereduction) { 403 ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 404 qn->m, PetscScalar, &qn->dFtdX, 405 qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 406 } 407 ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 408 409 /* set up the line search */ 410 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 411 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 412 } 413 PetscFunctionReturn(0); 414 } 415 416 #undef __FUNCT__ 417 #define __FUNCT__ "SNESReset_QN" 418 static PetscErrorCode SNESReset_QN(SNES snes) 419 { 420 PetscErrorCode ierr; 421 SNES_QN *qn; 422 PetscFunctionBegin; 423 if (snes->data) { 424 qn = (SNES_QN*)snes->data; 425 if (qn->dX) { 426 ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 427 } 428 if (qn->dF) { 429 ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr); 430 } 431 if (qn->singlereduction) { 432 ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 433 } 434 ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr); 435 } 436 PetscFunctionReturn(0); 437 } 438 439 #undef __FUNCT__ 440 #define __FUNCT__ "SNESDestroy_QN" 441 static PetscErrorCode SNESDestroy_QN(SNES snes) 442 { 443 PetscErrorCode ierr; 444 PetscFunctionBegin; 445 ierr = SNESReset_QN(snes);CHKERRQ(ierr); 446 ierr = PetscFree(snes->data);CHKERRQ(ierr); 447 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr); 448 PetscFunctionReturn(0); 449 } 450 451 #undef __FUNCT__ 452 #define __FUNCT__ "SNESSetFromOptions_QN" 453 static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 454 { 455 456 PetscErrorCode ierr; 457 SNES_QN *qn; 458 PetscBool monflg = PETSC_FALSE; 459 SNESLineSearch linesearch; 460 PetscFunctionBegin; 461 462 qn = (SNES_QN*)snes->data; 463 464 ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 465 ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,PETSC_NULL);CHKERRQ(ierr); 466 ierr = PetscOptionsInt("-snes_qn_restart","Maximum number of iterations between restarts","SNESQN",qn->restart_periodic,&qn->restart_periodic, PETSC_NULL);CHKERRQ(ierr); 467 ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr); 468 ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr); 469 ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr); 470 ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, PETSC_NULL);CHKERRQ(ierr); 471 ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)qn->scale_type,(PetscEnum*)&qn->scale_type,PETSC_NULL);CHKERRQ(ierr); 472 ierr = PetscOptionsEnum("-snes_qn_composition_type","Composition type","SNESQNSetCompositionType",SNESQNCompositionTypes, 473 (PetscEnum)qn->composition_type,(PetscEnum*)&qn->composition_type,PETSC_NULL);CHKERRQ(ierr); 474 ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes, 475 (PetscEnum)qn->restart_type,(PetscEnum*)&qn->restart_type,PETSC_NULL);CHKERRQ(ierr); 476 ierr = PetscOptionsTail();CHKERRQ(ierr); 477 if (!snes->linesearch) { 478 ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 479 if (!snes->pc || qn->composition_type == SNES_QN_SEQUENTIAL) { 480 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 481 } else { 482 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 483 } 484 } 485 if (monflg) { 486 qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 487 } 488 PetscFunctionReturn(0); 489 } 490 491 492 #undef __FUNCT__ 493 #define __FUNCT__ "SNESQNSetRestartType" 494 /*@ 495 SNESQNSetRestartType - Sets the restart type for SNESQN. 496 497 Logically Collective on SNES 498 499 Input Parameters: 500 + snes - the iterative context 501 - rtype - restart type 502 503 Options Database: 504 + -snes_qn_restart_type<powell,periodic,none> - set the restart type 505 - -snes_qn_restart[30] - sets the number of iterations before restart for periodic 506 507 Level: intermediate 508 509 SNESQNRestartTypes: 510 + SNES_QN_RESTART_NONE - never restart 511 . SNES_QN_RESTART_POWELL - restart based upon descent criteria 512 - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 513 514 Notes: 515 The default line search used is the L2 line search and it requires two additional function evaluations. 516 517 .keywords: SNES, SNESQN, restart, type, set SNESLineSearch 518 @*/ 519 PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) { 520 PetscErrorCode ierr; 521 PetscFunctionBegin; 522 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 523 ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 524 PetscFunctionReturn(0); 525 } 526 527 #undef __FUNCT__ 528 #define __FUNCT__ "SNESQNSetScaleType" 529 /*@ 530 SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN. 531 532 Logically Collective on SNES 533 534 Input Parameters: 535 + snes - the iterative context 536 - stype - scale type 537 538 Options Database: 539 . -snes_qn_scale_type<shanno,none,linesearch,jacobian> 540 541 Level: intermediate 542 543 SNESQNSelectTypes: 544 + SNES_QN_SCALE_NONE - don't scale the problem 545 . SNES_QN_SCALE_SHANNO - use shanno scaling 546 . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 547 - SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian. 548 549 .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 550 @*/ 551 552 PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) { 553 PetscErrorCode ierr; 554 PetscFunctionBegin; 555 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 556 ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 557 PetscFunctionReturn(0); 558 } 559 560 #undef __FUNCT__ 561 #define __FUNCT__ "SNESQNSetCompositionType" 562 /*@ 563 SNESQNSetCompositionType - Sets the composition type 564 565 Logically Collective on SNES 566 567 Input Parameters: 568 + snes - the iterative context 569 - stype - composition type 570 571 Options Database: 572 . -snes_qn_composition_type<sequential, composed> 573 574 Level: intermediate 575 576 SNESQNSelectTypes: 577 + SNES_QN_COMPOSITION_SEQUENTIAL - Solve the system with X = PC(X) and D = F(PC(X)) 578 - SNES_QN_COMPOSITION_COMPOSED - solve the system with X = X and D = PC(X) - X 579 580 .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 581 @*/ 582 583 PetscErrorCode SNESQNSetCompositionType(SNES snes, SNESQNCompositionType ctype) { 584 PetscErrorCode ierr; 585 PetscFunctionBegin; 586 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 587 ierr = PetscTryMethod(snes,"SNESQNSetCompositionType_C",(SNES,SNESQNCompositionType),(snes,ctype));CHKERRQ(ierr); 588 PetscFunctionReturn(0); 589 } 590 591 EXTERN_C_BEGIN 592 #undef __FUNCT__ 593 #define __FUNCT__ "SNESQNSetScaleType_QN" 594 PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) { 595 SNES_QN *qn = (SNES_QN *)snes->data; 596 PetscFunctionBegin; 597 qn->scale_type = stype; 598 PetscFunctionReturn(0); 599 } 600 601 #undef __FUNCT__ 602 #define __FUNCT__ "SNESQNSetRestartType_QN" 603 PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) { 604 SNES_QN *qn = (SNES_QN *)snes->data; 605 PetscFunctionBegin; 606 qn->restart_type = rtype; 607 PetscFunctionReturn(0); 608 } 609 610 #undef __FUNCT__ 611 #define __FUNCT__ "SNESQNSetCompositionType_QN" 612 613 PetscErrorCode SNESQNSetCompositionType_QN(SNES snes, SNESQNCompositionType ctype) { 614 SNES_QN *qn = (SNES_QN *)snes->data; 615 PetscFunctionBegin; 616 qn->composition_type = ctype; 617 PetscFunctionReturn(0); 618 } 619 EXTERN_C_END 620 621 /* -------------------------------------------------------------------------- */ 622 /*MC 623 SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 624 625 Options Database: 626 627 + -snes_qn_m - Number of past states saved for the L-Broyden methods. 628 . -snes_qn_powell_angle - Angle condition for restart. 629 . -snes_qn_powell_descent - Descent condition for restart. 630 . -snes_qn_composition <sequential, composed>- Type of composition. 631 . -snes_linesearch_type <cp, l2, basic> - Type of line search. 632 - -snes_qn_monitor - Monitors the quasi-newton jacobian. 633 634 Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to 635 form the approximate inverse Jacobian using a series of multiplicative rank-one updates. This will eventually be 636 generalized to implement several limited-memory Broyden methods. 637 638 When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 639 these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 640 iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 641 perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 642 643 References: 644 645 L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 646 International Journal of Computer Mathematics, vol. 86, 2009. 647 648 649 Level: beginner 650 651 .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 652 653 M*/ 654 EXTERN_C_BEGIN 655 #undef __FUNCT__ 656 #define __FUNCT__ "SNESCreate_QN" 657 PetscErrorCode SNESCreate_QN(SNES snes) 658 { 659 660 PetscErrorCode ierr; 661 SNES_QN *qn; 662 663 PetscFunctionBegin; 664 snes->ops->setup = SNESSetUp_QN; 665 snes->ops->solve = SNESSolve_QN; 666 snes->ops->destroy = SNESDestroy_QN; 667 snes->ops->setfromoptions = SNESSetFromOptions_QN; 668 snes->ops->view = 0; 669 snes->ops->reset = SNESReset_QN; 670 671 snes->usespc = PETSC_TRUE; 672 snes->usesksp = PETSC_FALSE; 673 674 if (!snes->tolerancesset) { 675 snes->max_funcs = 30000; 676 snes->max_its = 10000; 677 } 678 679 ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 680 snes->data = (void *) qn; 681 qn->m = 10; 682 qn->scaling = 1.0; 683 qn->dX = PETSC_NULL; 684 qn->dF = PETSC_NULL; 685 qn->dXtdF = PETSC_NULL; 686 qn->dFtdX = PETSC_NULL; 687 qn->dXdFmat = PETSC_NULL; 688 qn->monitor = PETSC_NULL; 689 qn->singlereduction = PETSC_FALSE; 690 qn->powell_gamma = 0.9; 691 qn->powell_downhill = 0.2; 692 qn->composition_type= SNES_QN_SEQUENTIAL; 693 qn->scale_type = SNES_QN_SCALE_SHANNO; 694 qn->restart_type = SNES_QN_RESTART_POWELL; 695 qn->restart_periodic= -1; 696 697 PetscFunctionReturn(0); 698 } 699 700 EXTERN_C_END 701