/*@ SNESGetNPCFunction - Gets the function from a preconditioner after SNESSolve() has been called. Collective on SNES Input Parameters: . snes - the SNES context Output Parameter: . F - function vector . fnorm - the norm of F Level: developer .keywords: SNES, nonlinear, function .seealso: SNESGetNPC(),SNESSetNPC(),SNESComputeFunction(),SNESApplyNPC(),SNESSolve() @*/ PetscErrorCode SNESGetNPCFunction(SNES snes,Vec F,PetscReal *fnorm) { PetscErrorCode ierr; PCSide npcside; SNESFunctionType functype; SNESNormSchedule normschedule; Vec FPC,XPC; PetscFunctionBegin; if (snes->pc) { ierr = SNESGetNPCSide(snes->pc,&npcside);CHKERRQ(ierr); ierr = SNESGetFunctionType(snes->pc,&functype);CHKERRQ(ierr); ierr = SNESGetNormSchedule(snes->pc,&normschedule);CHKERRQ(ierr); /* check if the function is valid based upon how the inner solver is preconditioned */ if (normschedule != SNES_NORM_NONE && normschedule != SNES_NORM_INITIAL_ONLY && (npcside == PC_RIGHT || functype == SNES_FUNCTION_UNPRECONDITIONED)) { ierr = SNESGetFunction(snes->pc,&FPC,NULL,NULL);CHKERRQ(ierr); if (FPC) { if (fnorm) {ierr = SNESGetFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr);} ierr = VecCopy(FPC,F);CHKERRQ(ierr); } else { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no function"); } } else { ierr = SNESGetSolution(snes->pc,&XPC);CHKERRQ(ierr); if (XPC) { ierr = SNESComputeFunction(snes->pc,XPC,F);CHKERRQ(ierr); if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Preconditioner has no solution"); } } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No preconditioner set"); PetscFunctionReturn(0); }
PetscErrorCode SNESSolve_Composite(SNES snes) { Vec F; Vec X; Vec B; PetscInt i; PetscReal fnorm = 0.0, xnorm = 0.0, snorm = 0.0; PetscErrorCode ierr; SNESNormSchedule normtype; SNES_Composite *comp = (SNES_Composite*)snes->data; PetscFunctionBegin; X = snes->vec_sol; F = snes->vec_func; B = snes->vec_rhs; ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); snes->iter = 0; snes->norm = 0.; comp->innerFailures = 0; ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); ierr = SNESSetWorkVecs(snes, 1);CHKERRQ(ierr); snes->reason = SNES_CONVERGED_ITERATING; ierr = SNESGetNormSchedule(snes, &normtype);CHKERRQ(ierr); if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { if (!snes->vec_func_init_set) { ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); } else snes->vec_func_init_set = PETSC_FALSE; if (snes->xl && snes->xu) { ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr); } else { ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ } SNESCheckFunctionNorm(snes,fnorm); ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); snes->iter = 0; snes->norm = fnorm; ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); /* test convergence */ ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); if (snes->reason) PetscFunctionReturn(0); } else { ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); } /* Call general purpose update function */ if (snes->ops->update) { ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); } for (i = 0; i < snes->max_its; i++) { /* Copy the state before modification by application of the composite solver; we will subtract the new state after application */ ierr = VecCopy(X, snes->work[0]);CHKERRQ(ierr); if (comp->type == SNES_COMPOSITE_ADDITIVE) { ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr); } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) { ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr); } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr); } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type"); if (snes->reason < 0) break; /* Compute the solution update for convergence testing */ ierr = VecAXPY(snes->work[0], -1.0, X);CHKERRQ(ierr); ierr = VecScale(snes->work[0], -1.0);CHKERRQ(ierr); if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) { ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); if (snes->xl && snes->xu) { ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr); ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); } else { ierr = VecNormBegin(F, NORM_2, &fnorm);CHKERRQ(ierr); ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); ierr = VecNormEnd(F, NORM_2, &fnorm);CHKERRQ(ierr); ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); } SNESCheckFunctionNorm(snes,fnorm); } else if (normtype == SNES_NORM_ALWAYS) { ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); } /* Monitor convergence */ ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); snes->iter = i+1; snes->norm = fnorm; ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); /* Test for convergence */ if (normtype == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,snorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} if (snes->reason) break; /* Call general purpose update function */ if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} } if (normtype == SNES_NORM_ALWAYS) { if (i == snes->max_its) { ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; } } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; PetscFunctionReturn(0); }
PetscErrorCode SNESSolve_NASM(SNES snes) { Vec F; Vec X; Vec B; Vec Y; PetscInt i; PetscReal fnorm = 0.0; PetscErrorCode ierr; SNESNormSchedule normschedule; SNES_NASM *nasm = (SNES_NASM*)snes->data; PetscFunctionBegin; X = snes->vec_sol; Y = snes->vec_sol_update; F = snes->vec_func; B = snes->vec_rhs; ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); snes->iter = 0; snes->norm = 0.; ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); snes->reason = SNES_CONVERGED_ITERATING; ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr); if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { /* compute the initial function and preconditioned update delX */ if (!snes->vec_func_init_set) { ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); if (snes->domainerror) { snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; PetscFunctionReturn(0); } } else snes->vec_func_init_set = PETSC_FALSE; ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ if (PetscIsInfOrNanReal(fnorm)) { snes->reason = SNES_DIVERGED_FNORM_NAN; PetscFunctionReturn(0); } ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); snes->iter = 0; snes->norm = fnorm; ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); /* test convergence */ ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); if (snes->reason) PetscFunctionReturn(0); } else { ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); } /* Call general purpose update function */ if (snes->ops->update) { ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); } /* copy the initial solution over for later */ if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);} for (i = 0; i < snes->max_its; i++) { ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); if (snes->domainerror) { snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; break; } ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ if (PetscIsInfOrNanReal(fnorm)) { snes->reason = SNES_DIVERGED_FNORM_NAN; break; } } /* Monitor convergence */ ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); snes->iter = i+1; snes->norm = fnorm; ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); /* Test for convergence */ if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} if (snes->reason) break; /* Call general purpose update function */ if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} } if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);} if (normschedule == SNES_NORM_ALWAYS) { if (i == snes->max_its) { ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; } } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ PetscFunctionReturn(0); }
PetscErrorCode SNESSolve_NGS(SNES snes) { Vec F; Vec X; Vec B; PetscInt i; PetscReal fnorm; PetscErrorCode ierr; SNESNormSchedule normschedule; PetscFunctionBegin; ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); X = snes->vec_sol; F = snes->vec_func; B = snes->vec_rhs; ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); snes->iter = 0; snes->norm = 0.; ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); snes->reason = SNES_CONVERGED_ITERATING; ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr); if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { /* compute the initial function and preconditioned update delX */ if (!snes->vec_func_init_set) { ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); if (snes->domainerror) { snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; PetscFunctionReturn(0); } } else snes->vec_func_init_set = PETSC_FALSE; ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ if (PetscIsInfOrNanReal(fnorm)) { snes->reason = SNES_DIVERGED_FNORM_NAN; PetscFunctionReturn(0); } ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); snes->iter = 0; snes->norm = fnorm; ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); /* test convergence */ ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); if (snes->reason) PetscFunctionReturn(0); } else { ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); } /* Call general purpose update function */ if (snes->ops->update) { ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); } for (i = 0; i < snes->max_its; i++) { ierr = SNESComputeNGS(snes, B, X);CHKERRQ(ierr); /* only compute norms if requested or about to exit due to maximum iterations */ if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); if (snes->domainerror) { snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; PetscFunctionReturn(0); } ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ if (PetscIsInfOrNanReal(fnorm)) { snes->reason = SNES_DIVERGED_FNORM_NAN; PetscFunctionReturn(0); } } /* Monitor convergence */ ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); snes->iter = i+1; snes->norm = fnorm; ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); /* Test for convergence */ if (normschedule == SNES_NORM_ALWAYS) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); if (snes->reason) PetscFunctionReturn(0); /* Call general purpose update function */ if (snes->ops->update) { ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); } } if (normschedule == SNES_NORM_ALWAYS) { if (i == snes->max_its) { ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; } } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */ PetscFunctionReturn(0); }