-
Notifications
You must be signed in to change notification settings - Fork 1
/
qs_collocate_density.F
2835 lines (2397 loc) · 120 KB
/
qs_collocate_density.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 Calculate the plane wave density by collocating the primitive Gaussian
!> functions (pgf).
!> \par History
!> - rewrote collocate for increased accuracy and speed
!> - introduced the PGI hack for increased speed with that compiler
!> (22.02.02)
!> - Added Multiple Grid feature
!> - new way to get over the grid (01.03.02)
!> - removed timing calls since they were getting expensive
!> - Updated with the new QS data structures (09.04.02,MK)
!> - introduction of the real space grid type ( prelim. version JVdV 05.02)
!> - parallel FFT (JGH 22.05.02)
!> - multigrid arrays independent from density (JGH 30.08.02)
!> - old density stored in g space (JGH 30.08.02)
!> - distributed real space code (JGH 17.07.03)
!> - refactoring and new loop ordering (JGH 23.11.03)
!> - OpenMP parallelization (JGH 03.12.03)
!> - Modified to compute tau (Joost 12.03)
!> - removed the incremental density rebuild (Joost 01.04)
!> - introduced realspace multigridding (Joost 02.04)
!> - introduced map_consistent (Joost 02.04)
!> - Addition of the subroutine calculate_atomic_charge_density (TdK, 08.05)
!> - rewrite of the collocate/integrate kernels (Joost VandeVondele, 03.07)
!> - Extended by the derivatives for DFPT [Sandra Luber, Edward Ditler, 2021]
!> \author Matthias Krack (03.04.2001)
!> 1) Joost VandeVondele (01.2002)
!> Thomas D. Kuehne (04.08.2005)
!> Ole Schuett (2020)
! **************************************************************************************************
MODULE qs_collocate_density
USE admm_types, ONLY: get_admm_env
USE ao_util, ONLY: exp_radius_very_extended
USE atomic_kind_types, ONLY: atomic_kind_type, &
get_atomic_kind, &
get_atomic_kind_set
USE basis_set_types, ONLY: get_gto_basis_set, &
gto_basis_set_type
USE cell_types, ONLY: cell_type, &
pbc
USE cp_control_types, ONLY: dft_control_type
USE cp_dbcsr_operations, ONLY: dbcsr_deallocate_matrix_set
USE cp_fm_types, ONLY: cp_fm_get_element, &
cp_fm_get_info, &
cp_fm_type
USE cp_dbcsr_api, ONLY: dbcsr_copy, &
dbcsr_get_block_p, &
dbcsr_p_type, &
dbcsr_type
USE external_potential_types, ONLY: get_potential, &
gth_potential_type
USE gaussian_gridlevels, ONLY: gaussian_gridlevel, &
gridlevel_info_type
USE grid_api, ONLY: &
GRID_FUNC_AB, GRID_FUNC_CORE_X, GRID_FUNC_CORE_Y, GRID_FUNC_CORE_Z, GRID_FUNC_DAB_X, &
GRID_FUNC_DAB_Y, GRID_FUNC_DAB_Z, GRID_FUNC_DABpADB_X, GRID_FUNC_DABpADB_Y, &
GRID_FUNC_DABpADB_Z, GRID_FUNC_DADB, GRID_FUNC_DX, GRID_FUNC_DXDX, GRID_FUNC_DXDY, &
GRID_FUNC_DY, GRID_FUNC_DYDY, GRID_FUNC_DYDZ, GRID_FUNC_DZ, GRID_FUNC_DZDX, &
GRID_FUNC_DZDZ, collocate_pgf_product, grid_collocate_task_list
USE input_constants, ONLY: &
orb_dx2, orb_dxy, orb_dy2, orb_dyz, orb_dz2, orb_dzx, orb_px, orb_py, orb_pz, orb_s
USE kinds, ONLY: default_string_length, &
dp
USE lri_environment_types, ONLY: lri_kind_type
USE memory_utilities, ONLY: reallocate
USE message_passing, ONLY: mp_comm_type
USE orbital_pointers, ONLY: coset, &
ncoset
USE particle_types, ONLY: particle_type
USE pw_env_types, ONLY: pw_env_get, &
pw_env_type
USE pw_methods, ONLY: pw_axpy, &
pw_integrate_function, &
pw_transfer, &
pw_zero
USE pw_pool_types, ONLY: pw_pool_p_type, &
pw_pool_type, &
pw_pools_create_pws, &
pw_pools_give_back_pws
USE pw_types, ONLY: pw_r3d_rs_type, &
pw_c1d_gs_type, &
pw_r3d_rs_type
USE qs_environment_types, ONLY: get_qs_env, &
qs_environment_type
USE qs_kind_types, ONLY: get_qs_kind, &
get_qs_kind_set, &
qs_kind_type
USE qs_ks_types, ONLY: get_ks_env, &
qs_ks_env_type
USE qs_neighbor_list_types, ONLY: neighbor_list_set_p_type
USE realspace_grid_types, ONLY: map_gaussian_here, &
realspace_grid_desc_p_type, &
realspace_grid_type, &
rs_grid_zero, &
transfer_rs2pw
USE rs_pw_interface, ONLY: density_rs2pw
USE task_list_methods, ONLY: rs_copy_to_buffer, &
rs_distribute_matrix, &
rs_scatter_matrices
USE task_list_types, ONLY: atom_pair_type, &
task_list_type, &
task_type
!$ USE OMP_LIB, ONLY: omp_get_max_threads, omp_get_thread_num, omp_get_num_threads
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_collocate_density'
! *** Public subroutines ***
PUBLIC :: calculate_ppl_grid, &
calculate_rho_core, &
calculate_lri_rho_elec, &
calculate_rho_single_gaussian, &
calculate_rho_metal, &
calculate_rho_resp_single, &
calculate_rho_resp_all, &
calculate_rho_elec, &
calculate_drho_elec, &
calculate_wavefunction, &
collocate_function, &
calculate_rho_nlcc, &
calculate_drho_elec_dR, &
calculate_drho_core, &
collocate_single_gaussian
INTERFACE calculate_rho_core
MODULE PROCEDURE calculate_rho_core_r3d_rs
MODULE PROCEDURE calculate_rho_core_c1d_gs
END INTERFACE
INTERFACE calculate_rho_resp_all
MODULE PROCEDURE calculate_rho_resp_all_r3d_rs, calculate_rho_resp_all_c1d_gs
END INTERFACE
CONTAINS
! **************************************************************************************************
!> \brief computes the density of the non-linear core correction on the grid
!> \param rho_nlcc ...
!> \param qs_env ...
! **************************************************************************************************
SUBROUTINE calculate_rho_nlcc(rho_nlcc, qs_env)
TYPE(pw_r3d_rs_type), INTENT(INOUT) :: rho_nlcc
TYPE(qs_environment_type), POINTER :: qs_env
CHARACTER(len=*), PARAMETER :: routineN = 'calculate_rho_nlcc'
INTEGER :: atom_a, handle, iatom, iexp_nlcc, ikind, &
ithread, j, n, natom, nc, nexp_nlcc, &
ni, npme, nthread, subpatch_pattern
INTEGER, DIMENSION(:), POINTER :: atom_list, cores, nct_nlcc
LOGICAL :: nlcc
REAL(KIND=dp) :: alpha, eps_rho_rspace, radius
REAL(KIND=dp), DIMENSION(3) :: ra
REAL(KIND=dp), DIMENSION(:), POINTER :: alpha_nlcc
REAL(KIND=dp), DIMENSION(:, :), POINTER :: cval_nlcc, pab
TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
TYPE(cell_type), POINTER :: cell
TYPE(dft_control_type), POINTER :: dft_control
TYPE(gth_potential_type), POINTER :: gth_potential
TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
TYPE(pw_env_type), POINTER :: pw_env
TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
TYPE(realspace_grid_type), POINTER :: rs_rho
CALL timeset(routineN, handle)
NULLIFY (cell, dft_control, pab, particle_set, atomic_kind_set, &
qs_kind_set, atom_list, pw_env, rs_rho, auxbas_pw_pool, cores)
CALL get_qs_env(qs_env=qs_env, &
atomic_kind_set=atomic_kind_set, &
qs_kind_set=qs_kind_set, &
cell=cell, &
dft_control=dft_control, &
particle_set=particle_set, &
pw_env=pw_env)
CALL pw_env_get(pw_env, auxbas_rs_grid=rs_rho, &
auxbas_pw_pool=auxbas_pw_pool)
! be careful in parallel nsmax is chosen with multigrid in mind!
CALL rs_grid_zero(rs_rho)
eps_rho_rspace = dft_control%qs_control%eps_rho_rspace
DO ikind = 1, SIZE(atomic_kind_set)
CALL get_atomic_kind(atomic_kind_set(ikind), natom=natom, atom_list=atom_list)
CALL get_qs_kind(qs_kind_set(ikind), gth_potential=gth_potential)
IF (.NOT. ASSOCIATED(gth_potential)) CYCLE
CALL get_potential(potential=gth_potential, nlcc_present=nlcc, nexp_nlcc=nexp_nlcc, &
alpha_nlcc=alpha_nlcc, nct_nlcc=nct_nlcc, cval_nlcc=cval_nlcc)
IF (.NOT. nlcc) CYCLE
DO iexp_nlcc = 1, nexp_nlcc
alpha = alpha_nlcc(iexp_nlcc)
nc = nct_nlcc(iexp_nlcc)
ni = ncoset(2*nc - 2)
ALLOCATE (pab(ni, 1))
pab = 0._dp
nthread = 1
ithread = 0
CALL reallocate(cores, 1, natom)
npme = 0
cores = 0
! prepare core function
DO j = 1, nc
SELECT CASE (j)
CASE (1)
pab(1, 1) = cval_nlcc(1, iexp_nlcc)
CASE (2)
n = coset(2, 0, 0)
pab(n, 1) = cval_nlcc(2, iexp_nlcc)/alpha**2
n = coset(0, 2, 0)
pab(n, 1) = cval_nlcc(2, iexp_nlcc)/alpha**2
n = coset(0, 0, 2)
pab(n, 1) = cval_nlcc(2, iexp_nlcc)/alpha**2
CASE (3)
n = coset(4, 0, 0)
pab(n, 1) = cval_nlcc(3, iexp_nlcc)/alpha**4
n = coset(0, 4, 0)
pab(n, 1) = cval_nlcc(3, iexp_nlcc)/alpha**4
n = coset(0, 0, 4)
pab(n, 1) = cval_nlcc(3, iexp_nlcc)/alpha**4
n = coset(2, 2, 0)
pab(n, 1) = 2._dp*cval_nlcc(3, iexp_nlcc)/alpha**4
n = coset(2, 0, 2)
pab(n, 1) = 2._dp*cval_nlcc(3, iexp_nlcc)/alpha**4
n = coset(0, 2, 2)
pab(n, 1) = 2._dp*cval_nlcc(3, iexp_nlcc)/alpha**4
CASE (4)
n = coset(6, 0, 0)
pab(n, 1) = cval_nlcc(4, iexp_nlcc)/alpha**6
n = coset(0, 6, 0)
pab(n, 1) = cval_nlcc(4, iexp_nlcc)/alpha**6
n = coset(0, 0, 6)
pab(n, 1) = cval_nlcc(4, iexp_nlcc)/alpha**6
n = coset(4, 2, 0)
pab(n, 1) = 3._dp*cval_nlcc(4, iexp_nlcc)/alpha**6
n = coset(4, 0, 2)
pab(n, 1) = 3._dp*cval_nlcc(4, iexp_nlcc)/alpha**6
n = coset(2, 4, 0)
pab(n, 1) = 3._dp*cval_nlcc(4, iexp_nlcc)/alpha**6
n = coset(2, 0, 4)
pab(n, 1) = 3._dp*cval_nlcc(4, iexp_nlcc)/alpha**6
n = coset(0, 4, 2)
pab(n, 1) = 3._dp*cval_nlcc(4, iexp_nlcc)/alpha**6
n = coset(0, 2, 4)
pab(n, 1) = 3._dp*cval_nlcc(4, iexp_nlcc)/alpha**6
n = coset(2, 2, 2)
pab(n, 1) = 6._dp*cval_nlcc(4, iexp_nlcc)/alpha**6
CASE DEFAULT
CPABORT("")
END SELECT
END DO
IF (dft_control%nspins == 2) pab = pab*0.5_dp
DO iatom = 1, natom
atom_a = atom_list(iatom)
ra(:) = pbc(particle_set(atom_a)%r, cell)
IF (rs_rho%desc%parallel .AND. .NOT. rs_rho%desc%distributed) THEN
! replicated realspace grid, split the atoms up between procs
IF (MODULO(iatom, rs_rho%desc%group_size) == rs_rho%desc%my_pos) THEN
npme = npme + 1
cores(npme) = iatom
END IF
ELSE
npme = npme + 1
cores(npme) = iatom
END IF
END DO
DO j = 1, npme
iatom = cores(j)
atom_a = atom_list(iatom)
ra(:) = pbc(particle_set(atom_a)%r, cell)
subpatch_pattern = 0
ni = 2*nc - 2
radius = exp_radius_very_extended(la_min=0, la_max=ni, lb_min=0, lb_max=0, &
ra=ra, rb=ra, rp=ra, &
zetp=1/(2*alpha**2), eps=eps_rho_rspace, &
pab=pab, o1=0, o2=0, & ! without map_consistent
prefactor=1.0_dp, cutoff=0.0_dp)
CALL collocate_pgf_product(ni, 1/(2*alpha**2), 0, 0, 0.0_dp, 0, ra, &
(/0.0_dp, 0.0_dp, 0.0_dp/), 1.0_dp, pab, 0, 0, rs_rho, &
ga_gb_function=GRID_FUNC_AB, radius=radius, &
use_subpatch=.TRUE., subpatch_pattern=subpatch_pattern)
END DO
DEALLOCATE (pab)
END DO
END DO
IF (ASSOCIATED(cores)) THEN
DEALLOCATE (cores)
END IF
CALL transfer_rs2pw(rs_rho, rho_nlcc)
CALL timestop(handle)
END SUBROUTINE calculate_rho_nlcc
! **************************************************************************************************
!> \brief computes the local pseudopotential (without erf term) on the grid
!> \param vppl ...
!> \param qs_env ...
! **************************************************************************************************
SUBROUTINE calculate_ppl_grid(vppl, qs_env)
TYPE(pw_r3d_rs_type), INTENT(INOUT) :: vppl
TYPE(qs_environment_type), POINTER :: qs_env
CHARACTER(len=*), PARAMETER :: routineN = 'calculate_ppl_grid'
INTEGER :: atom_a, handle, iatom, ikind, ithread, &
j, lppl, n, natom, ni, npme, nthread, &
subpatch_pattern
INTEGER, DIMENSION(:), POINTER :: atom_list, cores
REAL(KIND=dp) :: alpha, eps_rho_rspace, radius
REAL(KIND=dp), DIMENSION(3) :: ra
REAL(KIND=dp), DIMENSION(:), POINTER :: cexp_ppl
REAL(KIND=dp), DIMENSION(:, :), POINTER :: pab
TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
TYPE(cell_type), POINTER :: cell
TYPE(dft_control_type), POINTER :: dft_control
TYPE(gth_potential_type), POINTER :: gth_potential
TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
TYPE(pw_env_type), POINTER :: pw_env
TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
TYPE(realspace_grid_type), POINTER :: rs_rho
CALL timeset(routineN, handle)
NULLIFY (cell, dft_control, pab, atomic_kind_set, qs_kind_set, particle_set, &
atom_list, pw_env, rs_rho, auxbas_pw_pool, cores)
CALL get_qs_env(qs_env=qs_env, &
atomic_kind_set=atomic_kind_set, &
qs_kind_set=qs_kind_set, &
cell=cell, &
dft_control=dft_control, &
particle_set=particle_set, &
pw_env=pw_env)
CALL pw_env_get(pw_env, auxbas_rs_grid=rs_rho, &
auxbas_pw_pool=auxbas_pw_pool)
! be careful in parallel nsmax is chosen with multigrid in mind!
CALL rs_grid_zero(rs_rho)
eps_rho_rspace = dft_control%qs_control%eps_rho_rspace
DO ikind = 1, SIZE(atomic_kind_set)
CALL get_atomic_kind(atomic_kind_set(ikind), natom=natom, atom_list=atom_list)
CALL get_qs_kind(qs_kind_set(ikind), gth_potential=gth_potential)
IF (.NOT. ASSOCIATED(gth_potential)) CYCLE
CALL get_potential(potential=gth_potential, alpha_ppl=alpha, nexp_ppl=lppl, cexp_ppl=cexp_ppl)
IF (lppl <= 0) CYCLE
ni = ncoset(2*lppl - 2)
ALLOCATE (pab(ni, 1))
pab = 0._dp
nthread = 1
ithread = 0
CALL reallocate(cores, 1, natom)
npme = 0
cores = 0
! prepare core function
DO j = 1, lppl
SELECT CASE (j)
CASE (1)
pab(1, 1) = cexp_ppl(1)
CASE (2)
n = coset(2, 0, 0)
pab(n, 1) = cexp_ppl(2)
n = coset(0, 2, 0)
pab(n, 1) = cexp_ppl(2)
n = coset(0, 0, 2)
pab(n, 1) = cexp_ppl(2)
CASE (3)
n = coset(4, 0, 0)
pab(n, 1) = cexp_ppl(3)
n = coset(0, 4, 0)
pab(n, 1) = cexp_ppl(3)
n = coset(0, 0, 4)
pab(n, 1) = cexp_ppl(3)
n = coset(2, 2, 0)
pab(n, 1) = 2._dp*cexp_ppl(3)
n = coset(2, 0, 2)
pab(n, 1) = 2._dp*cexp_ppl(3)
n = coset(0, 2, 2)
pab(n, 1) = 2._dp*cexp_ppl(3)
CASE (4)
n = coset(6, 0, 0)
pab(n, 1) = cexp_ppl(4)
n = coset(0, 6, 0)
pab(n, 1) = cexp_ppl(4)
n = coset(0, 0, 6)
pab(n, 1) = cexp_ppl(4)
n = coset(4, 2, 0)
pab(n, 1) = 3._dp*cexp_ppl(4)
n = coset(4, 0, 2)
pab(n, 1) = 3._dp*cexp_ppl(4)
n = coset(2, 4, 0)
pab(n, 1) = 3._dp*cexp_ppl(4)
n = coset(2, 0, 4)
pab(n, 1) = 3._dp*cexp_ppl(4)
n = coset(0, 4, 2)
pab(n, 1) = 3._dp*cexp_ppl(4)
n = coset(0, 2, 4)
pab(n, 1) = 3._dp*cexp_ppl(4)
n = coset(2, 2, 2)
pab(n, 1) = 6._dp*cexp_ppl(4)
CASE DEFAULT
CPABORT("")
END SELECT
END DO
DO iatom = 1, natom
atom_a = atom_list(iatom)
ra(:) = pbc(particle_set(atom_a)%r, cell)
IF (rs_rho%desc%parallel .AND. .NOT. rs_rho%desc%distributed) THEN
! replicated realspace grid, split the atoms up between procs
IF (MODULO(iatom, rs_rho%desc%group_size) == rs_rho%desc%my_pos) THEN
npme = npme + 1
cores(npme) = iatom
END IF
ELSE
npme = npme + 1
cores(npme) = iatom
END IF
END DO
IF (npme .GT. 0) THEN
DO j = 1, npme
iatom = cores(j)
atom_a = atom_list(iatom)
ra(:) = pbc(particle_set(atom_a)%r, cell)
subpatch_pattern = 0
ni = 2*lppl - 2
radius = exp_radius_very_extended(la_min=0, la_max=ni, &
lb_min=0, lb_max=0, &
ra=ra, rb=ra, rp=ra, &
zetp=alpha, eps=eps_rho_rspace, &
pab=pab, o1=0, o2=0, & ! without map_consistent
prefactor=1.0_dp, cutoff=0.0_dp)
CALL collocate_pgf_product(ni, alpha, 0, 0, 0.0_dp, 0, ra, &
(/0.0_dp, 0.0_dp, 0.0_dp/), 1.0_dp, pab, 0, 0, rs_rho, &
radius=radius, ga_gb_function=GRID_FUNC_AB, &
use_subpatch=.TRUE., subpatch_pattern=subpatch_pattern)
END DO
END IF
DEALLOCATE (pab)
END DO
IF (ASSOCIATED(cores)) THEN
DEALLOCATE (cores)
END IF
CALL transfer_rs2pw(rs_rho, vppl)
CALL timestop(handle)
END SUBROUTINE calculate_ppl_grid
! **************************************************************************************************
!> \brief Collocates the fitted lri density on a grid.
!> \param lri_rho_g ...
!> \param lri_rho_r ...
!> \param qs_env ...
!> \param lri_coef ...
!> \param total_rho ...
!> \param basis_type ...
!> \param exact_1c_terms ...
!> \param pmat replicated block diagonal density matrix (optional)
!> \param atomlist list of atoms to be included (optional)
!> \par History
!> 04.2013
!> \author Dorothea Golze
! **************************************************************************************************
SUBROUTINE calculate_lri_rho_elec(lri_rho_g, lri_rho_r, qs_env, &
lri_coef, total_rho, basis_type, exact_1c_terms, pmat, atomlist)
TYPE(pw_c1d_gs_type), INTENT(INOUT) :: lri_rho_g
TYPE(pw_r3d_rs_type), INTENT(INOUT) :: lri_rho_r
TYPE(qs_environment_type), POINTER :: qs_env
TYPE(lri_kind_type), DIMENSION(:), POINTER :: lri_coef
REAL(KIND=dp), INTENT(OUT) :: total_rho
CHARACTER(len=*), INTENT(IN) :: basis_type
LOGICAL, INTENT(IN) :: exact_1c_terms
TYPE(dbcsr_type), OPTIONAL :: pmat
INTEGER, DIMENSION(:), OPTIONAL :: atomlist
CHARACTER(len=*), PARAMETER :: routineN = 'calculate_lri_rho_elec'
INTEGER :: atom_a, group_size, handle, iatom, igrid_level, ikind, ipgf, iset, jpgf, jset, &
m1, maxco, maxsgf_set, my_pos, na1, natom, nb1, ncoa, ncob, nseta, offset, sgfa, sgfb
INTEGER, DIMENSION(:), POINTER :: atom_list, la_max, la_min, npgfa, nsgfa
INTEGER, DIMENSION(:, :), POINTER :: first_sgfa
LOGICAL :: found
LOGICAL, ALLOCATABLE, DIMENSION(:) :: map_it
LOGICAL, ALLOCATABLE, DIMENSION(:, :) :: map_it2
REAL(KIND=dp) :: eps_rho_rspace, radius, zetp
REAL(KIND=dp), DIMENSION(3) :: ra
REAL(KIND=dp), DIMENSION(:), POINTER :: aci
REAL(KIND=dp), DIMENSION(:, :), POINTER :: p_block, pab, sphi_a, work, zeta
TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
TYPE(cell_type), POINTER :: cell
TYPE(dft_control_type), POINTER :: dft_control
TYPE(gridlevel_info_type), POINTER :: gridlevel_info
TYPE(gto_basis_set_type), POINTER :: lri_basis_set, orb_basis_set
TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
TYPE(pw_env_type), POINTER :: pw_env
TYPE(pw_pool_p_type), DIMENSION(:), POINTER :: pw_pools
TYPE(pw_c1d_gs_type), ALLOCATABLE, DIMENSION(:) :: mgrid_gspace
TYPE(pw_r3d_rs_type), ALLOCATABLE, DIMENSION(:) :: mgrid_rspace
TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
TYPE(realspace_grid_type), DIMENSION(:), POINTER :: rs_rho
TYPE(realspace_grid_type), POINTER :: rs_grid
NULLIFY (aci, atomic_kind_set, qs_kind_set, atom_list, cell, &
dft_control, first_sgfa, gridlevel_info, la_max, &
la_min, lri_basis_set, npgfa, nsgfa, &
pab, particle_set, pw_env, pw_pools, rs_grid, rs_rho, sphi_a, &
work, zeta)
CALL timeset(routineN, handle)
IF (exact_1c_terms) THEN
CPASSERT(PRESENT(pmat))
END IF
CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set, &
atomic_kind_set=atomic_kind_set, &
cell=cell, particle_set=particle_set, &
pw_env=pw_env, &
dft_control=dft_control)
eps_rho_rspace = dft_control%qs_control%eps_rho_rspace
gridlevel_info => pw_env%gridlevel_info
! *** set up the pw multi-grids *** !
CPASSERT(ASSOCIATED(pw_env))
CALL pw_env_get(pw_env=pw_env, rs_grids=rs_rho, pw_pools=pw_pools)
CALL pw_pools_create_pws(pw_pools, mgrid_rspace)
CALL pw_pools_create_pws(pw_pools, mgrid_gspace)
! *** set up the rs multi-grids *** !
DO igrid_level = 1, gridlevel_info%ngrid_levels
CALL rs_grid_zero(rs_rho(igrid_level))
END DO
!take maxco from the LRI basis set!
CALL get_qs_kind_set(qs_kind_set=qs_kind_set, &
maxco=maxco, basis_type=basis_type)
ALLOCATE (pab(maxco, 1))
offset = 0
my_pos = mgrid_rspace(1)%pw_grid%para%group%mepos
group_size = mgrid_rspace(1)%pw_grid%para%group%num_pe
DO ikind = 1, SIZE(atomic_kind_set)
CALL get_atomic_kind(atomic_kind_set(ikind), natom=natom, atom_list=atom_list)
CALL get_qs_kind(qs_kind_set(ikind), basis_set=lri_basis_set, basis_type=basis_type)
!Take the lri basis set here!
CALL get_gto_basis_set(gto_basis_set=lri_basis_set, lmax=la_max, &
lmin=la_min, zet=zeta, nset=nseta, npgf=npgfa, &
sphi=sphi_a, first_sgf=first_sgfa, nsgf_set=nsgfa)
DO iatom = 1, natom
atom_a = atom_list(iatom)
IF (PRESENT(ATOMLIST)) THEN
IF (atomlist(atom_a) == 0) CYCLE
END IF
ra(:) = pbc(particle_set(atom_a)%r, cell)
aci => lri_coef(ikind)%acoef(iatom, :)
m1 = MAXVAL(npgfa(1:nseta))
ALLOCATE (map_it(m1))
DO iset = 1, nseta
! collocate this set locally?
map_it = .FALSE.
DO ipgf = 1, npgfa(iset)
igrid_level = gaussian_gridlevel(gridlevel_info, zeta(ipgf, iset))
rs_grid => rs_rho(igrid_level)
map_it(ipgf) = map_gaussian_here(rs_grid, cell%h_inv, ra, offset, group_size, my_pos)
END DO
offset = offset + 1
IF (ANY(map_it(1:npgfa(iset)))) THEN
sgfa = first_sgfa(1, iset)
ncoa = npgfa(iset)*ncoset(la_max(iset))
m1 = sgfa + nsgfa(iset) - 1
ALLOCATE (work(nsgfa(iset), 1))
work(1:nsgfa(iset), 1) = aci(sgfa:m1)
pab = 0._dp
CALL dgemm("N", "N", ncoa, 1, nsgfa(iset), 1.0_dp, lri_basis_set%sphi(1, sgfa), &
SIZE(lri_basis_set%sphi, 1), work(1, 1), SIZE(work, 1), 0.0_dp, pab(1, 1), &
SIZE(pab, 1))
DO ipgf = 1, npgfa(iset)
na1 = (ipgf - 1)*ncoset(la_max(iset))
igrid_level = gaussian_gridlevel(gridlevel_info, zeta(ipgf, iset))
rs_grid => rs_rho(igrid_level)
IF (map_it(ipgf)) THEN
radius = exp_radius_very_extended(la_min=la_min(iset), la_max=la_max(iset), &
lb_min=0, lb_max=0, &
ra=ra, rb=ra, rp=ra, &
zetp=zeta(ipgf, iset), eps=eps_rho_rspace, &
prefactor=1.0_dp, cutoff=1.0_dp)
CALL collocate_pgf_product(la_max=la_max(iset), &
zeta=zeta(ipgf, iset), &
la_min=la_min(iset), &
lb_max=0, zetb=0.0_dp, lb_min=0, &
ra=ra, rab=(/0.0_dp, 0.0_dp, 0.0_dp/), &
scale=1._dp, &
pab=pab, o1=na1, o2=0, &
rsgrid=rs_grid, &
radius=radius, &
ga_gb_function=GRID_FUNC_AB)
END IF
END DO
DEALLOCATE (work)
END IF
END DO
DEALLOCATE (map_it)
END DO
END DO
DEALLOCATE (pab)
! process the one-center terms
IF (exact_1c_terms) THEN
! find maximum numbers
offset = 0
CALL get_qs_kind_set(qs_kind_set=qs_kind_set, &
maxco=maxco, &
maxsgf_set=maxsgf_set, &
basis_type="ORB")
ALLOCATE (pab(maxco, maxco), work(maxco, maxsgf_set))
DO ikind = 1, SIZE(atomic_kind_set)
CALL get_atomic_kind(atomic_kind_set(ikind), natom=natom, atom_list=atom_list)
CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set, basis_type="ORB")
CALL get_gto_basis_set(gto_basis_set=orb_basis_set, lmax=la_max, &
lmin=la_min, zet=zeta, nset=nseta, npgf=npgfa, &
sphi=sphi_a, first_sgf=first_sgfa, nsgf_set=nsgfa)
DO iatom = 1, natom
atom_a = atom_list(iatom)
ra(:) = pbc(particle_set(atom_a)%r, cell)
CALL dbcsr_get_block_p(matrix=pmat, row=atom_a, col=atom_a, BLOCK=p_block, found=found)
m1 = MAXVAL(npgfa(1:nseta))
ALLOCATE (map_it2(m1, m1))
DO iset = 1, nseta
DO jset = 1, nseta
! processor mappint
map_it2 = .FALSE.
DO ipgf = 1, npgfa(iset)
DO jpgf = 1, npgfa(jset)
zetp = zeta(ipgf, iset) + zeta(jpgf, jset)
igrid_level = gaussian_gridlevel(gridlevel_info, zetp)
rs_grid => rs_rho(igrid_level)
map_it2(ipgf, jpgf) = map_gaussian_here(rs_grid, cell%h_inv, ra, offset, group_size, my_pos)
END DO
END DO
offset = offset + 1
!
IF (ANY(map_it2(1:npgfa(iset), 1:npgfa(jset)))) THEN
ncoa = npgfa(iset)*ncoset(la_max(iset))
sgfa = first_sgfa(1, iset)
ncob = npgfa(jset)*ncoset(la_max(jset))
sgfb = first_sgfa(1, jset)
! decontract density block
CALL dgemm("N", "N", ncoa, nsgfa(jset), nsgfa(iset), &
1.0_dp, sphi_a(1, sgfa), SIZE(sphi_a, 1), &
p_block(sgfa, sgfb), SIZE(p_block, 1), &
0.0_dp, work(1, 1), maxco)
CALL dgemm("N", "T", ncoa, ncob, nsgfa(jset), &
1.0_dp, work(1, 1), maxco, &
sphi_a(1, sgfb), SIZE(sphi_a, 1), &
0.0_dp, pab(1, 1), maxco)
DO ipgf = 1, npgfa(iset)
DO jpgf = 1, npgfa(jset)
zetp = zeta(ipgf, iset) + zeta(jpgf, jset)
igrid_level = gaussian_gridlevel(gridlevel_info, zetp)
rs_grid => rs_rho(igrid_level)
na1 = (ipgf - 1)*ncoset(la_max(iset))
nb1 = (jpgf - 1)*ncoset(la_max(jset))
IF (map_it2(ipgf, jpgf)) THEN
radius = exp_radius_very_extended(la_min=la_min(iset), &
la_max=la_max(iset), &
lb_min=la_min(jset), &
lb_max=la_max(jset), &
ra=ra, rb=ra, rp=ra, &
zetp=zetp, eps=eps_rho_rspace, &
prefactor=1.0_dp, cutoff=1.0_dp)
CALL collocate_pgf_product( &
la_max(iset), zeta(ipgf, iset), la_min(iset), &
la_max(jset), zeta(jpgf, jset), la_min(jset), &
ra, (/0.0_dp, 0.0_dp, 0.0_dp/), 1.0_dp, pab, na1, nb1, &
rs_grid, &
radius=radius, ga_gb_function=GRID_FUNC_AB)
END IF
END DO
END DO
END IF
END DO
END DO
DEALLOCATE (map_it2)
!
END DO
END DO
DEALLOCATE (pab, work)
END IF
CALL pw_zero(lri_rho_g)
CALL pw_zero(lri_rho_r)
DO igrid_level = 1, gridlevel_info%ngrid_levels
CALL pw_zero(mgrid_rspace(igrid_level))
CALL transfer_rs2pw(rs=rs_rho(igrid_level), &
pw=mgrid_rspace(igrid_level))
END DO
DO igrid_level = 1, gridlevel_info%ngrid_levels
CALL pw_zero(mgrid_gspace(igrid_level))
CALL pw_transfer(mgrid_rspace(igrid_level), &
mgrid_gspace(igrid_level))
CALL pw_axpy(mgrid_gspace(igrid_level), lri_rho_g)
END DO
CALL pw_transfer(lri_rho_g, lri_rho_r)
total_rho = pw_integrate_function(lri_rho_r, isign=-1)
! *** give back the multi-grids *** !
CALL pw_pools_give_back_pws(pw_pools, mgrid_gspace)
CALL pw_pools_give_back_pws(pw_pools, mgrid_rspace)
CALL timestop(handle)
END SUBROUTINE calculate_lri_rho_elec
#:for kind in ["r3d_rs", "c1d_gs"]
! **************************************************************************************************
!> \brief computes the density of the core charges on the grid
!> \param rho_core ...
!> \param total_rho ...
!> \param qs_env ...
!> \param calpha ...
!> \param ccore ...
!> \param only_nopaw ...
! **************************************************************************************************
SUBROUTINE calculate_rho_core_${kind}$ (rho_core, total_rho, qs_env, calpha, ccore, only_nopaw)
TYPE(pw_${kind}$_type), INTENT(INOUT) :: rho_core
REAL(KIND=dp), INTENT(OUT) :: total_rho
TYPE(qs_environment_type), POINTER :: qs_env
REAL(KIND=dp), DIMENSION(:), OPTIONAL :: calpha, ccore
LOGICAL, INTENT(IN), OPTIONAL :: only_nopaw
CHARACTER(len=*), PARAMETER :: routineN = 'calculate_rho_core'
INTEGER :: atom_a, handle, iatom, ikind, ithread, &
j, natom, npme, nthread, &
subpatch_pattern
INTEGER, DIMENSION(:), POINTER :: atom_list, cores
LOGICAL :: my_only_nopaw, paw_atom
REAL(KIND=dp) :: alpha, eps_rho_rspace, radius
REAL(KIND=dp), DIMENSION(3) :: ra
REAL(KIND=dp), DIMENSION(:, :), POINTER :: pab
TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
TYPE(cell_type), POINTER :: cell
TYPE(dft_control_type), POINTER :: dft_control
TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
TYPE(pw_env_type), POINTER :: pw_env
TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
TYPE(pw_r3d_rs_type) :: rhoc_r
TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
TYPE(realspace_grid_type), POINTER :: rs_rho
CALL timeset(routineN, handle)
NULLIFY (cell, dft_control, pab, atomic_kind_set, qs_kind_set, particle_set, &
atom_list, pw_env, rs_rho, auxbas_pw_pool, cores)
ALLOCATE (pab(1, 1))
my_only_nopaw = .FALSE.
IF (PRESENT(only_nopaw)) my_only_nopaw = only_nopaw
IF (PRESENT(calpha)) THEN
CPASSERT(PRESENT(ccore))
END IF
CALL get_qs_env(qs_env=qs_env, &
atomic_kind_set=atomic_kind_set, &
qs_kind_set=qs_kind_set, &
cell=cell, &
dft_control=dft_control, &
particle_set=particle_set, &
pw_env=pw_env)
CALL pw_env_get(pw_env, auxbas_rs_grid=rs_rho, &
auxbas_pw_pool=auxbas_pw_pool)
! be careful in parallel nsmax is chosen with multigrid in mind!
CALL rs_grid_zero(rs_rho)
eps_rho_rspace = dft_control%qs_control%eps_rho_rspace
DO ikind = 1, SIZE(atomic_kind_set)
CALL get_atomic_kind(atomic_kind_set(ikind), natom=natom, atom_list=atom_list)
IF (PRESENT(calpha)) THEN
alpha = calpha(ikind)
pab(1, 1) = ccore(ikind)
ELSE
CALL get_qs_kind(qs_kind_set(ikind), paw_atom=paw_atom, &
alpha_core_charge=alpha, ccore_charge=pab(1, 1))
END IF
IF (my_only_nopaw .AND. paw_atom) CYCLE
IF (alpha == 0.0_dp .OR. pab(1, 1) == 0.0_dp) CYCLE
nthread = 1
ithread = 0
CALL reallocate(cores, 1, natom)
npme = 0
cores = 0
DO iatom = 1, natom
IF (rs_rho%desc%parallel .AND. .NOT. rs_rho%desc%distributed) THEN
! replicated realspace grid, split the atoms up between procs
IF (MODULO(iatom, rs_rho%desc%group_size) == rs_rho%desc%my_pos) THEN
npme = npme + 1
cores(npme) = iatom
END IF
ELSE
npme = npme + 1
cores(npme) = iatom
END IF
END DO
IF (npme .GT. 0) THEN
DO j = 1, npme
iatom = cores(j)
atom_a = atom_list(iatom)
ra(:) = pbc(particle_set(atom_a)%r, cell)
subpatch_pattern = 0
radius = exp_radius_very_extended(la_min=0, la_max=0, &
lb_min=0, lb_max=0, &
ra=ra, rb=ra, rp=ra, &
zetp=alpha, eps=eps_rho_rspace, &
pab=pab, o1=0, o2=0, & ! without map_consistent
prefactor=-1.0_dp, cutoff=0.0_dp)
CALL collocate_pgf_product(0, alpha, 0, 0, 0.0_dp, 0, ra, &
(/0.0_dp, 0.0_dp, 0.0_dp/), -1.0_dp, pab, 0, 0, rs_rho, &
radius=radius, ga_gb_function=GRID_FUNC_AB, &
use_subpatch=.TRUE., subpatch_pattern=subpatch_pattern)
END DO
END IF
END DO
IF (ASSOCIATED(cores)) THEN
DEALLOCATE (cores)
END IF
DEALLOCATE (pab)
CALL auxbas_pw_pool%create_pw(rhoc_r)
CALL transfer_rs2pw(rs_rho, rhoc_r)
total_rho = pw_integrate_function(rhoc_r, isign=-1)
CALL pw_transfer(rhoc_r, rho_core)
CALL auxbas_pw_pool%give_back_pw(rhoc_r)
CALL timestop(handle)
END SUBROUTINE calculate_rho_core_${kind}$
#:endfor
! *****************************************************************************
!> \brief Computes the derivative of the density of the core charges with
!> respect to the nuclear coordinates on the grid.
!> \param drho_core The resulting density derivative
!> \param qs_env ...
!> \param beta Derivative direction
!> \param lambda Atom index
!> \note SL November 2014, ED 2021
! **************************************************************************************************
SUBROUTINE calculate_drho_core(drho_core, qs_env, beta, lambda)
TYPE(pw_c1d_gs_type), INTENT(INOUT) :: drho_core
TYPE(qs_environment_type), POINTER :: qs_env
INTEGER, INTENT(IN) :: beta, lambda
CHARACTER(len=*), PARAMETER :: routineN = 'calculate_drho_core'
INTEGER :: atom_a, dabqadb_func, handle, iatom, &
ikind, ithread, j, natom, npme, &
nthread, subpatch_pattern
INTEGER, DIMENSION(:), POINTER :: atom_list, cores
REAL(KIND=dp) :: alpha, eps_rho_rspace, radius
REAL(KIND=dp), DIMENSION(3) :: ra
REAL(KIND=dp), DIMENSION(:, :), POINTER :: pab
TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
TYPE(cell_type), POINTER :: cell
TYPE(dft_control_type), POINTER :: dft_control
TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
TYPE(pw_env_type), POINTER :: pw_env
TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
TYPE(pw_r3d_rs_type) :: rhoc_r
TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
TYPE(realspace_grid_type), POINTER :: rs_rho
CALL timeset(routineN, handle)
NULLIFY (cell, dft_control, pab, atomic_kind_set, qs_kind_set, particle_set, &
atom_list, pw_env, rs_rho, auxbas_pw_pool, cores)
ALLOCATE (pab(1, 1))
CALL get_qs_env(qs_env=qs_env, &
atomic_kind_set=atomic_kind_set, &
qs_kind_set=qs_kind_set, &
cell=cell, &
dft_control=dft_control, &
particle_set=particle_set, &
pw_env=pw_env)
CALL pw_env_get(pw_env, auxbas_rs_grid=rs_rho, &
auxbas_pw_pool=auxbas_pw_pool)
! be careful in parallel nsmax is chosen with multigrid in mind!
CALL rs_grid_zero(rs_rho)
eps_rho_rspace = dft_control%qs_control%eps_rho_rspace
SELECT CASE (beta)
CASE (1)
dabqadb_func = GRID_FUNC_CORE_X
CASE (2)
dabqadb_func = GRID_FUNC_CORE_Y
CASE (3)
dabqadb_func = GRID_FUNC_CORE_Z
CASE DEFAULT
CPABORT("invalid beta")
END SELECT
DO ikind = 1, SIZE(atomic_kind_set)
CALL get_atomic_kind(atomic_kind_set(ikind), natom=natom, atom_list=atom_list)
CALL get_qs_kind(qs_kind_set(ikind), &
alpha_core_charge=alpha, ccore_charge=pab(1, 1))
IF (alpha == 0.0_dp .OR. pab(1, 1) == 0.0_dp) CYCLE
nthread = 1
ithread = 0
CALL reallocate(cores, 1, natom)
npme = 0