3#include "../afivo/src/cpp_macros.h"
11 integer,
parameter :: scalar_voltage = 1
12 integer,
parameter :: tabulated_voltage = 2
15 integer :: field_given_by = -1
18 real(dp),
allocatable :: field_table_times(:)
21 real(dp),
allocatable :: field_table_values(:)
24 real(dp) :: field_rise_time = 0.0_dp
30 integer :: field_num_pulses = 1
42 logical :: field_electrode_grounded = .false.
45 logical :: field_electrode2_grounded = .false.
48 real(dp) :: rod_r0(ndim) = -1.0e100_dp
51 real(dp) :: rod_r1(ndim) = -1.0e100_dp
54 real(dp) :: rod2_r0(ndim) = -1.0e100_dp
57 real(dp) :: rod2_r1(ndim) = -1.0e100_dp
60 real(dp) :: rod_radius = -1.0e100_dp
63 real(dp) :: rod2_radius = -1.0e100_dp
66 real(dp) :: cone_length_frac = -1.0e100_dp
69 real(dp) :: cone2_length_frac = -1.0e100_dp
72 real(dp) :: cone_tip_radius = -1.0e100_dp
75 real(dp) :: cone2_tip_radius = -1.0e100_dp
80 real(dp) :: cone_tip_center(ndim)
82 real(dp) :: cone2_tip_center(ndim)
85 real(dp) :: cone_tip_r_curvature
87 real(dp) :: cone2_tip_r_curvature
92 character(string_len) :: field_bc_type =
"homogeneous"
111 type(af_t),
intent(inout) :: tree
112 type(cfg_t),
intent(inout) :: cfg
113 type(mg_t),
intent(inout) ::
mg
114 character(len=string_len) :: given_by, user_value
115 character(len=string_len) :: electrode_type
116 integer :: first_blank
119 call cfg_add_get(cfg,
"field_amplitude", field_amplitude, &
120 "The (initial) vertical applied electric field (V/m)")
123 call cfg_add_get(cfg,
"field_given_by", given_by, &
124 "How the electric field or voltage is specified")
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)
132 select case (given_by)
134 field_given_by = scalar_voltage
135 read(user_value, *) field_voltage
137 field_given_by = scalar_voltage
138 read(user_value, *) field_voltage
141 case (
"voltage_table")
142 field_given_by = tabulated_voltage
146 field_table_times, field_table_values)
148 field_given_by = tabulated_voltage
152 field_table_times, field_table_values)
154 field_table_values = -
st_domain_len(ndim) * field_table_values
157 error stop
"field_amplitude not specified"
159 print *,
"Warning: field_amplitude is deprecated, use field_given_by"
160 field_given_by = scalar_voltage
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"
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)")
179 "Time of one complete voltage pulse (s)")
181 if (field_pulse_width <
huge_real .and. field_rise_time <= 0) &
182 error stop
"Set field_rise_time when using field_pulse_width"
184 if (field_num_pulses > 1)
then
186 error stop
"field_num_pulses > 1 requires field_pulse_period"
188 error stop
"field_num_pulses > 1 requires field_pulse_width"
190 error stop
"field_pulse_period shorter than one full pulse"
193 call cfg_add_get(cfg,
"field_bc_type", field_bc_type, &
194 "Boundary condition for electric potential")
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)")
227 electrode_type =
"rod"
228 call cfg_add_get(cfg,
"field_electrode_type", electrode_type, &
229 "Type of electrode (sphere, rod, rod_cone_top, rod_rod, user)")
236 select case (field_bc_type)
240 mg%sides_bc => field_bc_neumann
242 mg%sides_bc => field_bc_all_neumann
244 error stop
"field_bc_select error: invalid condition"
256 select case (electrode_type)
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"
266 call check_general_electrode_parameters()
268 case (
"rod_cone_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"
276 call get_conical_rod_properties(rod_r0, rod_r1, rod_radius, &
277 cone_tip_radius, cone_tip_center, cone_tip_r_curvature)
279 mg%lsf => conical_rod_lsf
282 call check_general_electrode_parameters()
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"
291 mg%lsf => rod_rod_lsf
294 mg%lsf_boundary_function => rod_rod_get_potential
295 case (
"two_rod_cone_electrodes")
298 call check_general_electrode_parameters()
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 if (cone_tip_radius <= 0 .or. cone_tip_radius > rod_radius) &
306 error stop
"cone tip radius should be smaller than rod radius"
307 if (cone2_tip_radius <= 0 .or. cone2_tip_radius > rod2_radius) &
308 error stop
"cone2 tip radius should be smaller than rod2 radius"
309 if (cone_length_frac < 0 .or. cone_length_frac > 1) &
310 error stop
"cone_length_frac not set correctly"
311 if (cone2_length_frac < 0 .or. cone2_length_frac > 1) &
312 error stop
"cone2_length_frac not set correctly"
314 call get_conical_rod_properties(rod_r0, rod_r1, rod_radius, &
315 cone_tip_radius, cone_tip_center, cone_tip_r_curvature)
316 call get_conical_rod_properties(rod2_r0, rod2_r1, rod2_radius, &
317 cone2_tip_radius, cone2_tip_center, cone2_tip_r_curvature)
319 mg%lsf => two_conical_rods_lsf
322 mg%lsf_boundary_function => two_conical_rods_get_potential
324 if (.not.
associated(
user_lsf))
then
325 error stop
"user_lsf not set"
334 print *,
"Electrode types: sphere, rod, rod_cone_top, rod_rod, user"
335 error stop
"Invalid electrode type"
338 call af_set_cc_methods(tree,
i_lsf, funcval=set_lsf_box)
339 tree%mg_i_lsf =
i_lsf
341 mg%lsf_dist => mg_lsf_dist_gss
343 if (rod_radius <= 0)
then
344 error stop
"set field_rod_radius to smallest length scale of electrode"
346 mg%lsf_length_scale = rod_radius
350 af_bc_neumann_zero, af_gc_interp)
354 subroutine check_general_electrode_parameters()
355 if (any(rod_r0 <= -1.0e10_dp)) &
356 error stop
"field_rod_r0 not set correctly"
357 if (any(rod_r1 <= -1.0e10_dp)) &
358 error stop
"field_rod_r1 not set correctly"
359 if (rod_radius <= 0) &
360 error stop
"field_rod_radius not set correctly"
361 end subroutine check_general_electrode_parameters
367 type(af_t),
intent(inout) :: tree
368 integer,
intent(in) :: s_in
371 integer :: lvl, i, id, nc, n, ix
377 do lvl = 1, tree%highest_lvl
379 do i = 1,
size(tree%lvls(lvl)%leaves)
380 id = tree%lvls(lvl)%leaves(i)
382 tree%boxes(id)%cc(dtimes(:),
i_rhs) = 0.0_dp
388 tree%boxes(id)%cc(dtimes(:),
i_rhs) = &
389 tree%boxes(id)%cc(dtimes(:),
i_rhs) + &
390 q * tree%boxes(id)%cc(dtimes(:), ix)
407 type(af_t),
intent(inout) :: tree
408 type(mg_t),
intent(inout) ::
mg
409 integer,
intent(in) :: s_in
410 real(dp),
intent(in) :: time
411 logical,
intent(in) :: have_guess
413 real(dp) :: max_rhs, residual_threshold, conv_fac
414 real(dp) :: residual_ratio
415 integer,
parameter :: max_initial_iterations = 100
416 real(dp),
parameter :: max_residual = 1e8_dp
417 real(dp),
parameter :: min_residual = 1e-6_dp
418 real(dp) :: residuals(max_initial_iterations)
423 call af_tree_maxabs_cc(tree,
mg%i_rhs, max_rhs)
435 residual_threshold = max(min_residual, &
440 if (field_electrode_grounded)
then
441 mg%lsf_boundary_value = 0.0_dp
448 if (.not. have_guess)
then
449 do i = 1, max_initial_iterations
450 call mg_fas_fmg(tree,
mg, .true., .true.)
451 call af_tree_maxabs_cc(tree,
mg%i_tmp, residuals(i))
453 if (residuals(i) < residual_threshold)
then
458 residual_ratio = minval(residuals(i-2:i)) / &
459 maxval(residuals(i-2:i))
460 if (residual_ratio < 2.0_dp .and. residual_ratio > 0.5_dp &
461 .and. residuals(i) < max_residual)
exit
466 if (i == max_initial_iterations + 1)
then
467 print *,
"Iteration residual"
468 do i = 1, max_initial_iterations
469 write(*,
"(I4,E18.2)") i, residuals(i)
471 print *,
"Maybe increase multigrid_max_rel_residual?"
472 error stop
"No convergence in initial field computation"
478 call mg_fas_vcycle(tree,
mg, .true.)
479 call af_tree_maxabs_cc(tree,
mg%i_tmp, residuals(i))
480 if (residuals(i) < residual_threshold)
exit
491 type(af_t),
intent(inout) :: tree
492 type(mg_t),
intent(in) ::
mg
512 type(af_t),
intent(in) :: tree
513 real(dp),
intent(in) :: time
522 select case (field_given_by)
523 case (scalar_voltage)
529 if (t < field_rise_time)
then
531 else if (t < field_pulse_width + field_rise_time)
then
534 tmp = t - (field_pulse_width + field_rise_time)
536 (1 - tmp/field_rise_time))
539 case (tabulated_voltage)
540 call lt_lin_interp_list(field_table_times, field_table_values, &
548 type(box_t),
intent(in) :: box
549 integer,
intent(in) :: nb
550 integer,
intent(in) :: iv
551 real(dp),
intent(in) :: coords(ndim, box%n_cell**(ndim-1))
552 real(dp),
intent(out) :: bc_val(box%n_cell**(ndim-1))
553 integer,
intent(out) :: bc_type
555 if (af_neighb_dim(nb) == ndim)
then
556 if (af_neighb_low(nb))
then
557 bc_type = af_bc_dirichlet
560 bc_type = af_bc_dirichlet
564 bc_type = af_bc_neumann
572 subroutine field_bc_neumann(box, nb, iv, coords, bc_val, bc_type)
573 type(box_t),
intent(in) :: box
574 integer,
intent(in) :: nb
575 integer,
intent(in) :: iv
576 real(dp),
intent(in) :: coords(ndim, box%n_cell**(ndim-1))
577 real(dp),
intent(out) :: bc_val(box%n_cell**(ndim-1))
578 integer,
intent(out) :: bc_type
580 if (af_neighb_dim(nb) == ndim)
then
581 if (af_neighb_low(nb))
then
582 bc_type = af_bc_dirichlet
585 bc_type = af_bc_neumann
589 bc_type = af_bc_neumann
592 end subroutine field_bc_neumann
595 subroutine field_bc_all_neumann(box, nb, iv, coords, bc_val, bc_type)
596 type(box_t),
intent(in) :: box
597 integer,
intent(in) :: nb
598 integer,
intent(in) :: iv
599 real(dp),
intent(in) :: coords(ndim, box%n_cell**(ndim-1))
600 real(dp),
intent(out) :: bc_val(box%n_cell**(ndim-1))
601 integer,
intent(out) :: bc_type
603 bc_type = af_bc_neumann
605 end subroutine field_bc_all_neumann
608 subroutine set_lsf_box(box, iv)
609 type(box_t),
intent(inout) :: box
610 integer,
intent(in) :: iv
616 rr = af_r_cc(box, [ijk])
617 box%cc(ijk, iv) =
mg%lsf(rr)
619 end subroutine set_lsf_box
621 real(dp) function sphere_lsf(r)
622 real(dp),
intent(in) :: r(ndim)
623 sphere_lsf = norm2(r - rod_r0) - rod_radius
624 end function sphere_lsf
626 real(dp) function rod_lsf(r)
628 real(dp),
intent(in) :: r(ndim)
629 rod_lsf =
gm_dist_line(r, rod_r0, rod_r1, ndim) - rod_radius
633 subroutine get_conical_rod_properties(r0, r1, &
634 rod_radius, tip_radius, cone_tip_center, cone_tip_r_curvature)
635 real(dp),
intent(in) :: r0(ndim)
636 real(dp),
intent(in) :: r1(ndim)
637 real(dp),
intent(in) :: rod_radius
638 real(dp),
intent(in) :: tip_radius
639 real(dp),
intent(out) :: cone_tip_center(ndim)
640 real(dp),
intent(out) :: cone_tip_r_curvature
641 real(dp) :: cone_angle, cone_length
645 cone_length = cone_length_frac * norm2(r1 - r0)
646 cone_angle = atan((rod_radius - tip_radius) / cone_length)
651 cone_tip_r_curvature = tip_radius/cos(cone_angle)
652 cone_tip_center = r1 - sin(cone_angle) * &
653 cone_tip_r_curvature * (r1 - r0)/norm2(r1 - r0)
654 end subroutine get_conical_rod_properties
657 pure function conical_rod_lsf_arg(r, r0, r1, cone_tip_center, &
658 rod_radius, tip_radius, cone_length_frac, r_curvature)
result (lsf)
660 real(dp),
intent(in) :: r(ndim)
661 real(dp),
intent(in) :: r0(ndim), r1(ndim), cone_tip_center(ndim)
662 real(dp),
intent(in) :: rod_radius, tip_radius, cone_length_frac
663 real(dp),
intent(in) :: r_curvature
665 real(dp) :: dist_vec(ndim), frac, radius_at_height, tmp
670 if (frac <= 1 - cone_length_frac)
then
672 lsf = norm2(dist_vec) - rod_radius
673 else if (frac < 1.0_dp)
then
675 tmp = (1 - frac) / cone_length_frac
676 radius_at_height = tip_radius + tmp * (rod_radius - tip_radius)
677 lsf = norm2(dist_vec) - radius_at_height
680 lsf = norm2(r - cone_tip_center) - r_curvature
682 end function conical_rod_lsf_arg
684 real(dp) function conical_rod_lsf(r)
685 real(dp),
intent(in) :: r(ndim)
686 conical_rod_lsf = conical_rod_lsf_arg(r, rod_r0, rod_r1, cone_tip_center, &
687 rod_radius, cone_tip_radius, cone_length_frac, cone_tip_r_curvature)
688 end function conical_rod_lsf
691 real(dp) function two_conical_rods_lsf(r)
692 real(dp),
intent(in) :: r(ndim)
693 real(dp) :: lsf_1, lsf_2
695 lsf_1 = conical_rod_lsf_arg(r, rod_r0, rod_r1, cone_tip_center, &
696 rod_radius, cone_tip_radius, cone_length_frac, cone_tip_r_curvature)
697 lsf_2 = conical_rod_lsf_arg(r, rod2_r0, rod2_r1, cone2_tip_center, &
698 rod2_radius, cone2_tip_radius, cone2_length_frac, cone2_tip_r_curvature)
699 two_conical_rods_lsf = min(lsf_1, lsf_2)
700 end function two_conical_rods_lsf
702 real(dp) function two_conical_rods_get_potential(r) result(phi)
703 real(dp),
intent(in) :: r(ndim)
704 real(dp) :: lsf_1, lsf_2
706 lsf_1 = conical_rod_lsf_arg(r, rod_r0, rod_r1, cone_tip_center, &
707 rod_radius, cone_tip_radius, cone_length_frac, cone_tip_r_curvature)
708 lsf_2 = conical_rod_lsf_arg(r, rod2_r0, rod2_r1, cone2_tip_center, &
709 rod2_radius, cone2_tip_radius, cone2_length_frac, cone2_tip_r_curvature)
711 if (lsf_1 < lsf_2)
then
713 if (field_electrode_grounded)
then
719 if (field_electrode2_grounded)
then
725 end function two_conical_rods_get_potential
728 real(dp) function rod_rod_lsf(r)
730 real(dp),
intent(in) :: r(ndim)
732 rod_rod_lsf = min(
gm_dist_line(r, rod_r0, rod_r1, ndim) - rod_radius, &
734 end function rod_rod_lsf
737 function rod_rod_get_potential(r)
result(phi)
739 real(dp),
intent(in) :: r(ndim)
740 real(dp) :: phi, lsf_1, lsf_2
743 lsf_1 =
gm_dist_line(r, rod_r0, rod_r1, ndim) - rod_radius
744 lsf_2 =
gm_dist_line(r, rod2_r0, rod2_r1, ndim) - rod2_radius
746 if (lsf_1 < lsf_2)
then
748 if (field_electrode_grounded)
then
754 if (field_electrode2_grounded)
then
760 end function rod_rod_get_potential
765 type(af_t),
intent(in) :: tree
766 real(dp),
intent(out) :: field_energy
768 call af_reduction(tree, field_energy_box, reduce_sum, 0.0_dp, field_energy)
772 real(dp) function field_energy_box(box)
774 type(box_t),
intent(in) :: box
777 real(dp),
parameter :: twopi = 2 * acos(-1.0_dp)
779 real(dp) :: w(dtimes(box%n_cell))
784 if (st_use_dielectric)
then
785 w = 0.5_dp *
uc_eps0 * box%cc(dtimes(1:nc), i_eps) * product(box%dr)
787 w = 0.5_dp *
uc_eps0 * product(box%dr)
791 if (box%coord_t == af_cyl)
then
794 w(i, :) = w(i, :) * twopi * af_cyl_radius_cc(box, i)
799 field_energy_box = sum(w * box%cc(dtimes(1:nc), i_electric_fld)**2)
800 end function field_energy_box
802 real(dp) function reduce_sum(a, b)
803 real(dp),
intent(in) :: a, b
805 end function reduce_sum
808 type(box_t),
intent(in) :: box
809 real(dp) :: e_vector(dtimes(1:box%n_cell), ndim)
815 e_vector(dtimes(:), 1) = 0.5_dp * (box%fc(1:nc, 1, electric_fld) + &
816 box%fc(2:nc+1, 1, electric_fld))
818 e_vector(dtimes(:), 1) = 0.5_dp * (box%fc(1:nc, 1:nc, 1, electric_fld) + &
819 box%fc(2:nc+1, 1:nc, 1, electric_fld))
820 e_vector(dtimes(:), 2) = 0.5_dp * (box%fc(1:nc, 1:nc, 2, electric_fld) + &
821 box%fc(1:nc, 2:nc+1, 2, electric_fld))
823 e_vector(dtimes(:), 1) = 0.5_dp * (box%fc(1:nc, 1:nc, 1:nc, 1, electric_fld) + &
824 box%fc(2:nc+1, 1:nc, 1:nc, 1, electric_fld))
825 e_vector(dtimes(:), 2) = 0.5_dp * (box%fc(1:nc, 1:nc, 1:nc, 2, electric_fld) + &
826 box%fc(1:nc, 2:nc+1, 1:nc, 2, electric_fld))
827 e_vector(dtimes(:), 3) = 0.5_dp * (box%fc(1:nc, 1:nc, 1:nc, 3, electric_fld) + &
828 box%fc(1:nc, 1:nc, 2:nc+1, 3, electric_fld))
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.
subroutine, public field_set_rhs(tree, s_in)
real(dp), public, protected field_pulse_period
Time of one complete voltage pulse.
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...
real(dp) function, dimension(dtimes(1:box%n_cell), ndim), public field_get_e_vector(box)
subroutine, public field_from_potential(tree, mg)
Compute field from potential.
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.
subroutine, public field_set_voltage(tree, time)
Set the current voltage.
real(dp), public, protected current_voltage
The current applied voltage.
subroutine, public field_initialize(tree, cfg, mg)
Initialize this module.
subroutine, public field_compute(tree, mg, s_in, time, have_guess)
Compute electric field on the tree. First perform multigrid to get electric potential,...
Module that provides routines for geometric operations and calculations. Methods and types have a pre...
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:
integer, public, protected i_eps
Index can be set to include a dielectric.
type(mg_t), public mg
Multigrid option structure.
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.
real(dp), parameter huge_real
Huge number.
character(len= *), parameter undefined_str
Undefined string.
real(dp), parameter undefined_real
Undefined number.
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.