WarpX
Loading...
Searching...
No Matches
JacobianFunctionMF.H
Go to the documentation of this file.
1/* Copyright 2024 Justin Angus
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7#ifndef JacobianFunctionMF_H_
8#define JacobianFunctionMF_H_
9
10#include "LinearFunction.H"
11#include "CurlCurlMLMGPC.H"
12#include "JacobiPC.H"
13#include "Utils/TextMsg.H"
14#include <AMReX_Config.H>
15#include <string>
16
23
24template <class T, class Ops>
26{
27 public:
28
29 using RT = typename T::value_type;
30
31 JacobianFunctionMF() = default;
33
34 // Default move and copy operations
38 JacobianFunctionMF& operator=(JacobianFunctionMF&&) noexcept = default;
39
40 void apply ( T& a_dF, const T& a_dU ) override;
41
42 inline
43 void precond ( T& a_U, const T& a_X ) override
44 {
45 if (m_usePreCond) {
46 a_U.zero();
47 m_preCond->Apply(a_U, a_X);
48 } else {
49 a_U.Copy(a_X);
50 }
51 }
52
53 inline
54 void updatePreCondMat ( const T& a_X ) override
55 {
56 if (m_usePreCond) { m_preCond->Update(a_X); }
57 }
58
59 T makeVecLHS () const override;
60 T makeVecRHS () const override;
61
62 inline
63 bool isDefined() const { return m_is_defined; }
64
65 [[nodiscard]] inline
66 bool usePreconditioner() const { return m_usePreCond; }
67
68 inline
69 void setBaseSolution ( const T& a_U )
70 {
71 m_Y0.Copy(a_U);
72 m_normY0 = this->norm2(m_Y0);
73 }
74
75 inline
76 void setBaseRHS ( const T& a_R )
77 {
78 m_R0.Copy(a_R);
79 }
80
81 inline
82 void setJFNKEps ( RT a_eps )
83 {
84 m_epsJFNK = a_eps;
85 }
86
87 inline
88 void setIsLinear ( bool a_isLinear )
89 {
90 m_is_linear = a_isLinear;
91 }
92
93 inline
94 void curTime ( RT a_time )
95 {
96 m_cur_time = a_time;
97 if (m_usePreCond) { m_preCond->CurTime(a_time); }
98 }
99
100 inline
101 void curTimeStep ( RT a_dt )
102 {
103 m_dt = a_dt;
104 if (m_usePreCond) { m_preCond->CurTimeStep(a_dt); }
105 }
106
107 inline
108 void printParams () const
109 {
111 m_preCond->printParameters();
112 }
113 }
114
115 void define( const T&, Ops*, const PreconditionerType& ) override;
116
117 inline
118 PreconditionerType pcType () const override { return m_pc_type; }
119
120 private:
121
122 bool m_is_defined = false;
123 bool m_is_linear = false;
124 bool m_usePreCond = false;
125 RT m_epsJFNK = RT(1.0e-6);
128
130
132 Ops* m_ops = nullptr;
133 std::unique_ptr<Preconditioner<T,Ops>> m_preCond = nullptr;
134};
135
136template <class T, class Ops>
138 Ops* a_ops,
139 const PreconditionerType& a_pc_type )
140{
141 BL_PROFILE("JacobianFunctionMF::::define()");
142 m_Z.Define(a_U);
143 m_Y0.Define(a_U);
144 m_R0.Define(a_U);
145 m_R.Define(a_U);
146
147 m_ops = a_ops;
148
149 m_usePreCond = (a_pc_type != PreconditionerType::none);
150 if (m_usePreCond) {
151 m_pc_type = a_pc_type;
153 m_preCond = std::make_unique<CurlCurlMLMGPC<T,Ops>>();
155 m_preCond = std::make_unique<JacobiPC<T,Ops>>();
157#ifdef AMREX_USE_PETSC
158 WARPX_ABORT_WITH_MESSAGE("JacobianFunctionMF::Define(): pc_petsc not yet implemented");
159#else
160 WARPX_ABORT_WITH_MESSAGE("JacobianFunctionMF::Define(): must compile with PETSc to use pc_petsc (AMREX_USE_PETSC must be defined)");
161#endif
162 } else {
163 std::stringstream convergenceMsg;
164 convergenceMsg << "JacobianFunctionMF::define(): " << amrex::getEnumNameString(m_pc_type)
165 << " is not a valid preconditioner type.";
166 WARPX_ABORT_WITH_MESSAGE(convergenceMsg.str());
167 }
168 m_preCond->Define(a_U, a_ops);
169 }
170
171 m_is_defined = true;
172}
173
174template <class T, class Ops>
176{
177 BL_PROFILE("JacobianFunctionMF::::makeVecRHS()");
178 T x;
179 x.Define(m_R);
180 return x;
181}
182
183template <class T, class Ops>
185{
186 BL_PROFILE("JacobianFunctionMF::::makeVecLHS()");
187 T x;
188 x.Define(m_R);
189 return x;
190}
191
192template <class T, class Ops>
193void JacobianFunctionMF<T,Ops>::apply (T& a_dF, const T& a_dU)
194{
195 BL_PROFILE("JacobianFunctionMF::apply()");
196 using namespace amrex::literals;
197 RT const normY = this->norm2(a_dU); // always 1 when called from GMRES
198
200 isDefined(),
201 "JacobianFunction::apply() called on undefined JacobianFunction");
202
203 if (normY < 1.0e-15) { a_dF.zero(); }
204 else {
205
206 RT eps;
207 if (m_is_linear) {
208 eps = 1.0_rt;
209 } else {
210 /* eps = error_rel * sqrt(1 + ||Y0||) / ||dU||
211 * M. Pernice and H. F. Walker, "NITSOL: A Newton Iterative Solver for
212 * Nonlinear Systems", SIAM J. Sci. Stat. Comput., 1998, vol 19,
213 * pp. 302--318. */
214 if (m_normY0==0.0) { eps = m_epsJFNK * this->norm2(m_R0) / normY; }
215 else {
216 // m_epsJFNK * sqrt(1.0 + m_normY0) / normY
217 // above commonly used form not recommend for poorly scaled Y0
218 eps = m_epsJFNK * m_normY0 / normY;
219 }
220 }
221 const RT eps_inv = 1.0_rt/eps;
222
223 m_Z.linComb( 1.0, m_Y0, eps, a_dU ); // Z = Y0 + eps*dU
224 m_ops->ComputeRHS(m_R, m_Z, m_cur_time, -1, true );
225
226 // F(Y) = Y - b - R(Y) ==> dF = dF/dY*dU = [1 - dR/dY]*dU
227 // = dU - (R(Z)-R(Y0))/eps
228 a_dF.linComb( 1.0, a_dU, eps_inv, m_R0 );
229 a_dF.increment(m_R,-eps_inv);
230
231 }
232
233}
234
235#endif
#define BL_PROFILE(a)
PreconditionerType
Types for preconditioners for field solvers.
Definition Preconditioner.H:20
@ none
Definition Preconditioner.H:20
@ pc_jacobi
Definition Preconditioner.H:20
@ pc_curl_curl_mlmg
Definition Preconditioner.H:20
@ pc_petsc
Definition Preconditioner.H:20
#define WARPX_ABORT_WITH_MESSAGE(MSG)
Definition TextMsg.H:15
#define WARPX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition TextMsg.H:13
std::unique_ptr< Preconditioner< T, Ops > > m_preCond
Definition JacobianFunctionMF.H:133
RT m_cur_time
Definition JacobianFunctionMF.H:127
bool m_is_defined
Definition JacobianFunctionMF.H:122
void updatePreCondMat(const T &a_X) override
update preconditioner
Definition JacobianFunctionMF.H:54
T m_R0
Definition JacobianFunctionMF.H:131
void setBaseRHS(const T &a_R)
Definition JacobianFunctionMF.H:76
void setBaseSolution(const T &a_U)
Definition JacobianFunctionMF.H:69
bool isDefined() const
Definition JacobianFunctionMF.H:63
~JacobianFunctionMF()=default
bool usePreconditioner() const
Definition JacobianFunctionMF.H:66
void setIsLinear(bool a_isLinear)
Definition JacobianFunctionMF.H:88
T m_Z
Definition JacobianFunctionMF.H:131
void precond(T &a_U, const T &a_X) override
apply the preconditioner on a given vector of type T
Definition JacobianFunctionMF.H:43
void define(const T &, Ops *, const PreconditionerType &) override
define the linear function object
Definition JacobianFunctionMF.H:137
JacobianFunctionMF(const JacobianFunctionMF &)=default
RT m_normY0
Definition JacobianFunctionMF.H:126
JacobianFunctionMF()=default
void curTimeStep(RT a_dt)
Definition JacobianFunctionMF.H:101
RT m_dt
Definition JacobianFunctionMF.H:127
void curTime(RT a_time)
Definition JacobianFunctionMF.H:94
T makeVecLHS() const override
make LHS vector
Definition JacobianFunctionMF.H:184
void apply(Vec &a_dF, const Vec &a_dU) override
Definition JacobianFunctionMF.H:193
PreconditionerType pcType() const override
returns the preconditioner type
Definition JacobianFunctionMF.H:118
bool m_is_linear
Definition JacobianFunctionMF.H:123
JacobianFunctionMF & operator=(const JacobianFunctionMF &)=default
T m_Y0
Definition JacobianFunctionMF.H:131
void printParams() const
Definition JacobianFunctionMF.H:108
void setJFNKEps(RT a_eps)
Definition JacobianFunctionMF.H:82
RT m_epsJFNK
Definition JacobianFunctionMF.H:125
typename T::value_type RT
Definition JacobianFunctionMF.H:29
Ops * m_ops
Definition JacobianFunctionMF.H:132
T makeVecRHS() const override
make RHS vector
Definition JacobianFunctionMF.H:175
bool m_usePreCond
Definition JacobianFunctionMF.H:124
JacobianFunctionMF(JacobianFunctionMF &&) noexcept=default
T m_R
Definition JacobianFunctionMF.H:131
PreconditionerType m_pc_type
Definition JacobianFunctionMF.H:129
RT norm2(const T &a_U)
compute 2-norm of a vector
Definition LinearFunction.H:100
LinearFunction()=default
default constructor
std::string getEnumNameString(T const &v)