WarpX
Loading...
Searching...
No Matches
SemiImplicitEM.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 SEMI_IMPLICIT_EM_H_
8#define SEMI_IMPLICIT_EM_H_
9
11
12#include <AMReX_Array.H>
13#include <AMReX_MultiFab.H>
14#include <AMReX_REAL.H>
15
16#include "ImplicitSolver.H"
17
41
43{
44public:
45
46 SemiImplicitEM() = default;
47
48 ~SemiImplicitEM() override = default;
49
50 // Prohibit Move and Copy operations
55
56 void Define ( WarpX* a_WarpX ) override;
57
58 void PrintParameters () const override;
59
60 void OneStep ( amrex::Real start_time,
61 amrex::Real a_dt,
62 int a_step ) override;
63
64 void ComputeRHS ( WarpXSolverVec& a_RHS,
65 const WarpXSolverVec& a_E,
66 amrex::Real start_time,
67 int a_nl_iter,
68 bool a_from_jacobian ) override;
69
70 // This parameter is used for the time-step fraction in the PC for implicit
71 // treatment of light waves in the curl-curl MLMG solver.
72 // This function should return zero if light waves are not treated implicitly
73 [[nodiscard]] amrex::Real GetThetaForPC () const override { return 0.; }
74
75private:
76
81
82};
83
84#endif
ImplicitSolver()=default
SemiImplicitEM()=default
void Define(WarpX *a_WarpX) override
Read user-provided parameters that control the implicit solver. Allocate internal arrays for intermed...
Definition SemiImplicitEM.cpp:15
WarpXSolverVec m_Eold
Definition SemiImplicitEM.H:80
SemiImplicitEM & operator=(SemiImplicitEM &&)=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 SemiImplicitEM.cpp:56
amrex::Real GetThetaForPC() const override
Definition SemiImplicitEM.H:73
SemiImplicitEM(SemiImplicitEM &&)=delete
~SemiImplicitEM() override=default
SemiImplicitEM(const SemiImplicitEM &)=delete
SemiImplicitEM & operator=(const SemiImplicitEM &)=delete
WarpXSolverVec m_E
Solver vectors for E and Eold.
Definition SemiImplicitEM.H:80
void PrintParameters() const override
Definition SemiImplicitEM.cpp:44
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 SemiImplicitEM.cpp:106
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