27 use,
intrinsic :: iso_fortran_env
51 integer :: coupling_mode
52 integer :: coupling_terms
53 integer :: dipole_field
54 real(real64),
allocatable :: e_field(:,:)
55 real(real64),
allocatable :: b_field(:,:)
56 real(real64),
allocatable :: vec_pot(:,:)
57 real(real64),
allocatable :: e_field_dip(:)
58 real(real64),
allocatable :: vec_pot_dip(:)
59 real(real64),
allocatable :: b_field_dip(:)
60 real(real64),
allocatable :: e_quadrupole_pot(:)
61 logical :: calc_field_dip = .false.
62 logical :: add_electric_dip = .false.
63 logical :: add_electric_quad = .false.
64 logical :: add_magnetic_dip = .false.
65 logical :: add_zeeman = .false.
66 logical :: add_non_uniform_vec_pot = .false.
67 real(real64) :: center_of_mass(1:3)
68 integer :: center_of_mass_ip
69 integer :: center_of_mass_rankmin
70 logical :: test_equad = .false.
71 type(derivatives_t),
pointer,
private :: der
72 real(real64),
private :: mass
76 integer,
public,
parameter :: &
77 NO_MAXWELL_COUPLING = 0, &
83 integer,
public,
parameter :: &
93 type(mxll_coupling_t),
intent(inout) :: this
94 type(states_elec_dim_t),
intent(in) :: d
95 type(grid_t),
target,
intent(in) :: gr
96 type(namespace_t),
intent(in) :: namespace
97 real(real64),
intent(in) :: mass
99 integer :: terms_default, fmc_default, multipolar_terms, pauli_terms
131 call parse_variable(namespace,
'MaxwellCouplingMode', no_maxwell_coupling, &
156 option__multipolarexpansionterms__electric_dipole + &
157 option__multipolarexpansionterms__electric_quadrupole + &
158 option__multipolarexpansionterms__magnetic_dipole
160 call parse_variable(namespace,
'MultipolarExpansionTerms', terms_default, multipolar_terms)
162 if (
bitand(multipolar_terms, option__multipolarexpansionterms__electric_dipole) /= 0)
then
165 if (
bitand(multipolar_terms, option__multipolarexpansionterms__electric_quadrupole) /=0)
then
166 this%add_electric_quad = .
true.
168 if (
bitand(multipolar_terms, option__multipolarexpansionterms__magnetic_dipole) /=0)
then
169 this%add_magnetic_dip = .
true.
188 if (d%nspin == 2)
then
189 fmc_default = option__paulihamiltonianterms__non_uniform_vector_potential + &
190 option__paulihamiltonianterms__zeeman_term
192 fmc_default = option__paulihamiltonianterms__non_uniform_vector_potential
195 call parse_variable(namespace,
'PauliHamiltonianTerms', fmc_default, pauli_terms)
196 if (
bitand(pauli_terms, option__paulihamiltonianterms__non_uniform_vector_potential) /= 0)
then
197 this%add_non_uniform_vec_pot = .
true.
199 if (
bitand(pauli_terms, option__paulihamiltonianterms__zeeman_term) /=0)
then
200 this%add_zeeman = .
true.
217 call parse_variable(namespace,
'MaxwellDipoleField', dipole_average, &
220 if (this%coupling_mode /= no_maxwell_coupling)
then
221 safe_allocate(this%e_field(1:gr%np, 1:gr%box%dim))
222 safe_allocate(this%b_field(1:gr%np, 1:gr%box%dim))
223 safe_allocate(this%vec_pot(1:gr%np, 1:gr%box%dim))
227 safe_allocate(this%e_field_dip(1:gr%box%dim))
229 safe_allocate(this%b_field_dip(1:gr%box%dim))
231 safe_allocate(this%vec_pot_dip(1:gr%box%dim))
236 safe_allocate(this%e_quadrupole_pot(1:gr%np))
251 type(
grid_t),
intent(in) :: gr
254 integer :: ip, rankmin
266 call parse_variable(namespace,
'MaxwellTestQuadrupole', .false., this%test_equad)
268 if (this%test_equad)
then
269 safe_allocate(this%e_field(1:gr%np, 1:gr%box%dim))
272 this%e_field(ip,1) = 0.02_real64 * gr%x(1, ip)
274 this%calc_field_dip = .
true.
275 this%center_of_mass(1:3) =
m_zero
277 this%center_of_mass_rankmin = rankmin
278 safe_allocate(this%e_quadrupole_pot(1:gr%np))
279 this%e_quadrupole_pot =
m_zero
291 type(
mesh_t),
intent(in) :: mesh
293 class(
space_t),
intent(in) :: space
295 integer :: ip, ispin, idir
300 select case (this%coupling_mode)
304 if (this%add_electric_dip)
then
305 do ispin = 1, d%spin_channels
308 hm_base%potential(ip, ispin) = hm_base%potential(ip, ispin) + &
309 sum(this%e_field_dip(1:space%dim) * (mesh%x(1:space%dim, ip) - this%center_of_mass(1:space%dim)))
314 if (this%add_electric_quad .and. this%calc_field_dip)
then
317 do ispin = 1, d%spin_channels
319 hm_base%potential(ip, ispin) = hm_base%potential(ip, ispin) + this%e_quadrupole_pot(ip)
329 hm_base%uniform_vector_potential(1:space%dim) = hm_base%uniform_vector_potential(1:space%dim) - &
330 this%vec_pot_dip(1:space%dim)
334 if (this%add_non_uniform_vec_pot)
then
335 do idir = 1, mesh%box%dim
336 hm_base%vector_potential(idir, :) = hm_base%vector_potential(idir, :) + this%vec_pot(:, idir)
339 if (this%add_zeeman)
then
341 hm_base%magnetic_field(:,:) = hm_base%magnetic_field(:,:) + this%b_field(:,:)
357 safe_deallocate_a(this%e_field)
358 safe_deallocate_a(this%b_field)
359 safe_deallocate_a(this%vec_pot)
360 safe_deallocate_a(this%e_field_dip)
361 safe_deallocate_a(this%b_field_dip)
362 safe_deallocate_a(this%vec_pot_dip)
363 safe_deallocate_a(this%e_quadrupole_pot)
377 this_out%coupling_mode = this_in%coupling_mode
378 this_out%coupling_terms = this_in%coupling_terms
379 this_out%dipole_field = this_in%dipole_field
380 this_out%calc_field_dip = this_in%calc_field_dip
381 this_out%add_electric_dip = this_in%add_electric_dip
382 this_out%add_electric_quad = this_in%add_electric_quad
383 this_out%add_magnetic_dip = this_in%add_magnetic_dip
384 this_out%add_zeeman = this_in%add_zeeman
385 this_out%add_non_uniform_vec_pot = this_in%add_non_uniform_vec_pot
386 this_out%center_of_mass = this_in%center_of_mass
387 this_out%center_of_mass_ip = this_in%center_of_mass_ip
388 this_out%center_of_mass_rankmin = this_in%center_of_mass_rankmin
389 this_out%test_equad = this_in%test_equad
391 this_out%der => der_target
392 this_out%mass = this_in%mass
394 if (
allocated(this_in%e_field))
then
395 safe_allocate(this_out%e_field(1:
size(this_in%e_field, dim=1), 1:
size(this_in%e_field, dim=2)))
396 this_out%e_field = this_in%e_field
399 if (
allocated(this_in%b_field))
then
400 safe_allocate(this_out%b_field(1:
size(this_in%b_field, dim=1), 1:
size(this_in%b_field, dim=2)))
401 this_out%b_field = this_in%b_field
404 if (
allocated(this_in%vec_pot))
then
405 safe_allocate(this_out%vec_pot(1:
size(this_in%vec_pot, dim=1), 1:
size(this_in%vec_pot, dim=2)))
406 this_out%vec_pot = this_in%vec_pot
409 if (
allocated(this_in%e_field_dip))
then
410 safe_allocate(this_out%e_field_dip(1:
size(this_in%e_field_dip)))
411 this_out%e_field_dip = this_in%e_field_dip
414 if (
allocated(this_in%vec_pot_dip))
then
415 safe_allocate(this_out%vec_pot_dip(1:
size(this_in%vec_pot_dip)))
416 this_out%vec_pot_dip = this_in%vec_pot_dip
419 if (
allocated(this_in%b_field_dip))
then
420 safe_allocate(this_out%b_field_dip(1:
size(this_in%b_field_dip)))
421 this_out%b_field_dip = this_in%b_field_dip
424 if (
allocated(this_in%e_quadrupole_pot))
then
425 safe_allocate(this_out%e_quadrupole_pot(1:
size(this_in%e_quadrupole_pot)))
426 this_out%e_quadrupole_pot = this_in%e_quadrupole_pot
439 type(
mesh_t),
intent(in) :: mesh
441 real(real64),
allocatable :: e_field_quadrupole_tensor(:,:), tmp_partial(:,:), this_e_field(:,:)
442 real(real64) :: r(3), tensor_dot_rr(3)
443 integer :: ip, i, dims
449 safe_allocate(e_field_quadrupole_tensor(1:mesh%box%dim, 1:mesh%box%dim))
450 safe_allocate(this_e_field(1:mesh%np_part, 1:mesh%box%dim))
451 safe_allocate(tmp_partial(1:mesh%np, 1:mesh%box%dim))
453 this_e_field(1:mesh%np,:) = this%e_field(1:mesh%np,:)
455 e_field_quadrupole_tensor =
m_zero
459 if (mesh%mpi_grp%rank == this%center_of_mass_rankmin)
then
460 e_field_quadrupole_tensor(i, 1:dims) = tmp_partial(this%center_of_mass_ip, 1:dims)
463 call mesh%allreduce(e_field_quadrupole_tensor)
466 r(:) = mesh%x(:, ip) - this%center_of_mass(:)
467 tensor_dot_rr = matmul(e_field_quadrupole_tensor, r)
469 this%e_quadrupole_pot(ip) =
m_half * dot_product(r, tensor_dot_rr)
472 safe_deallocate_a(e_field_quadrupole_tensor)
473 safe_deallocate_a(tmp_partial)
474 safe_deallocate_a(this_e_field)
493 real(real64) :: rtmp(3)
494 real(real64),
allocatable :: bxr(:,:)
500 safe_allocate(bxr(1:this%der%mesh%np,1:3))
503 do ip = 1, this%der%mesh%np
504 rtmp(:) = this%der%mesh%x(:, ip) - this%center_of_mass
510 call hpsib%copy_to(bxr_dot_gradb)
511 do idir = 1, this%der%dim
512 call batch_mul(this%der%mesh%np, bxr(:,idir), gradb(idir), bxr_dot_gradb)
516 call bxr_dot_gradb%end()
517 do idir = 1, this%der%dim
518 call gradb(idir)%end()
520 safe_deallocate_a(bxr)
batchified version of the BLAS axpy routine:
This module implements common operations on batches of mesh functions.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
subroutine, public zderivatives_batch_grad(der, ffb, opffb, ghost_update, set_bc, to_cartesian, factor)
apply the gradient to a batch of mesh functions
real(real64), parameter, public m_zero
complex(real64), parameter, public m_zi
real(real64), parameter, public m_half
This module implements the underlying real-space grid.
integer, parameter, public field_uniform_vector_potential
integer, parameter, public field_vector_potential
integer, parameter, public field_magnetic_field
This module is intended to contain "only mathematical" functions and procedures.
pure real(real64) function, dimension(1:3), public dcross_product(a, b)
This module defines the meshes, which are used in Octopus.
integer function, public mesh_nearest_point(mesh, pos, dmin, rankmin)
Returns the index of the point which is nearest to a given vector position pos.
subroutine, public mxll_coupling_init(this, d, gr, namespace, mass)
Parse variables and initialize Maxwell coupling.
subroutine, public set_electric_quadrupole_pot(this, mesh)
Computes the electric quadrupole potential where .
subroutine mxll_quadrupole_test_init(this, gr, namespace)
Initializes quadrupole test when requested. The test applies an electric field defined as E=(0....
integer, parameter, public length_gauge_dipole
subroutine, public mxll_coupling_end(this)
Finalize and deallocate Maxwell coupling arrays.
subroutine, public mxll_coupling_calc(this, hm_base, mesh, d, space)
Add the Maxwell coupling to the electronic Hamiltonian.
subroutine, public mxll_coupling_copy(this_out, this_in, der_target)
Deep-copy Maxwell-electron coupling data and rebind derivatives.
integer, parameter, public velocity_gauge_dipole
integer, parameter, public multipolar_expansion
subroutine, public magnetic_dipole_coupling(this, psib, hpsib)
Computes the magnetic dipole term of the Hamiltonian. This routine acually implements the equivalent...
integer, parameter, public dipole_at_com
integer, parameter, public full_minimal_coupling
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.
This module handles spin dimensions of the states and the k-point distribution.
class representing derivatives
Description of the grid, containing information on derivatives, stencil, and symmetries.
The basic Hamiltonian for electronic system.
Describes mesh distribution to nodes.
class for organizing spins and k-points
batches of electronic states