54 integer,
parameter,
public :: &
55 EXTERNAL_CURRENT_PARSER = 0, &
62 type(states_mxll_t),
intent(inout) :: st
63 class(space_t),
intent(in) :: space
64 class(mesh_t),
intent(in) :: mesh
65 real(real64),
intent(in) :: time
66 complex(real64),
optional,
intent(inout) :: rs_current_density_ext(:,:)
68 real(real64),
allocatable :: current(:,:,:)
72 safe_allocate(current(1:mesh%np, 1:space%dim, 1))
77 safe_deallocate_a(current)
85 type(states_mxll_t),
intent(inout) :: st
86 class(space_t),
intent(in) :: space
87 class(mesh_t),
intent(in) :: mesh
88 type(namespace_t),
intent(in) :: namespace
91 integer :: ip, il, nlines, ncols, idir, ierr
92 real(real64) :: j_vector(space%dim), dummy(space%dim), xx(1:space%dim), rr, omega
93 character(len=1024) :: tdf_expression, phase_expression
121 if (
parse_block(namespace,
'UserDefinedMaxwellExternalCurrent', blk) == 0)
then
126 st%external_current_number = nlines
127 safe_allocate(st%external_current_modus(1:nlines))
128 safe_allocate(st%external_current_string(1:space%dim, 1:nlines))
129 safe_allocate(st%external_current_amplitude(1:mesh%np, 1:space%dim, 1:nlines))
130 safe_allocate(st%external_current_td_function(1:nlines))
131 safe_allocate(st%external_current_omega(1:nlines))
132 safe_allocate(st%external_current_td_phase(1:nlines))
138 if ((ncols /= 4) .and. (ncols /= 5) .and. (ncols /= 6) .and. (ncols /= 7))
then
139 message(1) =
'Each line in the MaxwellExternalCurrent block must have'
140 message(2) =
'four, five, six or or seven columns.'
146 if (st%external_current_modus(il) == external_current_parser)
then
148 do idir = 1, space%dim
153 do idir = 1, space%dim
157 call mesh_r(mesh, ip, rr, coords = xx)
159 st%external_current_string(idir, il))
161 st%external_current_amplitude(ip, idir, il) = j_vector(idir)
165 st%external_current_omega(il) = omega
167 call tdf_read(st%external_current_td_function(il), namespace, trim(tdf_expression), ierr)
170 call tdf_read(st%external_current_td_phase(il), namespace, trim(phase_expression), ierr)
172 write(
message(1),
'(3A)')
'Error in the "', trim(tdf_expression),
'" field defined in the TDExternalFields block:'
173 write(
message(2),
'(3A)')
'Time-dependent phase function "', trim(phase_expression),
'" not found.'
192 class(
space_t),
intent(in) :: space
193 class(
mesh_t),
intent(in) :: mesh
194 real(real64),
intent(in) :: time
195 real(real64),
intent(inout) :: current(:,:)
197 integer :: ip, jn, idir
198 real(real64) :: xx(space%dim), rr, tt, j_vector(space%dim), dummy(space%dim), amp(space%dim)
199 complex(real64) :: exp_arg
200 real(real64) :: tmp_amp, phase
207 do jn = 1, st%external_current_number
208 if (st%external_current_modus(jn) == external_current_parser)
then
210 call mesh_r(mesh, ip, rr, coords = xx)
211 do idir = 1, space%dim
214 trim(st%external_current_string(idir,jn)))
217 current(ip, :) = current(ip, :) + j_vector
221 exp_arg = st%external_current_omega(jn) * time +
tdf(st%external_current_td_phase(jn),time)
222 phase = real(
exp(-
m_zi*exp_arg), real64)
223 tmp_amp =
tdf(st%external_current_td_function(jn), time)
225 amp = st%external_current_amplitude(ip, :, jn) * tmp_amp
226 current(ip, :) = current(ip, :) + amp(:) * phase
double exp(double __x) __attribute__((__nothrow__
subroutine, public get_rs_density_ext(st, space, mesh, time, rs_current_density_ext)
integer, parameter, public external_current_td_function
subroutine, public external_current_calculation(st, space, mesh, time, current)
subroutine, public external_current_init(st, space, namespace, mesh)
real(real64), parameter, public m_zero
complex(real64), parameter, public m_zi
This module is intended to contain "only mathematical" functions and procedures.
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
subroutine, public messages_warning(no_lines, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
integer function, public parse_block(namespace, name, blk, check_varinfo_)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
subroutine, public build_rs_current_state(current_state, mesh, rs_current_state, ep_field, np)
subroutine, public conv_to_c_string(str)
converts to c string
subroutine, public tdf_init(f)
subroutine, public tdf_read(f, namespace, function_name, ierr)
This function initializes "f" from the TDFunctions block.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
type(unit_t), public unit_one
some special units required for particular quantities
Describes mesh distribution to nodes.