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 <AMReX_Array.H>
11 #include <AMReX_IntVect.H>
12 #include <AMReX_REAL.H>
13 #include <AMReX_RealVect.H>
14 #include <AMReX_Vector.H>
15 
21 
22 public:
23 
53  void Init(
54  amrex::Real dt,
56  bool do_subcycling,
57  bool do_fdtd_nci_corr,
58  short grid_type,
59  bool do_moving_window,
60  int moving_window_dir,
61  int nox,
62  int nox_fft, int noy_fft, int noz_fft,
63  int nci_corr_stencil,
64  int electromagnetic_solver_id,
65  int max_level,
66  const amrex::Vector<amrex::Real>& v_galilean,
67  const amrex::Vector<amrex::Real>& v_comoving,
68  bool safe_guard_cells,
69  int do_multi_J,
70  bool fft_do_time_averaging,
71  bool do_pml,
72  int do_pml_in_domain,
73  int pml_ncell,
74  const amrex::Vector<amrex::IntVect>& ref_ratios,
75  bool use_filter,
76  const amrex::IntVect& bilinear_filter_stencil_length);
77 
78  // Guard cells allocated for MultiFabs E and B
80  // Guard cells allocated for MultiFab J
82  // Guard cells allocated for MultiFab Rho
84  // Guard cells allocated for MultiFab F
86  // Guard cells allocated for MultiFab G
88 
89  // Guard cells exchanged for specific parts of the PIC loop
90 
91  // Number of guard cells of E and B that must exchanged before Field Solver
93  // Number of guard cells of F that must exchanged before Field Solver
95  // Number of guard cells of G that must be exchanged before Field Solver
97  // Number of guard cells of E and B that must exchanged before Field Gather
99  // Number of guard cells of E and B that must exchanged before updating the Aux grid
101  // Number of guard cells of all MultiFabs that must exchanged before moving window
103  // Number of guard cells of E and B that are exchanged immediately after the main PSATD push
105 
106  // Number of guard cells for local deposition of J and rho
109 };
110 
111 #endif // WARPX_GUARDCELLMANAGER_H_
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE IntVect TheZeroVector() noexcept
This class computes and stores the number of guard cells needed for the allocation of the MultiFabs a...
Definition: GuardCellManager.H:20
amrex::IntVect ng_alloc_G
Definition: GuardCellManager.H:87
amrex::IntVect ng_afterPushPSATD
Definition: GuardCellManager.H:104
amrex::IntVect ng_MovingWindow
Definition: GuardCellManager.H:102
amrex::IntVect ng_FieldSolver
Definition: GuardCellManager.H:92
amrex::IntVect ng_depos_rho
Definition: GuardCellManager.H:108
amrex::IntVect ng_depos_J
Definition: GuardCellManager.H:107
amrex::IntVect ng_FieldSolverF
Definition: GuardCellManager.H:94
amrex::IntVect ng_alloc_F
Definition: GuardCellManager.H:85
amrex::IntVect ng_alloc_EB
Definition: GuardCellManager.H:79
amrex::IntVect ng_FieldSolverG
Definition: GuardCellManager.H:96
amrex::IntVect ng_alloc_J
Definition: GuardCellManager.H:81
void Init(amrex::Real dt, amrex::RealVect dx, bool do_subcycling, bool do_fdtd_nci_corr, short 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
amrex::IntVect ng_UpdateAux
Definition: GuardCellManager.H:100
amrex::IntVect ng_FieldGather
Definition: GuardCellManager.H:98
amrex::IntVect ng_alloc_Rho
Definition: GuardCellManager.H:83
float dt
Definition: stencil.py:442
tuple dx
lab frame
Definition: stencil.py:429