aboutsummaryrefslogtreecommitdiff
path: root/Source/Particles/Collision/BackgroundStopping/BackgroundStopping.cpp
blob: e832020683dad0679afe49fbf664fa95ec91d3e6 (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
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
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
/* Copyright 2022 David Grote
 *
 * This file is part of WarpX.
 *
 * License: BSD-3-Clause-LBNL
 */
#include "BackgroundStopping.H"

#include "Utils/Parser/ParserUtils.H"
#include "Utils/ParticleUtils.H"
#include "Utils/WarpXProfilerWrapper.H"
#include "WarpX.H"

#include <AMReX_ParmParse.H>
#include <AMReX_REAL.H>

#include <string>

BackgroundStopping::BackgroundStopping (std::string const collision_name)
    : CollisionBase(collision_name)
{
    using namespace amrex::literals;

    AMREX_ALWAYS_ASSERT_WITH_MESSAGE(m_species_names.size() == 1,
                                     "Background stopping must have exactly one species.");

    const amrex::ParmParse pp_collision_name(collision_name);

    std::string background_type_str;
    pp_collision_name.get("background_type", background_type_str);
    if (background_type_str == "electrons") {
        m_background_type = BackgroundStoppingType::ELECTRONS;
    } else if (background_type_str == "ions") {
        m_background_type = BackgroundStoppingType::IONS;
    } else {
        AMREX_ALWAYS_ASSERT_WITH_MESSAGE(false, "background_type must be either electrons or ions");
    }

    amrex::ParticleReal background_density;
    std::string background_density_str;
    if (utils::parser::queryWithParser(pp_collision_name, "background_density", background_density)) {
        AMREX_ALWAYS_ASSERT_WITH_MESSAGE(background_density > 0_prt,
                 "For background stopping, the background density must be greater than 0");
        m_background_density_parser =
            utils::parser::makeParser(std::to_string(background_density), {"x", "y", "z", "t"});
    } else if (pp_collision_name.query("background_density(x,y,z,t)", background_density_str)) {
        m_background_density_parser =
            utils::parser::makeParser(background_density_str, {"x", "y", "z", "t"});
    } else {
        AMREX_ALWAYS_ASSERT_WITH_MESSAGE(false,
                 "For background stopping, the background density must be specified.");
    }

    amrex::ParticleReal background_temperature;
    std::string background_temperature_str;
    if (utils::parser::queryWithParser(pp_collision_name, "background_temperature", background_temperature)) {
        AMREX_ALWAYS_ASSERT_WITH_MESSAGE(background_temperature > 0_prt,
                 "For background stopping, the background temperature must be greater than 0");
        m_background_temperature_parser =
            utils::parser::makeParser(std::to_string(background_temperature), {"x", "y", "z", "t"});
    } else if (pp_collision_name.query("background_temperature(x,y,z,t)", background_temperature_str)) {
        m_background_temperature_parser =
            utils::parser::makeParser(background_temperature_str, {"x", "y", "z", "t"});
    } else {
        AMREX_ALWAYS_ASSERT_WITH_MESSAGE(false,
                 "For background stopping, the background temperature must be specified.");
    }

    constexpr auto num_parser_args = 4;
    m_background_density_func = m_background_density_parser.compile<num_parser_args>();
    m_background_temperature_func = m_background_temperature_parser.compile<num_parser_args>();

    if (m_background_type == BackgroundStoppingType::ELECTRONS) {
        m_background_mass = PhysConst::m_e;
        utils::parser::queryWithParser(
            pp_collision_name, "background_mass", m_background_mass);
    } else if (m_background_type == BackgroundStoppingType::IONS) {
        utils::parser::getWithParser(
            pp_collision_name, "background_mass", m_background_mass);
        utils::parser::getWithParser(
            pp_collision_name, "background_charge_state", m_background_charge_state);
    }
    AMREX_ALWAYS_ASSERT_WITH_MESSAGE(m_background_mass > 0_prt,
             "For background stopping, the background mass must be greater than 0");

}

void
BackgroundStopping::doCollisions (amrex::Real cur_time, amrex::Real dt, MultiParticleContainer* mypc)
{
    WARPX_PROFILE("BackgroundStopping::doCollisions()");
    using namespace amrex::literals;

    auto& species = mypc->GetParticleContainerFromName(m_species_names[0]);
    amrex::ParticleReal const species_mass = species.getMass();
    amrex::ParticleReal const species_charge = species.getCharge();

    WARPX_ALWAYS_ASSERT_WITH_MESSAGE(species_mass > 0_prt, "Error: With background stopping, the species mass must be > 0");
    WARPX_ALWAYS_ASSERT_WITH_MESSAGE(species_charge != 0_prt, "Error: With background stopping, the species charge must be nonzero");

    const BackgroundStoppingType background_type = m_background_type;

    // Loop over refinement levels
    auto const flvl = species.finestLevel();
    for (int lev = 0; lev <= flvl; ++lev) {

        auto cost = WarpX::getCosts(lev);

        // loop over particles box by box
#ifdef _OPENMP
#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
#endif
        for (WarpXParIter pti(species, lev); pti.isValid(); ++pti) {
            if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers)
            {
                amrex::Gpu::synchronize();
            }
            amrex::Real wt = amrex::second();

            if (background_type == BackgroundStoppingType::ELECTRONS) {
                doBackgroundStoppingOnElectronsWithinTile(pti, dt, cur_time, species_mass, species_charge);
            } else if (background_type == BackgroundStoppingType::IONS) {
                doBackgroundStoppingOnIonsWithinTile(pti, dt, cur_time, species_mass, species_charge);
            }

            if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers)
            {
                amrex::Gpu::synchronize();
                wt = amrex::second() - wt;
                amrex::HostDevice::Atomic::Add(&(*cost)[pti.index()], wt);
            }
        }

    }
}

