WarpX
WarpXUtil.H
Go to the documentation of this file.
1 /* Copyright 2019-2020 Andrew Myers, Luca Fedeli, Maxence Thevenet
2  * Revathi Jambunathan
3  *
4  * This file is part of WarpX.
5  *
6  * License: BSD-3-Clause-LBNL
7  */
8 #ifndef WARPX_UTILS_H_
9 #define WARPX_UTILS_H_
10 
11 #include <AMReX_BoxArray.H>
13 #include <AMReX_Extension.H>
14 #include <AMReX_GpuQualifiers.H>
15 #include <AMReX_LayoutData.H>
16 #include <AMReX_ParmParse.H>
17 #include <AMReX_Parser.H>
18 #include <AMReX_REAL.H>
19 #include <AMReX_Utility.H>
20 #include <AMReX_Vector.H>
21 
22 #include <AMReX_BaseFwd.H>
23 
24 #include <cstddef>
25 #include <cstdint>
26 #include <string>
27 #include <vector>
28 
29 void ParseGeometryInput();
30 
31 void ReadBoostedFrameParameters(amrex::Real& gamma_boost, amrex::Real& beta_boost,
32  amrex::Vector<int>& boost_direction);
33 
35 
39 void ReadBCParams ();
40 
43 void CheckDims ();
44 
55 
56 void NullifyMF(amrex::MultiFab& mf, int lev, amrex::Real zmin,
57  amrex::Real zmax);
58 
59 
60 namespace WarpXUtilIO{
67 bool WriteBinaryDataOnFile(std::string filename, const amrex::Vector<char>& data);
68 
69 }
70 
71 namespace WarpXUtilAlgo{
72 
88 void getCellCoordinates (int i, int j, int k,
89  amrex::GpuArray<int, 3> const mf_type,
92  amrex::Real &x, amrex::Real &y, amrex::Real &z)
93 {
94  using namespace amrex::literals;
95  x = domain_lo[0] + i*dx[0] + (1._rt - mf_type[0]) * dx[0]*0.5_rt;
96 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
98  y = 0._rt;
99  z = domain_lo[1] + k*dx[1] + (1._rt - mf_type[1]) * dx[1]*0.5_rt;
100 #else
101  y = domain_lo[1] + j*dx[1] + (1._rt - mf_type[1]) * dx[1]*0.5_rt;
102  z = domain_lo[2] + k*dx[2] + (1._rt - mf_type[2]) * dx[2]*0.5_rt;
103 #endif
104 }
105 
106 }
107 
108 
109 namespace WarpXUtilLoadBalance
110 {
119  bool doCosts (const amrex::LayoutData<amrex::Real>* cost, const amrex::BoxArray ba,
120  const amrex::DistributionMapping& dm);
121 }
122 
123 #endif //WARPX_UTILS_H_
bool doCosts(const amrex::LayoutData< amrex::Real > *costs, const amrex::BoxArray ba, const amrex::DistributionMapping &dm)
We only want to update the cost data if the grids we are working on are the main grids, i.e. not the PML grids. This function returns whether this is the case or not.
Definition: WarpXUtil.cpp:470
void ReadBCParams()
Definition: WarpXUtil.cpp:385
data
Definition: run_alltests_1node.py:325
int gamma_boost
Definition: compute_domain.py:41
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void getCellCoordinates(int i, int j, int k, amrex::GpuArray< int, 3 > const mf_type, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const domain_lo, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const dx, amrex::Real &x, amrex::Real &y, amrex::Real &z)
Compute physical coordinates (x,y,z) that correspond to a given (i,j,k) and the corresponding stagger...
Definition: WarpXUtil.H:88
void ParseGeometryInput()
Definition: WarpXUtil.cpp:56
int dx
Definition: compute_domain.py:35
void ConvertLabParamsToBoost()
Definition: WarpXUtil.cpp:139
void ReadBoostedFrameParameters(amrex::Real &gamma_boost, amrex::Real &beta_boost, amrex::Vector< int > &boost_direction)
void CheckGriddingForRZSpectral()
Definition: WarpXUtil.cpp:298
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
filename
Definition: write_atomic_data_cpp.py:21
Definition: WarpXUtil.H:71
#define AMREX_FORCE_INLINE
Definition: WarpXUtil.cpp:468
#define AMREX_GPU_HOST_DEVICE
i
Definition: check_interp_points_and_weights.py:174
Definition: WarpXUtil.cpp:266
void CheckDims()
Definition: WarpXUtil.cpp:276
void NullifyMF(amrex::MultiFab &mf, int lev, amrex::Real zmin, amrex::Real zmax)
Definition: WarpXUtil.cpp:224
bool WriteBinaryDataOnFile(std::string filename, const amrex::Vector< char > &data)
Definition: WarpXUtil.cpp:267