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
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
|
/* Copyright 2019-2020 Andrew Myers, Burlen Loring, Luca Fedeli
* Maxence Thevenet, Remi Lehe, Revathi Jambunathan
* Revathi Jambunathan
*
* This file is part of WarpX.
*
* License: BSD-3-Clause-LBNL
*/
#include <WarpXUtil.H>
#include <WarpXConst.H>
#include <AMReX_ParmParse.H>
#include <WarpX.H>
#include <cmath>
#include <fstream>
using namespace amrex;
void ReadBoostedFrameParameters(Real& gamma_boost, Real& beta_boost,
Vector<int>& boost_direction)
{
ParmParse pp("warpx");
pp.query("gamma_boost", gamma_boost);
if( gamma_boost > 1. ) {
beta_boost = std::sqrt(1.-1./pow(gamma_boost,2));
std::string s;
pp.get("boost_direction", s);
if (s == "x" || s == "X") {
boost_direction[0] = 1;
}
#if (AMREX_SPACEDIM == 3)
else if (s == "y" || s == "Y") {
boost_direction[1] = 1;
}
#endif
else if (s == "z" || s == "Z") {
boost_direction[2] = 1;
}
else {
const std::string msg = "Unknown boost_dir: "+s;
Abort(msg.c_str());
}
AMREX_ALWAYS_ASSERT_WITH_MESSAGE( s == "z" || s == "Z" ,
"The boost must be in the z direction.");
}
}
void ConvertLabParamsToBoost()
{
Real gamma_boost = 1., beta_boost = 0.;
int max_level = 0;
Vector<int> boost_direction {0,0,0};
ReadBoostedFrameParameters(gamma_boost, beta_boost, boost_direction);
if (gamma_boost <= 1.) return;
Vector<Real> prob_lo(AMREX_SPACEDIM);
Vector<Real> prob_hi(AMREX_SPACEDIM);
Vector<Real> fine_tag_lo(AMREX_SPACEDIM);
Vector<Real> fine_tag_hi(AMREX_SPACEDIM);
Vector<Real> slice_lo(AMREX_SPACEDIM);
Vector<Real> slice_hi(AMREX_SPACEDIM);
ParmParse pp_geom("geometry");
ParmParse pp_wpx("warpx");
ParmParse pp_amr("amr");
ParmParse pp_slice("slice");
pp_geom.getarr("prob_lo",prob_lo,0,AMREX_SPACEDIM);
BL_ASSERT(prob_lo.size() == AMREX_SPACEDIM);
pp_geom.getarr("prob_hi",prob_hi,0,AMREX_SPACEDIM);
BL_ASSERT(prob_hi.size() == AMREX_SPACEDIM);
pp_slice.queryarr("dom_lo",slice_lo,0,AMREX_SPACEDIM);
BL_ASSERT(slice_lo.size() == AMREX_SPACEDIM);
pp_slice.queryarr("dom_hi",slice_hi,0,AMREX_SPACEDIM);
BL_ASSERT(slice_hi.size() == AMREX_SPACEDIM);
pp_amr.query("max_level", max_level);
if (max_level > 0){
pp_wpx.getarr("fine_tag_lo", fine_tag_lo);
pp_wpx.getarr("fine_tag_hi", fine_tag_hi);
}
#if (AMREX_SPACEDIM == 3)
Vector<int> dim_map {0, 1, 2};
#else
Vector<int> dim_map {0, 2};
#endif
for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
{
if (boost_direction[dim_map[idim]]) {
amrex::Real convert_factor;
// Assume that the window travels with speed +c
convert_factor = 1./( gamma_boost * ( 1 - beta_boost ) );
prob_lo[idim] *= convert_factor;
prob_hi[idim] *= convert_factor;
if (max_level > 0){
fine_tag_lo[idim] *= convert_factor;
fine_tag_hi[idim] *= convert_factor;
}
slice_lo[idim] *= convert_factor;
slice_hi[idim] *= convert_factor;
break;
}
}
pp_geom.addarr("prob_lo", prob_lo);
pp_geom.addarr("prob_hi", prob_hi);
if (max_level > 0){
pp_wpx.addarr("fine_tag_lo", fine_tag_lo);
pp_wpx.addarr("fine_tag_hi", fine_tag_hi);
}
pp_slice.addarr("dom_lo",slice_lo);
pp_slice.addarr("dom_hi",slice_hi);
}
/* \brief Function that sets the value of MultiFab MF to zero for z between
* zmin and zmax.
*/
void NullifyMF(amrex::MultiFab& mf, int lev, amrex::Real zmin, amrex::Real zmax){
BL_PROFILE("WarpX::NullifyMF()");
#ifdef _OPENMP
#pragma omp parallel if (Gpu::notInLaunchRegion())
#endif
for(amrex::MFIter mfi(mf, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi){
const amrex::Box& bx = mfi.tilebox();
// Get box lower and upper physical z bound, and dz
const amrex::Real zmin_box = WarpX::LowerCorner(bx, lev)[2];
const amrex::Real zmax_box = WarpX::UpperCorner(bx, lev)[2];
amrex::Real dz = WarpX::CellSize(lev)[2];
// Get box lower index in the z direction
#if (AMREX_SPACEDIM==3)
const int lo_ind = bx.loVect()[2];
#else
const int lo_ind = bx.loVect()[1];
#endif
// Check if box intersect with [zmin, zmax]
if ( (zmax>zmin_box && zmin<=zmax_box) ){
Array4<Real> arr = mf[mfi].array();
// Set field to 0 between zmin and zmax
ParallelFor(bx,
[=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept{
#if (AMREX_SPACEDIM==3)
const Real z_gridpoint = zmin_box+(k-lo_ind)*dz;
#else
const Real z_gridpoint = zmin_box+(j-lo_ind)*dz;
#endif
if ( (z_gridpoint >= zmin) && (z_gridpoint < zmax) ) {
arr(i,j,k) = 0.;
};
}
);
}
}
}
namespace WarpXUtilIO{
bool WriteBinaryDataOnFile(std::string filename, const amrex::Vector<char>& data)
{
std::ofstream of{filename, std::ios::binary};
of.write(data.data(), data.size());
of.close();
return of.good();
}
}
void Store_parserString(amrex::ParmParse& pp, std::string query_string,
std::string& stored_string)
{
char cstr[query_string.size()+1];
strcpy(cstr, query_string.c_str());
std::vector<std::string> f;
pp.getarr(cstr, f);
stored_string.clear();
for (auto const& s : f) {
stored_string += s;
}
f.clear();
}
WarpXParser makeParser (std::string const& parse_function)
{
WarpXParser parser(parse_function);
parser.registerVariables({"x","y","z","t"});
ParmParse pp("my_constants");
std::set<std::string> symbols = parser.symbols();
symbols.erase("x");
symbols.erase("y");
symbols.erase("z");
symbols.erase("t");
for (auto it = symbols.begin(); it != symbols.end(); ) {
Real v;
if (pp.query(it->c_str(), v)) {
parser.setConstant(*it, v);
it = symbols.erase(it);
} else {
++it;
}
}
for (auto const& s : symbols) {
amrex::Abort("makeParser::Unknown symbol "+s);
}
return parser;
}
|