afivo-streamer 1.1
1D/2D/3D streamer simulations with AMR
Loading...
Searching...
No Matches
m_field.f90
Go to the documentation of this file.
1!> Module to compute electric fields
2module m_field
3#include "../afivo/src/cpp_macros.h"
4 use m_af_all
5 use m_types
6 use m_streamer
7
8 implicit none
9 private
10
11 integer, parameter :: scalar_voltage = 1
12 integer, parameter :: tabulated_voltage = 2
13
14 !> How the electric field or voltage is specified
15 integer :: field_given_by = -1
16
17 !> List of times
18 real(dp), allocatable :: field_table_times(:)
19
20 !> List of voltages
21 real(dp), allocatable :: field_table_values(:)
22
23 !> Linear rise time of field (s)
24 real(dp) :: field_rise_time = 0.0_dp
25
26 !> Pulse width excluding rise and fall time
27 real(dp) :: field_pulse_width = huge_real
28
29 !> Number of voltage pulses
30 integer :: field_num_pulses = 1
31
32 !> Time of one complete voltage pulse
33 real(dp), public, protected :: field_pulse_period = huge_real
34
35 !> The (initial) vertical applied electric field
36 real(dp) :: field_amplitude = undefined_real
37
38 !> The applied voltage (vertical direction)
39 real(dp) :: field_voltage = undefined_real
40
41 !> Whether electrode 1 is grounded or at the applied voltage
42 logical :: field_electrode_grounded = .false.
43
44 !> Whether electrode 2 is grounded or at the applied voltage
45 logical :: field_electrode2_grounded = .false.
46
47 !> Electrode 1: first relative coordinate
48 real(dp) :: rod_r0(ndim) = -1.0e100_dp
49
50 !> Electrode 1: second relative coordinate
51 real(dp) :: rod_r1(ndim) = -1.0e100_dp
52
53 !> Electrode 2: first relative coordinate
54 real(dp) :: rod2_r0(ndim) = -1.0e100_dp
55
56 !> Electrode 2: second relative coordinate
57 real(dp) :: rod2_r1(ndim) = -1.0e100_dp
58
59 !> Electrode 1 radius (in m)
60 real(dp) :: rod_radius = -1.0e100_dp
61
62 !> Electrode 2 radius (in m)
63 real(dp) :: rod2_radius = -1.0e100_dp
64
65 !> Electrode 1: fraction of conical part (if conical)
66 real(dp) :: cone_length_frac = -1.0e100_dp
67
68 !> Electrode 2: fraction of conical part (if conical)
69 real(dp) :: cone2_length_frac = -1.0e100_dp
70
71 !> Electrode 1: tip radius (if conical)
72 real(dp) :: cone_tip_radius = -1.0e100_dp
73
74 !> Electrode 2: tip radius (if conical)
75 real(dp) :: cone2_tip_radius = -1.0e100_dp
76
77 ! Internal variables
78
79 !> Electrode 1: 'origin' of spherical tip (if conical)
80 real(dp) :: cone_tip_center(ndim)
81 !> Electrode 2: 'origin' of spherical tip (if conical)
82 real(dp) :: cone2_tip_center(ndim)
83
84 !> Electrode 1: radius of curvature of spherical tip (if conical)
85 real(dp) :: cone_tip_r_curvature
86 !> Electrode 2: radius of curvature of spherical tip (if conical)
87 real(dp) :: cone2_tip_r_curvature
88
89 !> The current applied voltage
90 real(dp), public, protected :: current_voltage = 0.0_dp
91
92 character(string_len) :: field_bc_type = "homogeneous"
93
94 public :: field_initialize
95 public :: field_compute
96 public :: field_set_rhs
97 public :: field_set_voltage
98
99 public :: field_bc_homogeneous
100 public :: field_from_potential
101 public :: field_compute_energy
102 public :: field_get_e_vector
103
104contains
105
106 !> Initialize this module
107 subroutine field_initialize(tree, cfg, mg)
108 use m_config
109 use m_table_data
111 type(af_t), intent(inout) :: tree
112 type(cfg_t), intent(inout) :: cfg !< Settings
113 type(mg_t), intent(inout) :: mg !< Multigrid option struct
114 character(len=string_len) :: given_by, user_value
115 character(len=string_len) :: electrode_type
116 integer :: first_blank
117
118 ! This is for backward compatibility
119 call cfg_add_get(cfg, "field_amplitude", field_amplitude, &
120 "The (initial) vertical applied electric field (V/m)")
121
122 given_by = undefined_str
123 call cfg_add_get(cfg, "field_given_by", given_by, &
124 "How the electric field or voltage is specified")
125
126 ! Split given_by string
127 given_by = adjustl(given_by)
128 first_blank = index(given_by, " ")
129 user_value = adjustl(given_by(first_blank:))
130 given_by = given_by(1:first_blank-1)
131
132 select case (given_by)
133 case ("voltage")
134 field_given_by = scalar_voltage
135 read(user_value, *) field_voltage
136 case ("field")
137 field_given_by = scalar_voltage
138 read(user_value, *) field_voltage
139 ! Convert to a voltage
140 field_voltage = -st_domain_len(ndim) * field_voltage
141 case ("voltage_table")
142 field_given_by = tabulated_voltage
143
144 ! Read in the table
145 call table_from_file(trim(user_value), "voltage_vs_time", &
146 field_table_times, field_table_values)
147 case ("field_table")
148 field_given_by = tabulated_voltage
149
150 ! Read in the table
151 call table_from_file(trim(user_value), "field_vs_time", &
152 field_table_times, field_table_values)
153 ! Convert to a voltage
154 field_table_values = -st_domain_len(ndim) * field_table_values
155 case (undefined_str)
156 if (field_amplitude <= undefined_real) then
157 error stop "field_amplitude not specified"
158 else
159 print *, "Warning: field_amplitude is deprecated, use field_given_by"
160 field_given_by = scalar_voltage
161 field_voltage = -st_domain_len(ndim) * field_amplitude
162 end if
163 case default
164 print *, "field_given_by value: ", trim(given_by), ", options are:"
165 print *, "1. voltage <value in V>"
166 print *, "2. field <value in V/m>"
167 print *, "3. voltage_table <filename>"
168 print *, "4. field_table <filename>"
169 error stop "Unknown field_given_by value"
170 end select
171
172 call cfg_add_get(cfg, "field_rise_time", field_rise_time, &
173 "Linear rise time of field (s)")
174 call cfg_add_get(cfg, "field_pulse_width", field_pulse_width, &
175 "Pulse width excluding rise and fall time (s)")
176 call cfg_add_get(cfg, "field_num_pulses", field_num_pulses, &
177 "Number of voltage pulses (default: 1)")
178 call cfg_add_get(cfg, "field_pulse_period", field_pulse_period, &
179 "Time of one complete voltage pulse (s)")
180
181 if (field_pulse_width < huge_real .and. field_rise_time <= 0) &
182 error stop "Set field_rise_time when using field_pulse_width"
183
184 if (field_num_pulses > 1) then
186 error stop "field_num_pulses > 1 requires field_pulse_period"
187 if (field_pulse_width >= huge_real) &
188 error stop "field_num_pulses > 1 requires field_pulse_width"
189 if (field_pulse_width + 2 * field_rise_time > field_pulse_period) &
190 error stop "field_pulse_period shorter than one full pulse"
191 end if
192
193 call cfg_add_get(cfg, "field_bc_type", field_bc_type, &
194 "Boundary condition for electric potential")
195
196 !< [electrode_settings]
197 call cfg_add_get(cfg, "field_electrode_grounded", field_electrode_grounded, &
198 "Whether electrode 1 is grounded or at the applied voltage")
199 call cfg_add_get(cfg, "field_electrode2_grounded", field_electrode2_grounded, &
200 "Whether electrode 2 is grounded or at the applied voltage")
201 call cfg_add_get(cfg, "field_rod_r0", rod_r0, &
202 "Electrode 1: first relative coordinate")
203 call cfg_add_get(cfg, "field_rod_r1", rod_r1, &
204 "Electrode 1: second relative coordinate")
205 call cfg_add_get(cfg, "field_rod2_r0", rod2_r0, &
206 "Electrode 2: first relative coordinate")
207 call cfg_add_get(cfg, "field_rod2_r1", rod2_r1, &
208 "Electrode 2: second relative coordinate")
209 call cfg_add_get(cfg, "field_rod_radius", rod_radius, &
210 "Electrode 1 radius (in m)")
211 call cfg_add_get(cfg, "field_rod2_radius", rod2_radius, &
212 "Electrode 2 radius (in m)")
213 call cfg_add_get(cfg, "cone_tip_radius", cone_tip_radius, &
214 "Electrode 1: tip radius (if conical)")
215 call cfg_add_get(cfg, "cone_length_frac", cone_length_frac, &
216 "Electrode 1: fraction of conical part (if conical)")
217 call cfg_add_get(cfg, "cone2_tip_radius", cone2_tip_radius, &
218 "Electrode 2: tip radius (if conical)")
219 call cfg_add_get(cfg, "cone2_length_frac", cone2_length_frac, &
220 "Electrode 2: fraction of conical part (if conical)")
221
222 rod_r0 = st_domain_origin + rod_r0 * st_domain_len
223 rod_r1 = st_domain_origin + rod_r1 * st_domain_len
224 rod2_r0 = st_domain_origin + rod2_r0 * st_domain_len
225 rod2_r1 = st_domain_origin + rod2_r1 * st_domain_len
226
227 electrode_type = "rod"
228 call cfg_add_get(cfg, "field_electrode_type", electrode_type, &
229 "Electrode: sphere, rod, rod_cone_top, rod_rod, sphere_rod, user")
230 !< [electrode_settings]
231
232 if (associated(user_potential_bc)) then
233 mg%sides_bc => user_potential_bc
234 else
235 ! Use one of the predefined boundary conditions
236 select case (field_bc_type)
237 case ("homogeneous")
238 mg%sides_bc => field_bc_homogeneous
239 case ("neumann")
240 mg%sides_bc => field_bc_neumann
241 case ("all_neumann")
242 mg%sides_bc => field_bc_all_neumann
243 case default
244 error stop "field_bc_select error: invalid condition"
245 end select
246 end if
247
248 ! Set the multigrid options. First define the variables to use
249 mg%i_phi = i_phi
250 mg%i_tmp = i_tmp
251 mg%i_rhs = i_rhs
252
253 if (st_use_dielectric) tree%mg_i_eps = i_eps
254
255 if (st_use_electrode) then
256 select case (electrode_type)
257 case ("sphere")
258 ! A single spherical electrode
259 if (any(rod_r0 <= -1.0e10_dp)) &
260 error stop "field_rod_r0 not set correctly"
261 if (rod_radius <= 0) &
262 error stop "field_rod_radius not set correctly"
263 mg%lsf => sphere_lsf
264 case ("rod")
265 ! A single rod electrode with a semi-spherical cap
266 call check_general_electrode_parameters()
267 mg%lsf => rod_lsf
268 case ("rod_cone_top")
269 ! A single rod-shaped electrode with a conical top
270 call check_general_electrode_parameters()
271 if (cone_tip_radius <= 0 .or. cone_tip_radius > rod_radius) &
272 error stop "cone_tip_radius should be smaller than rod radius"
273 if (cone_length_frac < 0 .or. cone_length_frac > 1) &
274 error stop "cone_length_frac not set correctly"
275
276 call get_conical_rod_properties(rod_r0, rod_r1, rod_radius, &
277 cone_tip_radius, cone_tip_center, cone_tip_r_curvature)
278
279 mg%lsf => conical_rod_lsf
280 case ("rod_rod")
281 ! Two rod electrodes with semi-spherical caps
282 call check_general_electrode_parameters()
283
284 if (any(rod2_r0 <= -1.0e10_dp)) &
285 error stop "field_rod2_r0 not set correctly"
286 if (any(rod2_r1 <= -1.0e10_dp)) &
287 error stop "field_rod2_r1 not set correctly"
288 if (rod2_radius <= 0) &
289 error stop "field_rod2_radius not set correctly"
290
291 mg%lsf => rod_rod_lsf
292
293 ! Provide a function to set the voltage on the electrodes
294 mg%lsf_boundary_function => rod_rod_get_potential
295 case ("sphere_rod")
296 ! Sphere and rod electrode with semi-spherical cap
297 call check_general_electrode_parameters()
298
299 if (any(rod2_r0 <= -1.0e10_dp)) &
300 error stop "field_rod2_r0 not set correctly"
301 if (any(rod2_r1 <= -1.0e10_dp)) &
302 error stop "field_rod2_r1 not set correctly"
303 if (rod2_radius <= 0) &
304 error stop "field_rod2_radius not set correctly"
305
306 mg%lsf => sphere_rod_lsf
307 mg%lsf_boundary_function => sphere_rod_get_potential
308 case ("two_rod_cone_electrodes")
309 ! Two rod-shaped electrodes with conical tops (for now assumed to have
310 ! the same shape)
311 call check_general_electrode_parameters()
312 if (any(rod2_r0 <= -1.0e10_dp)) &
313 error stop "field_rod2_r0 not set correctly"
314 if (any(rod2_r1 <= -1.0e10_dp)) &
315 error stop "field_rod2_r1 not set correctly"
316 if (rod2_radius <= 0) &
317 error stop "field_rod2_radius not set correctly"
318 if (cone_tip_radius <= 0 .or. cone_tip_radius > rod_radius) &
319 error stop "cone tip radius should be smaller than rod radius"
320 if (cone2_tip_radius <= 0 .or. cone2_tip_radius > rod2_radius) &
321 error stop "cone2 tip radius should be smaller than rod2 radius"
322 if (cone_length_frac < 0 .or. cone_length_frac > 1) &
323 error stop "cone_length_frac not set correctly"
324 if (cone2_length_frac < 0 .or. cone2_length_frac > 1) &
325 error stop "cone2_length_frac not set correctly"
326
327 call get_conical_rod_properties(rod_r0, rod_r1, rod_radius, &
328 cone_tip_radius, cone_tip_center, cone_tip_r_curvature)
329 call get_conical_rod_properties(rod2_r0, rod2_r1, rod2_radius, &
330 cone2_tip_radius, cone2_tip_center, cone2_tip_r_curvature)
331
332 mg%lsf => two_conical_rods_lsf
333
334 ! Provide a function to set the voltage on the electrodes
335 mg%lsf_boundary_function => two_conical_rods_get_potential
336 case ("coaxial")
337 ! A coaxial geometry, with an inner and outer conductor
338#if NDIM == 1
339 error stop "coaxial not supported in 1D"
340#else
341 if (rod_radius <= 0) &
342 error stop "field_rod_radius not set correctly"
343 if (rod2_radius <= 0) &
344 error stop "field_rod2_radius not set correctly"
345
346 ! Check if outer radius is correctly set
347 if (st_cylindrical) then
348 if (rod2_radius > st_domain_len(1)) &
349 error stop "field_rod2_radius too large"
350 else if (ndim == 2 .or. ndim == 3) then
351 if (rod2_radius > 0.5_dp * minval(st_domain_len(1:2))) &
352 error stop "field_rod2_radius too large"
353 end if
354
355 if (ndim == 3 .and. .not. st_periodic(ndim)) &
356 error stop "Coaxial requires periodic = F F T (last dim. periodic)"
357
358 mg%lsf => coaxial_lsf
359 mg%lsf_boundary_function => coaxial_get_potential
360
361 ! Ensure that outer electrode is grounded
362 field_electrode_grounded = .false.
363 field_electrode2_grounded = .true.
364 mg%sides_bc => field_bc_all_dirichlet
365#endif
366 case ("user")
367 if (.not. associated(user_lsf)) then
368 error stop "user_lsf not set"
369 else
370 mg%lsf => user_lsf
371 end if
372
373 if (associated(user_lsf_bc)) then
374 mg%lsf_boundary_function => user_lsf_bc
375 end if
376 case default
377 print *, "Electrode types: sphere, rod, rod_cone_top, rod_rod, user"
378 error stop "Invalid electrode type"
379 end select
380
381 call af_set_cc_methods(tree, i_lsf, funcval=set_lsf_box)
382 tree%mg_i_lsf = i_lsf
383
384 mg%lsf_dist => mg_lsf_dist_gss
385
386 if (rod_radius <= 0) then
387 error stop "set field_rod_radius to smallest length scale of electrode"
388 end if
389 mg%lsf_length_scale = rod_radius
390 end if
391
392 call af_set_cc_methods(tree, i_electric_fld, &
393 af_bc_neumann_zero, af_gc_interp)
394
395 end subroutine field_initialize
396
397 subroutine check_general_electrode_parameters()
398 if (any(rod_r0 <= -1.0e10_dp)) &
399 error stop "field_rod_r0 not set correctly"
400 if (any(rod_r1 <= -1.0e10_dp)) &
401 error stop "field_rod_r1 not set correctly"
402 if (rod_radius <= 0) &
403 error stop "field_rod_radius not set correctly"
404 end subroutine check_general_electrode_parameters
405
406 subroutine field_set_rhs(tree, s_in)
408 use m_chemistry
409 use m_dielectric
410 type(af_t), intent(inout) :: tree
411 integer, intent(in) :: s_in
412 real(dp), parameter :: fac = -uc_elem_charge / uc_eps0
413 real(dp) :: q
414 integer :: lvl, i, id, nc, n, ix
415
416 nc = tree%n_cell
417
418 ! Set the source term (rhs)
419 !$omp parallel private(lvl, i, id, n, ix, q)
420 do lvl = 1, tree%highest_lvl
421 !$omp do
422 do i = 1, size(tree%lvls(lvl)%leaves)
423 id = tree%lvls(lvl)%leaves(i)
424
425 tree%boxes(id)%cc(dtimes(:), i_rhs) = 0.0_dp
426
427 do n = 1, size(charged_species_itree)
428 ix = charged_species_itree(n) + s_in
429 q = charged_species_charge(n) * fac
430
431 tree%boxes(id)%cc(dtimes(:), i_rhs) = &
432 tree%boxes(id)%cc(dtimes(:), i_rhs) + &
433 q * tree%boxes(id)%cc(dtimes(:), ix)
434 end do
435 end do
436 !$omp end do nowait
437 end do
438 !$omp end parallel
439
440 if (st_use_dielectric) then
441 call surface_surface_charge_to_rhs(tree, diel, i_surf_dens, i_rhs, fac)
442 end if
443
444 end subroutine field_set_rhs
445
446 !> Compute electric field on the tree. First perform multigrid to get electric
447 !> potential, then take numerical gradient to geld field.
448 subroutine field_compute(tree, mg, s_in, time, have_guess)
449 use m_chemistry
450 type(af_t), intent(inout) :: tree
451 type(mg_t), intent(inout) :: mg ! Multigrid option struct
452 integer, intent(in) :: s_in
453 real(dp), intent(in) :: time
454 logical, intent(in) :: have_guess
455 integer :: i
456 real(dp) :: max_rhs, residual_threshold, conv_fac
457 real(dp) :: residual_ratio
458 integer, parameter :: max_initial_iterations = 100
459 real(dp), parameter :: max_residual = 1e8_dp
460 real(dp), parameter :: min_residual = 1e-6_dp
461 real(dp) :: residuals(max_initial_iterations)
462
463 call field_set_rhs(tree, s_in)
464 call field_set_voltage(tree, time)
465
466 call af_tree_maxabs_cc(tree, mg%i_rhs, max_rhs)
467
468 ! With an electrode, the initial convergence testing should be less strict
469 if (st_use_electrode) then
470 conv_fac = 1e-8_dp
471 else
472 conv_fac = 1e-10_dp
473 end if
474
475 ! Set threshold based on rhs and on estimate of round-off error, given by
476 ! delta phi / dx^2 = (phi/L * dx)/dx^2
477 ! Note that we use min_residual in case max_rhs and current_voltage are zero
478 residual_threshold = max(min_residual, &
480 conv_fac * abs(current_voltage)/(st_domain_len(ndim) * af_min_dr(tree)))
481
482 if (st_use_electrode) then
483 if (field_electrode_grounded) then
484 mg%lsf_boundary_value = 0.0_dp
485 else
486 mg%lsf_boundary_value = current_voltage
487 end if
488 end if
489
490 ! Perform a FMG cycle when we have no guess
491 if (.not. have_guess) then
492 do i = 1, max_initial_iterations
493 call mg_fas_fmg(tree, mg, .true., .true.)
494 call af_tree_maxabs_cc(tree, mg%i_tmp, residuals(i))
495
496 if (residuals(i) < residual_threshold) then
497 exit
498 else if (i > 2) then
499 ! Check if the residual is not changing much anymore, and if it is
500 ! small enough, in which case we assume convergence
501 residual_ratio = minval(residuals(i-2:i)) / &
502 maxval(residuals(i-2:i))
503 if (residual_ratio < 2.0_dp .and. residual_ratio > 0.5_dp &
504 .and. residuals(i) < max_residual) exit
505 end if
506 end do
507
508 ! Check for convergence
509 if (i == max_initial_iterations + 1) then
510 print *, "Iteration residual"
511 do i = 1, max_initial_iterations
512 write(*, "(I4,E18.2)") i, residuals(i)
513 end do
514 print *, "Maybe increase multigrid_max_rel_residual?"
515 error stop "No convergence in initial field computation"
516 end if
517 end if
518
519 ! Perform cheaper V-cycles
521 call mg_fas_vcycle(tree, mg, .true.)
522 call af_tree_maxabs_cc(tree, mg%i_tmp, residuals(i))
523 if (residuals(i) < residual_threshold) exit
524 end do
525
526 call field_from_potential(tree, mg)
527
528 end subroutine field_compute
529
530 !> Compute field from potential
531 subroutine field_from_potential(tree, mg)
533 use m_dielectric
534 type(af_t), intent(inout) :: tree
535 type(mg_t), intent(in) :: mg ! Multigrid option struct
536
537 if (st_use_dielectric) then
538 call mg_compute_phi_gradient(tree, mg, electric_fld, -1.0_dp)
539 call surface_correct_field_fc(tree, diel, i_surf_dens, &
541 call mg_compute_field_norm(tree, electric_fld, i_electric_fld)
542 else
543 call mg_compute_phi_gradient(tree, mg, electric_fld, -1.0_dp, i_electric_fld)
544 end if
545
546 ! Set the field norm also in ghost cells
547 call af_gc_tree(tree, [i_electric_fld])
548 end subroutine field_from_potential
549
550 !> Set the current voltage
551 subroutine field_set_voltage(tree, time)
553 use m_lookup_table
555 type(af_t), intent(in) :: tree
556 real(dp), intent(in) :: time
557 real(dp) :: electric_fld, t, tmp
558
559 if (associated(user_field_amplitude)) then
562 return
563 end if
564
565 select case (field_given_by)
566 case (scalar_voltage)
567 current_voltage = 0.0_dp
568
569 if (time < field_pulse_period * field_num_pulses) then
570 t = modulo(time, field_pulse_period)
571
572 if (t < field_rise_time) then
573 current_voltage = field_voltage * (t/field_rise_time)
574 else if (t < field_pulse_width + field_rise_time) then
575 current_voltage = field_voltage
576 else
577 tmp = t - (field_pulse_width + field_rise_time)
578 current_voltage = field_voltage * max(0.0_dp, &
579 (1 - tmp/field_rise_time))
580 end if
581 end if
582 case (tabulated_voltage)
583 call lt_lin_interp_list(field_table_times, field_table_values, &
584 time, current_voltage)
585 end select
586 end subroutine field_set_voltage
587
588 !> Dirichlet boundary conditions for the potential in the last dimension,
589 !> Neumann zero boundary conditions in the other directions
590 subroutine field_bc_homogeneous(box, nb, iv, coords, bc_val, bc_type)
591 type(box_t), intent(in) :: box
592 integer, intent(in) :: nb
593 integer, intent(in) :: iv
594 real(dp), intent(in) :: coords(ndim, box%n_cell**(ndim-1))
595 real(dp), intent(out) :: bc_val(box%n_cell**(ndim-1))
596 integer, intent(out) :: bc_type
597
598 if (af_neighb_dim(nb) == ndim) then
599 if (af_neighb_low(nb)) then
600 bc_type = af_bc_dirichlet
601 bc_val = 0.0_dp
602 else
603 bc_type = af_bc_dirichlet
604 bc_val = current_voltage
605 end if
606 else
607 bc_type = af_bc_neumann
608 bc_val = 0.0_dp
609 end if
610 end subroutine field_bc_homogeneous
611
612 !> A Dirichlet zero and non-zero Neumann boundary condition for the potential
613 !> in the last dimension, Neumann zero boundary conditions in the other
614 !> directions
615 subroutine field_bc_neumann(box, nb, iv, coords, bc_val, bc_type)
616 type(box_t), intent(in) :: box
617 integer, intent(in) :: nb
618 integer, intent(in) :: iv
619 real(dp), intent(in) :: coords(ndim, box%n_cell**(ndim-1))
620 real(dp), intent(out) :: bc_val(box%n_cell**(ndim-1))
621 integer, intent(out) :: bc_type
622
623 if (af_neighb_dim(nb) == ndim) then
624 if (af_neighb_low(nb)) then
625 bc_type = af_bc_dirichlet
626 bc_val = 0.0_dp
627 else
628 bc_type = af_bc_neumann
629 bc_val = current_voltage/st_domain_len(ndim)
630 end if
631 else
632 bc_type = af_bc_neumann
633 bc_val = 0.0_dp
634 end if
635 end subroutine field_bc_neumann
636
637 !> All Neumann zero boundary conditions for the potential
638 subroutine field_bc_all_neumann(box, nb, iv, coords, bc_val, bc_type)
639 type(box_t), intent(in) :: box
640 integer, intent(in) :: nb
641 integer, intent(in) :: iv
642 real(dp), intent(in) :: coords(ndim, box%n_cell**(ndim-1))
643 real(dp), intent(out) :: bc_val(box%n_cell**(ndim-1))
644 integer, intent(out) :: bc_type
645
646 bc_type = af_bc_neumann
647 bc_val = 0.0_dp
648 end subroutine field_bc_all_neumann
649
650 !> All Dirichlet zero boundary conditions for the potential
651 subroutine field_bc_all_dirichlet(box, nb, iv, coords, bc_val, bc_type)
652 type(box_t), intent(in) :: box
653 integer, intent(in) :: nb
654 integer, intent(in) :: iv
655 real(dp), intent(in) :: coords(ndim, box%n_cell**(ndim-1))
656 real(dp), intent(out) :: bc_val(box%n_cell**(ndim-1))
657 integer, intent(out) :: bc_type
658
659 if (st_cylindrical) then
660 if (nb == af_neighb_lowx) then
661 bc_type = af_bc_neumann
662 else
663 bc_type = af_bc_dirichlet
664 end if
665 bc_val = 0.0_dp
666 else
667 bc_type = af_bc_dirichlet
668 bc_val = 0.0_dp
669 end if
670 end subroutine field_bc_all_dirichlet
671
672 ! This routine sets the level set function for a simple rod
673 subroutine set_lsf_box(box, iv)
674 type(box_t), intent(inout) :: box
675 integer, intent(in) :: iv
676 integer :: ijk, nc
677 real(dp) :: rr(ndim)
678
679 nc = box%n_cell
680 do kji_do(0,nc+1)
681 rr = af_r_cc(box, [ijk])
682 box%cc(ijk, iv) = mg%lsf(rr)
683 end do; close_do
684 end subroutine set_lsf_box
685
686 real(dp) function sphere_lsf(r)
687 real(dp), intent(in) :: r(ndim)
688 sphere_lsf = norm2(r - rod_r0) - rod_radius
689 end function sphere_lsf
690
691 real(dp) function rod_lsf(r)
692 use m_geometry
693 real(dp), intent(in) :: r(ndim)
694 rod_lsf = gm_dist_line(r, rod_r0, rod_r1, ndim) - rod_radius
695 end function rod_lsf
696
697 !> Compute several parameters for a conical rod
698 subroutine get_conical_rod_properties(r0, r1, &
699 rod_radius, tip_radius, cone_tip_center, cone_tip_r_curvature)
700 real(dp), intent(in) :: r0(ndim) !< Beginning of rod
701 real(dp), intent(in) :: r1(ndim) !< End of conical part
702 real(dp), intent(in) :: rod_radius !< Radius of rod
703 real(dp), intent(in) :: tip_radius !< Radius of curvature of tip
704 real(dp), intent(out) :: cone_tip_center(ndim)
705 real(dp), intent(out) :: cone_tip_r_curvature
706 real(dp) :: cone_angle, cone_length
707
708 ! Determine (half) the opening angle of the top cone, which goes from
709 ! rod_radius to tip_radius over cone_length
710 cone_length = cone_length_frac * norm2(r1 - r0)
711 cone_angle = atan((rod_radius - tip_radius) / cone_length)
712
713 ! We have a point on a sphere with coordinates of the form (R*cos(a),
714 ! R*sin(a)) = (tip_radius, y), so we can get R and subtract R sin(a) to
715 ! obtain the center of the sphere
716 cone_tip_r_curvature = tip_radius/cos(cone_angle)
717 cone_tip_center = r1 - sin(cone_angle) * &
718 cone_tip_r_curvature * (r1 - r0)/norm2(r1 - r0)
719 end subroutine get_conical_rod_properties
720
721 !> Helper function to compute lsf for a rod with a conical tip
722 pure function conical_rod_lsf_arg(r, r0, r1, cone_tip_center, &
723 rod_radius, tip_radius, cone_length_frac, r_curvature) result (lsf)
724 use m_geometry
725 real(dp), intent(in) :: r(ndim)
726 real(dp), intent(in) :: r0(ndim), r1(ndim), cone_tip_center(ndim)
727 real(dp), intent(in) :: rod_radius, tip_radius, cone_length_frac
728 real(dp), intent(in) :: r_curvature
729 real(dp) :: lsf
730 real(dp) :: dist_vec(ndim), frac, radius_at_height, tmp
731
732 ! Project onto line from r0 to r1
733 call gm_dist_vec_line(r, r0, r1, ndim, dist_vec, frac)
734
735 if (frac <= 1 - cone_length_frac) then
736 ! Cylindrical part
737 lsf = norm2(dist_vec) - rod_radius
738 else if (frac < 1.0_dp) then
739 ! Conical part
740 tmp = (1 - frac) / cone_length_frac ! between 0 and 1
741 radius_at_height = tip_radius + tmp * (rod_radius - tip_radius)
742 lsf = norm2(dist_vec) - radius_at_height
743 else
744 ! Spherical tip
745 lsf = norm2(r - cone_tip_center) - r_curvature
746 end if
747 end function conical_rod_lsf_arg
748
749 real(dp) function conical_rod_lsf(r)
750 real(dp), intent(in) :: r(ndim)
751 conical_rod_lsf = conical_rod_lsf_arg(r, rod_r0, rod_r1, cone_tip_center, &
752 rod_radius, cone_tip_radius, cone_length_frac, cone_tip_r_curvature)
753 end function conical_rod_lsf
754
755 !> Get lsf for two conical rods
756 real(dp) function two_conical_rods_lsf(r)
757 real(dp), intent(in) :: r(ndim)
758 real(dp) :: lsf_1, lsf_2
759
760 lsf_1 = conical_rod_lsf_arg(r, rod_r0, rod_r1, cone_tip_center, &
761 rod_radius, cone_tip_radius, cone_length_frac, cone_tip_r_curvature)
762 lsf_2 = conical_rod_lsf_arg(r, rod2_r0, rod2_r1, cone2_tip_center, &
763 rod2_radius, cone2_tip_radius, cone2_length_frac, cone2_tip_r_curvature)
764 two_conical_rods_lsf = min(lsf_1, lsf_2)
765 end function two_conical_rods_lsf
766
767 real(dp) function two_conical_rods_get_potential(r) result(phi)
768 real(dp), intent(in) :: r(ndim)
769 real(dp) :: lsf_1, lsf_2
770
771 lsf_1 = conical_rod_lsf_arg(r, rod_r0, rod_r1, cone_tip_center, &
772 rod_radius, cone_tip_radius, cone_length_frac, cone_tip_r_curvature)
773 lsf_2 = conical_rod_lsf_arg(r, rod2_r0, rod2_r1, cone2_tip_center, &
774 rod2_radius, cone2_tip_radius, cone2_length_frac, cone2_tip_r_curvature)
775
776 if (lsf_1 < lsf_2) then
777 ! Closer to electrode 1
778 if (field_electrode_grounded) then
779 phi = 0.0_dp
780 else
781 phi = current_voltage
782 end if
783 else
784 if (field_electrode2_grounded) then
785 phi = 0.0_dp
786 else
787 phi = current_voltage
788 end if
789 end if
790 end function two_conical_rods_get_potential
791
792 !> Get level set function for case of two rods
793 real(dp) function rod_rod_lsf(r)
794 use m_geometry
795 real(dp), intent(in) :: r(ndim)
796
797 rod_rod_lsf = min(gm_dist_line(r, rod_r0, rod_r1, ndim) - rod_radius, &
798 gm_dist_line(r, rod2_r0, rod2_r1, ndim) - rod2_radius)
799 end function rod_rod_lsf
800
801 !> Get potential to apply at electrode when there are two rods
802 function rod_rod_get_potential(r) result(phi)
803 use m_geometry
804 real(dp), intent(in) :: r(ndim)
805 real(dp) :: phi, lsf_1, lsf_2
806
807 ! Determine distance to electrodes
808 lsf_1 = gm_dist_line(r, rod_r0, rod_r1, ndim) - rod_radius
809 lsf_2 = gm_dist_line(r, rod2_r0, rod2_r1, ndim) - rod2_radius
810
811 if (lsf_1 < lsf_2) then
812 ! Closer to electrode 1
813 if (field_electrode_grounded) then
814 phi = 0.0_dp
815 else
816 phi = current_voltage
817 end if
818 else
819 if (field_electrode2_grounded) then
820 phi = 0.0_dp
821 else
822 phi = current_voltage
823 end if
824 end if
825 end function rod_rod_get_potential
826
827 !> Get level set function for case of sphere and rod
828 real(dp) function sphere_rod_lsf(r)
829 use m_geometry
830 real(dp), intent(in) :: r(ndim)
831
832 sphere_rod_lsf = min(sphere_lsf(r), &
833 gm_dist_line(r, rod2_r0, rod2_r1, ndim) - rod2_radius)
834 end function sphere_rod_lsf
835
836 !> Get potential to apply at electrode for sphere-rod case
837 function sphere_rod_get_potential(r) result(phi)
838 use m_geometry
839 real(dp), intent(in) :: r(ndim)
840 real(dp) :: phi, lsf_1, lsf_2
841
842 ! Determine distance to electrodes
843 lsf_1 = sphere_lsf(r)
844 lsf_2 = gm_dist_line(r, rod2_r0, rod2_r1, ndim) - rod2_radius
845
846 if (lsf_1 < lsf_2) then
847 ! Closer to electrode 1
848 if (field_electrode_grounded) then
849 phi = 0.0_dp
850 else
851 phi = current_voltage
852 end if
853 else
854 if (field_electrode2_grounded) then
855 phi = 0.0_dp
856 else
857 phi = current_voltage
858 end if
859 end if
860 end function sphere_rod_get_potential
861
862 subroutine coaxial_get_lsf(r, lsf_1, lsf_2)
863 real(dp), intent(in) :: r(ndim)
864 real(dp), intent(out) :: lsf_1, lsf_2
865 real(dp) :: domain_center(ndim)
866
867#if NDIM == 1
868 lsf_1 = 0; lsf_2 = 0; domain_center = 0 ! prevent warnings
869 error stop "coaxial not supported in 1D"
870#else
871 if (st_cylindrical) then
872 domain_center(1) = 0
873 domain_center(2) = st_domain_origin(2) + 0.5_dp * st_domain_len(2)
874 lsf_1 = norm2(r(1:2) - domain_center(1:2)) - rod_radius
875 lsf_2 = rod2_radius - norm2(r(1:2) - domain_center(1:2))
876 else if (ndim == 2 .or. ndim == 3) then
877 domain_center = st_domain_origin + 0.5_dp * st_domain_len
878 ! Axis is parallel to last coordinate
879 lsf_1 = norm2(r(1:2) - domain_center(1:2)) - rod_radius
880 lsf_2 = rod2_radius - norm2(r(1:2) - domain_center(1:2))
881 end if
882#endif
883 end subroutine coaxial_get_lsf
884
885 real(dp) function coaxial_lsf(r)
886 real(dp), intent(in) :: r(ndim)
887 real(dp) :: lsf_1, lsf_2
888 call coaxial_get_lsf(r, lsf_1, lsf_2)
889 coaxial_lsf = min(lsf_1, lsf_2)
890 end function coaxial_lsf
891
892 !> Get potential to apply at electrode for sphere-rod case
893 function coaxial_get_potential(r) result(phi)
894 real(dp), intent(in) :: r(ndim)
895 real(dp) :: phi, lsf_1, lsf_2
896
897 call coaxial_get_lsf(r, lsf_1, lsf_2)
898
899 if (lsf_1 < lsf_2) then
900 phi = current_voltage
901 else
902 phi = 0.0_dp
903 end if
904 end function coaxial_get_potential
905
906 !> Compute total field energy in Joule, defined as the volume integral over
907 !> 1/2 * epsilon * E^2
908 subroutine field_compute_energy(tree, field_energy)
909 type(af_t), intent(in) :: tree
910 real(dp), intent(out) :: field_energy
911
912 call af_reduction(tree, field_energy_box, reduce_sum, 0.0_dp, field_energy)
913 end subroutine field_compute_energy
914
915 !> Get the electrostatic field energy in a box
916 real(dp) function field_energy_box(box)
918 type(box_t), intent(in) :: box
919#if NDIM == 2
920 integer :: i
921 real(dp), parameter :: twopi = 2 * acos(-1.0_dp)
922#endif
923 real(dp) :: w(dtimes(box%n_cell))
924 integer :: nc
925
926 nc = box%n_cell
927
928 if (st_use_dielectric) then
929 w = 0.5_dp * uc_eps0 * box%cc(dtimes(1:nc), i_eps) * product(box%dr)
930 else
931 w = 0.5_dp * uc_eps0 * product(box%dr)
932 end if
933
934#if NDIM == 2
935 if (box%coord_t == af_cyl) then
936 ! Weight by 2 * pi * r
937 do i = 1, nc
938 w(i, :) = w(i, :) * twopi * af_cyl_radius_cc(box, i)
939 end do
940 end if
941#endif
942
943 field_energy_box = sum(w * box%cc(dtimes(1:nc), i_electric_fld)**2)
944 end function field_energy_box
945
946 real(dp) function reduce_sum(a, b)
947 real(dp), intent(in) :: a, b
948 reduce_sum = a + b
949 end function reduce_sum
950
951 function field_get_e_vector(box) result(E_vector)
952 type(box_t), intent(in) :: box
953 real(dp) :: e_vector(dtimes(1:box%n_cell), ndim)
954 integer :: nc
955
956 nc = box%n_cell
957
958#if NDIM == 1
959 e_vector(dtimes(:), 1) = 0.5_dp * (box%fc(1:nc, 1, electric_fld) + &
960 box%fc(2:nc+1, 1, electric_fld))
961#elif NDIM == 2
962 e_vector(dtimes(:), 1) = 0.5_dp * (box%fc(1:nc, 1:nc, 1, electric_fld) + &
963 box%fc(2:nc+1, 1:nc, 1, electric_fld))
964 e_vector(dtimes(:), 2) = 0.5_dp * (box%fc(1:nc, 1:nc, 2, electric_fld) + &
965 box%fc(1:nc, 2:nc+1, 2, electric_fld))
966#elif NDIM == 3
967 e_vector(dtimes(:), 1) = 0.5_dp * (box%fc(1:nc, 1:nc, 1:nc, 1, electric_fld) + &
968 box%fc(2:nc+1, 1:nc, 1:nc, 1, electric_fld))
969 e_vector(dtimes(:), 2) = 0.5_dp * (box%fc(1:nc, 1:nc, 1:nc, 2, electric_fld) + &
970 box%fc(1:nc, 2:nc+1, 1:nc, 2, electric_fld))
971 e_vector(dtimes(:), 3) = 0.5_dp * (box%fc(1:nc, 1:nc, 1:nc, 3, electric_fld) + &
972 box%fc(1:nc, 1:nc, 2:nc+1, 3, electric_fld))
973#endif
974 end function field_get_e_vector
975
976end module m_field
Module for handling chemical reactions.
integer, dimension(:), allocatable, public, protected charged_species_itree
List with indices of charged species.
integer, dimension(:), allocatable, public, protected charged_species_charge
List with charges of charged species.
Module with settings and routines to handle dielectrics.
integer, parameter, public i_surf_dens
type(surfaces_t), public diel
To store dielectric surface.
Module to compute electric fields.
Definition m_field.f90:2
subroutine, public field_set_rhs(tree, s_in)
Definition m_field.f90:407
real(dp), public, protected field_pulse_period
Time of one complete voltage pulse.
Definition m_field.f90:33
subroutine, public field_bc_homogeneous(box, nb, iv, coords, bc_val, bc_type)
Dirichlet boundary conditions for the potential in the last dimension, Neumann zero boundary conditio...
Definition m_field.f90:591
real(dp) function, dimension(dtimes(1:box%n_cell), ndim), public field_get_e_vector(box)
Definition m_field.f90:952
subroutine, public field_from_potential(tree, mg)
Compute field from potential.
Definition m_field.f90:532
subroutine, public field_compute_energy(tree, field_energy)
Compute total field energy in Joule, defined as the volume integral over 1/2 * epsilon * E^2.
Definition m_field.f90:909
subroutine, public field_set_voltage(tree, time)
Set the current voltage.
Definition m_field.f90:552
real(dp), public, protected current_voltage
The current applied voltage.
Definition m_field.f90:90
subroutine, public field_initialize(tree, cfg, mg)
Initialize this module.
Definition m_field.f90:108
subroutine, public field_compute(tree, mg, s_in, time, have_guess)
Compute electric field on the tree. First perform multigrid to get electric potential,...
Definition m_field.f90:449
Module that provides routines for geometric operations and calculations. Methods and types have a pre...
Definition m_geometry.f90:3
real(dp) function, public gm_dist_line(r, r0, r1, n_dim)
pure subroutine, public gm_dist_vec_line(r, r0, r1, n_dim, dist_vec, frac)
Compute distance vector between point and its projection onto a line between r0 and r1.
This module contains several pre-defined variables like:
Definition m_streamer.f90:6
logical, public, protected st_cylindrical
Whether cylindrical coordinates are used.
integer, public, protected i_eps
Index can be set to include a dielectric.
type(mg_t), public mg
Multigrid option structure.
logical, dimension(ndim), public, protected st_periodic
Whether the domain is periodic (per dimension)
logical, public, protected st_use_electrode
Whether to include an electrode.
integer, public, protected i_lsf
Index can be set to include an electrode.
logical, public, protected st_use_dielectric
Whether a dielectric is used.
integer, public, protected electric_fld
Index of electric field vector.
integer, public, protected i_rhs
Index of source term Poisson.
integer, public, protected i_phi
Index of electrical potential.
integer, public, protected st_multigrid_num_vcycles
Number of V-cycles to perform per time step.
integer, public, protected i_electric_fld
Index of electric field norm.
real(dp), dimension(ndim), public, protected st_domain_len
Domain length per dimension.
real(dp), dimension(ndim), public, protected st_domain_origin
Origin of domain.
integer, public, protected i_tmp
Index of temporary variable.
real(dp), public, protected st_multigrid_max_rel_residual
Module with settings and routines for tabulated data.
subroutine, public table_from_file(file_name, data_name, x_data, y_data)
Routine to read in tabulated data from a file.
Module with basic types.
Definition m_types.f90:2
real(dp), parameter huge_real
Huge number.
Definition m_types.f90:16
character(len= *), parameter undefined_str
Undefined string.
Definition m_types.f90:10
real(dp), parameter undefined_real
Undefined number.
Definition m_types.f90:13
Module that contains physical and numerical constants.
real(dp), parameter uc_elem_charge
real(dp), parameter uc_eps0
This module contains all the methods that users can customize.
procedure(mg_func_lsf), pointer user_lsf_bc
Function to get boundary value for level set function.
procedure(field_func), pointer user_field_amplitude
To set the field amplitude manually.
procedure(af_subr_bc), pointer user_potential_bc
To set custom boundary conditions for the electric potential.
procedure(mg_func_lsf), pointer user_lsf
Custom level-set function to define an electrode.