7 #ifndef WARPX_SPECTRAL_SOLVER_H_
8 #define WARPX_SPECTRAL_SOLVER_H_
24 #ifdef WARPX_USE_PSATD
75 const int norder_x,
const int norder_y,
76 const int norder_z,
const short grid_type,
82 const bool periodic_single_box,
83 const bool update_with_rho,
84 const bool fft_do_time_averaging,
85 const int psatd_solution_type,
87 const int rho_in_time,
88 const bool dive_cleaning,
89 const bool divb_cleaning);
102 const int field_index,
103 const int i_comp = 0);
111 const int field_index,
113 const int i_comp=0 );
125 const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
Class that stores the fields in spectral space, and performs the Fourier transforms between real spac...
Definition: SpectralFieldData.H:122
SpectralField fields
Definition: SpectralFieldData.H:143
Definition: SpectralFieldData.H:33
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:33
void ComputeSpectralDivE(const int lev, const std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, amrex::MultiFab &divE)
Public interface to call the member function ComputeSpectralDivE of the base class SpectralBaseAlgori...
Definition: SpectralSolver.H:124
SpectralFieldIndex m_spectral_index
Definition: SpectralSolver.H:190
void VayDeposition()
Public interface to call the virtual function VayDeposition, declared in the base class SpectralBaseA...
Definition: SpectralSolver.H:147
void ZeroOutDataComp(const int icomp)
Set to zero the data on component icomp of field_data.fields.
Definition: SpectralSolver.H:171
void ForwardTransform(const int lev, const amrex::MultiFab &mf, const int field_index, const int i_comp=0)
Transform the component i_comp of the MultiFab mf to Fourier space, and store the result internally (...
Definition: SpectralSolver.cpp:118
void ScaleDataComp(const int icomp, const amrex::Real scale_factor)
Scale the data on component icomp of field_data.fields by a given scale factor.
Definition: SpectralSolver.H:184
SpectralSolver(const int lev, const amrex::BoxArray &realspace_ba, const amrex::DistributionMapping &dm, const int norder_x, const int norder_y, const int norder_z, const short grid_type, const amrex::Vector< amrex::Real > &v_galilean, const amrex::Vector< amrex::Real > &v_comoving, const amrex::RealVect dx, const amrex::Real dt, const bool pml, const bool periodic_single_box, const bool update_with_rho, const bool fft_do_time_averaging, const int psatd_solution_type, const int J_in_time, const int rho_in_time, const bool dive_cleaning, const bool divb_cleaning)
Constructor of the class SpectralSolver.
Definition: SpectralSolver.cpp:23
void pushSpectralFields()
Update the fields in spectral space, over one timestep.
Definition: SpectralSolver.cpp:139
void CurrentCorrection()
Public interface to call the virtual function CurrentCorrection, defined in the base class SpectralBa...
Definition: SpectralSolver.H:136
SpectralFieldData field_data
Definition: SpectralSolver.H:201
amrex::IntVect m_fill_guards
Definition: SpectralSolver.H:194
void CopySpectralDataComp(const int src_comp, const int dest_comp)
Copy spectral data from component src_comp to component dest_comp of field_data.fields.
Definition: SpectralSolver.H:159
std::unique_ptr< SpectralBaseAlgorithm > algorithm
Definition: SpectralSolver.H:206
void BackwardTransform(const int lev, amrex::MultiFab &mf, const int field_index, const amrex::IntVect &fill_guards, const int i_comp=0)
Transform spectral field specified by field_index back to real space, and store it in the component i...
Definition: SpectralSolver.cpp:128
void setVal(value_type val)
void mult(value_type val, int comp, int num_comp, int nghost=0)
int dx
Definition: stencil.py:436
int dt
Definition: stencil.py:440