void BackgroundStopping::doBackgroundStoppingOnElectronsWithinTile (WarpXParIter& pti, amrex::Real dt, amrex::Real t,
                                                                    amrex::ParticleReal species_mass, amrex::ParticleReal species_charge)
{
    using namespace amrex::literals;

    // So that CUDA code gets its intrinsic, not the host-only C++ library version
    using std::sqrt, std::abs, std::log, std::exp;

    // get particle count
    long const np = pti.numParticles();

    // get background particle mass
    amrex::ParticleReal const mass_e = m_background_mass;

    // setup parsers for the background density and temperature
    auto const n_e_func = m_background_density_func;
    auto const T_e_func = m_background_temperature_func;

    // get Struct-Of-Array particle data, also called attribs
    auto& attribs = pti.GetAttribs();
    amrex::ParticleReal* const AMREX_RESTRICT ux = attribs[PIdx::ux].dataPtr();
    amrex::ParticleReal* const AMREX_RESTRICT uy = attribs[PIdx::uy].dataPtr();
    amrex::ParticleReal* const AMREX_RESTRICT uz = attribs[PIdx::uz].dataPtr();

    // May be needed to evaluate the density and/or temperature functions
    auto const GetPosition = GetParticlePosition(pti);

    amrex::ParallelFor(np,
        [=] AMREX_GPU_HOST_DEVICE (long ip)
        {

            amrex::ParticleReal x, y, z;
            GetPosition.AsStored(ip, x, y, z);
            amrex::ParticleReal const n_e = n_e_func(x, y, z, t);
            amrex::ParticleReal const T_e = T_e_func(x, y, z, t)*PhysConst::kb;

            AMREX_ASSERT(n_e > 0_prt);
            AMREX_ASSERT(T_e > 0_prt);

            // This implements the equation 14.12 from Introduction to Plasma Physics,
            // Goldston and Rutherford, the slowing down of beam ions due to collisions with electrons.
            // The equation is written as dV/dt = -alpha*V, and integrated to
            // give V(t+dt) = V(t)*exp(-alpha*dt)

            amrex::ParticleReal constexpr pi = MathConst::pi;
            amrex::ParticleReal constexpr ep0 = PhysConst::ep0;
            amrex::ParticleReal constexpr q_e = PhysConst::q_e;
            amrex::ParticleReal constexpr q_e2 = q_e*q_e;
            amrex::ParticleReal constexpr ep02 = ep0*ep0;

            amrex::ParticleReal const Zb = abs(species_charge/q_e);

            amrex::ParticleReal const vth = sqrt(3_prt*T_e/mass_e);
            amrex::ParticleReal const wp = sqrt(n_e*q_e2/(ep0*mass_e));
            amrex::ParticleReal const lambdadb = vth/wp;
            amrex::ParticleReal const lambdadb3 = lambdadb*lambdadb*lambdadb;
            amrex::ParticleReal const loglambda = log((12_prt*pi/Zb)*(n_e*lambdadb3));

            AMREX_ASSERT(loglambda > 0_prt);

            amrex::ParticleReal const pi32 = pi*sqrt(pi);
            amrex::ParticleReal const q2 = species_charge*species_charge;
            amrex::ParticleReal const T32 = T_e*sqrt(T_e);

            amrex::ParticleReal const alpha = sqrt(2_prt)*n_e*q2*q_e2*sqrt(mass_e)*loglambda/(12_prt*pi32*ep02*species_mass*T32);

            ux[ip] *= exp(-alpha*dt);
            uy[ip] *= exp(-alpha*dt);
            uz[ip] *= exp(-alpha*dt);

        }
        );
}

