Octopus
mxll_elec_coupling.F90
Go to the documentation of this file.
1!! Copyright (C) 2023 F. Bonafé
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18#include "global.h"
19
22 use debug_oct_m
24 use global_oct_m
25 use grid_oct_m
27 use, intrinsic :: iso_fortran_env
28 use math_oct_m
29 use mesh_oct_m
32 use parser_oct_m
34 use space_oct_m
37
38 implicit none
39
40 private
41 public :: &
49
50 type :: mxll_coupling_t
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
73 end type mxll_coupling_t
74
75 ! See explanation in the MaxwellCouplingMode variable description.
76 integer, public, parameter :: &
77 NO_MAXWELL_COUPLING = 0, &
82
83 integer, public, parameter :: &
84 DIPOLE_AVERAGE = 0, &
86
87
88contains
89
90 ! ------------------------------------------------------------------------
92 subroutine mxll_coupling_init(this, d, gr, namespace, mass)
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
98
99 integer :: terms_default, fmc_default, multipolar_terms, pauli_terms
100
101 push_sub(mxll_coupling_init)
102
103 this%der => gr%der
104 this%mass = mass
105
106 !%Variable MaxwellCouplingMode
107 !%Type integer
108 !%Default none
109 !%Section Hamiltonian
110 !%Description
111 !% This variable selects the level of light-matter coupling with electromagnetic fields from the Maxwell solver.
112 !% Each option defines a coupling Hamiltonian <math> H_{int}</math>.
113 !%Option none 0
114 !% No coupling.
115 !%Option length_gauge_dipole 1
116 !% Dipole level in length gauge with transverse electric field E(t):
117 !% <math> H_{int} = -r.E </math>
118 !% The option MaxwellDipoleField should be set to define if <math>E</math> is calculated as an average at evaluated
119 !% at the center of mass (by default, dipole coupling with the average electric field will be used).
120 !%Option multipolar_expansion 2
121 !% The option MultipoleExpansionTerms selects the terms to be included (electric dipole, electric quadrupole
122 !% and/or magnetic dipole).
123 !%Option velocity_gauge_dipole 3
124 !% Coupling at dipole level in velocity gauge with transverse vector potential A(t).
125 !% <math> H_{int} = 1/2m (2 A(t).p + A^2(t)) </math>
126 !% The option MaxwellDipoleField should be set.
127 !%Option full_minimal_coupling 4
128 !% Full vector potential A(r,t) coupled to the momentum.
129 !% <math> H_{int} = 1/2m (2 A(r,t).p + A^2(r,t)) </math>
130 !%End
131 call parse_variable(namespace, 'MaxwellCouplingMode', no_maxwell_coupling, &
132 this%coupling_mode)
133 this%add_electric_dip = (this%coupling_mode == length_gauge_dipole)
134
135 if (this%coupling_mode == multipolar_expansion) then
136
137 !%Variable MultipolarExpansionTerms
138 !%Type flag
139 !%Default electric_dipole + electric_quadrupole + magnetic_dipole
140 !%Section Hamiltonian
141 !%Description
142 !% Terms to be included in the multipolar expansion.
143 !% For this type of coupling to make sense, the E field has to be calculated at the center of mass (not averaged).
144 !%Option electric_dipole bit(1)
145 !% Adds electric dipole term in length gauge, uses transverse electric field E(t):
146 !% <math> H_{int} = -e (r.E) </math>
147 !%Option electric_quadrupole bit(2)
148 !% Adds electric quadrupole term:
149 !% <math> H_{int} = \frac{1}{2} e (r . Q . r )] <\math>
150 !% where Q is the outer product of gradient and electric field: <math> Q_{ij} = \partial_i E_j |_{r=r_0} </math>
151 !%Option magnetic_dipole bit(3)
152 !% Adds magnetic dipole term:
153 !% <math> H_{int} = - i (e \hbar /2m) \sum_i (\vec{B}(r_0).(\vec{r} x \nabla)) </math>
154 !%End
155 terms_default = &
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
163 this%add_electric_dip = .true.
164 end if
165 if (bitand(multipolar_terms, option__multipolarexpansionterms__electric_quadrupole) /=0) then
166 this%add_electric_quad = .true.
167 end if
168 if (bitand(multipolar_terms, option__multipolarexpansionterms__magnetic_dipole) /=0) then
169 this%add_magnetic_dip = .true.
170 end if
172 end if
173
174 if (this%coupling_mode == full_minimal_coupling) then
175
176 !%Variable PauliHamiltonianTerms
177 !%Type flag
178 !%Default non_uniform_vector_potential + zeeman_term
179 !%Section Hamiltonian
180 !%Description
181 !% Terms to be included in the Pauli Hamiltonnian.
182 !%Option non_uniform_vector_potential bit(1)
183 !% Adds non-uniform vector potential to the canonical momentum.
184 !%Option zeeman_term bit(2)
185 !% Adds the non-uniform Zeeman potential.
186 !%End
188 if (d%nspin == 2) then
189 fmc_default = option__paulihamiltonianterms__non_uniform_vector_potential + &
190 option__paulihamiltonianterms__zeeman_term
191 else
192 fmc_default = option__paulihamiltonianterms__non_uniform_vector_potential
193 end if
194
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.
198 end if
199 if (bitand(pauli_terms, option__paulihamiltonianterms__zeeman_term) /=0) then
200 this%add_zeeman = .true.
201 end if
202
203 end if
204
205 !%Variable MaxwellDipoleField
206 !%Type integer
207 !%Default average
208 !%Section Hamiltonian
209 !%Description
210 !% This variable selects the method to get the E field vector at the position of the system
211 !% if Maxwell-matter coupling at dipole level within length gauge is done.
212 !%Option average 0
213 !% Transverse E field is averaged in the volume occupied by the matter system.
214 !%Option center_of_mass 1
215 !% Tranverse E field is evaluated at the center of mass of the matter system (at a single point).
216 !%End
217 call parse_variable(namespace, 'MaxwellDipoleField', dipole_average, &
218 this%dipole_field)
219
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))
224 this%e_field = m_zero
225 this%b_field = m_zero
226 this%vec_pot = m_zero
227 safe_allocate(this%e_field_dip(1:gr%box%dim))
228 this%e_field_dip = m_zero
229 safe_allocate(this%b_field_dip(1:gr%box%dim))
230 this%b_field_dip = m_zero
231 safe_allocate(this%vec_pot_dip(1:gr%box%dim))
232 this%vec_pot_dip = m_zero
233 end if
234
235 if (this%coupling_mode == multipolar_expansion .and. this%add_electric_quad) then
236 safe_allocate(this%e_quadrupole_pot(1:gr%np))
237 end if
238
239 call mxll_quadrupole_test_init(this, gr, namespace)
240
241 pop_sub(mxll_coupling_init)
242
243 end subroutine mxll_coupling_init
244
245 ! ------------------------------------------------------------------------
249 subroutine mxll_quadrupole_test_init(this, gr, namespace)
250 type(mxll_coupling_t), intent(inout) :: this
251 type(grid_t), intent(in) :: gr
252 type(namespace_t), intent(in) :: namespace
253
254 integer :: ip, rankmin
255 real(real64) :: dmin
256
258
259 !%Variable MaxwellTestQuadrupole
260 !%Type logical
261 !%Default no
262 !%Section Maxwell
263 !%Description
264 !% Override Maxwell field with linear E field for testing.
265 !%End
266 call parse_variable(namespace, 'MaxwellTestQuadrupole', .false., this%test_equad)
267
268 if (this%test_equad) then
269 safe_allocate(this%e_field(1:gr%np, 1:gr%box%dim))
270 this%e_field = m_zero
271 do ip = 1, gr%np
272 this%e_field(ip,1) = 0.02_real64 * gr%x(1, ip)
273 end do
274 this%calc_field_dip = .true.
275 this%center_of_mass(1:3) = m_zero
276 this%center_of_mass_ip = mesh_nearest_point(gr, this%center_of_mass, dmin, rankmin)
277 this%center_of_mass_rankmin = rankmin
278 safe_allocate(this%e_quadrupole_pot(1:gr%np))
279 this%e_quadrupole_pot = m_zero
280 end if
281
283
284 end subroutine mxll_quadrupole_test_init
285
286 ! ------------------------------------------------------------------------
288 subroutine mxll_coupling_calc(this, hm_base, mesh, d, space)
289 type(mxll_coupling_t), intent(inout) :: this
290 type(hamiltonian_elec_base_t), intent(inout) :: hm_base
291 type(mesh_t), intent(in) :: mesh
292 type(states_elec_dim_t), intent(in) :: d
293 class(space_t), intent(in) :: space
294
295 integer :: ip, ispin, idir
296
297 push_sub(mxll_coupling_calc)
298
299 ! coupling to Maxwell field
300 select case (this%coupling_mode)
301
303
304 if (this%add_electric_dip) then
305 do ispin = 1, d%spin_channels
306 do ip = 1, mesh%np
307 ! The -1 sign is missing here. Check epot.F90 for the explanation.
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)))
310 end do
311 end do
312 end if
313
314 if (this%add_electric_quad .and. this%calc_field_dip) then
315 ! the calc_field_dip condition prevents calling the routine before initializing the interactions
316 call set_electric_quadrupole_pot(this, mesh)
317 do ispin = 1, d%spin_channels
318 do ip = 1, mesh%np
319 hm_base%potential(ip, ispin) = hm_base%potential(ip, ispin) + this%e_quadrupole_pot(ip)
320 end do
321 end do
322 end if
323
324 ! magnetic_dipole term is added in hamiltonian_elec_inc.F90
325
327 call hm_base%allocate_field(mesh, field_uniform_vector_potential, .false.)
328 ! The minus sign is here is for the wrong convention of Octopus for the phase
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)
331
333 call hm_base%allocate_field(mesh, field_vector_potential, .false.)
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)
337 end do
338 end if
339 if (this%add_zeeman) then
340 call hm_base%allocate_field(mesh, field_magnetic_field, .false.)
341 hm_base%magnetic_field(:,:) = hm_base%magnetic_field(:,:) + this%b_field(:,:)
342 end if
343
344 end select
345
346 pop_sub(mxll_coupling_calc)
347
348 end subroutine mxll_coupling_calc
349
350 ! ------------------------------------------------------------------------
352 subroutine mxll_coupling_end(this)
353 type(mxll_coupling_t), intent(inout) :: this
354
355 push_sub(mxll_coupling_end)
356
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)
364
365 pop_sub(mxll_coupling_end)
366 end subroutine mxll_coupling_end
367
368 ! ------------------------------------------------------------------------
370 subroutine mxll_coupling_copy(this_out, this_in, der_target)
371 type(mxll_coupling_t), intent(out) :: this_out
372 type(mxll_coupling_t), intent(in) :: this_in
373 type(derivatives_t), target, intent(in) :: der_target
374
375 push_sub(mxll_coupling_copy)
376
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
390
391 this_out%der => der_target
392 this_out%mass = this_in%mass
393
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
397 end if
398
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
402 end if
403
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
407 end if
408
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
412 end if
413
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
417 end if
418
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
422 end if
423
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
427 end if
428
429 pop_sub(mxll_coupling_copy)
430 end subroutine mxll_coupling_copy
431
432
433 ! ------------------------------------------------------------------------
437 subroutine set_electric_quadrupole_pot(this, mesh)
438 type(mxll_coupling_t), intent(inout) :: this
439 type(mesh_t), intent(in) :: mesh
440
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
444
446
447 call profiling_in('SET_ELECTRIC_QUADRUPOLE')
448
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))
452
453 this_e_field(1:mesh%np,:) = this%e_field(1:mesh%np,:)
454
455 e_field_quadrupole_tensor = m_zero
456 dims = this%der%dim
457 do i = 1, dims
458 call dderivatives_grad(this%der, this_e_field(:, i), tmp_partial)
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)
461 end if
462 end do
463 call mesh%allreduce(e_field_quadrupole_tensor)
464
465 do ip = 1, mesh%np
466 r(:) = mesh%x(:, ip) - this%center_of_mass(:)
467 tensor_dot_rr = matmul(e_field_quadrupole_tensor, r)
468 ! This term is +1/2 r.Q.r and not -1/2 r.Q.r (as e=+1 in Octopus)
469 this%e_quadrupole_pot(ip) = m_half * dot_product(r, tensor_dot_rr)
470 end do
471
472 safe_deallocate_a(e_field_quadrupole_tensor)
473 safe_deallocate_a(tmp_partial)
474 safe_deallocate_a(this_e_field)
475
476 call profiling_out('SET_ELECTRIC_QUADRUPOLE')
478
479 end subroutine set_electric_quadrupole_pot
480
481
482 ! ------------------------------------------------------------------------
487 subroutine magnetic_dipole_coupling(this, psib, hpsib)
488 type(mxll_coupling_t), intent(in) :: this
489 type(wfs_elec_t), intent(inout) :: psib
490 type(wfs_elec_t), intent(inout) :: hpsib
491
492 integer :: ip, idir
493 real(real64) :: rtmp(3)
494 real(real64), allocatable :: bxr(:,:)
495 type(wfs_elec_t) :: gradb(this%der%dim)
496 type(wfs_elec_t) :: bxr_dot_gradb
497
499
500 safe_allocate(bxr(1:this%der%mesh%np,1:3))
501 call profiling_in('MAGNETIC_DIPOLE')
502
503 do ip = 1, this%der%mesh%np
504 rtmp(:) = this%der%mesh%x(:, ip) - this%center_of_mass
505 bxr(ip,:) = dcross_product(this%b_field_dip, rtmp)
506 end do
507
508 call zderivatives_batch_grad(this%der, psib, gradb)
509
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)
513 call batch_axpy(this%der%mesh%np, m_zi*m_half/this%mass, bxr_dot_gradb, hpsib)
514 end do
515
516 call bxr_dot_gradb%end()
517 do idir = 1, this%der%dim
518 call gradb(idir)%end()
519 end do
520 safe_deallocate_a(bxr)
521
522 call profiling_out('MAGNETIC_DIPOLE')
524 end subroutine magnetic_dipole_coupling
525
527
528!! Local Variables:
529!! mode: f90
530!! coding: utf-8
531!! End:
batchified version of the BLAS axpy routine:
Definition: batch_ops.F90:156
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
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
Definition: global.F90:191
complex(real64), parameter, public m_zi
Definition: global.F90:205
real(real64), parameter, public m_half
Definition: global.F90:197
This module implements the underlying real-space grid.
Definition: grid.F90:119
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.
Definition: math.F90:117
pure real(real64) function, dimension(1:3), public dcross_product(a, b)
Definition: math.F90:1941
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
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.
Definition: mesh.F90:386
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.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
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.
Definition: grid.F90:171
The basic Hamiltonian for electronic system.
Describes mesh distribution to nodes.
Definition: mesh.F90:187
class for organizing spins and k-points
batches of electronic states
Definition: wfs_elec.F90:141
int true(void)