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