void BackgroundStopping::doBackgroundStoppingOnIonsWithinTile (WarpXParIter& pti, amrex::Real dt, amrex::Real t,
                                                               amrex::ParticleReal species_mass, amrex::ParticleReal species_charge)
{
    using namespace amrex::literals;

    // So that CUDA code gets its intrinsic, not the host-only C++ library version
    using std::sqrt, std::abs, std::log, std::exp, std::pow;

    // get particle count
    long const np = pti.numParticles();

    // get background particle mass
    amrex::ParticleReal const mass_i = m_background_mass;
    amrex::ParticleReal const charge_state_i = m_background_charge_state;

    // setup parsers for the background density and temperature
    auto const n_i_func = m_background_density_func;
    auto const T_i_func = m_background_temperature_func;

    // get Struct-Of-Array particle data, also called attribs
    auto& attribs = pti.GetAttribs();
    amrex::ParticleReal* const AMREX_RESTRICT ux = attribs[PIdx::ux].dataPtr();
    amrex::ParticleReal* const AMREX_RESTRICT uy = attribs[PIdx::uy].dataPtr();
    amrex::ParticleReal* const AMREX_RESTRICT uz = attribs[PIdx::uz].dataPtr();

    // May be needed to evaluate the density function
    auto const GetPosition = GetParticlePosition(pti);

    amrex::ParallelFor(np,
        [=] AMREX_GPU_HOST_DEVICE (long ip)
        {

            amrex::ParticleReal x, y, z;
            GetPosition.AsStored(ip, x, y, z);
            amrex::ParticleReal const n_i = n_i_func(x, y, z, t);
            amrex::ParticleReal const T_i = T_i_func(x, y, z, t)*PhysConst::kb;

            AMREX_ASSERT(n_i > 0_prt);
            AMREX_ASSERT(T_i > 0_prt);

            // This implements the equation 14.20 from Introduction to Plasma Physics,
            // Goldston and Rutherford, the slowing down of beam ions due to collisions with electrons.
            // The equation is written with energy, W, as dW/dt = -alpha/W**0.5, and integrated to
            // give W(t+dt) = (W(t)**1.5 - 3./2.*alpha*dt)**(2/3)

            amrex::ParticleReal constexpr pi = MathConst::pi;
            amrex::ParticleReal constexpr q_e = PhysConst::q_e;
            amrex::ParticleReal constexpr q_e2 = q_e*q_e;
            amrex::ParticleReal constexpr ep0 = PhysConst::ep0;
            amrex::ParticleReal constexpr ep02 = ep0*ep0;

            amrex::ParticleReal const qi2 = charge_state_i*charge_state_i*q_e2;
            amrex::ParticleReal const qb2 = species_charge*species_charge;
            amrex::ParticleReal const Zb = abs(species_charge/q_e);

            amrex::ParticleReal const vth = sqrt(3_prt*T_i/mass_i);
            amrex::ParticleReal const wp = sqrt(n_i*q_e2/(ep0*mass_i));
            amrex::ParticleReal const lambdadb = vth/wp;
            amrex::ParticleReal const lambdadb3 = lambdadb*lambdadb*lambdadb;
            amrex::ParticleReal const loglambda = log((12_prt*pi/Zb)*(n_i*lambdadb3));

            AMREX_ASSERT(loglambda > 0_prt);

            amrex::ParticleReal const alpha = sqrt(2_prt)*n_i*qi2*qb2*sqrt(species_mass)*loglambda/(8_prt*pi*ep02*mass_i);

            amrex::ParticleReal const W0 = 0.5_prt*species_mass*(ux[ip]*ux[ip] + uy[ip]*uy[ip] + uz[ip]*uz[ip]);
            amrex::ParticleReal const f1 = pow(W0, 1.5_prt) - 1.5_prt*alpha*dt;
            // If f1 goes negative, the particle has fully stopped, so set W1 to 0.
            amrex::ParticleReal const W1 = pow((f1 > 0_prt ? f1 : 0_prt), 2_prt/3_prt);
            amrex::ParticleReal const vscale = (W0 > 0_prt ? std::sqrt(W1/W0) : 0_prt);

            ux[ip] *= vscale;
            uy[ip] *= vscale;
            uz[ip] *= vscale;

        }
        );
}