aboutsummaryrefslogtreecommitdiff
path: root/Examples/Modules/ParticleBoundaryProcess/PICMI_inputs_reflection.py
blob: c9903ad22d466fc93e8f4b0f7a71a5238a258df7 (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
#!/usr/bin/env python3
#
# --- Input file to test particle reflection off an absorbing boundary

import pywarpx
from pywarpx import picmi

constants = picmi.constants

##########################
# numerics parameters
##########################

dt = 7.5e-12

# --- Nb time steps

max_steps = 10

# --- grid

nx = 64
nz = 64

xmin = -125e-6
zmin = -149e-6
xmax = 125e-6
zmax = 1e-6


##########################
# numerics components
##########################

grid = picmi.Cartesian2DGrid(
    number_of_cells = [nx, nz],
    lower_bound = [xmin, zmin],
    upper_bound = [xmax, zmax],
    lower_boundary_conditions = ['dirichlet', 'dirichlet'],
    upper_boundary_conditions = ['dirichlet', 'dirichlet'],
    lower_boundary_conditions_particles = ['open', 'absorbing'],
    upper_boundary_conditions_particles = ['open', 'absorbing'],
    warpx_max_grid_size = 32
)

solver = picmi.ElectrostaticSolver(
    grid=grid, method='Multigrid', required_precision=1e-6,
    warpx_self_fields_verbosity=0
)

#embedded_boundary = picmi.EmbeddedBoundary(
#    implicit_function="-max(max(x-12.5e-6,-12.5e-6-x),max(z+6.15e-5,-8.65e-5-z))"
#)

##########################
# physics components
##########################

uniform_plasma_elec = picmi.UniformDistribution(
    density = 1e15, # number of electrons per m^3
    lower_bound = [-1e-5, -1e-5, -125e-6],
    upper_bound = [1e-5, 1e-5, -120e-6],
    directed_velocity = [0., 0., 5e6] # uth the std of the (unitless) momentum
)

electrons = picmi.Species(
    particle_type='electron', name='electrons',
    initial_distribution=uniform_plasma_elec,
    warpx_save_particles_at_zhi=1,
    warpx_save_particles_at_zlo=1,
    warpx_reflection_model_zhi="0.5"
)

##########################
# diagnostics
##########################

field_diag = picmi.ParticleDiagnostic(
    species=electrons,
    name = 'diag1',
    data_list=['previous_positions'],
    period = 10,
    write_dir = '.',
    warpx_file_prefix = 'Python_particle_reflection_plt'
)

##########################
# simulation setup
##########################

sim = picmi.Simulation(
    solver = solver,
    time_step_size = dt,
    max_steps = max_steps,
    # warpx_embedded_boundary=embedded_boundary,
    verbose = 1
)

sim.add_species(
    electrons,
    layout = picmi.GriddedLayout(
        n_macroparticle_per_cell=[5, 2], grid=grid
    )
)
sim.add_diagnostic(field_diag)

##########################
# simulation run
##########################

sim.step(max_steps)

################################################
# check that the wrappers to access the particle
# buffer functions as intended
################################################

n = pywarpx.get_particle_boundary_buffer_size("electrons", 'z_hi')
print("Number of electrons in upper buffer:", n)
assert n == 63

n = pywarpx.get_particle_boundary_buffer_size("electrons", 'z_lo')
print("Number of electrons in lower buffer:", n)
assert n == 67

scraped_steps = pywarpx.get_particle_boundary_buffer("electrons", 'z_hi', 'step_scraped', 0)
for arr in scraped_steps:
    # print(arr)
    assert all(arr == 4)

scraped_steps = pywarpx.get_particle_boundary_buffer("electrons", 'z_lo', 'step_scraped', 0)
for arr in scraped_steps:
    # print(arr)
    assert all(arr == 8)