WarpX
Loading...
Searching...
No Matches
MultiFluidContainer.H
Go to the documentation of this file.
1/* Copyright 2023 Grant Johnson, Remi Lehe
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7#ifndef WARPX_MultiFluidContainer_H_
8#define WARPX_MultiFluidContainer_H_
9
11
13
14#include <AMReX_Geometry.H>
15#include <AMReX_MultiFab.H>
16#include <AMReX_Vector.H>
17#include <AMReX_REAL.H>
18
19#include <string>
20
36{
37
38public:
39
41
43
48
49 [[nodiscard]] WarpXFluidContainer&
50 GetFluidContainer (int ispecies) const {return *allcontainers[ispecies];}
51
52#ifdef WARPX_USE_OPENPMD
53 std::unique_ptr<WarpXFluidContainer>& GetUniqueContainer(int ispecies) {
54 return allcontainers[ispecies];
55 }
56#endif
57
59
60 void InitData (
61 ablastr::fields::MultiFabRegister& m_fields, amrex::Box init_box, amrex::Real cur_time, int lev,
62 const amrex::Geometry& geom_lev, amrex::Real gamma_boost, amrex::Real beta_boost);
63
69 int lev,
70 std::string const& current_fp_string,
71 amrex::Real cur_time,
72 bool skip_deposition=false);
73
74 [[nodiscard]] int nSpecies() const {return static_cast<int>(species_names.size());}
75
79 int lev);
80
81private:
82
83 std::vector<std::string> species_names;
84
85 // Vector of fluid species
87
88};
89#endif /*WARPX_MultiFluidContainer_H_*/
MultiFluidContainer & operator=(MultiFluidContainer const &)=delete
void Evolve(ablastr::fields::MultiFabRegister &fields, int lev, std::string const &current_fp_string, amrex::Real cur_time, bool skip_deposition=false)
Definition MultiFluidContainer.cpp:66
MultiFluidContainer & operator=(MultiFluidContainer &&)=default
int nSpecies() const
Definition MultiFluidContainer.H:74
WarpXFluidContainer & GetFluidContainer(int ispecies) const
Definition MultiFluidContainer.H:50
MultiFluidContainer(MultiFluidContainer const &)=delete
void DepositCharge(ablastr::fields::MultiFabRegister &m_fields, amrex::MultiFab &rho, int lev)
Definition MultiFluidContainer.cpp:49
void DepositCurrent(ablastr::fields::MultiFabRegister &m_fields, amrex::MultiFab &jx, amrex::MultiFab &jy, amrex::MultiFab &jz, int lev)
Definition MultiFluidContainer.cpp:57
std::vector< std::string > species_names
Definition MultiFluidContainer.H:83
MultiFluidContainer(MultiFluidContainer &&)=default
void AllocateLevelMFs(ablastr::fields::MultiFabRegister &m_fields, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, int lev)
Definition MultiFluidContainer.cpp:30
std::unique_ptr< WarpXFluidContainer > & GetUniqueContainer(int ispecies)
Definition MultiFluidContainer.H:53
MultiFluidContainer()
Definition MultiFluidContainer.cpp:16
amrex::Vector< std::unique_ptr< WarpXFluidContainer > > allcontainers
Definition MultiFluidContainer.H:86
~MultiFluidContainer()=default
void InitData(ablastr::fields::MultiFabRegister &m_fields, amrex::Box init_box, amrex::Real cur_time, int lev, const amrex::Geometry &geom_lev, amrex::Real gamma_boost, amrex::Real beta_boost)
Definition MultiFluidContainer.cpp:38
Definition WarpXFluidContainer.H:30
BoxND< 3 > Box
Definition MultiFabRegister.H:262