WarpX
Loading...
Searching...
No Matches
StrangImplicitSpectralEM.H
Go to the documentation of this file.
1/* Copyright 2024 David Grote
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7#ifndef STRANG_IMPLICIT_SPECTRALEM_H_
8#define STRANG_IMPLICIT_SPECTRALEM_H_
9
11
12#include <AMReX_Array.H>
13#include <AMReX_MultiFab.H>
14#include <AMReX_REAL.H>
15
16#include "ImplicitSolver.H"
17
42
44{
45public:
46
48
49 ~StrangImplicitSpectralEM() override = default;
50
51 // Prohibit Move and Copy operations
56
57 void Define ( WarpX* a_WarpX ) override;
58
59 void PrintParameters () const override;
60
61 void OneStep ( amrex::Real start_time,
62 amrex::Real a_dt,
63 int a_step ) override;
64
65 void ComputeRHS ( WarpXSolverVec& a_RHS,
66 const WarpXSolverVec& a_E,
67 amrex::Real start_time,
68 int a_nl_iter,
69 bool a_from_jacobian ) override;
70
71 // This parameter is used for the time-step fraction in the PC for implicit
72 // treatment of light waves in the curl-curl MLMG solver.
73 // This function should return zero if light waves are not treated implicitly
74 [[nodiscard]] amrex::Real GetThetaForPC () const override { return 0.; }
75
76private:
77
87
95
99 void UpdateWarpXFields ( WarpXSolverVec const& a_E,
100 amrex::Real half_time );
101
106 void FinishFieldUpdate ( amrex::Real a_new_time );
107
108};
109
110#endif
ImplicitSolver()=default
amrex::Real GetThetaForPC() const override
Definition StrangImplicitSpectralEM.H:74
StrangImplicitSpectralEM()=default
void ComputeRHS(WarpXSolverVec &a_RHS, const WarpXSolverVec &a_E, amrex::Real start_time, int a_nl_iter, bool a_from_jacobian) override
Computes the RHS of the equation corresponding to the specified implicit algorithm....
Definition StrangImplicitSpectralEM.cpp:99
void PrintParameters() const override
Definition StrangImplicitSpectralEM.cpp:43
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > m_Bold
B is a derived variable from E. Need to save Bold to update B during the iterative nonlinear solve fo...
Definition StrangImplicitSpectralEM.H:94
void FinishFieldUpdate(amrex::Real a_new_time)
Nonlinear solver is for the time-centered values of E. After the solver, need to use m_E and m_Eold t...
Definition StrangImplicitSpectralEM.cpp:131
StrangImplicitSpectralEM(const StrangImplicitSpectralEM &)=delete
WarpXSolverVec m_E
Solver vectors to be used in the nonlinear solver to solve for the electric field E....
Definition StrangImplicitSpectralEM.H:86
StrangImplicitSpectralEM & operator=(StrangImplicitSpectralEM &&)=delete
void UpdateWarpXFields(WarpXSolverVec const &a_E, amrex::Real half_time)
Update the E and B fields owned by WarpX.
Definition StrangImplicitSpectralEM.cpp:122
~StrangImplicitSpectralEM() override=default
StrangImplicitSpectralEM(StrangImplicitSpectralEM &&)=delete
void OneStep(amrex::Real start_time, amrex::Real a_dt, int a_step) override
Advance fields and particles by one time step using the specified implicit algorithm.
Definition StrangImplicitSpectralEM.cpp:55
void Define(WarpX *a_WarpX) override
Read user-provided parameters that control the implicit solver. Allocate internal arrays for intermed...
Definition StrangImplicitSpectralEM.cpp:15
WarpXSolverVec m_Eold
Definition StrangImplicitSpectralEM.H:86
StrangImplicitSpectralEM & operator=(const StrangImplicitSpectralEM &)=delete
Definition WarpX.H:85
This is a wrapper class around a Vector of pointers to MultiFabs that contains basic math operators a...
Definition WarpXSolverVec.H:58