Octopus
pcm_eom.F90
Go to the documentation of this file.
1!! Copyright (C) 2017 Gabriel Gil, Stefano Corni, Silvio Pipolo, Carlo Andrea Rozzi
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
20module pcm_eom_oct_m
21 use debug_oct_m
22 use global_oct_m
23 use io_oct_m
24 use, intrinsic :: iso_fortran_env
25 use math_oct_m
29 implicit none
30
31 private
32 public :: &
41 pcm_nuclei, &
46 save
47
49 type :: pcm_tessera_t
50 real(real64) :: point(1:3)
51 real(real64) :: area
52 real(real64) :: normal(1:3)
53 real(real64) :: r_sphere
54 end type pcm_tessera_t
55
56 type(pcm_tessera_t), allocatable :: cts_act(:)
57 integer :: nts_act
58
60 type :: debye_param_t
61 real(real64) :: eps_0
62 real(real64) :: eps_d
63 real(real64) :: tau
64 end type debye_param_t
65
69 type :: drude_param_t
70 real(real64) :: aa
71 real(real64) :: gm
72 real(real64) :: w0
73 end type drude_param_t
74
75 integer, parameter :: &
76 PCM_DEBYE_MODEL = 1, &
78
79 integer, parameter :: &
80 PCM_ELECTRONS = 0, &
81 pcm_nuclei = 1, &
83 pcm_kick = 3, &
85
86
87 type(debye_param_t) :: deb
88 type(drude_param_t) :: drl
89
90
91 integer :: which_eom
93
94 integer :: which_eps
95
96 real(real64) :: dt
97
98 ! polarization charges and variation of it in the previous iteration:
99 real(real64), allocatable :: q_tp(:), dq_tp(:)
100 real(real64), allocatable :: qext_tp(:), dqext_tp(:)
101 real(real64), allocatable :: qkick_tp(:), dqkick_tp(:)
102
103 real(real64), allocatable :: pot_tp(:)
104 real(real64), allocatable :: potext_tp(:)
105
106 ! See Chem.Phys.Lett. 429 (2006) 310-316 for Velocity-Verlet (VV) algorithm...
107 real(real64) :: f1, f2, f3, f4, f5
108 ! analogous to force in the equation of motion for the pol.charges, prev. iter.
109 real(real64), allocatable :: force_tp(:)
110 real(real64), allocatable :: force_qext_tp(:)
111 real(real64), allocatable :: force_qkick_tp(:)
112
113 ! In Ref.1 - J.Phys.Chem.A 2015, 119, 5405-5416... - Debye dielectric func. (eps)
114 ! In Ref.2 - In J. Phys. Chem. C 2016, 120, 28774−28781... - Drude-Lorentz eps
115 real(real64), allocatable :: cals(:,:), cald(:,:)
116 real(real64), allocatable :: eigv(:), eigt(:,:)
117 real(real64), allocatable :: sm12(:,:), sp12(:,:)
119 real(real64), allocatable :: matq0(:,:), matqd(:,:)
120 real(real64), allocatable :: matq0_lf(:,:), matqd_lf(:,:)
121 ! !! \f$ \tilde{Q} \f$ and R matrices from Eq.(38)-(39) (Ref.1), respectively
122 ! !! Q_f and Q_{\omega} from Eq.(17)-(16) (Ref.2), respectively
123 real(real64), allocatable :: matqv(:,:), matqq(:,:)
124 real(real64), allocatable :: matqv_lf(:,:)
125 ! !! Q^{IEF(d)}_d, \f$ \tilde{Q} \f$ and R matrices enter the EOM in eq.(37), Ref.1
126 ! !! Q_f and Q_{\omega} matrices enter the EOM in eq.(15), Ref.2
127 ! !! N.B.: matrices R (in case of Debye) and Q_{\omega} (in case of Drude-Lorentz),
128 ! !! i.e., matqq in this implementation,
129 ! !! are the same in the EOMs for polarization charges due to the solute or ext.pot.
130
131 logical :: enough_initial = .false.
132
133contains
134
135 !------------------------------------------------------------------------------------------------------------------------------
138 subroutine pcm_charges_propagation(q_t, pot_t, this_dt, this_cts_act, input_asc, this_eom, &
139 this_eps, namespace, this_deb, this_drl)
140 real(real64), intent(out) :: q_t(:)
141 real(real64), intent(in) :: pot_t(:)
142 real(real64), intent(in) :: this_dt
143 type(pcm_tessera_t), intent(in) :: this_cts_act(:)
144 logical, intent(in) :: input_asc
145 integer, intent(in) :: this_eom
146 ! !! or due to external potential ('external')
147 integer, intent(in) :: this_eps
148 ! !!(PCM_DEBYE_MODEL) or Drude-Lorentz (PCM_DRUDE_MODEL)
149 type(namespace_t), intent(in) :: namespace
150 type(debye_param_t), optional, intent(in) :: this_deb
151 type(drude_param_t), optional, intent(in) :: this_drl
153 logical, save :: firsttime = .true.
154 logical, save :: initial_electron = .true.
155 logical, save :: initial_external = .true.
156 logical, save :: initial_kick = .true.
159
160 which_eom = this_eom
161 if (which_eom /= pcm_electrons .and. which_eom /= pcm_external_potential .and. &
162 which_eom /= pcm_external_plus_kick .and. which_eom /= pcm_kick) then
163 message(1) = "pcm_charges_propagation: EOM evolution of PCM charges can only be due to solute electrons"
164 message(2) = "or external potential (including the kick)."
165 call messages_fatal(2, namespace=namespace)
166 end if
168 if (firsttime) then
169 dt = this_dt
170 nts_act = size(this_cts_act)
171 if (size(q_t) /= nts_act) then
172 message(1) = "pcm_charges_propagation: Number of tesserae do not coincide with size of PCM charges array."
173 call messages_fatal(1, namespace=namespace)
174 end if
175
176 safe_allocate(cts_act(1:nts_act))
177 cts_act = this_cts_act
178
179 which_eps = this_eps
180 select case (which_eps)
181 case (pcm_debye_model)
182 if (present(this_deb)) then
183 deb = this_deb
184 else
185 message(1) = "pcm_charges_propagation: EOM-PCM error. Debye dielectric function requires three parameters."
186 call messages_fatal(1, namespace=namespace)
187 end if
188 case (pcm_drude_model)
189 if (present(this_drl)) then
190 drl = this_drl
191 else
192 message(1) = "pcm_charges_propagation: EOM-PCM error. Drude-Lorentz dielectric function requires three parameters."
193 call messages_fatal(1, namespace=namespace)
194 end if
195 case default
196 message(1) = "pcm_charges_propagation: EOM-PCM error. Only Debye or Drude-Lorent dielectric models are allowed."
197 call messages_fatal(1, namespace=namespace)
198 end select
200 if (abs(deb%tau) <= m_epsilon) then
201 message(1) = "pcm_charges_propagation: EOM-PCM error. Debye EOM-PCM require a non-null Debye relaxation time."
202 call messages_fatal(1, namespace=namespace)
203 end if
204 firsttime = .false.
205 end if
207
208 if (input_asc) then
209 select case (which_eps)
210 case (pcm_debye_model)
211 ! initialize pcm charges due to electrons, external potential or kick
212 call pcm_charges_from_input_file(q_t, pot_t, namespace)
213
215 return
216 case default
217 message(1) = "pcm_charges_propagation: EOM-PCM error. Only Debye EOM-PCM can startup from input charges."
218 call messages_fatal(1, namespace=namespace)
219 end select
220 end if
221
222 if ((initial_electron .and. which_eom == pcm_electrons) .or. &
223 (initial_external .and. (which_eom == pcm_external_potential .or. which_eom == pcm_external_plus_kick)) .or. &
224 (initial_kick .and. which_eom == pcm_kick)) then
225
226 ! initialize pcm matrices
227 call pcm_bem_init(namespace)
228
229 ! initialize pcm charges due to electrons, external potential or kick
230 call init_charges(q_t, pot_t, namespace)
231
232 if (initial_electron .and. which_eom == pcm_electrons) initial_electron = .false.
233 if (initial_external .and. (which_eom == pcm_external_potential .or. &
234 which_eom == pcm_external_plus_kick)) initial_external = .false.
235 if (initial_kick .and. which_eom == pcm_kick) initial_kick = .false.
236
237 else
238
239 ! propagate pcm charges due to electrons or external potential (including possible kick)
240 if (which_eps == pcm_debye_model) then
241 call pcm_ief_prop_deb(q_t, pot_t)
242 else if (which_eps == pcm_drude_model) then
243 call pcm_ief_prop_vv_ief_drl(q_t, pot_t)
244 end if
245
246 end if
247
249 end subroutine pcm_charges_propagation
250
251 !------------------------------------------------------------------------------------------------------------------------------
253 subroutine pcm_charges_from_input_file(q_t, pot_t, namespace)
254 real(real64), intent(out) :: q_t(:)
255 real(real64), intent(in) :: pot_t(:)
256 type(namespace_t), intent(in) :: namespace
257
258 real(real64) :: aux1(3)
259 integer :: aux2
260 integer :: asc_unit
261 integer :: ia
262
264
265 if (which_eom == pcm_electrons) then
266 safe_allocate(q_tp(1:nts_act))
267 asc_unit = io_open(pcm_dir//'ASC_e.dat', namespace, action='read')
268 do ia = 1, nts_act
269 read(asc_unit,*) aux1, q_t(ia), aux2
270 end do
271 q_tp = q_t
272
273 else if (which_eom == pcm_external_potential .or. which_eom == pcm_external_plus_kick) then
274 safe_allocate(qext_tp(1:nts_act))
275 asc_unit = io_open(pcm_dir//'ASC_ext.dat', namespace, action='read')
276 do ia = 1, nts_act
277 read(asc_unit,*) aux1, q_t(ia), aux2
278 end do
279 qext_tp = q_t
280
281 else if (which_eom == pcm_kick) then
282 safe_allocate(qkick_tp(1:nts_act))
283 asc_unit = io_open(pcm_dir//'ASC_kick.dat', namespace, action='read')
284 do ia = 1, nts_act
285 read(asc_unit,*) aux1, q_t(ia), aux2
286 end do
287 qkick_tp = q_t
288 end if
289 call io_close(asc_unit)
290
291 safe_allocate(pot_tp(1:nts_act))
292 pot_tp = pot_t
293
295 end subroutine pcm_charges_from_input_file
296
297 !------------------------------------------------------------------------------------------------------------------------------
299 subroutine init_charges(q_t,pot_t, namespace)
300 real(real64), intent(out) :: q_t(:)
301 real(real64), intent(in) :: pot_t(:)
302 type(namespace_t), intent(in) :: namespace
303
304 push_sub(init_charges)
305
306 if (which_eom == pcm_electrons) then
307 ! Here we consider the potential at any earlier time equal to the initial potential.
308 ! Therefore, we can suppose that the solvent is already in equilibrium with the initial solute potential.
309 ! The latter is valid when starting the propagation from the ground state but does not hold in general.
310 message(1) = 'EOM-PCM for solvent polarization due to solute electrons considers that you start from a ground state run.'
311 call messages_warning(1, namespace=namespace)
312
313 safe_allocate(pot_tp(1:nts_act))
314 pot_tp = pot_t
315
316 ! applying the static IEF-PCM response matrix (corresponging to epsilon_0) to the initial potential
317 q_t = matmul(matq0, pot_t)
318
319 safe_allocate(q_tp(1:nts_act))
320 q_tp = q_t
321
322 if (which_eps == pcm_drude_model) then
323 safe_allocate(dq_tp(1:nts_act))
324 safe_allocate(force_tp(1:nts_act))
325 dq_tp = m_zero
326 force_tp = m_zero
327 end if
328
329 else if (which_eom == pcm_external_potential .or. which_eom == pcm_external_plus_kick) then
330 ! Here (instead) we consider zero the potential at any earlier time.
331 ! Therefore, the solvent is not initially in equilibrium with the external potential unless its initial value is zero.
332
333 safe_allocate(potext_tp(1:nts_act))
334 potext_tp = pot_t
335
336 ! applying the dynamic IEF-PCM response matrix (for epsilon_d) to the initial potential
337 q_t = matmul(matqd_lf, pot_t)
338
339 safe_allocate(qext_tp(1:nts_act))
340 qext_tp = q_t
341
342 if (which_eps == pcm_drude_model) then
343 safe_allocate(dqext_tp(1:nts_act))
344 safe_allocate(force_qext_tp(1:nts_act))
345 dqext_tp = m_zero
346 force_qext_tp = m_zero
347 end if
349 else if (which_eom == pcm_kick) then
350
351 if (which_eps == pcm_drude_model) then
352 safe_allocate(dqkick_tp(1:nts_act))
353 safe_allocate(force_qkick_tp(1:nts_act))
354 dqkick_tp = m_zero
355 force_qkick_tp = m_zero
356 end if
357
358 if (which_eps == pcm_drude_model) then
359 dqkick_tp = matmul(matqv_lf, pot_t)*dt
360 else
361 q_t = matmul(matqv_lf - matmul(matqq, matqd_lf), pot_t)
362 end if
363
364 safe_allocate(qkick_tp(1:nts_act))
365 qkick_tp = q_t
366
367 end if
368
369 ! initializing Velocity-Verlet algorithm for the integration of EOM-PCM for Drude-Lorentz
370 if (which_eps == pcm_drude_model) call init_vv_propagator
371
372 pop_sub(init_charges)
373 end subroutine init_charges
374
375 !------------------------------------------------------------------------------------------------------------------------------
378 subroutine pcm_ief_prop_deb(q_t, pot_t)
379 real(real64), intent(out) :: q_t(:)
380 real(real64), intent(in) :: pot_t(:)
381
382 push_sub(pcm_ief_prop_deb)
383
384 if (which_eom == pcm_electrons) then
385 ! Eq.(47) in S. Corni, S. Pipolo and R. Cammi, J.Phys.Chem.A 2015, 119, 5405-5416.
386 q_t(:) = q_tp(:) - dt*matmul(matqq, q_tp) + &
387 dt*matmul(matqv, pot_tp) + &
388 matmul(matqd, pot_t - pot_tp)
389
390 q_tp = q_t
391
392 pot_tp = pot_t
393
394 else if (which_eom == pcm_external_potential .or. which_eom == pcm_external_plus_kick) then
395 ! analogous to Eq.(47) ibid. with different matrices
396 q_t(:) = qext_tp(:) - dt*matmul(matqq, qext_tp) + &
397 dt*matmul(matqv_lf, potext_tp) + &
398 matmul(matqd_lf, pot_t - potext_tp) ! N.B. matqq
399
400 qext_tp = q_t
401
402 potext_tp = pot_t
403
404 else if (which_eom == pcm_kick) then
405 ! simplifying for kick
406 q_t(:) = qkick_tp(:) - dt*matmul(matqq, qkick_tp) ! N.B. matqq
407
408 qkick_tp = q_t
409
410 end if
411
412 pop_sub(pcm_ief_prop_deb)
413 end subroutine pcm_ief_prop_deb
414
415 !------------------------------------------------------------------------------------------------------------------------------
421 subroutine init_vv_propagator
422 push_sub(init_vv_propagator)
423
424 f1 = dt*(m_one - dt*m_half*drl%gm)
425 f2 = dt*dt*m_half
426 f3 = m_one - dt*drl%gm*(m_one - dt*m_half*drl%gm)
427 f4 = m_half*dt
428 f5 = drl%gm*f2
429
430 pop_sub(init_vv_propagator)
431 end subroutine init_vv_propagator
432
433 !------------------------------------------------------------------------------------------------------------------------------
436 subroutine pcm_ief_prop_vv_ief_drl(q_t, pot_t)
437 real(real64), intent(out) :: q_t(:)
438 real(real64), intent(in) :: pot_t(:)
439
440 real(real64) :: force(nts_act)
441
443
444 if (which_eom == pcm_electrons) then
445 ! From Eq.(15) S. Pipolo and S. Corni, J.Phys.Chem.C 2016, 120, 28774-28781.
446 ! Using the scheme in Eq.(21) and (17) of E. Vanden-Eijnden, G. Ciccotti, Chem.Phys.Lett. 429 (2006) 310-316.
447 q_t = q_tp + f1*dq_tp + f2*force_tp
448 force = -matmul(matqq, q_t) + matmul(matqv, pot_t)
449 dq_tp = f3*dq_tp + f4*(force+force_tp) -f5*force_tp
450 force_tp = force
451 q_tp = q_t
452
453 else if (which_eom == pcm_external_potential .or. which_eom == pcm_external_plus_kick) then
454 ! analagous to Eq.(15) ibid. with different matrices
455 q_t = qext_tp + f1*dqext_tp + f2*force_qext_tp
456 force = -matmul(matqq, q_t) + matmul(matqv_lf, pot_t) ! N.B. matqq
457 dqext_tp = f3*dqext_tp + f4*(force + force_qext_tp) -f5*force_qext_tp
458 force_qext_tp = force
459 q_tp = q_t
460
461 else if (which_eom == pcm_kick) then
462 ! simplifying for kick
463 q_t = qkick_tp + f1*dqkick_tp + f2*force_qkick_tp
464 force = -matmul(matqq, q_t) ! N.B. matqq
465 dqkick_tp = f3*dqkick_tp + f4*(force + force_qkick_tp) -f5*force_qkick_tp
466 force_qkick_tp = force
467 q_tp = q_t
468
469 end if
470
471 pot_tp = pot_t
472
474 end subroutine pcm_ief_prop_vv_ief_drl
475
476 !------------------------------------------------------------------------------------------------------------------------------
478 subroutine pcm_bem_init(namespace)
479 type(namespace_t), intent(in) :: namespace
480 integer :: itess, jtess
481 integer :: pcmmat0_unit, pcmmatd_unit
482
483 push_sub(pcm_bem_init)
484
485 if (which_eom == pcm_electrons) then
486 safe_allocate(matq0(1:nts_act, 1:nts_act))
487 safe_allocate(matqd(1:nts_act, 1:nts_act))
488 safe_allocate(matqv(1:nts_act, 1:nts_act))
489 if (.not. allocated(matqq)) then
490 safe_allocate(matqq(1:nts_act, 1:nts_act))
491 end if
492 else if (which_eom == pcm_external_potential .or. which_eom == pcm_external_plus_kick .or. which_eom == pcm_kick) then
493 safe_allocate(matq0_lf(1:nts_act, 1:nts_act)) ! not used yet
494 safe_allocate(matqd_lf(1:nts_act, 1:nts_act))
495 safe_allocate(matqv_lf(1:nts_act, 1:nts_act))
496 if (.not. allocated(matqq)) then
497 safe_allocate(matqq(1:nts_act, 1:nts_act))
498 end if
499 end if
500 call do_pcm_propmat(namespace)
501
502 if (which_eom == pcm_electrons) then
503 pcmmat0_unit = io_open(pcm_dir//'pcm_matrix_static_from_eom.out', namespace, action='write')
504 pcmmatd_unit = io_open(pcm_dir//'pcm_matrix_dynamic_from_eom.out', namespace, action='write')
505 do jtess = 1, nts_act
506 do itess = 1, nts_act
507 write(pcmmat0_unit,*) matq0(itess, jtess)
508 write(pcmmatd_unit,*) matqd(itess, jtess)
509 end do
510 end do
511 call io_close(pcmmat0_unit)
512 call io_close(pcmmatd_unit)
513 else if (which_eom == pcm_external_potential .or. which_eom == pcm_external_plus_kick) then
514 pcmmat0_unit = io_open(pcm_dir//'pcm_matrix_static_lf_from_eom.out', namespace, action='write')
515 pcmmatd_unit = io_open(pcm_dir//'pcm_matrix_dynamic_lf_from_eom.out', namespace, action='write')
516 do jtess = 1, nts_act
517 do itess = 1, nts_act
518 write(pcmmat0_unit,*) matq0_lf(itess, jtess)
519 write(pcmmatd_unit,*) matqd_lf(itess, jtess)
520 end do
521 end do
522 call io_close(pcmmat0_unit)
523 call io_close(pcmmatd_unit)
524 end if
525
526 pop_sub(pcm_bem_init)
527 end subroutine pcm_bem_init
528
529 !------------------------------------------------------------------------------------------------------------------------------
530
531 subroutine pcm_eom_end()
532 push_sub(pcm_eom_end)
533
534 ! pcm charges
535 safe_deallocate_a(q_tp)
536 safe_deallocate_a(qext_tp)
537 safe_deallocate_a(qkick_tp)
538
539 ! increment pcm charges
540 safe_deallocate_a(dq_tp)
541 safe_deallocate_a(dqext_tp)
542 safe_deallocate_a(dqkick_tp)
543
544 ! force-like term in pcm-eom equation
545 safe_deallocate_a(force_tp)
546 safe_deallocate_a(force_qext_tp)
547 safe_deallocate_a(force_qkick_tp)
548
549 ! pcm-eom bem matrices
550 safe_deallocate_a(matq0)
551 safe_deallocate_a(matqd)
552 safe_deallocate_a(matq0_lf)
553 safe_deallocate_a(matqd_lf)
554 safe_deallocate_a(matqv)
555 safe_deallocate_a(matqv_lf)
556 safe_deallocate_a(matqq)
557
558 safe_deallocate_a(pot_tp)
559 safe_deallocate_a(cts_act)
560
562
563 pop_sub(pcm_eom_end)
564 end subroutine pcm_eom_end
565
566 !------------------------------------------------------------------------------------------------------------------------------
573 subroutine do_pcm_propmat(namespace)
574 type(namespace_t), intent(in) :: namespace
575 integer, save :: i, j
576 real(real64), allocatable, save :: scr4(:,:), scr1(:,:)
577 real(real64), allocatable, save :: scr2(:,:), scr3(:,:)
578 real(real64), allocatable, save :: fact1(:), fact2(:)
581 real(real64), allocatable, save :: Kdiag0(:), Kdiagd(:)
586 real(real64), save :: sgn,sgn_lf, fac_eps0, fac_epsd
587 real(real64), save :: temp
588
589 logical, save :: firsttime = .true.
590
591 push_sub(do_pcm_propmat)
592
593 sgn = m_one ! In the case of NP this is -1 because the normal to the cavity is always pointing outward
594 ! and the field acreated by a positive unit charge outside the cavity is directed inward.
595
596 sgn_lf = m_one
597 ! ''local field'' differ from ''standard'' PCM response matrix in some sign changes
598 if (which_eom == pcm_external_potential .or. which_eom == pcm_external_plus_kick .or. which_eom == pcm_kick) sgn_lf = -m_one
599
600 safe_allocate(cals(1:nts_act, 1:nts_act))
601 safe_allocate(cald(1:nts_act, 1:nts_act))
602 safe_allocate(kdiag0(1:nts_act))
603 safe_allocate(kdiagd(1:nts_act))
604 safe_allocate(fact1(1:nts_act))
605 safe_allocate(fact2(1:nts_act))
606
607 ! generate Calderon S and D matrices
608 do i = 1, nts_act
609 do j = 1, nts_act
610 call green_s(i, j, temp)
611 cals(i, j) = temp
612 call green_d(i, j, temp)
613 cald(i, j) = temp
614 end do
615 end do
616 if (firsttime) then
617 call allocate_ts_matrix()
618 call do_ts_matrix(namespace)
619 firsttime = .false.
620 end if
621 safe_deallocate_a(cals)
622 safe_deallocate_a(cald)
623
624 safe_allocate(scr4(1:nts_act, 1:nts_act))
625 safe_allocate(scr1(1:nts_act, 1:nts_act))
626 safe_allocate(scr2(1:nts_act, 1:nts_act))
627 safe_allocate(scr3(1:nts_act, 1:nts_act))
628
629 if (which_eps == pcm_debye_model) then
630 if (.not. is_close(deb%eps_0, m_one)) then
631 fac_eps0 = (deb%eps_0 + m_one)/(deb%eps_0 - m_one)
632 kdiag0(:) = sgn_lf*(m_two*m_pi - sgn*sgn_lf*eigv(:))/(m_two*m_pi*fac_eps0 - sgn*eigv(:)) ! Eq.(14) with eps_0 in Ref.1
633 else
634 kdiag0(:) = m_zero
635 end if
636 if (.not. is_close(deb%eps_d, m_one)) then
637 fac_epsd = (deb%eps_d + m_one)/(deb%eps_d - m_one)
638 kdiagd(:) = sgn_lf*(m_two*m_pi - sgn*sgn_lf*eigv(:))/(m_two*m_pi*fac_epsd - sgn*eigv(:)) ! Eq.(14) with eps_d, ibid.
639 else
640 kdiagd(:) = m_zero
641 end if
642 fact1(:) = ((m_two*m_pi - sgn*eigv(:))*deb%eps_0 + m_two*m_pi + eigv(:))/ & ! inverse of Eq.(32), ibid.
643 ((m_two*m_pi - sgn*eigv(:))*deb%eps_d + m_two*m_pi + eigv(:))/deb%tau
644 fact2(:) = kdiag0(:)*fact1(:) ! tau^{-1}K_0 in Eq.(38), ibid.
645
646 else if (which_eps == pcm_drude_model) then
647 kdiagd(:) = m_zero ! from Eq.(10) up in Ref.2
648 fact2(:) = (m_two*m_pi - sgn*eigv(:))*drl%aa/(m_four*m_pi) ! Eq.(10) down
649 do i = 1, nts_act
650 if (fact2(i) < m_zero)fact2(i) = m_zero ! check out
651 end do
652 if (abs(drl%w0) <= m_epsilon) drl%w0 = 1.0e-8_real64 ! check out
653 fact1(:) = fact2(:) + drl%w0*drl%w0 ! Eq.(19), ibid.
654 fact2(:) = sgn_lf*(m_two*m_pi - sgn*sgn_lf*eigv(:))*drl%aa/(m_four*m_pi) ! Eq.(10) down, local field analogous
655 kdiag0(:) = fact2(:)/fact1(:) ! from Eq.(10) up, ibid.
656 end if
657
658 scr3 = matmul(sm12, eigt)
659 scr2 = matmul(transpose(eigt), sp12)
660 scr4 = transpose(scr3)
661 do i = 1, nts_act
662 scr1(:,i) = scr3(:, i)*fact1(i)
663 end do
664 matqq = matmul(scr1, scr2) ! Eq.(39) in Ref.1 and Eq.(16) in Ref.2
665 do i = 1,nts_act
666 scr1(:,i) = scr3(:, i)*fact2(i)
667 end do
668 if (which_eom == pcm_electrons) then
669 matqv = -matmul(scr1, scr4) ! Eq.(38) in Ref.1 and Eq.(17) in Ref.2
670 else if (which_eom == pcm_external_potential .or. which_eom == pcm_external_plus_kick .or. which_eom == pcm_kick) then
671 matqv_lf = -matmul(scr1, scr4) ! local field analogous
672 end if
673 do i = 1, nts_act
674 scr1(:,i) = scr3(:,i)*kdiag0(i)
675 end do
676 if (which_eom == pcm_electrons) then
677 matq0 = -matmul(scr1, scr4) ! from Eq.(14) and (18) for eps_0 in Ref.1
678 else if (which_eom == pcm_external_potential .or. which_eom == pcm_external_plus_kick .or. which_eom == pcm_kick) then
679 matq0_lf = -matmul(scr1, scr4) ! local field analogous (not used yet)
680 end if
681 do i = 1, nts_act
682 scr1(:,i) = scr3(:, i)*kdiagd(i)
683 end do
684 if (which_eom == pcm_electrons) then
685 matqd = -matmul(scr1, scr4) ! from Eq.(14) and (18) for eps_d, ibid.
686 else if (which_eom == pcm_external_potential .or. which_eom == pcm_external_plus_kick .or. which_eom == pcm_kick) then
687 matqd_lf = -matmul(scr1, scr4) ! local field analogous
688 end if
689
690 safe_deallocate_a(scr4)
691 safe_deallocate_a(scr1)
692 safe_deallocate_a(scr2)
693 safe_deallocate_a(scr3)
694 safe_deallocate_a(fact1)
695 safe_deallocate_a(fact2)
696 safe_deallocate_a(kdiag0)
697 safe_deallocate_a(kdiagd)
698
699 if (enough_initial) then
701 end if
702
703 pop_sub(do_pcm_propmat)
704 end subroutine do_pcm_propmat
705
706 !------------------------------------------------------------------------------------------------------------------------------
711 subroutine green_d(i, j, value)
712 integer, intent(in) :: i, j
713 real(real64), intent(out) :: value
714
715 real(real64) :: dist,diff(3)
716
717 push_sub(green_d)
718
719 if (i /= j) then
720 diff = cts_act(i)%point - cts_act(j)%point
721 dist = norm2(diff)
722 value = dot_product(cts_act(j)%normal, diff)/dist**3 ! Eq.(5) in Refs.1-2
723 else
724 value = -1.0694_real64*sqrt(m_four*m_pi*cts_act(i)%area)
725 value = value/(m_two*cts_act(i)%r_sphere)/cts_act(i)%area ! diagonal part is a bit cumbersome
726 end if
727
728 pop_sub(green_d)
729 end subroutine green_d
730
731 !------------------------------------------------------------------------------------------------------------------------------
736 subroutine green_s(i, j, value)
737 integer, intent(in):: i,j
738 real(real64), intent(out) :: value
739
740 real(real64) :: dist,diff(3)
741
742 push_sub(green_s)
743
744 if (i /= j) then
745 diff = cts_act(i)%point - cts_act(j)%point
746 dist = norm2(diff)
747 value = m_one/dist ! Eq.(5) in Refs.1-2
748 else
749 value = 1.0694_real64*sqrt(m_four*m_pi/cts_act(i)%area) ! diagonal part is a bit cumbersome
750 end if
751
752 pop_sub(green_s)
753 end subroutine green_s
754
755 !------------------------------------------------------------------------------------------------------------------------------
756 subroutine allocate_ts_matrix()
757 push_sub(allocate_ts_matrix)
758
759 safe_allocate(eigv(1:nts_act))
760 safe_allocate(eigt(1:nts_act, 1:nts_act))
761 safe_allocate(sm12(1:nts_act, 1:nts_act))
762 safe_allocate(sp12(1:nts_act, 1:nts_act))
763
764 pop_sub(allocate_ts_matrix)
765 end subroutine allocate_ts_matrix
766
767 !------------------------------------------------------------------------------------------------------------------------------
768 subroutine deallocate_ts_matrix()
769 push_sub(deallocate_ts_matrix)
770
771 safe_deallocate_a(eigv)
772 safe_deallocate_a(eigt)
773 safe_deallocate_a(sm12)
774 safe_deallocate_a(sp12)
775
776 pop_sub(deallocate_ts_matrix)
777 end subroutine deallocate_ts_matrix
778
779 !------------------------------------------------------------------------------------------------------------------------------
781 subroutine do_ts_matrix(namespace)
782 type(namespace_t), intent(in) :: namespace
783 integer :: i, j
784 integer :: info, lwork, liwork
785 real(real64), allocatable :: scr1(:,:), scr2(:,:), eigt_t(:,:)
786 real(real64) :: sgn
787 character jobz, uplo
788 integer, allocatable :: iwork(:)
789 real(real64), allocatable :: work(:)
790
791 push_sub(do_ts_matrix)
792
793 safe_allocate(scr1(1:nts_act, 1:nts_act))
794 safe_allocate(scr2(1:nts_act, 1:nts_act))
795 safe_allocate(eigt_t(1:nts_act, 1:nts_act))
796 safe_allocate(work(1:1 + 6*nts_act + 2*nts_act*nts_act))
797 safe_allocate(iwork(1:3 + 5*nts_act))
798
799 sgn = m_one
800
801 jobz = 'V'
802 uplo = 'U'
803 lwork = 1 + 6*nts_act + 2*nts_act*nts_act
804 liwork = 3 + 5*nts_act
805 eigt = cals
806 call dsyevd(jobz, uplo, nts_act, eigt, nts_act, eigv, work, lwork, iwork, liwork, info)
807 do i = 1, nts_act
808 if (eigv(i) < m_zero) then
809 write(message(1),*) "Eigenvalue ", i, " of S when constructing the TS matrix is negative!"
810 write(message(2),*) "I put it to 1e-8"
811 call messages_warning(2, namespace=namespace)
812 eigv(i) = 1.0e-8_real64
813 end if
814 scr1(:,i) = eigt(:,i)*sqrt(eigv(i))
815 end do
816 eigt_t = transpose(eigt)
817 sp12 = matmul(scr1, eigt_t) ! building S^{1/2} to be used in R (Ref.1) and Q_{\omega} (Ref.2)
818 do i = 1, nts_act
819 scr1(:, i) = eigt(:, i)/sqrt(eigv(i))
820 end do
821 sm12 = matmul(scr1, eigt_t) ! building S^{-1/2} to be used in R and \tilde{Q} (Ref.1) and Q_{\omega} and Q_f (Ref.2)
822 do i = 1, nts_act
823 scr1(:,i) = cald(:, i)*cts_act(i)%area
824 end do
825 scr2 = matmul(matmul(sm12, scr1), sp12) ! Eq.(10) in Ref.1 (paragraph after Eq.(9), Ref.2)
826 do j = 1, nts_act
827 do i = 1, nts_act
828 eigt(i, j) = m_half*(scr2(i, j) + scr2(j, i)) ! re-symmetrizing for numerical reasons
829 end do
830 end do
831 call dsyevd(jobz,uplo,nts_act,eigt,nts_act,eigv,work,lwork, &
832 iwork,liwork,info) ! obtaining \Lambda (eigv) and T (eigt), Eq.(10), ibid.
833
834 safe_deallocate_a(work)
835 safe_deallocate_a(iwork)
836 safe_deallocate_a(scr1)
837 safe_deallocate_a(scr2)
838 safe_deallocate_a(eigt_t)
839
840 pop_sub(do_ts_matrix)
841 end subroutine do_ts_matrix
842
843 ! -----------------------------------------------------------------------------
844 subroutine pcm_eom_enough_initial(not_yet_called)
845 logical, intent(out) :: not_yet_called
846
847 push_sub(pcm_eom_enough_initial)
848
849 enough_initial = .true.
850 not_yet_called = .false.
853 end subroutine pcm_eom_enough_initial
854
855end module pcm_eom_oct_m
856
857!! Local Variables:
858!! mode: f90
859!! coding: utf-8
860!! End:
double sqrt(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_two
Definition: global.F90:192
real(real64), parameter, public m_zero
Definition: global.F90:190
real(real64), parameter, public m_four
Definition: global.F90:194
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:188
character(len= *), parameter, public pcm_dir
Definition: global.F90:264
real(real64), parameter, public m_epsilon
Definition: global.F90:206
real(real64), parameter, public m_half
Definition: global.F90:196
real(real64), parameter, public m_one
Definition: global.F90:191
Definition: io.F90:116
subroutine, public io_close(iunit, grp)
Definition: io.F90:466
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:401
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:531
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:416
integer, parameter, public pcm_external_plus_kick
Definition: pcm_eom.F90:174
subroutine pcm_bem_init(namespace)
Boundary Element Method (BEM) EOM-IEF-PCM matrices initialization.
Definition: pcm_eom.F90:574
integer, parameter, public pcm_electrons
Definition: pcm_eom.F90:174
integer, parameter, public pcm_kick
Definition: pcm_eom.F90:174
subroutine do_ts_matrix(namespace)
Subroutine to build matrices , , and (notation of Refs.1-2)
Definition: pcm_eom.F90:877
integer, parameter, public pcm_nuclei
Definition: pcm_eom.F90:174
subroutine do_pcm_propmat(namespace)
Subroutine to build the required BEM matrices for the EOM-IEF-PCM for Debye and Drude-Lorentz cases.
Definition: pcm_eom.F90:669
subroutine init_vv_propagator
Subroutine to initialize numerical constants required by the Velocity-Verlet (VV) algorithm.
Definition: pcm_eom.F90:517
subroutine, public pcm_charges_propagation(q_t, pot_t, this_dt, this_cts_act, input_asc, this_eom, this_eps, namespace, this_deb, this_drl)
Driving subroutine for the Equation of Motion (EOM) propagation of the polarization charges within th...
Definition: pcm_eom.F90:235
subroutine pcm_ief_prop_deb(q_t, pot_t)
Euler method for integrating first order EOM for the polarization charges within IEF-PCM in the case ...
Definition: pcm_eom.F90:474
subroutine green_d(i, j, value)
Subroutine to build BEM matrix corresponding to the Calderon operator D, using the Green function of ...
Definition: pcm_eom.F90:807
integer, parameter, public pcm_external_potential
Definition: pcm_eom.F90:174
subroutine deallocate_ts_matrix()
Definition: pcm_eom.F90:864
subroutine init_charges(q_t, pot_t, namespace)
Polarization charges initialization (in equilibrium with the initial potential for electrons)
Definition: pcm_eom.F90:395
integer, parameter, public pcm_drude_model
Definition: pcm_eom.F90:170
subroutine, public pcm_eom_end()
Definition: pcm_eom.F90:627
subroutine pcm_charges_from_input_file(q_t, pot_t, namespace)
Polarization charges initialization from input file.
Definition: pcm_eom.F90:349
subroutine allocate_ts_matrix()
Definition: pcm_eom.F90:852
subroutine green_s(i, j, value)
Subroutine to build BEM matrix corresponding to the Calderon operator S, using the Green function of ...
Definition: pcm_eom.F90:832
subroutine, public pcm_eom_enough_initial(not_yet_called)
Definition: pcm_eom.F90:940
subroutine pcm_ief_prop_vv_ief_drl(q_t, pot_t)
VV algorithm for integrating second order EOM for the polarization charges within IEF-PCM in the case...
Definition: pcm_eom.F90:532
integer, parameter, public pcm_debye_model
Definition: pcm_eom.F90:170
set of parameters for Debye dielectric model
Definition: pcm_eom.F90:155
set of parameters for Drude-Lorentz dielectric model
Definition: pcm_eom.F90:164
tesselation derived type
Definition: pcm_eom.F90:144
int true(void)