aboutsummaryrefslogtreecommitdiff
path: root/Source/Utils/WarpXUtil.H
blob: c306a5a648aa7f28cfe5bfa2ad797708ea1de96a (plain) (blame)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
/* Copyright 2019-2020 Andrew Myers, Luca Fedeli, Maxence Thevenet
 * Revathi Jambunathan
 *
 * This file is part of WarpX.
 *
 * License: BSD-3-Clause-LBNL
 */
#ifndef WARPX_UTILS_H_
#define WARPX_UTILS_H_

#include <AMReX_BoxArray.H>
#include <AMReX_DistributionMapping.H>
#include <AMReX_Extension.H>
#include <AMReX_GpuQualifiers.H>
#include <AMReX_LayoutData.H>
#include <AMReX_ParmParse.H>
#include <AMReX_Parser.H>
#include <AMReX_REAL.H>
#include <AMReX_Utility.H>
#include <AMReX_Vector.H>

#include <AMReX_BaseFwd.H>

#include <cstddef>
#include <cstdint>
#include <string>
#include <vector>

void ParseGeometryInput();

void ReadBoostedFrameParameters(amrex::Real& gamma_boost, amrex::Real& beta_boost,
                                amrex::Vector<int>& boost_direction);

void ConvertLabParamsToBoost();

/**
 * Reads the user-defined field and particle boundary condition parameters
 */
void ReadBCParams ();

/** Check the warpx.dims matches the binary name
 */
void CheckDims ();

/** Check the warpx.dims matches the binary name & set up RZ gridding
 *
 * Ensures that the blocks are setup correctly for the RZ spectral solver
 * When using the RZ spectral solver, the Hankel transform cannot be
 * divided among multiple blocks. Each block must extend over the
 * entire radial extent.
 * The grid can be divided up along z, but the number of blocks
 * must be >= the number of processors.
 */
void CheckGriddingForRZSpectral ();

void NullifyMF(amrex::MultiFab& mf, int lev, amrex::Real zmin,
               amrex::Real zmax);


namespace WarpXUtilIO{
/**
 * A helper function to write binary data on disk.
 * @param[in] filename where to write
 * @param[in] data Vector containing binary data to write on disk
 * return true if it succeeds, false otherwise
 */
bool WriteBinaryDataOnFile(std::string filename, const amrex::Vector<char>& data);

}

namespace WarpXUtilAlgo{

/** \brief Compute physical coordinates (x,y,z) that correspond to a given (i,j,k) and
 *  the corresponding staggering, mf_type.
 *
 * \param[in] i    index along x
 * \param[in] j    index along y
 * \param[in] k    index along z
 * \param[in] mf_type GpuArray containing the staggering type to convert (i,j,k) to (x,y,z)
 * \param[in] domain_lo Physical coordinates of the lowest corner of the simulation domain
 * \param[in] dx   Cell size of the simulation domain
 *
 * \param[out] x   physical coordinate along x
 * \param[out] y   physical coordinate along y
 * \param[out] z   physical coordinate along z
 */
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)
{
    using namespace amrex::literals;
    x = domain_lo[0] + i*dx[0] + (1._rt - mf_type[0]) * dx[0]*0.5_rt;
#if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
    amrex::ignore_unused(j);
    y = 0._rt;
    z = domain_lo[1] + k*dx[1] + (1._rt - mf_type[1]) * dx[1]*0.5_rt;
#else
    y = domain_lo[1] + j*dx[1] + (1._rt - mf_type[1]) * dx[1]*0.5_rt;
    z = domain_lo[2] + k*dx[2] + (1._rt - mf_type[2]) * dx[2]*0.5_rt;
#endif
}

}


namespace WarpXUtilLoadBalance
{
    /** \brief 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.
     * @param[in] cost pointer to the cost data
     * @param[in] ba the grids to check
     * @param[in] dm the dmap to check
     * @return consistent whether the grids are consistent or not.
     */
    bool doCosts (const amrex::LayoutData<amrex::Real>* cost, const amrex::BoxArray ba,
                  const amrex::DistributionMapping& dm);
}

#endif //WARPX_UTILS_H_