Nektar++
NekNonlinSysNewton.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: NekNonlinSysNewton.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: NekNonlinSysNewton definition
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
38using namespace std;
39
41{
42/**
43 * @class NekNonlinSysNewton
44 *
45 * Solves a nonlinear system using iterative methods.
46 */
49 "Newton", NekNonlinSysNewton::create, "NekNonlinSysNewton solver.");
50
53 const LibUtilities::CommSharedPtr &vRowComm, const int nscale,
54 const NekSysKey &pKey)
55 : NekNonlinSys(pSession, vRowComm, nscale, pKey)
56{
57}
58
60{
64}
65
67{
69 m_linsol->SetSysOperators(in);
70}
71
72/**
73 *
74 **/
76 const int nGlobal, const Array<OneD, const NekDouble> &pInput,
77 Array<OneD, NekDouble> &pOutput, const int nDir, const NekDouble tol,
78 [[maybe_unused]] const NekDouble factor)
79{
80 ASSERTL0(nDir == 0, "nDir != 0 not tested");
81 ASSERTL0(nGlobal == m_SysDimen, "nGlobal != m_SysDimen");
82
83 int ntotal = nGlobal - nDir;
84
85 m_SourceVec = pInput;
86 m_Solution = pOutput;
88 m_converged = false;
89
90 Vmath::Vcopy(ntotal, pInput, 1, m_Solution, 1);
91
92 NekDouble resnormOld = 0.0;
93 int NttlNonlinIte = 0;
94 for (; NttlNonlinIte < m_maxiter; ++NttlNonlinIte)
95 {
97
98 m_converged = v_ConvergenceCheck(NttlNonlinIte, m_Residual, tol);
99 if (m_converged)
100 {
101 break;
102 }
103
104 NekDouble LinSysRelativeIteTol =
105 CalcInexactNewtonForcing(NttlNonlinIte, resnormOld, m_SysResNorm);
106 resnormOld = m_SysResNorm;
107 m_linsol->setRhsMagnitude(m_SysResNorm);
108 int ntmpLinSysIts = m_linsol->SolveSystem(
109 ntotal, m_Residual, m_DeltSltn, 0, LinSysRelativeIteTol);
110 m_NtotLinSysIts += ntmpLinSysIts;
111
112 Vmath::Vsub(ntotal, m_Solution, 1, m_DeltSltn, 1, m_Solution, 1);
113 }
114
115 if (((!m_converged) || m_verbose) && m_root && m_FlagWarnings)
116 {
117 int nwidthcolm = 11;
118
120 " # Nonlinear solver not converge in DoImplicitSolve");
121 cout << right << scientific << setw(nwidthcolm)
122 << setprecision(nwidthcolm - 6)
123 << " * Newton-Its converged (RES=" << sqrt(m_SysResNorm)
124 << " Res/(DtRHS): " << sqrt(m_SysResNorm / m_SysResNorm0)
125 << " with " << setw(3) << NttlNonlinIte << " Non-Its)" << endl;
126 }
127
128 return NttlNonlinIte;
129}
130
132 const int nIteration, const Array<OneD, const NekDouble> &Residual,
133 const NekDouble tol)
134{
135 m_SysResNorm = Vmath::Dot(Residual.size(), Residual, Residual);
137
138 if (nIteration == 0)
139 {
141 }
142
145
146 return resratio < restol * restol || m_SysResNorm < tol * tol;
147}
148
150 const int &k, const NekDouble &resnormOld, const NekDouble &resnorm)
151{
152 if (k == 0 || !m_InexactNewtonForcing)
153 {
155 }
156 else
157 {
158 NekDouble tmpForc =
159 m_forcingGamma * pow((resnorm / resnormOld), m_forcingAlpha);
160 return max(min(m_LinSysRelativeTolInNonlin, tmpForc), 1.0E-6);
161 }
162}
163
164} // namespace Nektar::LibUtilities
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:215
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
Array< OneD, NekDouble > m_Residual
Definition: NekNonlinSys.h:124
Array< OneD, NekDouble > m_DeltSltn
Definition: NekNonlinSys.h:125
Array< OneD, NekDouble > m_SourceVec
Definition: NekNonlinSys.h:126
NekLinSysIterSharedPtr m_linsol
Definition: NekNonlinSys.h:113
Array< OneD, NekDouble > m_Solution
Definition: NekNonlinSys.h:123
static NekNonlinSysSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
NekDouble CalcInexactNewtonForcing(const int &k, const NekDouble &resnormOld, const NekDouble &resnorm)
NekNonlinSysNewton(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nscale, const NekSysKey &pKey)
int v_SolveSystem(const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir, const NekDouble tol, const NekDouble factor) override
bool v_ConvergenceCheck(const int nIteration, const Array< OneD, const NekDouble > &Residual, const NekDouble tol) override
void v_SetSysOperators(const NekSysOperators &in) override
bool m_root
Root if parallel.
Definition: NekSys.h:297
LibUtilities::CommSharedPtr m_rowComm
Communicate.
Definition: NekSys.h:293
virtual void v_InitObject()
Definition: NekSys.h:306
bool m_verbose
Verbose.
Definition: NekSys.h:299
NekSysOperators m_operator
Operators.
Definition: NekSys.h:302
int m_SysDimen
The dimension of the system.
Definition: NekSys.h:304
virtual void v_SetSysOperators(const NekSysOperators &in)
Definition: NekSys.h:310
bool m_converged
Whether the iteration has been converged.
Definition: NekSys.h:295
int m_maxiter
Maximum iterations.
Definition: NekSys.h:289
void DoNekSysResEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:134
NekNonlinSysFactory & GetNekNonlinSysFactory()
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
double NekDouble
T Dot(int n, const T *w, const T *x)
dot product
Definition: Vmath.hpp:761
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294