-
Notifications
You must be signed in to change notification settings - Fork 1
/
ewalds_multipole.F
2159 lines (2007 loc) · 102 KB
/
ewalds_multipole.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!--------------------------------------------------------------------------------------------------!
! CP2K: A general program to perform molecular dynamics simulations !
! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
! !
! SPDX-License-Identifier: GPL-2.0-or-later !
!--------------------------------------------------------------------------------------------------!
! **************************************************************************************************
!> \brief Treats the electrostatic for multipoles (up to quadrupoles)
!> \author Teodoro Laino [tlaino] - 12.2007 - University of Zurich
!> inclusion of optional electric field damping for the polarizable atoms
!> Rodolphe Vuilleumier and Mathieu Salanne - 12.2009
! **************************************************************************************************
MODULE ewalds_multipole
USE atomic_kind_types, ONLY: atomic_kind_type
USE bibliography, ONLY: Aguado2003, &
Laino2008, &
cite_reference
USE cell_types, ONLY: cell_type, &
pbc
USE cp_log_handling, ONLY: cp_get_default_logger, &
cp_logger_type
USE damping_dipole_types, ONLY: damping_type, &
no_damping, &
tang_toennies
USE dg_rho0_types, ONLY: dg_rho0_type
USE dg_types, ONLY: dg_get, &
dg_type
USE distribution_1d_types, ONLY: distribution_1d_type
USE ewald_environment_types, ONLY: ewald_env_get, &
ewald_environment_type
USE ewald_pw_types, ONLY: ewald_pw_get, &
ewald_pw_type
USE fist_neighbor_list_control, ONLY: list_control
USE fist_neighbor_list_types, ONLY: fist_neighbor_type, &
neighbor_kind_pairs_type
USE fist_nonbond_env_types, ONLY: fist_nonbond_env_get, &
fist_nonbond_env_type, &
pos_type
USE input_section_types, ONLY: section_vals_type
USE kinds, ONLY: dp
USE mathconstants, ONLY: fourpi, &
oorootpi, &
pi, &
sqrthalf
USE message_passing, ONLY: mp_comm_type
USE parallel_rng_types, ONLY: UNIFORM, &
rng_stream_type
USE particle_types, ONLY: particle_type
USE pw_grid_types, ONLY: pw_grid_type
USE pw_pool_types, ONLY: pw_pool_type
USE structure_factor_types, ONLY: structure_factor_type
USE structure_factors, ONLY: structure_factor_allocate, &
structure_factor_deallocate, &
structure_factor_evaluate
#include "./base/base_uses.f90"
#:include "ewalds_multipole_sr.fypp"
IMPLICIT NONE
PRIVATE
LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .FALSE.
LOGICAL, PRIVATE, PARAMETER :: debug_r_space = .FALSE.
LOGICAL, PRIVATE, PARAMETER :: debug_g_space = .FALSE.
LOGICAL, PRIVATE, PARAMETER :: debug_e_field = .FALSE.
LOGICAL, PRIVATE, PARAMETER :: debug_e_field_en = .FALSE.
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'ewalds_multipole'
PUBLIC :: ewald_multipole_evaluate
CONTAINS
! **************************************************************************************************
!> \brief Computes the potential and the force for a lattice sum of multipoles (up to quadrupole)
!> \param ewald_env ...
!> \param ewald_pw ...
!> \param nonbond_env ...
!> \param cell ...
!> \param particle_set ...
!> \param local_particles ...
!> \param energy_local ...
!> \param energy_glob ...
!> \param e_neut ...
!> \param e_self ...
!> \param task ...
!> \param do_correction_bonded ...
!> \param do_forces ...
!> \param do_stress ...
!> \param do_efield ...
!> \param radii ...
!> \param charges ...
!> \param dipoles ...
!> \param quadrupoles ...
!> \param forces_local ...
!> \param forces_glob ...
!> \param pv_local ...
!> \param pv_glob ...
!> \param efield0 ...
!> \param efield1 ...
!> \param efield2 ...
!> \param iw ...
!> \param do_debug ...
!> \param atomic_kind_set ...
!> \param mm_section ...
!> \par Note
!> atomic_kind_set and mm_section are between the arguments only
!> for debug purpose (therefore optional) and can be avoided when this
!> function is called in other part of the program
!> \par Note
!> When a gaussian multipole is used instead of point multipole, i.e.
!> when radii(i)>0, the electrostatic fields (efield0, efield1, efield2)
!> become derivatives of the electrostatic potential energy towards
!> these gaussian multipoles.
!> \author Teodoro Laino [tlaino] - 12.2007 - University of Zurich
! **************************************************************************************************
RECURSIVE SUBROUTINE ewald_multipole_evaluate(ewald_env, ewald_pw, nonbond_env, &
cell, particle_set, local_particles, energy_local, energy_glob, e_neut, e_self, &
task, do_correction_bonded, do_forces, do_stress, &
do_efield, radii, charges, dipoles, &
quadrupoles, forces_local, forces_glob, pv_local, pv_glob, efield0, efield1, &
efield2, iw, do_debug, atomic_kind_set, mm_section)
TYPE(ewald_environment_type), POINTER :: ewald_env
TYPE(ewald_pw_type), POINTER :: ewald_pw
TYPE(fist_nonbond_env_type), POINTER :: nonbond_env
TYPE(cell_type), POINTER :: cell
TYPE(particle_type), POINTER :: particle_set(:)
TYPE(distribution_1d_type), POINTER :: local_particles
REAL(KIND=dp), INTENT(INOUT) :: energy_local, energy_glob
REAL(KIND=dp), INTENT(OUT) :: e_neut, e_self
LOGICAL, DIMENSION(3), INTENT(IN) :: task
LOGICAL, INTENT(IN) :: do_correction_bonded, do_forces, &
do_stress, do_efield
REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: radii, charges
REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: dipoles
REAL(KIND=dp), DIMENSION(:, :, :), OPTIONAL, &
POINTER :: quadrupoles
REAL(KIND=dp), DIMENSION(:, :), INTENT(INOUT), &
OPTIONAL :: forces_local, forces_glob, pv_local, &
pv_glob
REAL(KIND=dp), DIMENSION(:), INTENT(OUT), OPTIONAL :: efield0
REAL(KIND=dp), DIMENSION(:, :), INTENT(OUT), &
OPTIONAL :: efield1, efield2
INTEGER, INTENT(IN) :: iw
LOGICAL, INTENT(IN) :: do_debug
TYPE(atomic_kind_type), DIMENSION(:), OPTIONAL, &
POINTER :: atomic_kind_set
TYPE(section_vals_type), OPTIONAL, POINTER :: mm_section
CHARACTER(len=*), PARAMETER :: routineN = 'ewald_multipole_evaluate'
INTEGER :: handle, i, j, size1, size2
LOGICAL :: check_debug, check_efield, check_forces, &
do_task(3)
LOGICAL, DIMENSION(3, 3) :: my_task
REAL(KIND=dp) :: e_bonded, e_bonded_t, e_rspace, &
e_rspace_t, energy_glob_t
REAL(KIND=dp), DIMENSION(:), POINTER :: efield0_lr, efield0_sr
REAL(KIND=dp), DIMENSION(:, :), POINTER :: efield1_lr, efield1_sr, efield2_lr, &
efield2_sr
TYPE(mp_comm_type) :: group
CALL cite_reference(Aguado2003)
CALL cite_reference(Laino2008)
CALL timeset(routineN, handle)
CPASSERT(ASSOCIATED(nonbond_env))
check_debug = (debug_this_module .OR. debug_r_space .OR. debug_g_space .OR. debug_e_field .OR. debug_e_field_en) &
.EQV. debug_this_module
CPASSERT(check_debug)
check_forces = do_forces .EQV. (PRESENT(forces_local) .AND. PRESENT(forces_glob))
CPASSERT(check_forces)
check_efield = do_efield .EQV. (PRESENT(efield0) .OR. PRESENT(efield1) .OR. PRESENT(efield2))
CPASSERT(check_efield)
! Debugging this module
IF (debug_this_module .AND. do_debug) THEN
! Debug specifically real space part
IF (debug_r_space) THEN
CALL debug_ewald_multipoles(ewald_env, ewald_pw, nonbond_env, cell, &
particle_set, local_particles, iw, debug_r_space)
CPABORT("Debug Multipole Requested: Real Part!")
END IF
! Debug electric fields and gradients as pure derivatives
IF (debug_e_field) THEN
CPASSERT(PRESENT(atomic_kind_set))
CPASSERT(PRESENT(mm_section))
CALL debug_ewald_multipoles_fields(ewald_env, ewald_pw, nonbond_env, &
cell, particle_set, local_particles, radii, charges, dipoles, &
quadrupoles, task, iw, atomic_kind_set, mm_section)
CPABORT("Debug Multipole Requested: POT+EFIELDS+GRAD!")
END IF
! Debug the potential, electric fields and electric fields gradient in oder
! to retrieve the correct energy
IF (debug_e_field_en) THEN
CALL debug_ewald_multipoles_fields2(ewald_env, ewald_pw, nonbond_env, &
cell, particle_set, local_particles, radii, charges, dipoles, &
quadrupoles, task, iw)
CPABORT("Debug Multipole Requested: POT+EFIELDS+GRAD to give the correct energy!")
END IF
END IF
! Setup the tasks (needed to skip useless parts in the real-space term)
do_task = task
DO i = 1, 3
IF (do_task(i)) THEN
SELECT CASE (i)
CASE (1)
do_task(1) = ANY(charges /= 0.0_dp)
CASE (2)
do_task(2) = ANY(dipoles /= 0.0_dp)
CASE (3)
do_task(3) = ANY(quadrupoles /= 0.0_dp)
END SELECT
END IF
END DO
DO i = 1, 3
DO j = i, 3
my_task(j, i) = do_task(i) .AND. do_task(j)
my_task(i, j) = my_task(j, i)
END DO
END DO
! Allocate arrays for the evaluation of the potential, fields and electrostatic field gradients
NULLIFY (efield0_sr, efield0_lr, efield1_sr, efield1_lr, efield2_sr, efield2_lr)
IF (do_efield) THEN
IF (PRESENT(efield0)) THEN
size1 = SIZE(efield0)
ALLOCATE (efield0_sr(size1))
ALLOCATE (efield0_lr(size1))
efield0_sr = 0.0_dp
efield0_lr = 0.0_dp
END IF
IF (PRESENT(efield1)) THEN
size1 = SIZE(efield1, 1)
size2 = SIZE(efield1, 2)
ALLOCATE (efield1_sr(size1, size2))
ALLOCATE (efield1_lr(size1, size2))
efield1_sr = 0.0_dp
efield1_lr = 0.0_dp
END IF
IF (PRESENT(efield2)) THEN
size1 = SIZE(efield2, 1)
size2 = SIZE(efield2, 2)
ALLOCATE (efield2_sr(size1, size2))
ALLOCATE (efield2_lr(size1, size2))
efield2_sr = 0.0_dp
efield2_lr = 0.0_dp
END IF
END IF
e_rspace = 0.0_dp
e_bonded = 0.0_dp
IF ((.NOT. debug_g_space) .AND. (nonbond_env%do_nonbonded)) THEN
! Compute the Real Space (Short-Range) part of the Ewald sum.
! This contribution is only added when the nonbonded flag in the input
! is set, because these contributions depend. the neighborlists.
CALL ewald_multipole_SR(nonbond_env, ewald_env, atomic_kind_set, &
particle_set, cell, e_rspace, my_task, &
do_forces, do_efield, do_stress, radii, charges, dipoles, quadrupoles, &
forces_glob, pv_glob, efield0_sr, efield1_sr, efield2_sr)
energy_glob = energy_glob + e_rspace
IF (do_correction_bonded) THEN
! The corrections for bonded interactions are stored in the Real Space
! (Short-Range) part of the fields array.
CALL ewald_multipole_bonded(nonbond_env, particle_set, ewald_env, &
cell, e_bonded, my_task, do_forces, do_efield, do_stress, &
charges, dipoles, quadrupoles, forces_glob, pv_glob, &
efield0_sr, efield1_sr, efield2_sr)
energy_glob = energy_glob + e_bonded
END IF
END IF
e_neut = 0.0_dp
e_self = 0.0_dp
energy_local = 0.0_dp
IF (.NOT. debug_r_space) THEN
! Compute the Reciprocal Space (Long-Range) part of the Ewald sum
CALL ewald_multipole_LR(ewald_env, ewald_pw, cell, particle_set, &
local_particles, energy_local, my_task, do_forces, do_efield, do_stress, &
charges, dipoles, quadrupoles, forces_local, pv_local, efield0_lr, efield1_lr, &
efield2_lr)
! Self-Interactions corrections
CALL ewald_multipole_self(ewald_env, cell, local_particles, e_self, &
e_neut, my_task, do_efield, radii, charges, dipoles, quadrupoles, &
efield0_lr, efield1_lr, efield2_lr)
END IF
! Sumup energy contributions for possible IO
CALL ewald_env_get(ewald_env, group=group)
energy_glob_t = energy_glob
e_rspace_t = e_rspace
e_bonded_t = e_bonded
CALL group%sum(energy_glob_t)
CALL group%sum(e_rspace_t)
CALL group%sum(e_bonded_t)
! Print some info about energetics
CALL ewald_multipole_print(iw, energy_local, e_rspace_t, e_bonded_t, e_self, e_neut)
! Gather the components of the potential, fields and electrostatic field gradients
IF (do_efield) THEN
IF (PRESENT(efield0)) THEN
efield0 = efield0_sr + efield0_lr
CALL group%sum(efield0)
DEALLOCATE (efield0_sr)
DEALLOCATE (efield0_lr)
END IF
IF (PRESENT(efield1)) THEN
efield1 = efield1_sr + efield1_lr
CALL group%sum(efield1)
DEALLOCATE (efield1_sr)
DEALLOCATE (efield1_lr)
END IF
IF (PRESENT(efield2)) THEN
efield2 = efield2_sr + efield2_lr
CALL group%sum(efield2)
DEALLOCATE (efield2_sr)
DEALLOCATE (efield2_lr)
END IF
END IF
CALL timestop(handle)
END SUBROUTINE ewald_multipole_evaluate
! **************************************************************************************************
!> \brief computes the potential and the force for a lattice sum of multipoles
!> up to quadrupole - Short Range (Real Space) Term
!> \param nonbond_env ...
!> \param ewald_env ...
!> \param atomic_kind_set ...
!> \param particle_set ...
!> \param cell ...
!> \param energy ...
!> \param task ...
!> \param do_forces ...
!> \param do_efield ...
!> \param do_stress ...
!> \param radii ...
!> \param charges ...
!> \param dipoles ...
!> \param quadrupoles ...
!> \param forces ...
!> \param pv ...
!> \param efield0 ...
!> \param efield1 ...
!> \param efield2 ...
!> \author Teodoro Laino [tlaino] - 12.2007 - University of Zurich
! **************************************************************************************************
SUBROUTINE ewald_multipole_SR(nonbond_env, ewald_env, atomic_kind_set, &
particle_set, cell, energy, task, &
do_forces, do_efield, do_stress, radii, charges, dipoles, quadrupoles, &
forces, pv, efield0, efield1, efield2)
TYPE(fist_nonbond_env_type), POINTER :: nonbond_env
TYPE(ewald_environment_type), POINTER :: ewald_env
TYPE(atomic_kind_type), DIMENSION(:), OPTIONAL, &
POINTER :: atomic_kind_set
TYPE(particle_type), POINTER :: particle_set(:)
TYPE(cell_type), POINTER :: cell
REAL(KIND=dp), INTENT(INOUT) :: energy
LOGICAL, DIMENSION(3, 3), INTENT(IN) :: task
LOGICAL, INTENT(IN) :: do_forces, do_efield, do_stress
REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: radii, charges
REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: dipoles
REAL(KIND=dp), DIMENSION(:, :, :), OPTIONAL, &
POINTER :: quadrupoles
REAL(KIND=dp), DIMENSION(:, :), INTENT(INOUT), &
OPTIONAL :: forces, pv
REAL(KIND=dp), DIMENSION(:), POINTER :: efield0
REAL(KIND=dp), DIMENSION(:, :), POINTER :: efield1, efield2
CHARACTER(len=*), PARAMETER :: routineN = 'ewald_multipole_SR'
INTEGER :: a, atom_a, atom_b, b, c, d, e, handle, i, iend, igrp, ikind, ilist, ipair, &
istart, itype_ij, itype_ji, jkind, k, kind_a, kind_b, kk, nkdamp_ij, nkdamp_ji, nkinds, &
npairs
INTEGER, DIMENSION(:, :), POINTER :: list
LOGICAL :: do_efield0, do_efield1, do_efield2, &
force_eval
REAL(KIND=dp) :: alpha, beta, ch_i, ch_j, dampa_ij, dampa_ji, dampaexpi, dampaexpj, &
dampfac_ij, dampfac_ji, dampfuncdiffi, dampfuncdiffj, dampfunci, dampfuncj, dampsumfi, &
dampsumfj, ef0_i, ef0_j, eloc, fac, fac_ij, factorial, ir, irab2, ptens11, ptens12, &
ptens13, ptens21, ptens22, ptens23, ptens31, ptens32, ptens33, r, rab2, rab2_max, radius, &
rcut, tij, tmp, tmp1, tmp11, tmp12, tmp13, tmp2, tmp21, tmp22, tmp23, tmp31, tmp32, &
tmp33, tmp_ij, tmp_ji, xf
REAL(KIND=dp), DIMENSION(0:5) :: f
REAL(KIND=dp), DIMENSION(3) :: cell_v, cvi, damptij_a, damptji_a, dp_i, &
dp_j, ef1_i, ef1_j, fr, rab, tij_a
REAL(KIND=dp), DIMENSION(3, 3) :: damptij_ab, damptji_ab, ef2_i, ef2_j, &
qp_i, qp_j, tij_ab
REAL(KIND=dp), DIMENSION(3, 3, 3) :: tij_abc
REAL(KIND=dp), DIMENSION(3, 3, 3, 3) :: tij_abcd
REAL(KIND=dp), DIMENSION(3, 3, 3, 3, 3) :: tij_abcde
TYPE(damping_type) :: damping_ij, damping_ji
TYPE(fist_neighbor_type), POINTER :: nonbonded
TYPE(neighbor_kind_pairs_type), POINTER :: neighbor_kind_pair
TYPE(pos_type), DIMENSION(:), POINTER :: r_last_update, r_last_update_pbc
CALL timeset(routineN, handle)
NULLIFY (nonbonded, r_last_update, r_last_update_pbc)
do_efield0 = do_efield .AND. ASSOCIATED(efield0)
do_efield1 = do_efield .AND. ASSOCIATED(efield1)
do_efield2 = do_efield .AND. ASSOCIATED(efield2)
IF (do_stress) THEN
ptens11 = 0.0_dp; ptens12 = 0.0_dp; ptens13 = 0.0_dp
ptens21 = 0.0_dp; ptens22 = 0.0_dp; ptens23 = 0.0_dp
ptens31 = 0.0_dp; ptens32 = 0.0_dp; ptens33 = 0.0_dp
END IF
! Get nonbond_env info
CALL fist_nonbond_env_get(nonbond_env, nonbonded=nonbonded, natom_types=nkinds, &
r_last_update=r_last_update, r_last_update_pbc=r_last_update_pbc)
CALL ewald_env_get(ewald_env, alpha=alpha, rcut=rcut)
rab2_max = rcut**2
IF (debug_r_space) THEN
rab2_max = HUGE(0.0_dp)
END IF
! Starting the force loop
Lists: DO ilist = 1, nonbonded%nlists
neighbor_kind_pair => nonbonded%neighbor_kind_pairs(ilist)
npairs = neighbor_kind_pair%npairs
IF (npairs == 0) CYCLE
list => neighbor_kind_pair%list
cvi = neighbor_kind_pair%cell_vector
cell_v = MATMUL(cell%hmat, cvi)
Kind_Group_Loop: DO igrp = 1, neighbor_kind_pair%ngrp_kind
istart = neighbor_kind_pair%grp_kind_start(igrp)
iend = neighbor_kind_pair%grp_kind_end(igrp)
ikind = neighbor_kind_pair%ij_kind(1, igrp)
jkind = neighbor_kind_pair%ij_kind(2, igrp)
itype_ij = no_damping
nkdamp_ij = 1
dampa_ij = 1.0_dp
dampfac_ij = 0.0_dp
itype_ji = no_damping
nkdamp_ji = 1
dampa_ji = 1.0_dp
dampfac_ji = 0.0_dp
IF (PRESENT(atomic_kind_set)) THEN
IF (ASSOCIATED(atomic_kind_set(jkind)%damping)) THEN
damping_ij = atomic_kind_set(jkind)%damping%damp(ikind)
itype_ij = damping_ij%itype
nkdamp_ij = damping_ij%order
dampa_ij = damping_ij%bij
dampfac_ij = damping_ij%cij
END IF
IF (ASSOCIATED(atomic_kind_set(ikind)%damping)) THEN
damping_ji = atomic_kind_set(ikind)%damping%damp(jkind)
itype_ji = damping_ji%itype
nkdamp_ji = damping_ji%order
dampa_ji = damping_ji%bij
dampfac_ji = damping_ji%cij
END IF
END IF
Pairs: DO ipair = istart, iend
IF (ipair <= neighbor_kind_pair%nscale) THEN
! scale the electrostatic interaction if needed
! (most often scaled to zero)
fac_ij = neighbor_kind_pair%ei_scale(ipair)
IF (fac_ij <= 0) CYCLE
ELSE
fac_ij = 1.0_dp
END IF
atom_a = list(1, ipair)
atom_b = list(2, ipair)
kind_a = particle_set(atom_a)%atomic_kind%kind_number
kind_b = particle_set(atom_b)%atomic_kind%kind_number
IF (atom_a == atom_b) fac_ij = 0.5_dp
rab = r_last_update_pbc(atom_b)%r - r_last_update_pbc(atom_a)%r
rab = rab + cell_v
rab2 = rab(1)**2 + rab(2)**2 + rab(3)**2
IF (rab2 <= rab2_max) THEN
IF (PRESENT(radii)) THEN
radius = SQRT(radii(atom_a)*radii(atom_a) + radii(atom_b)*radii(atom_b))
ELSE
radius = 0.0_dp
END IF
IF (radius > 0.0_dp) THEN
beta = sqrthalf/radius
$: ewalds_multipole_sr_macro(mode="SCREENED_COULOMB_GAUSS", damping=False, store_energy=True, store_forces=True)
ELSE
$: ewalds_multipole_sr_macro(mode="SCREENED_COULOMB_ERFC", damping=True, store_energy=True, store_forces=True )
END IF
END IF
END DO Pairs
END DO Kind_Group_Loop
END DO Lists
IF (do_stress) THEN
pv(1, 1) = pv(1, 1) + ptens11
pv(1, 2) = pv(1, 2) + (ptens12 + ptens21)*0.5_dp
pv(1, 3) = pv(1, 3) + (ptens13 + ptens31)*0.5_dp
pv(2, 1) = pv(1, 2)
pv(2, 2) = pv(2, 2) + ptens22
pv(2, 3) = pv(2, 3) + (ptens23 + ptens32)*0.5_dp
pv(3, 1) = pv(1, 3)
pv(3, 2) = pv(2, 3)
pv(3, 3) = pv(3, 3) + ptens33
END IF
CALL timestop(handle)
END SUBROUTINE ewald_multipole_SR
! **************************************************************************************************
!> \brief computes the bonded correction for the potential and the force for a
!> lattice sum of multipoles up to quadrupole
!> \param nonbond_env ...
!> \param particle_set ...
!> \param ewald_env ...
!> \param cell ...
!> \param energy ...
!> \param task ...
!> \param do_forces ...
!> \param do_efield ...
!> \param do_stress ...
!> \param charges ...
!> \param dipoles ...
!> \param quadrupoles ...
!> \param forces ...
!> \param pv ...
!> \param efield0 ...
!> \param efield1 ...
!> \param efield2 ...
!> \author Teodoro Laino [tlaino] - 05.2009
! **************************************************************************************************
SUBROUTINE ewald_multipole_bonded(nonbond_env, particle_set, ewald_env, &
cell, energy, task, do_forces, do_efield, do_stress, charges, &
dipoles, quadrupoles, forces, pv, efield0, efield1, efield2)
TYPE(fist_nonbond_env_type), POINTER :: nonbond_env
TYPE(particle_type), POINTER :: particle_set(:)
TYPE(ewald_environment_type), POINTER :: ewald_env
TYPE(cell_type), POINTER :: cell
REAL(KIND=dp), INTENT(INOUT) :: energy
LOGICAL, DIMENSION(3, 3), INTENT(IN) :: task
LOGICAL, INTENT(IN) :: do_forces, do_efield, do_stress
REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: charges
REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: dipoles
REAL(KIND=dp), DIMENSION(:, :, :), OPTIONAL, &
POINTER :: quadrupoles
REAL(KIND=dp), DIMENSION(:, :), INTENT(INOUT), &
OPTIONAL :: forces, pv
REAL(KIND=dp), DIMENSION(:), POINTER :: efield0
REAL(KIND=dp), DIMENSION(:, :), POINTER :: efield1, efield2
CHARACTER(len=*), PARAMETER :: routineN = 'ewald_multipole_bonded'
INTEGER :: a, atom_a, atom_b, b, c, d, e, handle, &
i, iend, igrp, ilist, ipair, istart, &
k, nscale
INTEGER, DIMENSION(:, :), POINTER :: list
LOGICAL :: do_efield0, do_efield1, do_efield2, &
force_eval
REAL(KIND=dp) :: alpha, ch_i, ch_j, ef0_i, ef0_j, eloc, fac, fac_ij, ir, irab2, ptens11, &
ptens12, ptens13, ptens21, ptens22, ptens23, ptens31, ptens32, ptens33, r, rab2, tij, &
tmp, tmp1, tmp11, tmp12, tmp13, tmp2, tmp21, tmp22, tmp23, tmp31, tmp32, tmp33, tmp_ij, &
tmp_ji
REAL(KIND=dp), DIMENSION(0:5) :: f
REAL(KIND=dp), DIMENSION(3) :: dp_i, dp_j, ef1_i, ef1_j, fr, rab, tij_a
REAL(KIND=dp), DIMENSION(3, 3) :: ef2_i, ef2_j, qp_i, qp_j, tij_ab
REAL(KIND=dp), DIMENSION(3, 3, 3) :: tij_abc
REAL(KIND=dp), DIMENSION(3, 3, 3, 3) :: tij_abcd
REAL(KIND=dp), DIMENSION(3, 3, 3, 3, 3) :: tij_abcde
TYPE(fist_neighbor_type), POINTER :: nonbonded
TYPE(neighbor_kind_pairs_type), POINTER :: neighbor_kind_pair
CALL timeset(routineN, handle)
do_efield0 = do_efield .AND. ASSOCIATED(efield0)
do_efield1 = do_efield .AND. ASSOCIATED(efield1)
do_efield2 = do_efield .AND. ASSOCIATED(efield2)
IF (do_stress) THEN
ptens11 = 0.0_dp; ptens12 = 0.0_dp; ptens13 = 0.0_dp
ptens21 = 0.0_dp; ptens22 = 0.0_dp; ptens23 = 0.0_dp
ptens31 = 0.0_dp; ptens32 = 0.0_dp; ptens33 = 0.0_dp
END IF
CALL ewald_env_get(ewald_env, alpha=alpha)
CALL fist_nonbond_env_get(nonbond_env, nonbonded=nonbonded)
! Starting the force loop
Lists: DO ilist = 1, nonbonded%nlists
neighbor_kind_pair => nonbonded%neighbor_kind_pairs(ilist)
nscale = neighbor_kind_pair%nscale
IF (nscale == 0) CYCLE
list => neighbor_kind_pair%list
Kind_Group_Loop: DO igrp = 1, neighbor_kind_pair%ngrp_kind
istart = neighbor_kind_pair%grp_kind_start(igrp)
IF (istart > nscale) CYCLE
iend = MIN(neighbor_kind_pair%grp_kind_end(igrp), nscale)
Pairs: DO ipair = istart, iend
! only use pairs that are (partially) excluded for electrostatics
fac_ij = -1.0_dp + neighbor_kind_pair%ei_scale(ipair)
IF (fac_ij >= 0) CYCLE
atom_a = list(1, ipair)
atom_b = list(2, ipair)
rab = particle_set(atom_b)%r - particle_set(atom_a)%r
rab = pbc(rab, cell)
rab2 = rab(1)**2 + rab(2)**2 + rab(3)**2
$: ewalds_multipole_sr_macro(mode="SCREENED_COULOMB_ERF", store_energy=True, store_forces=True)
END DO Pairs
END DO Kind_Group_Loop
END DO Lists
IF (do_stress) THEN
pv(1, 1) = pv(1, 1) + ptens11
pv(1, 2) = pv(1, 2) + (ptens12 + ptens21)*0.5_dp
pv(1, 3) = pv(1, 3) + (ptens13 + ptens31)*0.5_dp
pv(2, 1) = pv(1, 2)
pv(2, 2) = pv(2, 2) + ptens22
pv(2, 3) = pv(2, 3) + (ptens23 + ptens32)*0.5_dp
pv(3, 1) = pv(1, 3)
pv(3, 2) = pv(2, 3)
pv(3, 3) = pv(3, 3) + ptens33
END IF
CALL timestop(handle)
END SUBROUTINE ewald_multipole_bonded
! **************************************************************************************************
!> \brief computes the potential and the force for a lattice sum of multipoles
!> up to quadrupole - Long Range (Reciprocal Space) Term
!> \param ewald_env ...
!> \param ewald_pw ...
!> \param cell ...
!> \param particle_set ...
!> \param local_particles ...
!> \param energy ...
!> \param task ...
!> \param do_forces ...
!> \param do_efield ...
!> \param do_stress ...
!> \param charges ...
!> \param dipoles ...
!> \param quadrupoles ...
!> \param forces ...
!> \param pv ...
!> \param efield0 ...
!> \param efield1 ...
!> \param efield2 ...
!> \author Teodoro Laino [tlaino] - 12.2007 - University of Zurich
! **************************************************************************************************
SUBROUTINE ewald_multipole_LR(ewald_env, ewald_pw, cell, particle_set, &
local_particles, energy, task, do_forces, do_efield, do_stress, &
charges, dipoles, quadrupoles, forces, pv, efield0, efield1, efield2)
TYPE(ewald_environment_type), POINTER :: ewald_env
TYPE(ewald_pw_type), POINTER :: ewald_pw
TYPE(cell_type), POINTER :: cell
TYPE(particle_type), POINTER :: particle_set(:)
TYPE(distribution_1d_type), POINTER :: local_particles
REAL(KIND=dp), INTENT(INOUT) :: energy
LOGICAL, DIMENSION(3, 3), INTENT(IN) :: task
LOGICAL, INTENT(IN) :: do_forces, do_efield, do_stress
REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: charges
REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: dipoles
REAL(KIND=dp), DIMENSION(:, :, :), OPTIONAL, &
POINTER :: quadrupoles
REAL(KIND=dp), DIMENSION(:, :), INTENT(INOUT), &
OPTIONAL :: forces, pv
REAL(KIND=dp), DIMENSION(:), POINTER :: efield0
REAL(KIND=dp), DIMENSION(:, :), POINTER :: efield1, efield2
CHARACTER(len=*), PARAMETER :: routineN = 'ewald_multipole_LR'
COMPLEX(KIND=dp) :: atm_factor, atm_factor_st(3), cnjg_fac, &
fac, summe_tmp
COMPLEX(KIND=dp), ALLOCATABLE, DIMENSION(:) :: summe_ef
COMPLEX(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: summe_st
INTEGER :: gpt, handle, iparticle, iparticle_kind, iparticle_local, lp, mp, nnodes, &
node, np, nparticle_kind, nparticle_local
INTEGER, DIMENSION(:, :), POINTER :: bds
LOGICAL :: do_efield0, do_efield1, do_efield2
REAL(KIND=dp) :: alpha, denom, dipole_t(3), f0, factor, &
four_alpha_sq, gauss, pref, q_t, tmp, &
trq_t
REAL(KIND=dp), DIMENSION(3) :: tmp_v, vec
REAL(KIND=dp), DIMENSION(3, 3) :: pv_tmp
REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: rho0
TYPE(dg_rho0_type), POINTER :: dg_rho0
TYPE(dg_type), POINTER :: dg
TYPE(pw_grid_type), POINTER :: pw_grid
TYPE(pw_pool_type), POINTER :: pw_pool
TYPE(structure_factor_type) :: exp_igr
TYPE(mp_comm_type) :: group
CALL timeset(routineN, handle)
do_efield0 = do_efield .AND. ASSOCIATED(efield0)
do_efield1 = do_efield .AND. ASSOCIATED(efield1)
do_efield2 = do_efield .AND. ASSOCIATED(efield2)
! Gathering data from the ewald environment
CALL ewald_env_get(ewald_env, alpha=alpha, group=group)
CALL ewald_pw_get(ewald_pw, pw_big_pool=pw_pool, dg=dg)
CALL dg_get(dg, dg_rho0=dg_rho0)
rho0 => dg_rho0%density%array
pw_grid => pw_pool%pw_grid
bds => pw_grid%bounds
! Allocation of working arrays
nparticle_kind = SIZE(local_particles%n_el)
nnodes = 0
DO iparticle_kind = 1, nparticle_kind
nnodes = nnodes + local_particles%n_el(iparticle_kind)
END DO
CALL structure_factor_allocate(pw_grid%bounds, nnodes, exp_igr)
ALLOCATE (summe_ef(1:pw_grid%ngpts_cut))
summe_ef = CMPLX(0.0_dp, 0.0_dp, KIND=dp)
! Stress Tensor
IF (do_stress) THEN
pv_tmp = 0.0_dp
ALLOCATE (summe_st(3, 1:pw_grid%ngpts_cut))
summe_st = CMPLX(0.0_dp, 0.0_dp, KIND=dp)
END IF
! Defining four_alpha_sq
four_alpha_sq = 4.0_dp*alpha**2
dipole_t = 0.0_dp
q_t = 0.0_dp
trq_t = 0.0_dp
! Zero node count
node = 0
DO iparticle_kind = 1, nparticle_kind
nparticle_local = local_particles%n_el(iparticle_kind)
DO iparticle_local = 1, nparticle_local
node = node + 1
iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
vec = MATMUL(cell%h_inv, particle_set(iparticle)%r)
CALL structure_factor_evaluate(vec, exp_igr%lb, &
exp_igr%ex(:, node), exp_igr%ey(:, node), exp_igr%ez(:, node))
! Computing the total charge, dipole and quadrupole trace (if any)
IF (ANY(task(1, :))) THEN
q_t = q_t + charges(iparticle)
END IF
IF (ANY(task(2, :))) THEN
dipole_t = dipole_t + dipoles(:, iparticle)
END IF
IF (ANY(task(3, :))) THEN
trq_t = trq_t + quadrupoles(1, 1, iparticle) + &
quadrupoles(2, 2, iparticle) + &
quadrupoles(3, 3, iparticle)
END IF
END DO
END DO
! Looping over the positive g-vectors
DO gpt = 1, pw_grid%ngpts_cut_local
lp = pw_grid%mapl%pos(pw_grid%g_hat(1, gpt))
mp = pw_grid%mapm%pos(pw_grid%g_hat(2, gpt))
np = pw_grid%mapn%pos(pw_grid%g_hat(3, gpt))
lp = lp + bds(1, 1)
mp = mp + bds(1, 2)
np = np + bds(1, 3)
! Initializing sum to be used in the energy and force
node = 0
DO iparticle_kind = 1, nparticle_kind
nparticle_local = local_particles%n_el(iparticle_kind)
DO iparticle_local = 1, nparticle_local
node = node + 1
iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
! Density for energy and forces
CALL get_atom_factor(atm_factor, pw_grid, gpt, iparticle, task, charges, &
dipoles, quadrupoles)
summe_tmp = exp_igr%ex(lp, node)*exp_igr%ey(mp, node)*exp_igr%ez(np, node)
summe_ef(gpt) = summe_ef(gpt) + atm_factor*summe_tmp
! Precompute pseudo-density for stress tensor calculation
IF (do_stress) THEN
CALL get_atom_factor_stress(atm_factor_st, pw_grid, gpt, iparticle, task, &
dipoles, quadrupoles)
summe_st(1:3, gpt) = summe_st(1:3, gpt) + atm_factor_st(1:3)*summe_tmp
END IF
END DO
END DO
END DO
CALL group%sum(q_t)
CALL group%sum(dipole_t)
CALL group%sum(trq_t)
CALL group%sum(summe_ef)
IF (do_stress) CALL group%sum(summe_st)
! Looping over the positive g-vectors
DO gpt = 1, pw_grid%ngpts_cut_local
! computing the potential energy
lp = pw_grid%mapl%pos(pw_grid%g_hat(1, gpt))
mp = pw_grid%mapm%pos(pw_grid%g_hat(2, gpt))
np = pw_grid%mapn%pos(pw_grid%g_hat(3, gpt))
lp = lp + bds(1, 1)
mp = mp + bds(1, 2)
np = np + bds(1, 3)
IF (pw_grid%gsq(gpt) == 0.0_dp) THEN
! G=0 vector for dipole-dipole and charge-quadrupole
energy = energy + (1.0_dp/6.0_dp)*DOT_PRODUCT(dipole_t, dipole_t) &
- (1.0_dp/9.0_dp)*q_t*trq_t
! Stress tensor
IF (do_stress) THEN
pv_tmp(1, 1) = pv_tmp(1, 1) + (1.0_dp/6.0_dp)*DOT_PRODUCT(dipole_t, dipole_t)
pv_tmp(2, 2) = pv_tmp(2, 2) + (1.0_dp/6.0_dp)*DOT_PRODUCT(dipole_t, dipole_t)
pv_tmp(3, 3) = pv_tmp(3, 3) + (1.0_dp/6.0_dp)*DOT_PRODUCT(dipole_t, dipole_t)
END IF
! Corrections for G=0 to potential, field and field gradient
IF (do_efield .AND. (debug_e_field_en .OR. (.NOT. debug_this_module))) THEN
! This term is important and may give problems if one is debugging
! VS finite differences since it comes from a residual integral in
! the complex plane (cannot be reproduced with finite differences)
node = 0
DO iparticle_kind = 1, nparticle_kind
nparticle_local = local_particles%n_el(iparticle_kind)
DO iparticle_local = 1, nparticle_local
node = node + 1
iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
! Potential
IF (do_efield0) THEN
efield0(iparticle) = efield0(iparticle)
END IF
! Electrostatic field
IF (do_efield1) THEN
efield1(1:3, iparticle) = efield1(1:3, iparticle) - (1.0_dp/6.0_dp)*dipole_t(1:3)
END IF
! Electrostatic field gradients
IF (do_efield2) THEN
efield2(1, iparticle) = efield2(1, iparticle) - (1.0_dp/(18.0_dp))*q_t
efield2(5, iparticle) = efield2(5, iparticle) - (1.0_dp/(18.0_dp))*q_t
efield2(9, iparticle) = efield2(9, iparticle) - (1.0_dp/(18.0_dp))*q_t
END IF
END DO
END DO
END IF
CYCLE
END IF
gauss = (rho0(lp, mp, np)*pw_grid%vol)**2/pw_grid%gsq(gpt)
factor = gauss*REAL(summe_ef(gpt)*CONJG(summe_ef(gpt)), KIND=dp)
energy = energy + factor
IF (do_forces .OR. do_efield) THEN
node = 0
DO iparticle_kind = 1, nparticle_kind
nparticle_local = local_particles%n_el(iparticle_kind)
DO iparticle_local = 1, nparticle_local
node = node + 1
iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
fac = exp_igr%ex(lp, node)*exp_igr%ey(mp, node)*exp_igr%ez(np, node)
cnjg_fac = CONJG(fac)
! Forces
IF (do_forces) THEN
CALL get_atom_factor(atm_factor, pw_grid, gpt, iparticle, task, charges, &
dipoles, quadrupoles)
tmp = gauss*AIMAG(summe_ef(gpt)*(cnjg_fac*CONJG(atm_factor)))
forces(1, node) = forces(1, node) + tmp*pw_grid%g(1, gpt)
forces(2, node) = forces(2, node) + tmp*pw_grid%g(2, gpt)
forces(3, node) = forces(3, node) + tmp*pw_grid%g(3, gpt)
END IF
! Electric field
IF (do_efield) THEN
! Potential
IF (do_efield0) THEN
efield0(iparticle) = efield0(iparticle) + gauss*REAL(fac*CONJG(summe_ef(gpt)), KIND=dp)
END IF
! Electric field
IF (do_efield1) THEN
tmp = AIMAG(fac*CONJG(summe_ef(gpt)))*gauss
efield1(1, iparticle) = efield1(1, iparticle) - tmp*pw_grid%g(1, gpt)
efield1(2, iparticle) = efield1(2, iparticle) - tmp*pw_grid%g(2, gpt)
efield1(3, iparticle) = efield1(3, iparticle) - tmp*pw_grid%g(3, gpt)
END IF
! Electric field gradient
IF (do_efield2) THEN
tmp_v(1) = REAL(fac*CONJG(summe_ef(gpt)), KIND=dp)*pw_grid%g(1, gpt)*gauss
tmp_v(2) = REAL(fac*CONJG(summe_ef(gpt)), KIND=dp)*pw_grid%g(2, gpt)*gauss
tmp_v(3) = REAL(fac*CONJG(summe_ef(gpt)), KIND=dp)*pw_grid%g(3, gpt)*gauss
efield2(1, iparticle) = efield2(1, iparticle) + tmp_v(1)*pw_grid%g(1, gpt)
efield2(2, iparticle) = efield2(2, iparticle) + tmp_v(1)*pw_grid%g(2, gpt)
efield2(3, iparticle) = efield2(3, iparticle) + tmp_v(1)*pw_grid%g(3, gpt)
efield2(4, iparticle) = efield2(4, iparticle) + tmp_v(2)*pw_grid%g(1, gpt)
efield2(5, iparticle) = efield2(5, iparticle) + tmp_v(2)*pw_grid%g(2, gpt)
efield2(6, iparticle) = efield2(6, iparticle) + tmp_v(2)*pw_grid%g(3, gpt)
efield2(7, iparticle) = efield2(7, iparticle) + tmp_v(3)*pw_grid%g(1, gpt)
efield2(8, iparticle) = efield2(8, iparticle) + tmp_v(3)*pw_grid%g(2, gpt)
efield2(9, iparticle) = efield2(9, iparticle) + tmp_v(3)*pw_grid%g(3, gpt)
END IF
END IF
END DO
END DO
END IF
! Compute the virial P*V
IF (do_stress) THEN
! The Stress Tensor can be decomposed in two main components.
! The first one is just a normal ewald component for reciprocal space
denom = 1.0_dp/four_alpha_sq + 1.0_dp/pw_grid%gsq(gpt)
pv_tmp(1, 1) = pv_tmp(1, 1) + factor*(1.0_dp - 2.0_dp*pw_grid%g(1, gpt)*pw_grid%g(1, gpt)*denom)
pv_tmp(1, 2) = pv_tmp(1, 2) - factor*(2.0_dp*pw_grid%g(1, gpt)*pw_grid%g(2, gpt)*denom)
pv_tmp(1, 3) = pv_tmp(1, 3) - factor*(2.0_dp*pw_grid%g(1, gpt)*pw_grid%g(3, gpt)*denom)
pv_tmp(2, 1) = pv_tmp(2, 1) - factor*(2.0_dp*pw_grid%g(2, gpt)*pw_grid%g(1, gpt)*denom)
pv_tmp(2, 2) = pv_tmp(2, 2) + factor*(1.0_dp - 2.0_dp*pw_grid%g(2, gpt)*pw_grid%g(2, gpt)*denom)
pv_tmp(2, 3) = pv_tmp(2, 3) - factor*(2.0_dp*pw_grid%g(2, gpt)*pw_grid%g(3, gpt)*denom)
pv_tmp(3, 1) = pv_tmp(3, 1) - factor*(2.0_dp*pw_grid%g(3, gpt)*pw_grid%g(1, gpt)*denom)
pv_tmp(3, 2) = pv_tmp(3, 2) - factor*(2.0_dp*pw_grid%g(3, gpt)*pw_grid%g(2, gpt)*denom)
pv_tmp(3, 3) = pv_tmp(3, 3) + factor*(1.0_dp - 2.0_dp*pw_grid%g(3, gpt)*pw_grid%g(3, gpt)*denom)
! The second one can be written in the following way
f0 = 2.0_dp*gauss
pv_tmp(1, 1) = pv_tmp(1, 1) + f0*pw_grid%g(1, gpt)*REAL(summe_st(1, gpt)*CONJG(summe_ef(gpt)), KIND=dp)
pv_tmp(1, 2) = pv_tmp(1, 2) + f0*pw_grid%g(1, gpt)*REAL(summe_st(2, gpt)*CONJG(summe_ef(gpt)), KIND=dp)
pv_tmp(1, 3) = pv_tmp(1, 3) + f0*pw_grid%g(1, gpt)*REAL(summe_st(3, gpt)*CONJG(summe_ef(gpt)), KIND=dp)
pv_tmp(2, 1) = pv_tmp(2, 1) + f0*pw_grid%g(2, gpt)*REAL(summe_st(1, gpt)*CONJG(summe_ef(gpt)), KIND=dp)
pv_tmp(2, 2) = pv_tmp(2, 2) + f0*pw_grid%g(2, gpt)*REAL(summe_st(2, gpt)*CONJG(summe_ef(gpt)), KIND=dp)
pv_tmp(2, 3) = pv_tmp(2, 3) + f0*pw_grid%g(2, gpt)*REAL(summe_st(3, gpt)*CONJG(summe_ef(gpt)), KIND=dp)
pv_tmp(3, 1) = pv_tmp(3, 1) + f0*pw_grid%g(3, gpt)*REAL(summe_st(1, gpt)*CONJG(summe_ef(gpt)), KIND=dp)
pv_tmp(3, 2) = pv_tmp(3, 2) + f0*pw_grid%g(3, gpt)*REAL(summe_st(2, gpt)*CONJG(summe_ef(gpt)), KIND=dp)
pv_tmp(3, 3) = pv_tmp(3, 3) + f0*pw_grid%g(3, gpt)*REAL(summe_st(3, gpt)*CONJG(summe_ef(gpt)), KIND=dp)
END IF
END DO
pref = fourpi/pw_grid%vol
energy = energy*pref
CALL structure_factor_deallocate(exp_igr)
DEALLOCATE (summe_ef)
IF (do_stress) THEN
pv_tmp = pv_tmp*pref
! Symmetrize the tensor
pv(1, 1) = pv(1, 1) + pv_tmp(1, 1)
pv(1, 2) = pv(1, 2) + (pv_tmp(1, 2) + pv_tmp(2, 1))*0.5_dp
pv(1, 3) = pv(1, 3) + (pv_tmp(1, 3) + pv_tmp(3, 1))*0.5_dp
pv(2, 1) = pv(1, 2)
pv(2, 2) = pv(2, 2) + pv_tmp(2, 2)
pv(2, 3) = pv(2, 3) + (pv_tmp(2, 3) + pv_tmp(3, 2))*0.5_dp
pv(3, 1) = pv(1, 3)
pv(3, 2) = pv(2, 3)
pv(3, 3) = pv(3, 3) + pv_tmp(3, 3)
DEALLOCATE (summe_st)
END IF
IF (do_forces) THEN
forces = 2.0_dp*forces*pref
END IF
IF (do_efield0) THEN
efield0 = 2.0_dp*efield0*pref
END IF
IF (do_efield1) THEN
efield1 = 2.0_dp*efield1*pref
END IF
IF (do_efield2) THEN
efield2 = 2.0_dp*efield2*pref
END IF
CALL timestop(handle)
END SUBROUTINE ewald_multipole_LR
! **************************************************************************************************
!> \brief Computes the atom factor including charge, dipole and quadrupole
!> \param atm_factor ...
!> \param pw_grid ...
!> \param gpt ...
!> \param iparticle ...
!> \param task ...
!> \param charges ...
!> \param dipoles ...
!> \param quadrupoles ...
!> \par History
!> none
!> \author Teodoro Laino [tlaino] - 12.2007 - University of Zurich
! **************************************************************************************************
SUBROUTINE get_atom_factor(atm_factor, pw_grid, gpt, iparticle, task, charges, &
dipoles, quadrupoles)
COMPLEX(KIND=dp), INTENT(OUT) :: atm_factor
TYPE(pw_grid_type), POINTER :: pw_grid
INTEGER, INTENT(IN) :: gpt
INTEGER :: iparticle
LOGICAL, DIMENSION(3, 3), INTENT(IN) :: task
REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: charges
REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: dipoles
REAL(KIND=dp), DIMENSION(:, :, :), OPTIONAL, &
POINTER :: quadrupoles
COMPLEX(KIND=dp) :: tmp
INTEGER :: i, j
atm_factor = CMPLX(0.0_dp, 0.0_dp, KIND=dp)
IF (task(1, 1)) THEN
! Charge
atm_factor = atm_factor + charges(iparticle)
END IF
IF (task(2, 2)) THEN
! Dipole
tmp = CMPLX(0.0_dp, 0.0_dp, KIND=dp)
DO i = 1, 3
tmp = tmp + dipoles(i, iparticle)*pw_grid%g(i, gpt)
END DO
atm_factor = atm_factor + tmp*CMPLX(0.0_dp, -1.0_dp, KIND=dp)
END IF
IF (task(3, 3)) THEN
! Quadrupole