WarpX
GuardCellManager.H
Go to the documentation of this file.
1 /* Copyright 2019-2020 Maxence Thevenet
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef WARPX_GUARDCELLMANAGER_H_
8 #define WARPX_GUARDCELLMANAGER_H_
9 
10 #include <ablastr/utils/Enums.H>
11 
12 #include <AMReX_Array.H>
13 #include <AMReX_IntVect.H>
14 #include <AMReX_REAL.H>
15 #include <AMReX_RealVect.H>
16 #include <AMReX_Vector.H>
17 
23 
24 public:
25 
55  void Init(
56  amrex::Real dt,
57  const amrex::Real * dx,
58  bool do_subcycling,
59  bool do_fdtd_nci_corr,
61  bool do_moving_window,
62  int moving_window_dir,
63  int nox,
64  int nox_fft, int noy_fft, int noz_fft,
65  int nci_corr_stencil,
66  int electromagnetic_solver_id,
67  int max_level,
68  const amrex::Vector<amrex::Real>& v_galilean,
69  const amrex::Vector<amrex::Real>& v_comoving,
70  bool safe_guard_cells,
71  int do_multi_J,
72  bool fft_do_time_averaging,
73  bool do_pml,
74  int do_pml_in_domain,
75  int pml_ncell,
76  const amrex::Vector<amrex::IntVect>& ref_ratios,
77  bool use_filter,
78  const amrex::IntVect& bilinear_filter_stencil_length);
79 
80  // Guard cells allocated for MultiFabs E and B
82  // Guard cells allocated for MultiFab J
84  // Guard cells allocated for MultiFab Rho
86  // Guard cells allocated for MultiFab F
88  // Guard cells allocated for MultiFab G
90 
91  // Guard cells exchanged for specific parts of the PIC loop
92 
93  // Number of guard cells of E and B that must exchanged before Field Solver
95  // Number of guard cells of F that must exchanged before Field Solver
97  // Number of guard cells of G that must be exchanged before Field Solver
99  // Number of guard cells of E and B that must exchanged before Field Gather
101  // Number of guard cells of E and B that must exchanged before updating the Aux grid
103  // Number of guard cells of all MultiFabs that must exchanged before moving window
105  // Number of guard cells of E and B that are exchanged immediately after the main PSATD push
107 
108  // Number of guard cells for local deposition of J and rho
111 };
112 
113 #endif // WARPX_GUARDCELLMANAGER_H_
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE IntVectND< dim > TheZeroVector() noexcept
This class computes and stores the number of guard cells needed for the allocation of the MultiFabs a...
Definition: GuardCellManager.H:22
amrex::IntVect ng_alloc_G
Definition: GuardCellManager.H:89
amrex::IntVect ng_afterPushPSATD
Definition: GuardCellManager.H:106
amrex::IntVect ng_MovingWindow
Definition: GuardCellManager.H:104
amrex::IntVect ng_FieldSolver
Definition: GuardCellManager.H:94
amrex::IntVect ng_depos_rho
Definition: GuardCellManager.H:110
amrex::IntVect ng_depos_J
Definition: GuardCellManager.H:109
amrex::IntVect ng_FieldSolverF
Definition: GuardCellManager.H:96
amrex::IntVect ng_alloc_F
Definition: GuardCellManager.H:87
amrex::IntVect ng_alloc_EB
Definition: GuardCellManager.H:81
amrex::IntVect ng_FieldSolverG
Definition: GuardCellManager.H:98
amrex::IntVect ng_alloc_J
Definition: GuardCellManager.H:83
amrex::IntVect ng_UpdateAux
Definition: GuardCellManager.H:102
amrex::IntVect ng_FieldGather
Definition: GuardCellManager.H:100
amrex::IntVect ng_alloc_Rho
Definition: GuardCellManager.H:85
void Init(amrex::Real dt, const amrex::Real *dx, bool do_subcycling, bool do_fdtd_nci_corr, ablastr::utils::enums::GridType grid_type, bool do_moving_window, int moving_window_dir, int nox, int nox_fft, int noy_fft, int noz_fft, int nci_corr_stencil, int electromagnetic_solver_id, int max_level, const amrex::Vector< amrex::Real > &v_galilean, const amrex::Vector< amrex::Real > &v_comoving, bool safe_guard_cells, int do_multi_J, bool fft_do_time_averaging, bool do_pml, int do_pml_in_domain, int pml_ncell, const amrex::Vector< amrex::IntVect > &ref_ratios, bool use_filter, const amrex::IntVect &bilinear_filter_stencil_length)
Initialize number of guard cells depending on the options used.
Definition: GuardCellManager.cpp:34
GridType
Definition: Enums.H:17
float dt
Definition: stencil.py:442
tuple dx
lab frame
Definition: stencil.py:429