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
|
#!/usr/bin/env python3
#
# --- Input file for binary Coulomb collision testing. This input script
# --- runs the same test as inputs_2d but via Python, therefore the input
# --- values where directly copied from inputs_2d.
from pywarpx import picmi
constants = picmi.constants
#################################
####### GENERAL PARAMETERS ######
#################################
max_steps = 150
max_grid_size = 8
max_level = 0
nx = max_grid_size
ny = max_grid_size
xmin = 0
xmax = 4.154046151855669e2
ymin = xmin
ymax = xmax
plasma_density = 1e21
elec_rms_velocity = 0.044237441120300*constants.c
ion_rms_velocity = 0.006256118919701*constants.c
number_per_cell = 200
#################################
############ NUMERICS ###########
#################################
serialize_initial_conditions = 1
verbose = 1
cfl = 1.0
# Order of particle shape factors
particle_shape = 1
#################################
############ PLASMA #############
#################################
elec_dist = picmi.UniformDistribution(
density=plasma_density,
rms_velocity=[elec_rms_velocity]*3,
directed_velocity=[elec_rms_velocity, 0., 0.]
)
ion_dist = picmi.UniformDistribution(
density=plasma_density,
rms_velocity=[ion_rms_velocity]*3,
)
electrons = picmi.Species(
particle_type='electron', name='electron',
warpx_do_not_deposit=1,
initial_distribution=elec_dist,
)
ions = picmi.Species(
particle_type='H',
name='ion', charge='q_e',
mass="5*m_e",
warpx_do_not_deposit=1,
initial_distribution=ion_dist
)
#################################
########## COLLISIONS ###########
#################################
collision1 = picmi.CoulombCollisions(
name='collisions1',
species=[electrons, ions],
CoulombLog=15.9
)
collision2 = picmi.CoulombCollisions(
name='collisions2',
species=[electrons, electrons],
CoulombLog=15.9
)
collision3 = picmi.CoulombCollisions(
name='collisions3',
species=[ions, ions],
CoulombLog=15.9
)
#################################
###### GRID AND SOLVER ##########
#################################
grid = picmi.Cartesian2DGrid(
number_of_cells=[nx, ny],
warpx_max_grid_size=max_grid_size,
warpx_blocking_factor=max_grid_size,
lower_bound=[xmin, ymin],
upper_bound=[xmax, ymax],
lower_boundary_conditions=['periodic', 'periodic'],
upper_boundary_conditions=['periodic', 'periodic'],
)
solver = picmi.ElectromagneticSolver(grid=grid, cfl=cfl)
#################################
######### DIAGNOSTICS ###########
#################################
field_diag = picmi.FieldDiagnostic(
name='diag1',
grid=grid,
period=10,
data_list=[],
write_dir='.',
warpx_file_prefix='Python_collisionXZ_plt'
)
#################################
####### SIMULATION SETUP ########
#################################
sim = picmi.Simulation(
solver=solver,
max_steps=max_steps,
verbose=verbose,
warpx_serialize_initial_conditions=serialize_initial_conditions,
warpx_collisions=[collision1, collision2, collision3]
)
sim.add_species(
electrons,
layout=picmi.PseudoRandomLayout(
n_macroparticles_per_cell=number_per_cell, grid=grid
)
)
sim.add_species(
ions,
layout=picmi.PseudoRandomLayout(
n_macroparticles_per_cell=number_per_cell, grid=grid
)
)
sim.add_diagnostic(field_diag)
#################################
##### SIMULATION EXECUTION ######
#################################
#sim.write_input_file('PICMI_inputs_2d')
sim.step(max_steps)
|