Octopus
kick.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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
19#include "global.h"
20
21module kick_oct_m
22 use debug_oct_m
23 use global_oct_m
25 use iso_c_binding
27 use ions_oct_m
28 use, intrinsic :: iso_fortran_env
31 use math_oct_m
32 use mesh_oct_m
35 use parser_oct_m
37 use pcm_oct_m
41 use space_oct_m
46 use unit_oct_m
49
50 implicit none
51
52 private
53 public :: &
54 kick_t, &
55 kick_init, &
56 kick_copy, &
57 kick_end, &
58 kick_read, &
59 kick_write, &
60 kick_apply, &
63
64
65 integer, public, parameter :: &
66 KICK_FUNCTION_DIPOLE = 0, &
69
70 integer, public, parameter :: &
71 KICK_DENSITY_MODE = 0, &
72 kick_spin_mode = 1, &
75
76 integer, public, parameter :: &
77 QKICKMODE_NONE = 0, &
78 qkickmode_exp = 1, &
79 qkickmode_cos = 2, &
80 qkickmode_sin = 3, &
82
83
84 type kick_t
85 ! Components are public by default
86
88 integer :: dim
90 real(real64) :: time
92 integer, private :: delta_strength_mode
93 real(real64) :: delta_strength
95 real(real64), allocatable :: pol(:,:)
96 integer :: pol_dir
97 integer :: pol_equiv_axes
98 real(real64), allocatable :: wprime(:)
99 real(real64) :: easy_axis(3)
111 integer :: n_multipoles
112 integer, allocatable :: l(:), m(:)
113 real(real64), allocatable :: weight(:)
114 integer :: nqmult(3)
115 integer :: nqvec
116 real(real64), allocatable :: qvector(:,:)
117 real(real64) :: trans_vec(3, 2)
118 real(real64) :: qlength
119 integer :: qkick_mode
120 integer :: qbessel_l, qbessel_m
122 integer :: function_mode
123 character(len=200), private:: user_defined_function
124 end type kick_t
125
126contains
127
128 ! ---------------------------------------------------------
129 subroutine kick_init(kick, namespace, space, kpoints, nspin)
130 type(kick_t), intent(inout) :: kick
131 type(namespace_t), intent(in) :: namespace
132 class(space_t), intent(in) :: space
133 type(kpoints_t), intent(in) :: kpoints
134 integer, intent(in) :: nspin
135
136 type(block_t) :: blk
137 integer :: n_rows, irow, idir, iop, iq, iqx, iqy, iqz
138 real(real64) :: norm, dot
139 real(real64) :: qtemp(space%dim)
140 integer :: periodic_dim
141
142 push_sub(kick_init)
143
144 kick%dim = space%dim
145 periodic_dim = space%periodic_dim
146 kick%nqvec = 0
147
148 safe_allocate(kick%pol(1:space%dim, 1:space%dim))
149 safe_allocate(kick%wprime(1:space%dim))
150
151 !%Variable TDDeltaKickTime
152 !%Type float
153 !%Default 0.0
154 !%Section Time-Dependent::Response
155 !%Description
156 !% The delta-perturbation that can be applied by making use of the <tt>TDDeltaStrength</tt> variable,
157 !% can be applied at the time given by this variable. Usually, this time is zero, since one wants
158 !% to apply the delta pertubation or "kick" at a system at equilibrium, and no other time-dependent
159 !% external potential is used. However, one may want to apply a kick on top of a laser field,
160 !% for example.
161 !%End
162 call parse_variable(namespace, 'TDDeltaKickTime', m_zero, kick%time, units_inp%time)
163 if (kick%time > m_zero) then
164 call messages_experimental('TDDeltaKickTime > 0', namespace=namespace)
165 end if
166
167 !%Variable TDDeltaStrength
168 !%Type float
169 !%Default 0
170 !%Section Time-Dependent::Response
171 !%Description
172 !% When no laser is applied, a delta (in time) perturbation with
173 !% strength <tt>TDDeltaStrength</tt> can be applied. This is used to
174 !% calculate, <i>e.g.</i>, the linear optical spectra. If the ions are
175 !% allowed to move, the kick will affect them also.
176 !% The electric field is <math>-(\hbar k / e) \delta(t)</math> for a dipole with
177 !% zero wavevector, where <i>k</i> = <tt>TDDeltaStrength</tt>, which causes
178 !% the wavefunctions instantaneously to acquire a phase <math>e^{ikx}</math>.
179 !% The unit is inverse length.
180 !%End
181 call parse_variable(namespace, 'TDDeltaStrength', m_zero, kick%delta_strength, units_inp%length**(-1))
182
183 kick%function_mode = kick_function_dipole
184
185 if (abs(kick%delta_strength) <= m_epsilon) then
186 kick%delta_strength_mode = 0
187 kick%pol_equiv_axes = 0
188 kick%pol = diagonal_matrix(kick%dim, m_one)
189 kick%pol_dir = 0
190 kick%wprime = m_zero
191 kick%n_multipoles = 0
192 kick%qkick_mode = qkickmode_none
193 kick%easy_axis = m_zero
194 pop_sub(kick_init)
195 return
196 end if
197
198 !%Variable TDDeltaStrengthMode
199 !%Type integer
200 !%Default kick_density
201 !%Section Time-Dependent::Response
202 !%Description
203 !% When calculating the density response via real-time propagation,
204 !% one needs to perform an initial kick on the KS system, at
205 !% time zero. Depending on what kind of response property one wants to obtain,
206 !% this kick may be done in several modes. For use to calculate triplet excitations,
207 !% see MJT Oliveira, A Castro, MAL Marques, and A Rubio, <i>J. Nanoscience and Nanotechnology</i> <b>8</b>, 3392 (2008).
208 !%Option kick_density 0
209 !% The total density of the system is perturbed. This mode is appropriate for
210 !% electric dipole response, as for optical absorption.
211 !%Option kick_spin 1
212 !% The individual spin densities are perturbed oppositely. Note that this mode
213 !% is only possible if the run is done in spin-polarized mode, or with spinors.
214 !% This mode is appropriate for the paramagnetic dipole response, which can couple
215 !% to triplet excitations.
216 !%Option kick_spin_and_density 2
217 !% A combination of the two above. Note that this mode
218 !% is only possible if the run is done in spin-polarized mode, or with spinors.
219 !% This mode is intended for use with symmetries to obtain both of the responses
220 !% at once, at described in the reference above.
221 !%Option kick_magnon 3
222 !% Rotates the magnetization. Only works for spinors.
223 !% Can be used in a supercell or my making use of the generalized Bloch theorem.
224 !% In the later case (see <tt>SpiralBoundaryConditions</tt>) spin-orbit coupling cannot be used.
225 !%End
226 call parse_variable(namespace, 'TDDeltaStrengthMode', kick_density_mode, kick%delta_strength_mode)
227 select case (kick%delta_strength_mode)
228 case (kick_density_mode)
230 case (kick_magnon_mode)
231 if (nspin /= spinors) then
232 call messages_input_error(namespace, 'TDDeltaStrengthMode', 'Magnon kick is incompatible with spinors')
233 end if
234 case default
235 call messages_input_error(namespace, 'TDDeltaStrengthMode', 'Unknown mode')
236 end select
237
238 if (parse_is_defined(namespace, 'TDDeltaUserDefined')) then
239
240 kick%function_mode = kick_function_user_defined
241 kick%n_multipoles = 0
242
243 !%Variable TDDeltaUserDefined
244 !%Type string
245 !%Section Time-Dependent::Response
246 !%Description
247 !% By default, the kick function will be a dipole. This will change if (1) the variable
248 !% <tt>TDDeltaUserDefined</tt> is present in the inp file, or (2) if the block <tt>TDKickFunction</tt>
249 !% is present in the <tt>inp</tt> file. If both are present in the <tt>inp</tt> file, the <tt>TDKickFunction</tt>
250 !% block will be ignored. The value of <tt>TDDeltaUserDefined</tt> should be a string describing
251 !% the function that is going to be used as delta perturbation.
252 !%End
253 call parse_variable(namespace, 'TDDeltaUserDefined', '0', kick%user_defined_function)
254
255 !%Variable TDKickFunction
256 !%Type block
257 !%Section Time-Dependent::Response
258 !%Description
259 !% If the block <tt>TDKickFunction</tt> is present in the input file, and the variable
260 !% <tt>TDDeltaUserDefined</tt> is not present in the input file, the kick function to
261 !% be applied at time zero of the time-propagation will not be a "dipole" function
262 !% (<i>i.e.</i> <math>\phi \rightarrow e^{ikx} \phi</math>, but a general multipole in the form <math>r^l Y_{lm}(r)</math>.
263 !%
264 !% Each line has three columns: integers <i>l</i> and <i>m</i> that defines the
265 !% multipole, and a weight. Any number of lines may be given, and the kick will be the sum of those
266 !% multipoles with the given weights.
267 !%
268 !% This feature allows calculation of quadrupole, octupole, etc., response functions.
269 !%End
270 else if (parse_block(namespace, 'TDKickFunction', blk) == 0) then
271
272 kick%function_mode = kick_function_multipole
273 n_rows = parse_block_n(blk)
274 kick%n_multipoles = n_rows
275 safe_allocate( kick%l(1:n_rows))
276 safe_allocate( kick%m(1:n_rows))
277 safe_allocate(kick%weight(1:n_rows))
278 do irow = 1, n_rows
279 call parse_block_integer(blk, irow - 1, 0, kick%l(irow))
280 call parse_block_integer(blk, irow - 1, 1, kick%m(irow))
281 call parse_block_float(blk, irow - 1, 2, kick%weight(irow))
282 if ((kick%l(irow) < 0) .or. (abs(kick%m(irow)) > abs(kick%l(irow)))) then
283 call messages_input_error(namespace, 'TDkickFunction')
284 end if
285 end do
286
287 else
288
289 kick%function_mode = kick_function_dipole
290 kick%n_multipoles = 0
291
292 ! Find out how many equivalent axes we have...
293 !%Variable TDPolarizationEquivAxes
294 !%Type integer
295 !%Default 0
296 !%Section Time-Dependent::Response::Dipole
297 !%Description
298 !% Defines how many of the <tt>TDPolarization</tt> axes are equivalent. This information is stored in a file and then
299 !% used by <tt>oct-propagation_spectrum</tt> to rebuild the full polarizability tensor from just the
300 !% first <tt>TDPolarizationEquivAxes</tt> directions. This variable is also used by <tt>CalculationMode = vdw</tt>.
301 !%End
302 call parse_variable(namespace, 'TDPolarizationEquivAxes', 0, kick%pol_equiv_axes)
303
304 !%Variable TDPolarizationDirection
305 !%Type integer
306 !%Section Time-Dependent::Response::Dipole
307 !%Description
308 !% When a delta potential is included in a time-dependent run, this
309 !% variable defines in which direction the field will be applied
310 !% by selecting one of the lines of <tt>TDPolarization</tt>. In a
311 !% typical run (without using symmetry), the <tt>TDPolarization</tt> block
312 !% would contain the three Cartesian unit vectors (the default
313 !% value), and one would make 3 runs varying
314 !% <tt>TDPolarization</tt> from 1 to 3.
315 !% If one is using symmetry, <tt>TDPolarization</tt> should run only from 1
316 !% to <tt>TDPolarizationEquivAxes</tt>.
317 !%End
318
319 call parse_variable(namespace, 'TDPolarizationDirection', 0, kick%pol_dir)
320
321 if (kick%delta_strength_mode /= kick_magnon_mode) then
322 if (kick%pol_dir < 1 .or. kick%pol_dir > kick%dim) call messages_input_error(namespace, 'TDPolarizationDirection')
323 end if
324
325 !%Variable TDPolarization
326 !%Type block
327 !%Section Time-Dependent::Response::Dipole
328 !%Description
329 !% The (real) polarization of the delta electric field. Normally
330 !% one needs three perpendicular polarization directions to calculate a
331 !% spectrum (unless symmetry is used).
332 !% The format of the block is:
333 !%
334 !% <tt>%TDPolarization
335 !% <br>&nbsp;&nbsp;pol1x | pol1y | pol1z
336 !% <br>&nbsp;&nbsp;pol2x | pol2y | pol2z
337 !% <br>&nbsp;&nbsp;pol3x | pol3y | pol3z
338 !% <br>%</tt>
339 !%
340 !% <tt>Octopus</tt> uses both this block and the variable
341 !% <tt>TDPolarizationDirection</tt> to determine the polarization
342 !% vector for the run. For example, if
343 !% <tt>TDPolarizationDirection=2</tt> the polarization <tt>(pol2x,
344 !% pol2y, pol2z)</tt> would be used.
345 !% These directions may not be in periodic directions.
346 !%
347 !% The default value for <tt>TDPolarization</tt> is the three
348 !% Cartesian unit vectors (1,0,0), (0,1,0), and (0,0,1).
349 !%
350 !% Note that the directions do not necessarily need to be perpendicular
351 !% when symmetries are used.
352 !%
353 !% WARNING: If you want to obtain the cross-section tensor, the
354 !% <tt>TDPolarization</tt> block must be exactly the same for the run in
355 !% each direction. The direction must be selected by the
356 !% <tt>TDPolarizationDirection</tt> variable.
357 !%
358 !%End
359
360 ! Default basis is the Cartesian unit vectors.
361 ! FIXME: Here the symmetry of the system should be analyzed, and the polarization
362 ! basis built accordingly.
363 kick%pol(1:kick%dim, 1:kick%dim) = diagonal_matrix(kick%dim, m_one)
364 if (parse_block(namespace, 'TDPolarization', blk) == 0) then
365 n_rows = parse_block_n(blk)
366
367 if (n_rows /= kick%dim) then
368 call messages_input_error(namespace, 'TDPolarization', 'There should be exactly one line per dimension')
369 end if
370
371 do irow = 1, n_rows
372 do idir = 1, kick%dim
373 call parse_block_float(blk, irow - 1, idir - 1, kick%pol(idir, irow))
374 end do
375 end do
376 call parse_block_end(blk)
377 end if
378
379 ! Normalize
380 do idir = 1, kick%dim
381 kick%pol(1:kick%dim, idir) = kick%pol(1:kick%dim, idir) / norm2(kick%pol(1:kick%dim, idir))
382 end do
383
384 if (kick%delta_strength_mode /= kick_magnon_mode) then
385 if (any(abs(kick%pol(1:periodic_dim, :)) > m_epsilon)) then
386 message(1) = "Kick cannot be applied in a periodic direction. Use GaugeVectorField instead."
387 call messages_fatal(1, namespace=namespace)
388 end if
389 end if
390
391 !%Variable TDPolarizationWprime
392 !%Type block
393 !%Section Time-Dependent::Response::Dipole
394 !%Description
395 !% This block is needed only when
396 !% <tt>TDPolarizationEquivAxes</tt> is set to 3. In such a case,
397 !% the three directions (<i>pol1</i>, <i>pol2</i>, and <i>pol3</i>) defined in
398 !% the <tt>TDPolarization</tt> block should be related by symmetry
399 !% operations. If <i>A</i> is the symmetry operation that takes you
400 !% from <i>pol1</i> to <i>pol2</i>, then <tt>TDPolarizationWprime</tt>
401 !% should be set to the direction defined by <i>A</i><math>^{-1}</math><i>pol3</i>.
402 !% For more information see MJT Oliveira
403 !% <i>et al.</i>, <i>J. Nanoscience and Nanotechnology</i> <b>8</b>,
404 !% 3392 (2008).
405 !%End
406 if (parse_block(namespace, 'TDPolarizationWprime', blk) == 0) then
407 do idir = 1, kick%dim
408 call parse_block_float(blk, 0, idir - 1, kick%wprime(idir))
409 end do
410 kick%wprime(1:kick%dim) = kick%wprime(1:kick%dim) / norm2(kick%wprime(1:kick%dim))
411 call parse_block_end(blk)
412 else
413 kick%wprime(1:kick%dim-1) = m_zero
414 kick%wprime(kick%dim) = m_one
415 end if
416 end if
417
418 ! for non-dipole, it is more complicated to check whether it is actually in the periodic direction
419 if (periodic_dim > 0 .and. kick%delta_strength_mode /= kick_magnon_mode) then
420 message(1) = "Kicks cannot be applied correctly in periodic directions."
421 call messages_warning(1, namespace=namespace)
422 end if
423
424
425 !%Variable TDReducedMomentumTransfer
426 !%Type block
427 !%Section Time-Dependent::Response
428 !%Description
429 !% The same as TDMomentumTransfer, but the momentum is specified in reduced coordinates.
430 !% Only available for magnon kicks at the moment, and only with an exponential kick.
431 !%End
432 if (periodic_dim > 0 .and. kick%delta_strength_mode == kick_magnon_mode .and. &
433 parse_block(namespace, 'TDReducedMomentumTransfer', blk) == 0) then
434
435 kick%nqvec = 1
436 safe_allocate(kick%qvector(1:kick%dim, 1))
437 do idir = 1, kick%dim
438 call parse_block_float(blk, 0, idir - 1, qtemp(idir))
439 end do
440 call kpoints_to_absolute(kpoints%latt, qtemp(1:kick%dim), kick%qvector(:, 1))
441 kick%qkick_mode = qkickmode_exp
442
443 call parse_block_end(blk)
444
445 if (kpoints%use_symmetries) then
446 do iop = 1, symmetries_number(kpoints%symm)
447 if (iop == symmetries_identity_index(kpoints%symm)) cycle
448 if (.not. symm_op_invariant_cart(kpoints%symm%ops(iop), kick%qvector(:,1), 1e-5_real64)) then
449 message(1) = "The TDMomentumTransfer breaks (at least) one of the symmetries used to reduce the k-points."
450 message(2) = "Set SymmetryBreakDir equal to TDMomemtumTransfer."
451 call messages_fatal(2, namespace=namespace)
452 end if
453 end do
454 end if
455
456 !%Variable TDMomentumTransfer
457 !%Type block
458 !%Section Time-Dependent::Response
459 !%Description
460 !% Momentum-transfer vector for the calculation of the dynamic structure factor.
461 !% When this variable is set, a non-dipole field is applied, and an output file
462 !% <tt>ftchd</tt> is created (it contains the Fourier transform of the charge density
463 !% at each time). The type of the applied external field can be set by
464 !% an optional last number. Possible options are <tt>qexp</tt> (default), <tt>qcos</tt>,
465 !% <tt>qsin</tt>, or <tt>qcos+qsin</tt>. In the formulae below,
466 !% <math>\vec{q}</math> is the momentum-transfer vector.
467 !%Option qexp 1
468 !% External field is <math>e^{i \vec{q} \cdot \vec{r}}</math>.
469 !%Option qcos 2
470 !% External field is <math>\cos \left( i \vec{q} \cdot \vec{r} \right)</math>.
471 !%Option qsin 3
472 !% External field is <math>\sin \left( i \vec{q} \cdot \vec{r} \right)</math>.
473 !%Option qbessel 4
474 !% External field is <math>j_l \left( \vec{q} \cdot \vec{r} \right) Y_{lm} \left(\vec{r} \right)</math>.
475 !% In this case, the block has to include two extra values (<i>l</i> and <i>m</i>).
476 !%End
477
478 else if (parse_block(namespace, 'TDMomentumTransfer', blk) == 0) then
479 kick%nqvec = 1
480 safe_allocate(kick%qvector(1:kick%dim, 1))
481 do idir = 1, kick%dim
482 call parse_block_float(blk, 0, idir - 1, kick%qvector(idir,1))
483 kick%qvector(idir,1) = units_to_atomic(unit_one / units_inp%length, kick%qvector(idir,1))
484 end do
485
486 ! Read the calculation mode (exp, cos, sin, or bessel)
487 if (parse_block_cols(blk, 0) > 3) then
488
489 call parse_block_integer(blk, 0, idir - 1, kick%qkick_mode)
490
491 ! Read l and m if bessel mode (j_l*Y_lm) is used
492 if (kick%qkick_mode == qkickmode_bessel .and. parse_block_cols(blk, 0) == 6) then
493 call parse_block_integer(blk, 0, idir + 0, kick%qbessel_l)
494 call parse_block_integer(blk, 0, idir + 1, kick%qbessel_m)
495 else
496 kick%qbessel_l = 0
497 kick%qbessel_m = 0
498 end if
499
500 else
501 kick%qkick_mode = qkickmode_exp
502 end if
503
504 call parse_block_end(blk)
505
506 if (kpoints%use_symmetries) then
507 do iop = 1, symmetries_number(kpoints%symm)
508 if (iop == symmetries_identity_index(kpoints%symm)) cycle
509 if (.not. symm_op_invariant_cart(kpoints%symm%ops(iop), kick%qvector(:,1), 1e-5_real64)) then
510 message(1) = "The TDMomentumTransfer breaks (at least) one of the symmetries used to reduce the k-points."
511 message(2) = "Set SymmetryBreakDir equal to TDMomemtumTransfer."
512 call messages_fatal(2, namespace=namespace)
513 end if
514 end do
515 end if
516
517 else
518 kick%qkick_mode = qkickmode_none
519 kick%nqvec = 1
520 safe_allocate(kick%qvector(1:kick%dim,1))
521 kick%qvector(:,1) = m_zero
522 end if
523
524 kick%qlength = norm2(kick%qvector(:,1))
525
526 if (kick%delta_strength_mode == kick_magnon_mode) then
527 !%Variable TDEasyAxis
528 !%Type block
529 !%Section Time-Dependent::Response::Dipole
530 !%Description
531 !% For magnon kicks only.
532 !% This variable defines the direction of the easy axis of the crystal.
533 !% The magnetization is kicked in the plane transverse to this vector
534 !%End
535 if (parse_block(namespace, 'TDEasyAxis', blk) == 0) then
536 n_rows = parse_block_n(blk)
537
538 do idir = 1, 3
539 call parse_block_float(blk, 0, idir - 1, kick%easy_axis(idir))
540 end do
541 norm = norm2(kick%easy_axis(1:3))
542 if (norm < 1e-9_real64) then
543 message(1) = "TDEasyAxis norm is too small."
544 call messages_fatal(1, namespace=namespace)
545 end if
546 kick%easy_axis(1:3) = kick%easy_axis(1:3)/norm
547 call parse_block_end(blk)
548 else
549 message(1) = "For magnons, the variable TDEasyAxis must be defined."
550 call messages_fatal(1, namespace=namespace)
551 end if
552
553 !We first two vectors defining a basis in the transverse plane
554 !For this we take two vectors not align with the first one
555 !and we perform a Gram-Schmidt orthogonalization
556 kick%trans_vec(1,1) = -kick%easy_axis(2)
557 kick%trans_vec(2,1) = m_two*kick%easy_axis(3)
558 kick%trans_vec(3,1) = m_three*kick%easy_axis(1)
559
560 dot = sum(kick%easy_axis(1:3)*kick%trans_vec(1:3,1))
561 kick%trans_vec(1:3,1) = kick%trans_vec(1:3,1) - dot*kick%easy_axis(1:3)
562 norm = sum(kick%trans_vec(1:3,1)**2)
563 kick%trans_vec(1:3,1) = kick%trans_vec(1:3,1)/sqrt(norm)
564
565 !To get a direct basis, the last vector is obtained by the cross product
566 kick%trans_vec(1,2) = kick%easy_axis(2) * kick%trans_vec(3,1) - kick%easy_axis(3) * kick%trans_vec(2,1)
567 kick%trans_vec(2,2) = kick%easy_axis(3) * kick%trans_vec(1,1) - kick%easy_axis(1) * kick%trans_vec(3,1)
568 kick%trans_vec(3,2) = kick%easy_axis(1) * kick%trans_vec(2,1) - kick%easy_axis(2) * kick%trans_vec(1,1)
569
570 !The perturbation direction is defined as
571 !cos(q.r)*uvec + sin(q.r)*vvec
572
573
574 if (parse_is_defined(namespace, 'TDMomentumTransfer') &
575 .and. parse_is_defined(namespace, 'TDMultipleMomentumTransfer')) then
576 message(1) = "TDMomentumTransfer and TDMultipleMomentumTransfer cannot be defined at the same time."
577 call messages_fatal(1, namespace=namespace)
578 end if
579
580 if (parse_is_defined(namespace, 'TDMultipleMomentumTransfer')) then
581
582 kick%qkick_mode = qkickmode_exp
583
584 !%Variable TDMultipleMomentumTransfer
585 !%Type block
586 !%Section Time-Dependent::Response
587 !%Description
588 !% For magnon kicks only.
589 !% A simple way to specify momentum-transfer vectors for the calculation of
590 !% the magnetization dynamics. This variable should be used for a supercell.
591 !% For each reciprocal lattice vectors, the code will kick the original magnetization
592 !% using all the multiples of it.
593 !% The syntax reads:
594 !%
595 !% <tt>%TDMultipleMomentumTransfer
596 !% <br>&nbsp;&nbsp;N_x | N_y | N_z
597 !% <br>%</tt>
598 !%
599 !% and will include the (2N_x+1)*(2N_y+1)*(2N_z+1) multiples vectors of the reciprocal
600 !% lattice vectors of the current cell.
601 !%End
602 if (parse_block(namespace, 'TDMultipleMomentumTransfer', blk) /= 0) then
603 write(message(1),'(a)') 'Internal error while reading TDMultipleMomentumTransfer.'
604 call messages_fatal(1, namespace=namespace)
605 end if
606
607 do idir = 1, 3
608 call parse_block_integer(blk, 0, idir-1, kick%nqmult(idir))
609 end do
610
611 call parse_block_end(blk)
612
613
614 kick%nqvec = (2*kick%nqmult(1)+1)*(2*kick%nqmult(2)+1)*(2*kick%nqmult(3)+1)
615 !qvector has been allocated by default to a null vector before
616 safe_deallocate_a(kick%qvector)
617 safe_allocate(kick%qvector(1:3, 1:kick%nqvec))
618 iq = 0
619 do iqx = -kick%nqmult(1), kick%nqmult(1)
620 do iqy = -kick%nqmult(2), kick%nqmult(2)
621 do iqz = -kick%nqmult(3), kick%nqmult(3)
622 iq = iq + 1
623 qtemp(1:3) = (/iqx, iqy, iqz/)
624 call kpoints_to_absolute(kpoints%latt, qtemp, kick%qvector(1:3, iq))
625
626 !Checking symmetries for all G vectors
627 if (kpoints%use_symmetries) then
628 do iop = 1, symmetries_number(kpoints%symm)
629 if (iop == symmetries_identity_index(kpoints%symm)) cycle
630 if (.not. symm_op_invariant_cart(kpoints%symm%ops(iop), kick%qvector(:,iq), 1e-5_real64)) then
631 message(1) = "The TDMultipleMomentumTransfer breaks (at least) one " &
632 // "of the symmetries used to reduce the k-points."
633 message(2) = "Set SymmetryBreakDir accordingly."
634 call messages_fatal(2, namespace=namespace)
635 end if
636 end do
637 end if
638 end do
639 end do
640 end do
641
642 end if
643
644 else
645 kick%easy_axis(:) = m_zero
646 end if
647
648 if (kick%delta_strength_mode == kick_magnon_mode .and. kick%qkick_mode /= qkickmode_exp) then
649 message(1) = "For magnons, the kick mode must be exponential."
650 call messages_fatal(1, namespace=namespace)
651 end if
652
653 pop_sub(kick_init)
654 end subroutine kick_init
655
656
657 ! ---------------------------------------------------------
658 subroutine kick_copy(kick_out, kick_in)
659 type(kick_t), intent(inout) :: kick_out
660 type(kick_t), intent(in) :: kick_in
661
662 push_sub(kick_copy)
663
664 kick_out%dim = kick_in%dim
665
666 kick_out%time = kick_in%time
667
668 kick_out%delta_strength_mode = kick_in%delta_strength_mode
669 kick_out%delta_strength = kick_in%delta_strength
670
672 safe_allocate_source_a(kick_out%pol, kick_in%pol)
673 kick_out%pol_dir = kick_in%pol_dir
674 kick_out%pol_equiv_axes = kick_in%pol_equiv_axes
675 safe_allocate_source_a(kick_out%wprime, kick_in%wprime)
676
678 kick_out%n_multipoles = kick_in%n_multipoles
679 if (kick_out%n_multipoles > 0) then
680 safe_allocate_source_a(kick_out%l, kick_in%l)
681 safe_allocate_source_a(kick_out%m, kick_in%m)
682 safe_allocate_source_a(kick_out%weight, kick_in%weight)
683 end if
684 kick_out%nqvec = kick_in%nqvec
685 safe_allocate_source_a(kick_out%qvector, kick_in%qvector)
686 kick_out%qlength = kick_in%qlength
687 kick_out%qkick_mode = kick_in%qkick_mode
688 kick_out%qbessel_l = kick_in%qbessel_l
689 kick_out%qbessel_m = kick_in%qbessel_m
690
692 kick_out%function_mode = kick_in%function_mode
693 kick_out%user_defined_function = kick_in%user_defined_function
694
695 kick_out%easy_axis = kick_in%easy_axis
696
697 pop_sub(kick_copy)
698 end subroutine kick_copy
699
700 ! ---------------------------------------------------------
701 subroutine kick_end(kick)
702 type(kick_t), intent(inout) :: kick
703
704 push_sub(kick_end)
705
706 kick%delta_strength_mode = 0
707 kick%dim = 0
708 kick%pol_equiv_axes = 0
709 safe_deallocate_a(kick%pol)
710 kick%pol_dir = 0
711 safe_deallocate_a(kick%wprime)
712 if (kick%n_multipoles > 0) then
713 safe_deallocate_a(kick%l)
714 safe_deallocate_a(kick%m)
715 safe_deallocate_a(kick%weight)
716 end if
717 kick%n_multipoles = 0
718 kick%qkick_mode = qkickmode_none
719 safe_deallocate_a(kick%qvector)
720
721 pop_sub(kick_end)
722 end subroutine kick_end
723
724 ! ---------------------------------------------------------
725 subroutine kick_read(kick, iunit, namespace)
726 type(kick_t), intent(inout) :: kick
727 integer, intent(in) :: iunit
728 type(namespace_t), intent(in) :: namespace
729
730 integer :: idir, im, ierr
731 character(len=100) :: line
732
733 push_sub(kick_read)
734
735 kick%function_mode = -1
736
737 read(iunit, '(15x,i2)') kick%delta_strength_mode
738 read(iunit, '(15x,f18.12)') kick%delta_strength
739 read(iunit, '(15x,i2)') kick%dim
740 read(iunit, '(a)') line
741 if (index(line,'defined') /= 0) then
742 kick%function_mode = kick_function_user_defined
743 ! "# User defined: "
744 read(line,'(16x,a)') kick%user_defined_function
745 elseif (index(line,'multipole') /= 0) then
746 kick%function_mode = kick_function_multipole
747 ! "# N multipoles "
748 read(line, '(15x,i3)') kick%n_multipoles
749 safe_allocate( kick%l(1:kick%n_multipoles))
750 safe_allocate( kick%m(1:kick%n_multipoles))
751 safe_allocate(kick%weight(1:kick%n_multipoles))
752 do im = 1, kick%n_multipoles
753 ! "# multipole "
754 read(iunit, '(15x,2i3,f18.12)') kick%l(im), kick%m(im), kick%weight(im)
755 end do
756 else
757 kick%function_mode = kick_function_dipole
758 kick%n_multipoles = 0
759 backspace(iunit)
760
761 safe_allocate(kick%pol(1:kick%dim, 1:kick%dim))
762 safe_allocate(kick%wprime(1:kick%dim))
763 do idir = 1, kick%dim
764 read(iunit, '(15x,99f18.12)') kick%pol(1:kick%dim, idir)
765 end do
766 read(iunit, '(15x,i2)') kick%pol_dir
767 read(iunit, '(15x,i2)') kick%pol_equiv_axes
768 read(iunit, '(15x,99f18.12)') kick%wprime(1:kick%dim)
769 end if
770 if (kick%delta_strength_mode == kick_magnon_mode) then
771 read(iunit, '(15x,i3)') kick%nqvec
772 safe_allocate(kick%qvector(1:3, 1:kick%nqvec))
773 do im = 1, kick%nqvec
774 read(iunit, '(15x,3f18.12)') kick%qvector(1:kick%dim, im)
775 end do
776 read(iunit, '(15x,3f18.12)') kick%easy_axis(1:3)
777 read(iunit, '(15x,3f18.12)') kick%trans_vec(1:3,1)
778 read(iunit, '(15x,3f18.12)') kick%trans_vec(1:3,2)
779 end if
780 read(iunit, '(15x,f18.12)', iostat = ierr) kick%time
781
782 if (ierr /= 0) then
783 kick%time = m_zero
784 backspace(iunit)
785 end if
786
787 if (kick%function_mode < 0) then
788 message(1) = "No kick could be read from file."
789 call messages_fatal(1, namespace=namespace)
790 end if
791
792 pop_sub(kick_read)
793 end subroutine kick_read
795
796 ! ---------------------------------------------------------
797 subroutine kick_write(kick, iunit, out)
798 type(kick_t), intent(in) :: kick
799 integer, optional, intent(in) :: iunit
800 type(c_ptr), optional, intent(inout) :: out
801
802 integer :: idir, im
803 character(len=120) :: aux
804
805 push_sub(kick_write)
806
807 if (present(iunit)) then
808 write(iunit, '(a15,i1)') '# kick mode ', kick%delta_strength_mode
809 write(iunit, '(a15,f18.12)') '# kick strength', kick%delta_strength
810 write(iunit, '(a15,i2)') '# dim ', kick%dim
811 ! if this were to be read by humans, we would want units_from_atomic(units_out%length**(-1))
812 if (kick%function_mode == kick_function_user_defined) then
813 write(iunit,'(a15,1x,a)') '# User defined:', trim(kick%user_defined_function)
814 elseif (kick%n_multipoles > 0) then
815 write(iunit, '(a15,i3)') '# N multipoles ', kick%n_multipoles
816 do im = 1, kick%n_multipoles
817 write(iunit, '(a15,2i3,f18.12)') '# multipole ', kick%l(im), kick%m(im), kick%weight(im)
818 end do
819 else
820 do idir = 1, kick%dim
821 write(iunit, '(a6,i1,a8,99f18.12)') '# pol(', idir, ') ', kick%pol(1:kick%dim, idir)
822 end do
823 write(iunit, '(a15,i1)') '# direction ', kick%pol_dir
824 write(iunit, '(a15,i1)') '# Equiv. axes ', kick%pol_equiv_axes
825 write(iunit, '(a15,99f18.12)') '# wprime ', kick%wprime(1:kick%dim)
826 end if
827 write(iunit, '(a15,f18.12)') "# kick time ", kick%time
828
829 else if (present(out)) then
830 write(aux, '(a15,i2)') '# kick mode ', kick%delta_strength_mode
831 call write_iter_string(out, aux)
832 call write_iter_nl(out)
833 write(aux, '(a15,f18.12)') '# kick strength', kick%delta_strength
834 call write_iter_string(out, aux)
835 call write_iter_nl(out)
836 write(aux, '(a15,i2)') '# dim ', kick%dim
837 call write_iter_string(out, aux)
838 call write_iter_nl(out)
839 if (kick%function_mode == kick_function_user_defined) then
840 write(aux,'(a15,1x,a)') '# User defined:', trim(kick%user_defined_function)
841 call write_iter_string(out, aux)
842 call write_iter_nl(out)
843 elseif (kick%n_multipoles > 0) then
844 write(aux, '(a15,i3)') '# N multipoles ', kick%n_multipoles
845 call write_iter_string(out, aux)
846 call write_iter_nl(out)
847 do im = 1, kick%n_multipoles
848 write(aux, '(a15,2i3,f18.12)') '# multipole ', kick%l(im), kick%m(im), kick%weight(im)
849 call write_iter_string(out, aux)
850 call write_iter_nl(out)
851 end do
852 else
853 do idir = 1, kick%dim
854 write(aux, '(a6,i1,a8,99f18.12)') '# pol(', idir, ') ', kick%pol(1:kick%dim, idir)
855 call write_iter_string(out, aux)
856 call write_iter_nl(out)
857 end do
858 write(aux, '(a15,i2)') '# direction ', kick%pol_dir
859 call write_iter_string(out, aux)
860 call write_iter_nl(out)
861 write(aux, '(a15,i2)') '# Equiv. axes ', kick%pol_equiv_axes
862 call write_iter_string(out, aux)
863 call write_iter_nl(out)
864 write(aux, '(a15,99f18.12)') '# wprime ', kick%wprime(1:kick%dim)
865 call write_iter_string(out, aux)
866 call write_iter_nl(out)
867 end if
868 if (present(out) .and. kick%delta_strength_mode == kick_magnon_mode) then
869 write(aux, '(a15,i3)') '# N q-vectors ', kick%nqvec
870 call write_iter_string(out, aux)
871 call write_iter_nl(out)
872 do im = 1, kick%nqvec
873 write(aux, '(a15,3f18.12)') '# q-vector ', kick%qvector(1:kick%dim, im)
874 call write_iter_string(out, aux)
875 call write_iter_nl(out)
876 end do
877 write(aux, '(a15,3f18.12)') '# Easy axis ', kick%easy_axis(1:3)
878 call write_iter_string(out, aux)
879 call write_iter_nl(out)
880 write(aux, '(a15,3f18.12)') '# Trans. dir 1 ', kick%trans_vec(1:3,1)
881 call write_iter_string(out, aux)
882 call write_iter_nl(out)
883 write(aux, '(a15,3f18.12)') '# Trans. dir 2 ', kick%trans_vec(1:3,2)
884 call write_iter_string(out, aux)
885 call write_iter_nl(out)
886 end if
887 write(aux, '(a15,f18.12)') "# kick time ", kick%time
888 call write_iter_string(out, aux)
889 call write_iter_nl(out)
891 end if
892
893 pop_sub(kick_write)
894 end subroutine kick_write
895
896
897 ! ---------------------------------------------------------
898 !
899 subroutine kick_function_get(space, mesh, kick, kick_function, iq, to_interpolate)
900 class(space_t), intent(in) :: space
901 class(mesh_t), intent(in) :: mesh
902 type(kick_t), intent(in) :: kick
903 complex(real64), intent(out) :: kick_function(:)
904 integer, intent(in) :: iq
905 logical, optional, intent(in) :: to_interpolate
906
907 integer :: ip, im
908 real(real64) :: xx(space%dim)
909 real(real64) :: rkick, ikick, rr, ylm
910
911 integer :: np
912
913 push_sub(kick_function_get)
914
915 np = mesh%np
916 if (present(to_interpolate)) then
917 if (to_interpolate) np = mesh%np_part
918 end if
919
920 if (abs(kick%qlength) > m_epsilon .or. kick%delta_strength_mode == kick_magnon_mode) then ! q-vector is set
921
922 select case (kick%qkick_mode)
923 case (qkickmode_cos)
924 write(message(1), '(a,3F9.5,a)') 'Info: Using cos(q.r) field with q = (', kick%qvector(1:3, iq), ')'
925 case (qkickmode_sin)
926 write(message(1), '(a,3F9.5,a)') 'Info: Using sin(q.r) field with q = (', kick%qvector(1:3, iq), ')'
928 write(message(1), '(a,3F9.5,a)') 'Info: Using sin(q.r)+cos(q.r) field with q = (', kick%qvector(1:3, iq), ')'
929 case (qkickmode_exp)
930 select case(kick%dim)
931 case (1)
932 write(message(1), '(a,1F9.5,a)') 'Info: Using exp(iq.r) field with q = (', kick%qvector(1:kick%dim, iq), ')'
933 case (2)
934 write(message(1), '(a,2F9.5,a)') 'Info: Using exp(iq.r) field with q = (', kick%qvector(1:kick%dim, iq), ')'
935 case (3)
936 write(message(1), '(a,3F9.5,a)') 'Info: Using exp(iq.r) field with q = (', kick%qvector(1:kick%dim, iq), ')'
937 end select
938 case (qkickmode_bessel)
939 write(message(1), '(a,I2,a,I2,a,F9.5)') 'Info: Using j_l(qr)*Y_lm(r) field with (l,m)= (', &
940 kick%qbessel_l, ",", kick%qbessel_m,') and q = ', kick%qlength
941 case default
942 write(message(1), '(a)') 'Info: Unknown field type!'
943 end select
944 call messages_info(1)
945
946 kick_function = m_z0
947 do ip = 1, np
948 xx = mesh%x(ip, :)
949 select case (kick%qkick_mode)
950 case (qkickmode_cos)
951 kick_function(ip) = kick_function(ip) + cos(dot_product(kick%qvector(1:kick%dim, iq), xx))
952 case (qkickmode_sin)
953 kick_function(ip) = kick_function(ip) + sin(dot_product(kick%qvector(1:kick%dim, iq), xx))
955 kick_function(ip) = kick_function(ip) + sin(dot_product(kick%qvector(1:kick%dim, iq), xx))
956 case (qkickmode_exp)
957 kick_function(ip) = kick_function(ip) + exp(m_zi * dot_product(kick%qvector(:, iq), xx))
958 case (qkickmode_bessel)
959 call ylmr_real(xx, kick%qbessel_l, kick%qbessel_m, ylm)
960 kick_function(ip) = kick_function(ip) + loct_sph_bessel(kick%qbessel_l, kick%qlength*norm2(xx))*ylm
961 end select
962 end do
963
964 else
965 if (kick%function_mode == kick_function_user_defined) then
966
967 kick_function = m_z0
968 do ip = 1, np
969 call mesh_r(mesh, ip, rr, coords = xx)
970 rkick = m_zero
971 ikick = m_zero
972 call parse_expression(rkick, ikick, space%dim, xx, rr, m_zero, trim(kick%user_defined_function))
973 kick_function(ip) = rkick
974 end do
975
976 elseif (kick%n_multipoles > 0) then
977
978 kick_function = m_z0
979 do im = 1, kick%n_multipoles
980 do ip = 1, np
981 call mesh_r(mesh, ip, rr, coords = xx)
982 call loct_ylm(1, xx(1), xx(2), xx(3), kick%l(im), kick%m(im), ylm)
983 kick_function(ip) = kick_function(ip) + kick%weight(im) * (rr**kick%l(im)) * ylm
984 end do
985 end do
986 else
987 do ip = 1, np
988 kick_function(ip) = sum(mesh%x(ip, 1:space%dim) * &
989 kick%pol(1:space%dim, kick%pol_dir))
990 end do
991 end if
992 end if
993
994 pop_sub(kick_function_get)
995 end subroutine kick_function_get
996
997
998 ! ---------------------------------------------------------
999 !
1000 subroutine kick_pcm_function_get(space, mesh, kick, psolver, pcm, kick_pcm_function)
1001 class(space_t), intent(in) :: space
1002 class(mesh_t), intent(in) :: mesh
1003 type(kick_t), intent(in) :: kick
1004 type(poisson_t), intent(in) :: psolver
1005 type(pcm_t), intent(inout) :: pcm
1006 complex(real64), intent(out) :: kick_pcm_function(:)
1007
1008 complex(real64), allocatable :: kick_function_interpolate(:)
1009 real(real64), allocatable :: kick_function_real(:)
1010
1011 push_sub(kick_pcm_function_get)
1012
1013 kick_pcm_function = m_zero
1014 if (pcm%localf) then
1015 safe_allocate(kick_function_interpolate(1:mesh%np_part))
1016 kick_function_interpolate = m_zero
1017 call kick_function_get(space, mesh, kick, kick_function_interpolate, 1, to_interpolate = .true.)
1018 safe_allocate(kick_function_real(1:mesh%np_part))
1019 kick_function_real = dreal(kick_function_interpolate)
1020 if (pcm%kick_like) then
1021 ! computing kick-like polarization due to kick
1022 call pcm_calc_pot_rs(pcm, mesh, psolver, kick = kick%delta_strength * kick_function_real, kick_time = .true.)
1023 else if (.not. pcm%kick_like .and. pcm%which_eps == pcm_debye_model) then
1024 ! computing the kick-like part of polarization due to kick for Debye dielectric model
1025 pcm%kick_like = .true.
1026 call pcm_calc_pot_rs(pcm, mesh, psolver, kick = kick%delta_strength * kick_function_real, kick_time = .true.)
1027 pcm%kick_like = .false.
1028 else if (.not. pcm%kick_like .and. pcm%which_eps == pcm_drude_model) then
1029 pop_sub(kick_pcm_function_get)
1030 return
1031 end if
1032 kick_pcm_function = pcm%v_kick_rs / kick%delta_strength
1033 end if
1034
1035 pop_sub(kick_pcm_function_get)
1036 end subroutine kick_pcm_function_get
1037
1038
1039 ! ---------------------------------------------------------
1042 subroutine kick_apply(space, mesh, st, ions_dyn, ions, kick, psolver, kpoints, pcm)
1043 class(space_t), intent(in) :: space
1044 class(mesh_t), intent(in) :: mesh
1045 type(states_elec_t), intent(inout) :: st
1046 type(ion_dynamics_t), intent(in) :: ions_dyn
1047 type(ions_t), intent(inout) :: ions
1048 type(kick_t), intent(in) :: kick
1049 type(poisson_t), intent(in) :: psolver
1050 type(kpoints_t), intent(in) :: kpoints
1051 type(pcm_t), optional, intent(inout) :: pcm
1052
1053 integer :: iqn, ist, idim, ip, ispin, iatom
1054 complex(real64) :: cc(2), kick_value
1055 complex(real64), allocatable :: kick_function(:), psi(:, :)
1056
1057 complex(real64), allocatable :: kick_pcm_function(:)
1058 integer :: ns, iq
1059 real(real64) :: uvec(3), vvec(3), gvec(3, 3)
1060 real(real64) :: xx(space%dim), rr
1061
1062 push_sub(kick_apply)
1063
1064 ! The wavefunctions at time delta t read
1065 ! psi(delta t) = psi(t) exp(i k x)
1066 delta_strength: if (abs(kick%delta_strength) > m_epsilon) then
1067
1068 safe_allocate(kick_function(1:mesh%np))
1069 if (kick%delta_strength_mode /= kick_magnon_mode .or. kick%nqvec == 1) then
1070 call kick_function_get(space, mesh, kick, kick_function, 1)
1071 end if
1072
1073 ! PCM - computing polarization due to kick
1074 if (present(pcm)) then
1075 safe_allocate(kick_pcm_function(1:mesh%np))
1076 call kick_pcm_function_get(space, mesh, kick, psolver, pcm, kick_pcm_function)
1077 kick_function = kick_function + kick_pcm_function
1078 end if
1079
1080 write(message(1),'(a,f11.6)') 'Info: Applying delta kick: k = ', kick%delta_strength
1081 select case (kick%function_mode)
1082 case (kick_function_dipole)
1083 message(2) = "Info: kick function: dipole."
1085 message(2) = "Info: kick function: multipoles."
1087 message(2) = "Info: kick function: user defined function."
1088 end select
1089 select case (kick%delta_strength_mode)
1090 case (kick_density_mode)
1091 message(3) = "Info: Delta kick mode: Density mode"
1092 case (kick_spin_mode)
1093 message(3) = "Info: Delta kick mode: Spin mode"
1095 message(3) = "Info: Delta kick mode: Density + Spin modes"
1096 end select
1097 call messages_info(3)
1098
1099 ns = 1
1100 if (st%d%nspin == 2) ns = 2
1101
1102 safe_allocate(psi(1:mesh%np, 1:st%d%dim))
1103
1104 if (kick%delta_strength_mode /= kick_magnon_mode) then
1105
1106 do iqn = st%d%kpt%start, st%d%kpt%end
1107 do ist = st%st_start, st%st_end
1108 call states_elec_get_state(st, mesh, ist, iqn, psi)
1109
1110 select case (kick%delta_strength_mode)
1111 case (kick_density_mode)
1112 do idim = 1, st%d%dim
1113 do ip = 1, mesh%np
1114 psi(ip, idim) = exp(m_zi*kick%delta_strength*kick_function(ip))*psi(ip, idim)
1115 end do
1116 end do
1117
1118 case (kick_spin_mode)
1119 ispin = st%d%get_spin_index(iqn)
1120 do ip = 1, mesh%np
1121 kick_value = m_zi*kick%delta_strength*kick_function(ip)
1122
1123 cc(1) = exp(kick_value)
1124 cc(2) = exp(-kick_value)
1125
1126 select case (st%d%ispin)
1127 case (spin_polarized)
1128 psi(ip, 1) = cc(ispin)*psi(ip, 1)
1129 case (spinors)
1130 psi(ip, 1) = cc(1)*psi(ip, 1)
1131 psi(ip, 2) = cc(2)*psi(ip, 2)
1132 end select
1133 end do
1134
1136 do ip = 1, mesh%np
1137 kick_value = m_zi*kick%delta_strength*kick_function(ip)
1138 cc(1) = exp(m_two*kick_value)
1139
1140 select case (st%d%ispin)
1141 case (spin_polarized)
1142 if (is_spin_up(iqn)) then
1143 psi(ip, 1) = cc(1)*psi(ip, 1)
1144 end if
1145 case (spinors)
1146 psi(ip, 1) = cc(1)*psi(ip, 1)
1147 end select
1148 end do
1149 end select
1150
1151 call states_elec_set_state(st, mesh, ist, iqn, psi)
1152
1153 end do
1154 end do
1155
1156 else
1157 assert(st%d%ispin == spinors)
1158
1159 if (kick%nqvec == 1) then
1160 !The perturbation direction is defined as
1161 !cos(q.r)*uvec + sin(q.r)*vvec
1162 uvec(1:3) = kick%trans_vec(1:3,1)
1163 vvec(1:3) = kick%trans_vec(1:3,2)
1164
1165 do iqn = st%d%kpt%start, st%d%kpt%end, ns
1166 do ist = st%st_start, st%st_end
1167
1168 call states_elec_get_state(st, mesh, ist, iqn, psi)
1169
1170 do ip = 1, mesh%np
1171
1172 cc(1) = psi(ip, 1)
1173 cc(2) = psi(ip, 2)
1174
1175 !First part: 1I*cos(\lambda)
1176 psi(ip, 1) = cos(kick%delta_strength)* cc(1)
1177 psi(ip, 2) = cos(kick%delta_strength)* cc(2)
1178
1179 !We now add -i sin(\lambda) u.\sigma
1180 ! (u_z u_x-i*u_y) (v_z v_x-i*v_y)
1181 ! =cos(q.r) ( ) + sin(q.r)( )
1182 ! (u_x+i*u_y -u_z ) (v_x+i*v_y -v_z )
1183 psi(ip, 1) = psi(ip, 1) -m_zi*sin(kick%delta_strength)*( real(kick_function(ip), real64) &
1184 * (uvec(3)*cc(1) + (uvec(1)-m_zi*uvec(2))*cc(2)) &
1185 + aimag(kick_function(ip)) * (vvec(3)*cc(1) + (vvec(1)-m_zi*vvec(2))*cc(2)))
1186 psi(ip, 2) = psi(ip, 2) -m_zi*sin(kick%delta_strength)*( real(kick_function(ip), real64) &
1187 * (-uvec(3)*cc(2) + (uvec(1)+m_zi*uvec(2))*cc(1)) &
1188 + aimag(kick_function(ip)) * (-vvec(3)*cc(2) + (vvec(1)+m_zi*vvec(2))*cc(1)))
1189
1190 end do
1191
1192 call states_elec_set_state(st, mesh, ist, iqn, psi)
1193
1194 end do
1195 end do
1196
1197 else ! Multi-q kick
1198
1199 call kpoints_to_absolute(kpoints%latt, (/m_one,m_zero,m_zero/), gvec(1:3, 1))
1200 call kpoints_to_absolute(kpoints%latt, (/m_zero,m_one,m_zero/), gvec(1:3, 2))
1201 call kpoints_to_absolute(kpoints%latt, (/m_zero,m_zero,m_one/), gvec(1:3, 3))
1202
1203 kick_function = m_one
1204 do ip = 1, mesh%np
1205 call mesh_r(mesh, ip, rr, coords = xx)
1206 do iq = 1, 3
1207 if (kick%nqmult(iq) == 0) cycle
1208 if (abs(sin(m_half*sum(gvec(1:3, iq) * xx(1:3)))) <= m_epsilon) cycle
1209
1210 kick_function(ip) = kick_function(ip)*sin(m_half*(2*kick%nqmult(iq)+1) &
1211 *sum(gvec(1:3, iq) * xx(1:3)))/sin(m_half*sum(gvec(1:3, iq) * xx(1:3)))
1212 end do
1213 kick_function(ip) = kick_function(ip)*kick%delta_strength
1214 end do
1215
1216 do iqn = st%d%kpt%start, st%d%kpt%end, ns
1217 do ist = st%st_start, st%st_end
1218
1219 call states_elec_get_state(st, mesh, ist, iqn, psi)
1220
1221 do ip = 1, mesh%np
1222
1223 cc(1) = psi(ip, 1)
1224 cc(2) = psi(ip, 2)
1225
1226 ! (cos(F) + in_x sin(F) sin(F)(u_y (u_x-iu_y)/(1+u_z) - iu_z))
1227 ! = ( )
1228 ! (-sin(F)(u_y (u_x+iu_y)/(1+u_z)+iu_z) cos(F) - in_x sin(F) )
1229
1230 psi(ip, 1) = (cos(kick_function(ip))+m_zi*kick%easy_axis(1)*sin(kick_function(ip)))*cc(1) &
1231 +sin(kick_function(ip))*(kick%easy_axis(2)*(kick%easy_axis(1) &
1232 -m_zi*kick%easy_axis(2))/(1+kick%easy_axis(3))-m_zi*kick%easy_axis(3))*cc(2)
1233 psi(ip, 2) =-sin(kick_function(ip))*(kick%easy_axis(2)*(kick%easy_axis(1) &
1234 +m_zi*kick%easy_axis(2))/(1+kick%easy_axis(3))+m_zi*kick%easy_axis(3))*cc(1) &
1235 + (cos(kick_function(ip))-m_zi*kick%easy_axis(1)*sin(kick_function(ip)))*cc(2)
1236
1237 end do
1238
1239 call states_elec_set_state(st, mesh, ist, iqn, psi)
1240
1241 end do
1242 end do
1243
1244 end if
1245 end if
1246
1247 safe_deallocate_a(psi)
1248
1249 ! The nuclear velocities will be changed by
1250 ! Delta v_z = ( Z*e*E_0 / M) = - ( Z*k*\hbar / M)
1251 ! where M and Z are the ionic mass and charge, respectively.
1252 if (ion_dynamics_ions_move(ions_dyn) .and. abs(kick%delta_strength) > m_epsilon) then
1253 if (kick%delta_strength_mode /= kick_magnon_mode) then
1254 do iatom = 1, ions%natoms
1255 ions%vel(:, iatom) = ions%vel(:, iatom) + &
1256 kick%delta_strength * kick%pol(1:ions%space%dim, kick%pol_dir) * &
1257 p_proton_charge * ions%charge(iatom) / ions%mass(iatom)
1258 end do
1259 end if
1260 end if
1261
1262 safe_deallocate_a(kick_function)
1263 end if delta_strength
1264
1265 pop_sub(kick_apply)
1266 end subroutine kick_apply
1267
1268 pure integer function kick_get_type(kick) result(kick_type)
1269 type(kick_t), intent(in) :: kick
1270
1271 kick_type = kick%delta_strength_mode
1272
1273 end function kick_get_type
1274
1275end module kick_oct_m
1276
1277!! Local Variables:
1278!! mode: f90
1279!! coding: utf-8
1280!! End:
double exp(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public p_proton_charge
Definition: global.F90:226
real(real64), parameter, public m_zero
Definition: global.F90:187
complex(real64), parameter, public m_z0
Definition: global.F90:197
complex(real64), parameter, public m_zi
Definition: global.F90:201
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
real(real64), parameter, public m_three
Definition: global.F90:190
logical pure function, public ion_dynamics_ions_move(this)
subroutine, public kick_copy(kick_out, kick_in)
Definition: kick.F90:752
integer, parameter, public kick_magnon_mode
Definition: kick.F90:163
integer, parameter, public kick_function_multipole
Definition: kick.F90:158
subroutine, public kick_read(kick, iunit, namespace)
Definition: kick.F90:819
subroutine kick_pcm_function_get(space, mesh, kick, psolver, pcm, kick_pcm_function)
Definition: kick.F90:1094
integer, parameter, public kick_spin_mode
Definition: kick.F90:163
integer, parameter, public qkickmode_exp
Definition: kick.F90:169
integer, parameter, public qkickmode_cos
Definition: kick.F90:169
subroutine, public kick_end(kick)
Definition: kick.F90:795
integer, parameter, public qkickmode_sin
Definition: kick.F90:169
subroutine, public kick_init(kick, namespace, space, kpoints, nspin)
Definition: kick.F90:223
subroutine, public kick_function_get(space, mesh, kick, kick_function, iq, to_interpolate)
Definition: kick.F90:993
pure integer function, public kick_get_type(kick)
Definition: kick.F90:1362
subroutine, public kick_apply(space, mesh, st, ions_dyn, ions, kick, psolver, kpoints, pcm)
Applies the delta-function electric field where k = kick%delta_strength.
Definition: kick.F90:1136
subroutine, public kick_write(kick, iunit, out)
Definition: kick.F90:891
integer, parameter, public kick_function_user_defined
Definition: kick.F90:158
integer, parameter, public qkickmode_bessel
Definition: kick.F90:169
integer, parameter, public kick_spin_density_mode
Definition: kick.F90:163
subroutine, public kpoints_to_absolute(latt, kin, kout)
Definition: kpoints.F90:1031
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
subroutine, public ylmr_real(xx, li, mi, ylm)
This is a Numerical Recipes-based subroutine computes real spherical harmonics ylm at position (x,...
Definition: math.F90:372
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:336
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
integer, parameter, public pcm_drude_model
Definition: pcm_eom.F90:168
integer, parameter, public pcm_debye_model
Definition: pcm_eom.F90:168
subroutine, public pcm_calc_pot_rs(pcm, mesh, psolver, ions, v_h, v_ext, kick, time_present, kick_time)
Definition: pcm.F90:1215
This module handles spin dimensions of the states and the k-point distribution.
logical pure function, public is_spin_up(ik)
Returns true if k-point ik denotes spin-up, in spin-polarized case.
integer pure function, public symmetries_identity_index(this)
Definition: symmetries.F90:602
integer pure function, public symmetries_number(this)
Definition: symmetries.F90:556
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
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
Explicit interfaces to C functions, defined in write_iter_low.cc.
Definition: write_iter.F90:114
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.
int true(void)