WarpX
Loading...
Searching...
No Matches
SpectralSolverRZ.H
Go to the documentation of this file.
1/* Copyright 2019-2020 David Grote
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7#ifndef WARPX_SPECTRAL_SOLVER_RZ_H_
8#define WARPX_SPECTRAL_SOLVER_RZ_H_
9
11
13#include "SpectralFieldDataRZ.H"
14
16#include <ablastr/utils/Enums.H>
17
18
19/* \brief Top-level class for the electromagnetic spectral solver
20 *
21 * Stores the field in spectral space, and has member functions
22 * to Fourier-transform the fields between real space and spectral space
23 * and to update fields in spectral space over one time step.
24 */
26{
27 public:
28 // Inline definition of the member functions of `SpectralSolverRZ`,
29 // except the constructor (see `SpectralSolverRZ.cpp`)
30 // The body of these functions is short, since the work is done in the
31 // underlying classes `SpectralFieldData` and `PsatdAlgorithm`
32
33 // Constructor
34 SpectralSolverRZ (int lev,
35 amrex::BoxArray const & realspace_ba,
37 int n_rz_azimuthal_modes,
38 int norder_z,
40 const amrex::Vector<amrex::Real>& v_galilean,
41 amrex::RealVect dx, amrex::Real dt,
42 bool with_pml,
43 bool update_with_rho,
44 bool fft_do_time_averaging,
45 TimeDependencyJ time_dependency_J,
46 TimeDependencyRho time_dependency_rho,
47 bool dive_cleaning,
48 bool divb_cleaning);
49
50 /* \brief Transform the component `i_comp` of MultiFab `field_mf`
51 * to spectral space, and store the corresponding result internally
52 * (in the spectral field specified by `field_index`) */
53 void ForwardTransform (int lev, amrex::MultiFab const & field_mf, int field_index,
54 int i_comp=0);
55
56 /* \brief Transform the two MultiFabs `field_mf1` and `field_mf2`
57 * to spectral space, and store the corresponding results internally
58 * (in the spectral field specified by `field_index1` and `field_index2`) */
59 void ForwardTransform (int lev, amrex::MultiFab const & field_mf1, int field_index1,
60 amrex::MultiFab const & field_mf2, int field_index2);
61
62 /* \brief Transform spectral field specified by `field_index` back to
63 * real space, and store it in the component `i_comp` of `field_mf` */
64 void BackwardTransform (int lev, amrex::MultiFab& field_mf, int field_index,
65 int i_comp=0);
66
67 /* \brief Transform spectral fields specified by `field_index1` and `field_index2`
68 * back to real space, and store it in `field_mf1` and `field_mf2`*/
69 void BackwardTransform (int lev, amrex::MultiFab& field_mf1, int field_index1,
70 amrex::MultiFab& field_mf2, int field_index2);
71
72 /* \brief Update the fields in spectral space, over one timestep */
73 void pushSpectralFields (bool doing_pml=false);
74
75 /* \brief Initialize K space filtering arrays */
76 void InitFilter (amrex::IntVect const & filter_npass_each_dir,
77 bool compensation);
78
79 /* \brief Apply K space filtering for a scalar */
80 void ApplyFilter (int lev, int field_index);
81
82 /* \brief Apply K space filtering for a vector */
83 void ApplyFilter (int lev, int field_index1,
84 int field_index2, int field_index3);
85
90 void ComputeSpectralDivE (int lev,
91 ablastr::fields::VectorField const & Efield,
92 amrex::MultiFab& divE);
93
100 void CurrentCorrection ();
101
108 void VayDeposition ();
109
117 void CopySpectralDataComp (int src_comp, int dest_comp);
118
124 void ZeroOutDataComp (int icomp);
125
133 void ScaleDataComp (int icomp, amrex::Real scale_factor);
134
136
137 // Solve time step size
138 amrex::Real m_dt;
139
140 private:
141
142 SpectralKSpaceRZ k_space; // Save the instance to initialize filtering
143 SpectralFieldDataRZ field_data; // Store field in spectral space
144 // and perform the Fourier transforms
145 std::unique_ptr<SpectralBaseAlgorithmRZ> algorithm;
146 std::unique_ptr<SpectralBaseAlgorithmRZ> PML_algorithm;
147 // Defines field update equation in spectral space,
148 // and the associated coefficients.
149 // SpectralBaseAlgorithmRZ is a base class ; this pointer is meant
150 // to point an instance of a *sub-class* defining a specific algorithm
151
152};
153
154#endif // WARPX_SPECTRAL_SOLVER_RZ_H_
TimeDependencyJ
Definition WarpXAlgorithmSelection.H:106
TimeDependencyRho
Definition WarpXAlgorithmSelection.H:112
Definition SpectralFieldDataRZ.H:23
Definition SpectralFieldData.H:35
Definition SpectralKSpaceRZ.H:21
SpectralFieldIndex m_spectral_index
Definition SpectralSolverRZ.H:135
void ApplyFilter(int lev, int field_index)
Definition SpectralSolverRZ.cpp:146
void CopySpectralDataComp(int src_comp, int dest_comp)
Copy spectral data from component src_comp to component dest_comp of field_data.fields.
Definition SpectralSolverRZ.cpp:191
SpectralKSpaceRZ k_space
Definition SpectralSolverRZ.H:142
void ScaleDataComp(int icomp, amrex::Real scale_factor)
Scale the data on component icomp of field_data.fields by a given scale factor.
Definition SpectralSolverRZ.cpp:203
SpectralFieldDataRZ field_data
Definition SpectralSolverRZ.H:143
void VayDeposition()
Public interface to call the virtual function VayDeposition, declared in the base class SpectralBaseA...
Definition SpectralSolverRZ.cpp:184
void BackwardTransform(int lev, amrex::MultiFab &field_mf, int field_index, int i_comp=0)
Definition SpectralSolverRZ.cpp:106
std::unique_ptr< SpectralBaseAlgorithmRZ > algorithm
Definition SpectralSolverRZ.H:145
SpectralSolverRZ(int lev, amrex::BoxArray const &realspace_ba, amrex::DistributionMapping const &dm, int n_rz_azimuthal_modes, int norder_z, ablastr::utils::enums::GridType grid_type, const amrex::Vector< amrex::Real > &v_galilean, amrex::RealVect dx, amrex::Real dt, bool with_pml, bool update_with_rho, bool fft_do_time_averaging, TimeDependencyJ time_dependency_J, TimeDependencyRho time_dependency_rho, bool dive_cleaning, bool divb_cleaning)
Definition SpectralSolverRZ.cpp:28
void pushSpectralFields(bool doing_pml=false)
Definition SpectralSolverRZ.cpp:127
void ForwardTransform(int lev, amrex::MultiFab const &field_mf, int field_index, int i_comp=0)
Definition SpectralSolverRZ.cpp:83
void ComputeSpectralDivE(int lev, ablastr::fields::VectorField const &Efield, amrex::MultiFab &divE)
Public interface to call the member function ComputeSpectralDivE of the base class SpectralBaseAlgori...
Definition SpectralSolverRZ.cpp:163
void ZeroOutDataComp(int icomp)
Set to zero the data on component icomp of field_data.fields.
Definition SpectralSolverRZ.cpp:197
std::unique_ptr< SpectralBaseAlgorithmRZ > PML_algorithm
Definition SpectralSolverRZ.H:146
void CurrentCorrection()
Public interface to call the virtual function CurrentCorrection, defined in the base class SpectralBa...
Definition SpectralSolverRZ.cpp:178
void InitFilter(amrex::IntVect const &filter_npass_each_dir, bool compensation)
Definition SpectralSolverRZ.cpp:139
amrex::Real m_dt
Definition SpectralSolverRZ.H:138
std::array< amrex::MultiFab *, 3 > VectorField
Definition MultiFabRegister.H:191
GridType
Definition Enums.H:23
IntVectND< 3 > IntVect
RealVectND< 3 > RealVect