WarpX
Loading...
Searching...
No Matches
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
11#include <ablastr/utils/Enums.H>
12
13#include <AMReX_Array.H>
14#include <AMReX_IntVect.H>
15#include <AMReX_REAL.H>
16#include <AMReX_RealVect.H>
17#include <AMReX_Vector.H>
18
24
25public:
26
58 void Init(
59 amrex::Real dt,
60 const amrex::Real * dx,
61 bool do_subcycling,
62 bool do_fdtd_nci_corr,
64 bool do_moving_window,
65 int moving_window_dir,
66 int particle_max_grid_crossings,
67 int nox,
68 int nox_fft, int noy_fft, int noz_fft,
69 int nci_corr_stencil,
70 ElectromagneticSolverAlgo electromagnetic_solver_id,
71 EvolveScheme evolve_scheme,
72 int max_level,
73 const amrex::Vector<amrex::Real>& v_galilean,
74 const amrex::Vector<amrex::Real>& v_comoving,
75 bool safe_guard_cells,
76 int do_psatd_JRhom,
77 bool fft_do_time_averaging,
78 bool do_pml,
79 int do_pml_in_domain,
80 int pml_ncell,
81 const amrex::Vector<amrex::IntVect>& ref_ratios,
82 bool use_filter,
83 const amrex::IntVect& bilinear_filter_stencil_length);
84
85 // Guard cells allocated for MultiFabs E and B
87 // Guard cells allocated for MultiFab J
89 // Guard cells allocated for MultiFab Rho
91 // Guard cells allocated for MultiFab F
93 // Guard cells allocated for MultiFab G
95
96 // Guard cells exchanged for specific parts of the PIC loop
97
98 // Number of guard cells of E and B that must exchanged before Field Solver
100 // Number of guard cells of F that must exchanged before Field Solver
102 // Number of guard cells of G that must be exchanged before Field Solver
104 // Number of guard cells of E and B that must exchanged before Field Gather
106 // Number of guard cells of E and B that must exchanged before updating the Aux grid
108 // Number of guard cells of all MultiFabs that must exchanged before moving window
110 // Number of guard cells of E and B that are exchanged immediately after the main PSATD push
112
113 // Number of guard cells for local deposition of J and rho
116};
117
118#endif // WARPX_GUARDCELLMANAGER_H_
EvolveScheme
struct to select the overall evolve scheme
Definition WarpXAlgorithmSelection.H:37
ElectromagneticSolverAlgo
Definition WarpXAlgorithmSelection.H:58
__host__ static __device__ constexpr 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:23
amrex::IntVect ng_alloc_G
Definition GuardCellManager.H:94
amrex::IntVect ng_afterPushPSATD
Definition GuardCellManager.H:111
amrex::IntVect ng_MovingWindow
Definition GuardCellManager.H:109
amrex::IntVect ng_FieldSolver
Definition GuardCellManager.H:99
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 particle_max_grid_crossings, int nox, int nox_fft, int noy_fft, int noz_fft, int nci_corr_stencil, ElectromagneticSolverAlgo electromagnetic_solver_id, EvolveScheme evolve_scheme, int max_level, const amrex::Vector< amrex::Real > &v_galilean, const amrex::Vector< amrex::Real > &v_comoving, bool safe_guard_cells, int do_psatd_JRhom, 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:36
amrex::IntVect ng_depos_rho
Definition GuardCellManager.H:115
amrex::IntVect ng_depos_J
Definition GuardCellManager.H:114
amrex::IntVect ng_FieldSolverF
Definition GuardCellManager.H:101
amrex::IntVect ng_alloc_F
Definition GuardCellManager.H:92
amrex::IntVect ng_alloc_EB
Definition GuardCellManager.H:86
amrex::IntVect ng_FieldSolverG
Definition GuardCellManager.H:103
amrex::IntVect ng_alloc_J
Definition GuardCellManager.H:88
amrex::IntVect ng_UpdateAux
Definition GuardCellManager.H:107
amrex::IntVect ng_FieldGather
Definition GuardCellManager.H:105
amrex::IntVect ng_alloc_Rho
Definition GuardCellManager.H:90
GridType
Definition Enums.H:23
IntVectND< 3 > IntVect