#ifndef WARPX_MACROSCOPICPROPERTIES_H_ #define WARPX_MACROSCOPICPROPERTIES_H_ #include "MacroscopicProperties_fwd.H" #include "Utils/WarpXConst.H" #include #include #include #include #include #include #include #include /** * \brief This class contains the macroscopic properties of the medium needed to * evaluate macroscopic Maxwell equation. */ class MacroscopicProperties { public: MacroscopicProperties (); // constructor /** Read user-defined macroscopic properties. Called in constructor. */ void ReadParameters (); /** Initialize multifabs storing macroscopic multifabs */ void InitData (); /** return MultiFab, sigma (conductivity) of the medium. */ amrex::MultiFab& getsigma_mf () {return (*m_sigma_mf);} /** return MultiFab, epsilon (permittivity) of the medium. */ amrex::MultiFab& getepsilon_mf () {return (*m_eps_mf);} /** return MultiFab, mu (permeability) of the medium. */ amrex::MultiFab& getmu_mf () {return (*m_mu_mf);} /** Initializes the Multifabs storing macroscopic properties * with user-defined functions(x,y,z). */ void InitializeMacroMultiFabUsingParser (amrex::MultiFab *macro_mf, amrex::ParserExecutor<3> const& macro_parser, int lev); /** Gpu Vector with index type of the conductivity multifab */ amrex::GpuArray sigma_IndexType; /** Gpu Vector with index type of the permittivity multifab */ amrex::GpuArray epsilon_IndexType; /** Gpu Vector with index type of the permeability multifab */ amrex::GpuArray mu_IndexType; /** Gpu Vector with index type of the Ex multifab */ amrex::GpuArray Ex_IndexType; /** Gpu Vector with index type of the Ey multifab */ amrex::GpuArray Ey_IndexType; /** Gpu Vector with index type of the Ez multifab */ amrex::GpuArray Ez_IndexType; /** Gpu Vector with index type of coarsening ratio with default value (1,1,1) */ amrex::GpuArray macro_cr_ratio; private: /** Conductivity, sigma, of the medium */ amrex::Real m_sigma = 0.0; /** Permittivity, epsilon, of the medium */ amrex::Real m_epsilon = PhysConst::ep0; /** Permeability, mu, of the medium */ amrex::Real m_mu = PhysConst::mu0; /** Multifab for m_sigma */ std::unique_ptr m_sigma_mf; /** Multifab for m_epsilon */ std::unique_ptr m_eps_mf; /** Multifab for m_mu */ std::unique_ptr m_mu_mf; /** Stores initialization type for conductivity : constant or parser */ std::string m_sigma_s = "constant"; /** Stores initialization type for permittivity : constant or parser */ std::string m_epsilon_s = "constant"; /** Stores initialization type for permeability : constant or parser */ std::string m_mu_s = "constant"; /** string for storing parser function */ std::string m_str_sigma_function; std::string m_str_epsilon_function; std::string m_str_mu_function; /** Parser Wrappers */ std::unique_ptr m_sigma_parser; std::unique_ptr m_epsilon_parser; std::unique_ptr m_mu_parser; }; /** * \brief * This struct contains only static functions to compute the co-efficients for the * Lax-Wendroff scheme of macroscopic Maxwells equations using * macroscopic properties, namely, conductivity (sigma), permittivity (epsilon). * Permeability of the material, mu, is used as (beta/mu) for the E-update * defined in MacroscopicEvolveECartesian(). */ struct LaxWendroffAlgo { AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE static amrex::Real alpha (amrex::Real const sigma, amrex::Real const epsilon, amrex::Real dt) { using namespace amrex; amrex::Real fac1 = 0.5_rt * sigma * dt / epsilon; amrex::Real alpha = (1._rt - fac1)/(1._rt + fac1); return alpha; } AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE static amrex::Real beta (amrex::Real const sigma, amrex::Real const epsilon, amrex::Real dt) { using namespace amrex; amrex::Real fac1 = 0.5_rt * sigma * dt / epsilon; amrex::Real beta = dt / ( epsilon * (1._rt + fac1) ); return beta; } }; /** * \brief * This struct contains only static functions to compute the co-efficients for the * BackwardEuler scheme of macroscopic Maxwells equations using * macroscopic properties, namely, conductivity (sigma) and permittivity (epsilon). * Permeability of the material, mu, is used as (beta/mu) for the E-update * defined in MacroscopicEvolveECartesian(). */ struct BackwardEulerAlgo { AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE static amrex::Real alpha (amrex::Real const sigma, amrex::Real const epsilon, amrex::Real dt) { using namespace amrex; amrex::Real fac1 = sigma * dt / epsilon; amrex::Real alpha = (1._rt)/(1._rt + fac1); return alpha; } AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE static amrex::Real beta (amrex::Real const sigma, amrex::Real const epsilon, amrex::Real dt) { using namespace amrex; amrex::Real fac1 = sigma * dt / epsilon; amrex::Real beta = dt / ( epsilon * (1._rt + fac1) ); return beta; } }; #endif // WARPX_MACROSCOPIC_PROPERTIES_H_