| OOFEM
    2.4
    OOFEM.org - Object Oriented Finite Element Solver | 
This class implements Newton-Raphson Method, derived from abstract NumericalMethod class for solving non-linear problems. More...
#include <nrsolver.h>
 Inheritance diagram for oofem::NRSolver:
 Inheritance diagram for oofem::NRSolver: Collaboration diagram for oofem::NRSolver:
 Collaboration diagram for oofem::NRSolver:| Public Member Functions | |
| NRSolver (Domain *d, EngngModel *m) | |
| virtual | ~NRSolver () | 
| virtual NM_Status | solve (SparseMtrx &k, FloatArray &R, FloatArray *R0, FloatArray &X, FloatArray &dX, FloatArray &F, const FloatArray &internalForcesEBENorm, double &l, referenceLoadInputModeType rlm, int &nite, TimeStep *) | 
| Solves the given sparse linear system of equations  .  More... | |
| virtual void | printState (FILE *outputStream) | 
| Prints status message of receiver to output stream.  More... | |
| virtual IRResultType | initializeFrom (InputRecord *ir) | 
| virtual const char * | giveClassName () const | 
| virtual const char * | giveInputRecordName () const | 
| virtual void | setDomain (Domain *d) | 
| virtual void | reinitialize () | 
| Reinitializes the receiver.  More... | |
| virtual SparseLinearSystemNM * | giveLinearSolver () | 
| Constructs (if necessary) and returns a linear solver.  More... | |
|  Public Member Functions inherited from oofem::SparseNonLinearSystemNM | |
| SparseNonLinearSystemNM (Domain *d, EngngModel *m) | |
| Constructor.  More... | |
| virtual | ~SparseNonLinearSystemNM () | 
| Destructor.  More... | |
| virtual double | giveCurrentStepLength () | 
| Returns step length.  More... | |
| virtual void | setStepLength (double s) | 
| Sets the step length.  More... | |
| virtual bool | referenceLoad () const | 
| Returns true if reference loads are used (i.e.  More... | |
| IRResultType | initializeFrom (InputRecord *ir) | 
| virtual void | convertPertMap () | 
| virtual void | applyPerturbation (FloatArray *displacement) | 
| std::string | errorInfo (const char *func) const | 
| Error printing helper.  More... | |
|  Public Member Functions inherited from oofem::NumericalMethod | |
| NumericalMethod (Domain *d, EngngModel *m) | |
| Constructor.  More... | |
| virtual | ~NumericalMethod () | 
| Destructor.  More... | |
| EngngModel * | giveEngngModel () | 
| virtual contextIOResultType | saveContext (DataStream &stream, ContextMode mode, void *obj=NULL) | 
| virtual contextIOResultType | restoreContext (DataStream &stream, ContextMode mode, void *obj=NULL) | 
| Protected Types | |
| enum | nrsolver_ModeType { nrsolverModifiedNRM, nrsolverFullNRM, nrsolverAccelNRM } | 
| Protected Member Functions | |
| LineSearchNM * | giveLineSearchSolver () | 
| Constructs and returns a line search solver.  More... | |
| void | initPrescribedEqs () | 
| Initiates prescribed equations.  More... | |
| void | applyConstraintsToStiffness (SparseMtrx &k) | 
| void | applyConstraintsToLoadIncrement (int nite, const SparseMtrx &k, FloatArray &R, referenceLoadInputModeType rlm, TimeStep *tStep) | 
| bool | checkConvergence (FloatArray &RT, FloatArray &F, FloatArray &rhs, FloatArray &ddX, FloatArray &X, double RRT, const FloatArray &internalForcesEBENorm, int nite, bool &errorOutOfRange) | 
| Determines whether or not the solution has reached convergence.  More... | |
| Protected Attributes | |
| int | nsmax | 
| int | minIterations | 
| double | minStepLength | 
| nrsolver_ModeType | NR_Mode | 
| nrsolver_ModeType | NR_OldMode | 
| int | NR_ModeTick | 
| int | MANRMSteps | 
| std::unique_ptr< SparseLinearSystemNM > | linSolver | 
| linear system solver  More... | |
| LinSystSolverType | solverType | 
| linear system solver ID  More... | |
| SparseMtrx::SparseMtrxVersionType | smConstraintVersion | 
| sparse matrix version, used to control constrains application to stiffness  More... | |
| int | numberOfPrescribedDofs | 
| number of prescribed displacement  More... | |
| bool | prescribedDofsFlag | 
| Flag indicating that some dofs are controlled under displacement control.  More... | |
| IntArray | prescribedDofs | 
| Array of pairs identifying prescribed dofs (node, dof)  More... | |
| FloatArray | prescribedDofsValues | 
| Array of prescribed values.  More... | |
| int | prescribedDisplacementTF | 
| Load Time Function of prescribed values.  More... | |
| IntArray | prescribedEqs | 
| Array of prescribed equations.  More... | |
| bool | prescribedEqsInitFlag | 
| Flag indicating that prescribedEqs were initialized.  More... | |
| FloatArray | lastReactions | 
| Computed reactions. They are stored in order to print them in printState method.  More... | |
| bool | lsFlag | 
| Flag indicating whether to use line-search.  More... | |
| std::unique_ptr< LineSearchNM > | linesearchSolver | 
| Line search solver.  More... | |
| bool | mCalcStiffBeforeRes | 
| Flag indicating if the stiffness should be evaluated before the residual in the first iteration.  More... | |
| bool | constrainedNRFlag | 
| Flag indicating whether to use constrained Newton.  More... | |
| double | constrainedNRalpha | 
| Scale factor for dX, dX_new = alpha * dX.  More... | |
| int | constrainedNRminiter | 
| Minimum number of iterations before constraint is activated.  More... | |
| FloatArray | rtolf | 
| Relative unbalanced force tolerance for each group.  More... | |
| FloatArray | rtold | 
| Relative iterative displacement change tolerance for each group.  More... | |
| FloatArray | forceErrVec | 
| FloatArray | forceErrVecOld | 
| std::map< int, double > | dg_forceScale | 
| Optional user supplied scale of forces used in convergence check.  More... | |
| double | maxIncAllowed | 
|  Protected Attributes inherited from oofem::SparseNonLinearSystemNM | |
| double | deltaL | 
| Load level.  More... | |
| double | randPertAmplitude | 
| Amplitude of a random perturbation applied on the solution before the iteration process.  More... | |
| int | randSeed | 
| bool | pert_init_needed | 
| IntArray | igp_PertDmanDofSrcArray | 
| FloatArray | igp_PertWeightArray | 
| IntArray | igp_Map | 
| FloatArray | igp_Weight | 
|  Protected Attributes inherited from oofem::NumericalMethod | |
| Domain * | domain | 
| Pointer to domain.  More... | |
| EngngModel * | engngModel | 
| Pointer to engineering model.  More... | |
| Additional Inherited Members | |
|  Public Types inherited from oofem::SparseNonLinearSystemNM | |
| enum | referenceLoadInputModeType { rlm_total =0, rlm_incremental =1 } | 
| The following parameter allows to specify how the reference load vector is obtained from given totalLoadVector and initialLoadVector.  More... | |
This class implements Newton-Raphson Method, derived from abstract NumericalMethod class for solving non-linear problems.
Except the traditional load control, it also provides direct displacement control without requiring BC applied, so that the equation renumbering is not required when combining arc-length and Newton-Raphson solvers within a simulation.
The direct displacement control is achieved by adding a large number alpha to the corresponding diagonal member of K and replacing the right-hand side by alpha*prescribed_value. If alpha is very much larger than other stiffness coefficients then this alteration effectively replaces the corresponding equation by alpha*unknown_value = alpha*prescribed_value that is, the required condition, but the whole system remains symmetric and minimal changes are necessary in the computational sequence. The above artifice has been introduced by Payne and Irons.
Definition at line 95 of file nrsolver.h.
| 
 | protected | 
| Enumerator | |
|---|---|
| nrsolverModifiedNRM | |
| nrsolverFullNRM | |
| nrsolverAccelNRM | |
Definition at line 98 of file nrsolver.h.
| oofem::NRSolver::NRSolver | ( | Domain * | d, | 
| EngngModel * | m | ||
| ) | 
| 
 | virtual | 
Definition at line 97 of file nrsolver.C.
| 
 | protected | 
Definition at line 498 of file nrsolver.C.
References oofem::IntArray::at(), oofem::FloatArray::at(), oofem::SparseMtrx::at(), oofem::PetscSparseMtrx::createVecGlobal(), oofem::NumericalMethod::engngModel, oofem::Function::evaluateAtTime(), oofem::EngngModel::giveDomain(), oofem::Domain::giveFunction(), oofem::IntArray::giveSize(), oofem::TimeStep::giveTargetTime(), oofem::TimeStep::giveTimeIncrement(), oofem::TimeStep::isTheFirstStep(), numberOfPrescribedDofs, prescribedDisplacementTF, prescribedDofsValues, prescribedEqs, oofem::SparseNonLinearSystemNM::rlm_total, solverType, and oofem::ST_Petsc.
Referenced by solve().
| 
 | protected | 
Definition at line 454 of file nrsolver.C.
References oofem::IntArray::at(), oofem::SparseMtrx::at(), oofem::PetscSparseMtrx::createVecGlobal(), oofem::PetscSparseMtrx::giveMtrx(), oofem::SparseMtrx::giveVersion(), numberOfPrescribedDofs, prescribedEqs, and smConstraintVersion.
Referenced by solve().
| 
 | protected | 
Determines whether or not the solution has reached convergence.
Definition at line 580 of file nrsolver.C.
References oofem::__DofIDItemToString(), oofem::ParallelContext::accumulate(), oofem::IntArray::at(), oofem::FloatArray::at(), constrainedNRFlag, dg_forceScale, oofem::NumericalMethod::domain, oofem::Element_local, oofem::NumericalMethod::engngModel, forceErrVec, forceErrVecOld, oofem::Domain::giveBcs(), oofem::Domain::giveDofManagers(), oofem::Domain::giveElements(), oofem::Domain::giveMaxDofID(), oofem::Domain::giveNumber(), oofem::EngngModel::giveParallelContext(), oofem::EngngModel::giveProblemScale(), oofem::FloatArray::giveSize(), oofem::ParallelContext::isLocal(), oofem::ParallelContext::localNorm(), oofem::macroScale, nrsolver_ERROR_NORM_SMALL_NUM, NRSOLVER_MAX_REL_ERROR_BOUND, OOFEM_LOG_INFO, oofem::FloatArray::resize(), rtold, rtolf, oofem::FloatArray::zero(), and oofem::IntArray::zero().
Referenced by oofem::DynamicRelaxationSolver::solve(), and solve().
| 
 | inlinevirtual | 
Reimplemented from oofem::SparseNonLinearSystemNM.
Reimplemented in oofem::StaggeredSolver, and oofem::DynamicRelaxationSolver.
Definition at line 172 of file nrsolver.h.
| 
 | inlinevirtual | 
Reimplemented in oofem::StaggeredSolver, and oofem::DynamicRelaxationSolver.
Definition at line 173 of file nrsolver.h.
References _IFT_NRSolver_Name.
| 
 | virtual | 
Constructs (if necessary) and returns a linear solver.
Public method because some problems require it for sensitivity analysis, etc. even for nonlinear problems (e.g. tangent relations in multiscale simulations).
Reimplemented from oofem::SparseNonLinearSystemNM.
Definition at line 390 of file nrsolver.C.
References oofem::classFactory, oofem::ClassFactory::createSparseLinSolver(), oofem::NumericalMethod::domain, oofem::NumericalMethod::engngModel, linSolver, OOFEM_ERROR, and solverType.
Referenced by initializeFrom(), and solve().
| 
 | protected | 
Constructs and returns a line search solver.
Definition at line 410 of file nrsolver.C.
References oofem::NumericalMethod::domain, oofem::NumericalMethod::engngModel, and linesearchSolver.
Referenced by initializeFrom(), and solve().
| 
 | virtual | 
Reimplemented from oofem::NumericalMethod.
Reimplemented in oofem::StaggeredSolver, and oofem::DynamicRelaxationSolver.
Definition at line 103 of file nrsolver.C.
References _IFT_NRSolver_calcstiffbeforeres, _IFT_NRSolver_constrainedNRminiter, _IFT_NRSolver_ddfunc, _IFT_NRSolver_ddm, _IFT_NRSolver_ddv, _IFT_NRSolver_forceScale, _IFT_NRSolver_forceScaleDofs, _IFT_NRSolver_linesearch, _IFT_NRSolver_lstype, _IFT_NRSolver_manrmsteps, _IFT_NRSolver_maxinc, _IFT_NRSolver_maxiter, _IFT_NRSolver_miniterations, _IFT_NRSolver_minsteplength, _IFT_NRSolver_rtold, _IFT_NRSolver_rtolf, _IFT_NRSolver_rtolv, oofem::FloatArray::at(), oofem::IntArray::clear(), oofem::FloatArray::clear(), constrainedNRFlag, constrainedNRminiter, dg_forceScale, giveLinearSolver(), giveLineSearchSolver(), oofem::IntArray::giveSize(), oofem::FloatArray::giveSize(), oofem::InputRecord::hasField(), oofem::LineSearchNM::initializeFrom(), oofem::NumericalMethod::initializeFrom(), oofem::SparseNonLinearSystemNM::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, lsFlag, MANRMSteps, maxIncAllowed, mCalcStiffBeforeRes, minIterations, minStepLength, NR_Mode, NR_OldMode, nrsolverAccelNRM, nrsolverModifiedNRM, nsmax, numberOfPrescribedDofs, OOFEM_ERROR, prescribedDisplacementTF, prescribedDofs, prescribedDofsFlag, prescribedDofsValues, oofem::FloatArray::resize(), rtold, rtolf, and solverType.
Referenced by oofem::DynamicRelaxationSolver::initializeFrom(), and oofem::StaggeredSolver::initializeFrom().
| 
 | protected | 
Initiates prescribed equations.
Definition at line 420 of file nrsolver.C.
References oofem::IntArray::at(), oofem::NumericalMethod::domain, oofem::NumericalMethod::engngModel, oofem::DofManager::giveDofWithID(), oofem::Dof::giveEquationNumber(), oofem::DofManager::giveGlobalNumber(), oofem::Domain::giveNode(), oofem::Domain::giveNumber(), oofem::Domain::giveNumberOfDofManagers(), oofem::EngngModel::giveParallelContext(), oofem::ParallelContext::isLocal(), numberOfPrescribedDofs, prescribedDofs, prescribedEqs, prescribedEqsInitFlag, and oofem::IntArray::resize().
Referenced by solve().
| 
 | virtual | 
Prints status message of receiver to output stream.
Prints the message corresponding to last solve.
| outputStream | Stream to print state to. | 
Reimplemented from oofem::SparseNonLinearSystemNM.
Definition at line 558 of file nrsolver.C.
References oofem::IntArray::at(), oofem::FloatArray::at(), oofem::FloatArray::giveSize(), lastReactions, numberOfPrescribedDofs, and prescribedDofs.
| 
 | inlinevirtual | 
Reinitializes the receiver.
This is used, when topology of problem has changed (for example after adaptive refinement or load transfer in parallel applications). This is necessary for numerical methods, that cache some data between solution steps and that may depend on domain or problem topology. The caching of data by receiver is intended only for speeding up the calculation, but numerical method must be always able to generate this data again. This method clears receiver cached data dependent on topology, when it changes.
Reimplemented from oofem::NumericalMethod.
Definition at line 184 of file nrsolver.h.
| 
 | inlinevirtual | 
Reimplemented from oofem::NumericalMethod.
Definition at line 175 of file nrsolver.h.
| 
 | virtual | 
Solves the given sparse linear system of equations  .
. 
Total load vector is not passed, it is defined as  , where s is scale factor.
, where s is scale factor. 
| K | Coefficient matrix (  ; stiffness matrix). | 
| R | Reference incremental RHS (incremental load). | 
| R0 | Initial RHS (initial load). | 
| X | Total solution (total displacement). | 
| dX | Increment of solution (incremental displacements). | 
| F | InternalRhs (real internal forces). | 
| internalForcesEBENorm | Norm of internal nodal forces (evaluated on element by element basis) (split into each DOF id). | 
| s | RHS scale factor (load level). | 
| rlm | Reference load mode. | 
| nite | Number of iterations needed. | 
| tStep | Time step to solve for. | 
Implements oofem::SparseNonLinearSystemNM.
Reimplemented in oofem::StaggeredSolver, and oofem::DynamicRelaxationSolver.
Definition at line 202 of file nrsolver.C.
References oofem::FloatArray::add(), applyConstraintsToLoadIncrement(), applyConstraintsToStiffness(), oofem::IntArray::at(), oofem::FloatArray::at(), oofem::FloatArray::beDifferenceOf(), checkConvergence(), oofem::FloatArray::computeSquaredNorm(), constrainedNRalpha, constrainedNRFlag, constrainedNRminiter, oofem::SparseNonLinearSystemNM::deltaL, oofem::NumericalMethod::domain, oofem::ExportModuleManager::doOutput(), oofem::NumericalMethod::engngModel, forceErrVec, forceErrVecOld, oofem::EngngModel::giveExportModuleManager(), giveLinearSolver(), giveLineSearchSolver(), oofem::Domain::giveNumber(), oofem::EngngModel::giveParallelContext(), oofem::EngngModel::giveProblemScale(), oofem::FloatArray::giveSize(), oofem::TimeStep::incrementStateCounter(), oofem::TimeStep::incrementSubStepNumber(), initPrescribedEqs(), oofem::InternalRhs, lastReactions, linSolver, oofem::ParallelContext::localNorm(), lsFlag, oofem::macroScale, MANRMSteps, maxIncAllowed, mCalcStiffBeforeRes, minIterations, NM_None, NM_NoSuccess, NM_Success, oofem::NonLinearLhs, NR_Mode, nrsolverAccelNRM, nrsolverFullNRM, nsmax, numberOfPrescribedDofs, OOFEM_LOG_DEBUG, OOFEM_LOG_INFO, OOFEM_WARNING, prescribedDofs, prescribedDofsFlag, prescribedEqs, prescribedEqsInitFlag, oofem::FloatArray::resize(), rtold, rtolf, oofem::LineSearchNM::solve(), oofem::FloatArray::times(), oofem::EngngModel::updateComponent(), and oofem::FloatArray::zero().
| 
 | protected | 
Scale factor for dX, dX_new = alpha * dX.
Definition at line 143 of file nrsolver.h.
Referenced by solve().
| 
 | protected | 
Flag indicating whether to use constrained Newton.
Definition at line 141 of file nrsolver.h.
Referenced by checkConvergence(), initializeFrom(), and solve().
| 
 | protected | 
Minimum number of iterations before constraint is activated.
Definition at line 145 of file nrsolver.h.
Referenced by initializeFrom(), and solve().
| 
 | protected | 
Optional user supplied scale of forces used in convergence check.
Definition at line 157 of file nrsolver.h.
Referenced by checkConvergence(), and initializeFrom().
| 
 | protected | 
Definition at line 153 of file nrsolver.h.
Referenced by checkConvergence(), and solve().
| 
 | protected | 
Definition at line 154 of file nrsolver.h.
Referenced by checkConvergence(), and solve().
| 
 | protected | 
Computed reactions. They are stored in order to print them in printState method.
Definition at line 133 of file nrsolver.h.
Referenced by printState(), and solve().
| 
 | protected | 
Line search solver.
Definition at line 137 of file nrsolver.h.
Referenced by giveLineSearchSolver().
| 
 | protected | 
linear system solver
Definition at line 107 of file nrsolver.h.
Referenced by giveLinearSolver(), and solve().
| 
 | protected | 
Flag indicating whether to use line-search.
Definition at line 135 of file nrsolver.h.
Referenced by initializeFrom(), and solve().
| 
 | protected | 
Definition at line 104 of file nrsolver.h.
Referenced by initializeFrom(), and solve().
| 
 | protected | 
Definition at line 159 of file nrsolver.h.
Referenced by initializeFrom(), and solve().
| 
 | protected | 
Flag indicating if the stiffness should be evaluated before the residual in the first iteration.
Definition at line 139 of file nrsolver.h.
Referenced by initializeFrom(), and solve().
| 
 | protected | 
Definition at line 100 of file nrsolver.h.
Referenced by initializeFrom(), oofem::DynamicRelaxationSolver::solve(), and solve().
| 
 | protected | 
Definition at line 101 of file nrsolver.h.
Referenced by initializeFrom().
| 
 | protected | 
Definition at line 102 of file nrsolver.h.
Referenced by initializeFrom(), and solve().
| 
 | protected | 
Definition at line 103 of file nrsolver.h.
| 
 | protected | 
Definition at line 102 of file nrsolver.h.
Referenced by initializeFrom().
| 
 | protected | 
Definition at line 100 of file nrsolver.h.
Referenced by initializeFrom(), oofem::DynamicRelaxationSolver::solve(), and solve().
| 
 | protected | 
number of prescribed displacement
Definition at line 113 of file nrsolver.h.
Referenced by applyConstraintsToLoadIncrement(), applyConstraintsToStiffness(), initializeFrom(), initPrescribedEqs(), printState(), and solve().
| 
 | protected | 
Load Time Function of prescribed values.
Definition at line 127 of file nrsolver.h.
Referenced by applyConstraintsToLoadIncrement(), and initializeFrom().
| 
 | protected | 
Array of pairs identifying prescribed dofs (node, dof)
Definition at line 123 of file nrsolver.h.
Referenced by initializeFrom(), initPrescribedEqs(), printState(), and solve().
| 
 | protected | 
Flag indicating that some dofs are controlled under displacement control.
In parallel mode, numberOfPrescribedDofs is local (related to specific partition) so its nonzero value does not mean that there are no prescribed dofs on other partitions.
Definition at line 120 of file nrsolver.h.
Referenced by initializeFrom(), and solve().
| 
 | protected | 
Array of prescribed values.
Definition at line 125 of file nrsolver.h.
Referenced by applyConstraintsToLoadIncrement(), and initializeFrom().
| 
 | protected | 
Array of prescribed equations.
Definition at line 129 of file nrsolver.h.
Referenced by applyConstraintsToLoadIncrement(), applyConstraintsToStiffness(), initPrescribedEqs(), and solve().
| 
 | protected | 
Flag indicating that prescribedEqs were initialized.
Definition at line 131 of file nrsolver.h.
Referenced by initPrescribedEqs(), and solve().
| 
 | protected | 
Relative iterative displacement change tolerance for each group.
Definition at line 150 of file nrsolver.h.
Referenced by checkConvergence(), initializeFrom(), oofem::DynamicRelaxationSolver::solve(), and solve().
| 
 | protected | 
Relative unbalanced force tolerance for each group.
Definition at line 148 of file nrsolver.h.
Referenced by checkConvergence(), initializeFrom(), oofem::DynamicRelaxationSolver::solve(), and solve().
| 
 | protected | 
sparse matrix version, used to control constrains application to stiffness
Definition at line 111 of file nrsolver.h.
Referenced by applyConstraintsToStiffness().
| 
 | protected | 
linear system solver ID
Definition at line 109 of file nrsolver.h.
Referenced by applyConstraintsToLoadIncrement(), giveLinearSolver(), and initializeFrom().