-
Notifications
You must be signed in to change notification settings - Fork 1
/
optimize_embedding_potential.F
3634 lines (3035 loc) · 157 KB
/
optimize_embedding_potential.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 !
!--------------------------------------------------------------------------------------------------!
MODULE optimize_embedding_potential
USE atomic_kind_types, ONLY: atomic_kind_type,&
get_atomic_kind,&
get_atomic_kind_set
USE cell_types, ONLY: cell_type
USE cp_blacs_env, ONLY: cp_blacs_env_create,&
cp_blacs_env_release,&
cp_blacs_env_type
USE cp_control_types, ONLY: dft_control_type
USE cp_dbcsr_api, ONLY: dbcsr_p_type
USE cp_dbcsr_operations, ONLY: copy_dbcsr_to_fm,&
dbcsr_deallocate_matrix_set
USE cp_files, ONLY: close_file,&
open_file
USE cp_fm_basic_linalg, ONLY: cp_fm_column_scale,&
cp_fm_scale,&
cp_fm_scale_and_add,&
cp_fm_trace
USE cp_fm_diag, ONLY: choose_eigv_solver
USE cp_fm_struct, ONLY: cp_fm_struct_create,&
cp_fm_struct_release,&
cp_fm_struct_type
USE cp_fm_types, ONLY: &
cp_fm_copy_general, cp_fm_create, cp_fm_get_element, cp_fm_get_info, cp_fm_release, &
cp_fm_set_all, cp_fm_to_fm, cp_fm_to_fm_submat, cp_fm_type, cp_fm_write_unformatted
USE cp_log_handling, ONLY: cp_get_default_logger,&
cp_logger_type
USE cp_output_handling, ONLY: cp_p_file,&
cp_print_key_finished_output,&
cp_print_key_should_output,&
cp_print_key_unit_nr
USE cp_realspace_grid_cube, ONLY: cp_cube_to_pw,&
cp_pw_to_cube,&
cp_pw_to_simple_volumetric
USE embed_types, ONLY: opt_embed_pot_type
USE force_env_types, ONLY: force_env_type
USE input_constants, ONLY: &
embed_diff, embed_fa, embed_grid_angstrom, embed_grid_bohr, embed_level_shift, embed_none, &
embed_quasi_newton, embed_resp, embed_steep_desc
USE input_section_types, ONLY: section_get_ival,&
section_get_ivals,&
section_get_rval,&
section_vals_get_subs_vals,&
section_vals_type,&
section_vals_val_get
USE kinds, ONLY: default_path_length,&
dp
USE lri_environment_types, ONLY: lri_kind_type
USE mathconstants, ONLY: pi
USE message_passing, ONLY: mp_para_env_type
USE mixed_environment_utils, ONLY: get_subsys_map_index
USE parallel_gemm_api, ONLY: parallel_gemm
USE particle_list_types, ONLY: particle_list_type
USE particle_types, ONLY: particle_type
USE pw_env_types, ONLY: pw_env_get,&
pw_env_type
USE pw_methods, ONLY: &
pw_axpy, pw_copy, pw_derive, pw_dr2, pw_integral_ab, pw_integrate_function, pw_scale, &
pw_transfer, pw_zero
USE pw_poisson_methods, ONLY: pw_poisson_solve
USE pw_poisson_types, ONLY: pw_poisson_type
USE pw_pool_types, ONLY: pw_pool_type
USE pw_types, ONLY: pw_c1d_gs_type,&
pw_r3d_rs_type
USE qs_collocate_density, ONLY: calculate_rho_resp_all,&
calculate_wavefunction,&
collocate_function
USE qs_environment_types, ONLY: get_qs_env,&
qs_environment_type,&
set_qs_env
USE qs_integrate_potential_single, ONLY: integrate_v_rspace_one_center
USE qs_kind_types, ONLY: get_qs_kind,&
qs_kind_type
USE qs_kinetic, ONLY: build_kinetic_matrix
USE qs_ks_types, ONLY: qs_ks_env_type
USE qs_mo_types, ONLY: get_mo_set,&
mo_set_type
USE qs_neighbor_list_types, ONLY: neighbor_list_set_p_type
USE qs_rho_types, ONLY: qs_rho_get,&
qs_rho_type
USE qs_subsys_types, ONLY: qs_subsys_get,&
qs_subsys_type
USE xc, ONLY: smooth_cutoff
USE xc_rho_cflags_types, ONLY: xc_rho_cflags_setall,&
xc_rho_cflags_type
USE xc_rho_set_types, ONLY: xc_rho_set_create,&
xc_rho_set_release,&
xc_rho_set_type,&
xc_rho_set_update
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'optimize_embedding_potential'
PUBLIC :: prepare_embed_opt, init_embed_pot, release_opt_embed, calculate_embed_pot_grad, &
opt_embed_step, print_rho_diff, step_control, max_dens_diff, print_emb_opt_info, &
conv_check_embed, make_subsys_embed_pot, print_embed_restart, find_aux_dimen, &
read_embed_pot, understand_spin_states, given_embed_pot, print_rho_spin_diff, &
print_pot_simple_grid, get_prev_density, get_max_subsys_diff, Coulomb_guess
CONTAINS
! **************************************************************************************************
!> \brief Find out whether we need to swap alpha- and beta- spind densities in the second subsystem
!> \brief It's only needed because by default alpha-spins go first in a subsystem.
!> \brief By swapping we impose the constraint:
!> \brief rho_1(alpha) + rho_2(alpha) = rho_total(alpha)
!> \brief rho_1(beta) + rho_2(beta) = rho_total(beta)
!> \param force_env ...
!> \param ref_subsys_number ...
!> \param change_spin ...
!> \param open_shell_embed ...
!> \param all_nspins ...
!> \return ...
!> \author Vladimir Rybkin
! **************************************************************************************************
SUBROUTINE understand_spin_states(force_env, ref_subsys_number, change_spin, open_shell_embed, all_nspins)
TYPE(force_env_type), POINTER :: force_env
INTEGER :: ref_subsys_number
LOGICAL :: change_spin, open_shell_embed
INTEGER, ALLOCATABLE, DIMENSION(:) :: all_nspins
INTEGER :: i_force_eval, nspins, sub_spin_1, &
sub_spin_2, total_spin
INTEGER, DIMENSION(2) :: nelectron_spin
INTEGER, DIMENSION(2, 3) :: all_spins
TYPE(dft_control_type), POINTER :: dft_control
change_spin = .FALSE.
open_shell_embed = .FALSE.
ALLOCATE (all_nspins(ref_subsys_number))
IF (ref_subsys_number .EQ. 3) THEN
all_spins = 0
DO i_force_eval = 1, ref_subsys_number
CALL get_qs_env(qs_env=force_env%sub_force_env(i_force_eval)%force_env%qs_env, &
nelectron_spin=nelectron_spin, dft_control=dft_control)
all_spins(:, i_force_eval) = nelectron_spin
nspins = dft_control%nspins
all_nspins(i_force_eval) = nspins
END DO
! Find out whether we need a spin-dependend embedding potential
IF (.NOT. ((all_nspins(1) .EQ. 1) .AND. (all_nspins(2) .EQ. 1) .AND. (all_nspins(3) .EQ. 1))) &
open_shell_embed = .TRUE.
! If it's open shell, we need to check spin states
IF (open_shell_embed) THEN
IF (all_nspins(3) .EQ. 1) THEN
total_spin = 0
ELSE
total_spin = all_spins(1, 3) - all_spins(2, 3)
END IF
IF (all_nspins(1) .EQ. 1) THEN
sub_spin_1 = 0
ELSE
sub_spin_1 = all_spins(1, 1) - all_spins(2, 1)
END IF
IF (all_nspins(2) .EQ. 1) THEN
sub_spin_2 = 0
ELSE
sub_spin_2 = all_spins(1, 2) - all_spins(2, 2)
END IF
IF ((sub_spin_1 + sub_spin_2) .EQ. total_spin) THEN
change_spin = .FALSE.
ELSE
IF (ABS(sub_spin_1 - sub_spin_2) .EQ. total_spin) THEN
change_spin = .TRUE.
ELSE
CPABORT("Spin states of subsystems are not compatible.")
END IF
END IF
END IF ! not open_shell
ELSE
CPABORT("Reference subsystem must be the third FORCE_EVAL.")
END IF
END SUBROUTINE understand_spin_states
! **************************************************************************************************
!> \brief ...
!> \param qs_env ...
!> \param embed_pot ...
!> \param add_const_pot ...
!> \param Fermi_Amaldi ...
!> \param const_pot ...
!> \param open_shell_embed ...
!> \param spin_embed_pot ...
!> \param pot_diff ...
!> \param Coulomb_guess ...
!> \param grid_opt ...
! **************************************************************************************************
SUBROUTINE init_embed_pot(qs_env, embed_pot, add_const_pot, Fermi_Amaldi, const_pot, open_shell_embed, &
spin_embed_pot, pot_diff, Coulomb_guess, grid_opt)
TYPE(qs_environment_type), POINTER :: qs_env
TYPE(pw_r3d_rs_type), POINTER :: embed_pot
LOGICAL :: add_const_pot, Fermi_Amaldi
TYPE(pw_r3d_rs_type), POINTER :: const_pot
LOGICAL :: open_shell_embed
TYPE(pw_r3d_rs_type), POINTER :: spin_embed_pot, pot_diff
LOGICAL :: Coulomb_guess, grid_opt
INTEGER :: nelectrons
INTEGER, DIMENSION(2) :: nelectron_spin
REAL(KIND=dp) :: factor
TYPE(pw_env_type), POINTER :: pw_env
TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
TYPE(pw_r3d_rs_type), POINTER :: v_hartree_r_space
! Extract plane waves environment
CALL get_qs_env(qs_env=qs_env, pw_env=pw_env, &
nelectron_spin=nelectron_spin, &
v_hartree_rspace=v_hartree_r_space)
! Prepare plane-waves pool
CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
! Create embedding potential and set to zero
NULLIFY (embed_pot)
ALLOCATE (embed_pot)
CALL auxbas_pw_pool%create_pw(embed_pot)
CALL pw_zero(embed_pot)
! Spin embedding potential if asked
IF (open_shell_embed) THEN
NULLIFY (spin_embed_pot)
ALLOCATE (spin_embed_pot)
CALL auxbas_pw_pool%create_pw(spin_embed_pot)
CALL pw_zero(spin_embed_pot)
END IF
! Coulomb guess/constant potential
IF (Coulomb_guess) THEN
NULLIFY (pot_diff)
ALLOCATE (pot_diff)
CALL auxbas_pw_pool%create_pw(pot_diff)
CALL pw_zero(pot_diff)
END IF
! Initialize constant part of the embedding potential
IF (add_const_pot .AND. (.NOT. grid_opt)) THEN
! Now the constant potential is the Coulomb one
NULLIFY (const_pot)
ALLOCATE (const_pot)
CALL auxbas_pw_pool%create_pw(const_pot)
CALL pw_zero(const_pot)
END IF
! Add Fermi-Amaldi potential if requested
IF (Fermi_Amaldi) THEN
! Extract Hartree potential
NULLIFY (v_hartree_r_space)
CALL get_qs_env(qs_env=qs_env, pw_env=pw_env, &
v_hartree_rspace=v_hartree_r_space)
CALL pw_copy(v_hartree_r_space, embed_pot)
! Calculate the number of electrons
nelectrons = nelectron_spin(1) + nelectron_spin(2)
factor = (REAL(nelectrons, dp) - 1.0_dp)/(REAL(nelectrons, dp))
! Scale the Hartree potential to get Fermi-Amaldi
CALL pw_scale(embed_pot, a=factor)
! Copy Fermi-Amaldi to embedding potential for basis-based optimization
IF (.NOT. grid_opt) CALL pw_copy(embed_pot, embed_pot)
END IF
END SUBROUTINE init_embed_pot
! **************************************************************************************************
!> \brief Creates and allocates objects for optimization of embedding potential
!> \param qs_env ...
!> \param opt_embed ...
!> \param opt_embed_section ...
!> \author Vladimir Rybkin
! **************************************************************************************************
SUBROUTINE prepare_embed_opt(qs_env, opt_embed, opt_embed_section)
TYPE(qs_environment_type), POINTER :: qs_env
TYPE(opt_embed_pot_type) :: opt_embed
TYPE(section_vals_type), POINTER :: opt_embed_section
INTEGER :: diff_size, i_dens, size_prev_dens
TYPE(cp_blacs_env_type), POINTER :: blacs_env
TYPE(cp_fm_struct_type), POINTER :: fm_struct
TYPE(mp_para_env_type), POINTER :: para_env
TYPE(pw_env_type), POINTER :: pw_env
TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
!TYPE(pw_env_type), POINTER :: pw_env
!TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
! First, read the input
CALL read_opt_embed_section(opt_embed, opt_embed_section)
! All these are needed for optimization in a finite Gaussian basis
IF (.NOT. opt_embed%grid_opt) THEN
! Create blacs environment
CALL get_qs_env(qs_env=qs_env, &
para_env=para_env)
NULLIFY (blacs_env)
CALL cp_blacs_env_create(blacs_env=blacs_env, para_env=para_env)
! Reveal the dimension of the RI basis
CALL find_aux_dimen(qs_env, opt_embed%dimen_aux)
! Prepare the object for integrals
CALL make_lri_object(qs_env, opt_embed%lri)
! In case if spin embedding potential has to be optimized,
! the dimension of variational space is two times larger
IF (opt_embed%open_shell_embed) THEN
opt_embed%dimen_var_aux = 2*opt_embed%dimen_aux
ELSE
opt_embed%dimen_var_aux = opt_embed%dimen_aux
END IF
! Allocate expansion coefficients and gradient
NULLIFY (opt_embed%embed_pot_grad, opt_embed%embed_pot_coef, opt_embed%step, fm_struct)
NULLIFY (opt_embed%prev_embed_pot_grad, opt_embed%prev_embed_pot_coef, opt_embed%prev_step)
CALL cp_fm_struct_create(fm_struct, para_env=para_env, context=blacs_env, &
nrow_global=opt_embed%dimen_var_aux, ncol_global=1)
ALLOCATE (opt_embed%embed_pot_grad, opt_embed%embed_pot_coef, &
opt_embed%prev_embed_pot_grad, opt_embed%prev_embed_pot_coef, &
opt_embed%step, opt_embed%prev_step)
CALL cp_fm_create(opt_embed%embed_pot_grad, fm_struct, name="pot_grad")
CALL cp_fm_create(opt_embed%embed_pot_coef, fm_struct, name="pot_coef")
CALL cp_fm_create(opt_embed%prev_embed_pot_grad, fm_struct, name="prev_pot_grad")
CALL cp_fm_create(opt_embed%prev_embed_pot_coef, fm_struct, name="prev_pot_coef")
CALL cp_fm_create(opt_embed%step, fm_struct, name="step")
CALL cp_fm_create(opt_embed%prev_step, fm_struct, name="prev_step")
CALL cp_fm_struct_release(fm_struct)
CALL cp_fm_set_all(opt_embed%embed_pot_grad, 0.0_dp)
CALL cp_fm_set_all(opt_embed%prev_embed_pot_grad, 0.0_dp)
CALL cp_fm_set_all(opt_embed%embed_pot_coef, 0.0_dp)
CALL cp_fm_set_all(opt_embed%prev_embed_pot_coef, 0.0_dp)
CALL cp_fm_set_all(opt_embed%step, 0.0_dp)
CALL cp_fm_set_all(opt_embed%prev_step, 0.0_dp)
! Allocate Hessian
NULLIFY (opt_embed%embed_pot_hess, opt_embed%prev_embed_pot_hess, fm_struct)
CALL cp_fm_struct_create(fm_struct, para_env=para_env, context=blacs_env, &
nrow_global=opt_embed%dimen_var_aux, ncol_global=opt_embed%dimen_var_aux)
ALLOCATE (opt_embed%embed_pot_hess, opt_embed%prev_embed_pot_hess)
CALL cp_fm_create(opt_embed%embed_pot_hess, fm_struct, name="pot_Hess")
CALL cp_fm_create(opt_embed%prev_embed_pot_hess, fm_struct, name="prev_pot_Hess")
CALL cp_fm_struct_release(fm_struct)
! Special structure for the kinetic energy matrix
NULLIFY (fm_struct, opt_embed%kinetic_mat)
CALL cp_fm_struct_create(fm_struct, para_env=para_env, context=blacs_env, &
nrow_global=opt_embed%dimen_aux, ncol_global=opt_embed%dimen_aux)
ALLOCATE (opt_embed%kinetic_mat)
CALL cp_fm_create(opt_embed%kinetic_mat, fm_struct, name="kinetic_mat")
CALL cp_fm_struct_release(fm_struct)
CALL cp_fm_set_all(opt_embed%kinetic_mat, 0.0_dp)
! Hessian is set as a unit matrix
CALL cp_fm_set_all(opt_embed%embed_pot_hess, 0.0_dp, -1.0_dp)
CALL cp_fm_set_all(opt_embed%prev_embed_pot_hess, 0.0_dp, -1.0_dp)
! Release blacs environment
CALL cp_blacs_env_release(blacs_env)
END IF
CALL get_qs_env(qs_env=qs_env, pw_env=pw_env)
CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
NULLIFY (opt_embed%prev_subsys_dens)
size_prev_dens = SUM(opt_embed%all_nspins(1:(SIZE(opt_embed%all_nspins) - 1)))
ALLOCATE (opt_embed%prev_subsys_dens(size_prev_dens))
DO i_dens = 1, size_prev_dens
CALL auxbas_pw_pool%create_pw(opt_embed%prev_subsys_dens(i_dens))
CALL pw_zero(opt_embed%prev_subsys_dens(i_dens))
END DO
ALLOCATE (opt_embed%max_subsys_dens_diff(size_prev_dens))
! Array to store functional values
ALLOCATE (opt_embed%w_func(opt_embed%n_iter))
opt_embed%w_func = 0.0_dp
! Allocate max_diff and int_diff
diff_size = 1
IF (opt_embed%open_shell_embed) diff_size = 2
ALLOCATE (opt_embed%max_diff(diff_size))
ALLOCATE (opt_embed%int_diff(diff_size))
ALLOCATE (opt_embed%int_diff_square(diff_size))
! FAB update
IF (opt_embed%fab) THEN
NULLIFY (opt_embed%prev_embed_pot)
ALLOCATE (opt_embed%prev_embed_pot)
CALL auxbas_pw_pool%create_pw(opt_embed%prev_embed_pot)
CALL pw_zero(opt_embed%prev_embed_pot)
IF (opt_embed%open_shell_embed) THEN
NULLIFY (opt_embed%prev_spin_embed_pot)
ALLOCATE (opt_embed%prev_spin_embed_pot)
CALL auxbas_pw_pool%create_pw(opt_embed%prev_spin_embed_pot)
CALL pw_zero(opt_embed%prev_spin_embed_pot)
END IF
END IF
! Set allowed energy decrease parameter
opt_embed%allowed_decrease = 0.0001_dp
! Regularization contribution is set to zero
opt_embed%reg_term = 0.0_dp
! Step is accepted in the beginning
opt_embed%accept_step = .TRUE.
opt_embed%newton_step = .FALSE.
opt_embed%last_accepted = 1
! Set maximum and minimum trust radii
opt_embed%max_trad = opt_embed%trust_rad*7.900_dp
opt_embed%min_trad = opt_embed%trust_rad*0.125*0.065_dp
END SUBROUTINE prepare_embed_opt
! **************************************************************************************************
!> \brief ...
!> \param opt_embed ...
!> \param opt_embed_section ...
! **************************************************************************************************
SUBROUTINE read_opt_embed_section(opt_embed, opt_embed_section)
TYPE(opt_embed_pot_type) :: opt_embed
TYPE(section_vals_type), POINTER :: opt_embed_section
INTEGER :: embed_guess, embed_optimizer
! Read keywords
CALL section_vals_val_get(opt_embed_section, "REG_LAMBDA", &
r_val=opt_embed%lambda)
CALL section_vals_val_get(opt_embed_section, "N_ITER", &
i_val=opt_embed%n_iter)
CALL section_vals_val_get(opt_embed_section, "TRUST_RAD", &
r_val=opt_embed%trust_rad)
CALL section_vals_val_get(opt_embed_section, "DENS_CONV_MAX", &
r_val=opt_embed%conv_max)
CALL section_vals_val_get(opt_embed_section, "DENS_CONV_INT", &
r_val=opt_embed%conv_int)
CALL section_vals_val_get(opt_embed_section, "SPIN_DENS_CONV_MAX", &
r_val=opt_embed%conv_max_spin)
CALL section_vals_val_get(opt_embed_section, "SPIN_DENS_CONV_INT", &
r_val=opt_embed%conv_int_spin)
CALL section_vals_val_get(opt_embed_section, "CHARGE_DISTR_WIDTH", &
r_val=opt_embed%eta)
CALL section_vals_val_get(opt_embed_section, "READ_EMBED_POT", &
l_val=opt_embed%read_embed_pot)
CALL section_vals_val_get(opt_embed_section, "READ_EMBED_POT_CUBE", &
l_val=opt_embed%read_embed_pot_cube)
CALL section_vals_val_get(opt_embed_section, "GRID_OPT", &
l_val=opt_embed%grid_opt)
CALL section_vals_val_get(opt_embed_section, "LEEUWEN-BAERENDS", &
l_val=opt_embed%leeuwen)
CALL section_vals_val_get(opt_embed_section, "FAB", &
l_val=opt_embed%fab)
CALL section_vals_val_get(opt_embed_section, "VW_CUTOFF", &
r_val=opt_embed%vw_cutoff)
CALL section_vals_val_get(opt_embed_section, "VW_SMOOTH_CUT_RANGE", &
r_val=opt_embed%vw_smooth_cutoff_range)
CALL section_vals_val_get(opt_embed_section, "OPTIMIZER", i_val=embed_optimizer)
SELECT CASE (embed_optimizer)
CASE (embed_steep_desc)
opt_embed%steep_desc = .TRUE.
CASE (embed_quasi_newton)
opt_embed%steep_desc = .FALSE.
opt_embed%level_shift = .FALSE.
CASE (embed_level_shift)
opt_embed%steep_desc = .FALSE.
opt_embed%level_shift = .TRUE.
CASE DEFAULT
opt_embed%steep_desc = .TRUE.
END SELECT
CALL section_vals_val_get(opt_embed_section, "POT_GUESS", i_val=embed_guess)
SELECT CASE (embed_guess)
CASE (embed_none)
opt_embed%add_const_pot = .FALSE.
opt_embed%Fermi_Amaldi = .FALSE.
opt_embed%Coulomb_guess = .FALSE.
opt_embed%diff_guess = .FALSE.
CASE (embed_diff)
opt_embed%add_const_pot = .TRUE.
opt_embed%Fermi_Amaldi = .FALSE.
opt_embed%Coulomb_guess = .FALSE.
opt_embed%diff_guess = .TRUE.
CASE (embed_fa)
opt_embed%add_const_pot = .TRUE.
opt_embed%Fermi_Amaldi = .TRUE.
opt_embed%Coulomb_guess = .FALSE.
opt_embed%diff_guess = .FALSE.
CASE (embed_resp)
opt_embed%add_const_pot = .TRUE.
opt_embed%Fermi_Amaldi = .TRUE.
opt_embed%Coulomb_guess = .TRUE.
opt_embed%diff_guess = .FALSE.
CASE DEFAULT
opt_embed%add_const_pot = .FALSE.
opt_embed%Fermi_Amaldi = .FALSE.
opt_embed%Coulomb_guess = .FALSE.
opt_embed%diff_guess = .FALSE.
END SELECT
END SUBROUTINE read_opt_embed_section
! **************************************************************************************************
!> \brief Find the dimension of the auxiliary basis for the expansion of the embedding potential
!> \param qs_env ...
!> \param dimen_aux ...
! **************************************************************************************************
SUBROUTINE find_aux_dimen(qs_env, dimen_aux)
TYPE(qs_environment_type), POINTER :: qs_env
INTEGER :: dimen_aux
INTEGER :: iatom, ikind, nsgf
INTEGER, ALLOCATABLE, DIMENSION(:) :: kind_of
TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
! First, reveal the dimension of the RI basis
CALL get_qs_env(qs_env=qs_env, &
particle_set=particle_set, &
qs_kind_set=qs_kind_set, &
atomic_kind_set=atomic_kind_set)
CALL get_atomic_kind_set(atomic_kind_set, kind_of=kind_of)
dimen_aux = 0
DO iatom = 1, SIZE(particle_set)
ikind = kind_of(iatom)
CALL get_qs_kind(qs_kind=qs_kind_set(ikind), nsgf=nsgf, basis_type="RI_AUX")
dimen_aux = dimen_aux + nsgf
END DO
END SUBROUTINE find_aux_dimen
! **************************************************************************************************
!> \brief Prepare the lri_kind_type object for integrals between density and aux. basis functions
!> \param qs_env ...
!> \param lri ...
! **************************************************************************************************
SUBROUTINE make_lri_object(qs_env, lri)
TYPE(qs_environment_type), POINTER :: qs_env
TYPE(lri_kind_type), DIMENSION(:), POINTER :: lri
INTEGER :: ikind, natom, nkind, nsgf
TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
TYPE(atomic_kind_type), POINTER :: atomic_kind
TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
NULLIFY (atomic_kind, lri)
CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set, &
qs_kind_set=qs_kind_set)
nkind = SIZE(atomic_kind_set)
ALLOCATE (lri(nkind))
! Here we need only v_int and acoef (the latter as dummies)
DO ikind = 1, nkind
NULLIFY (lri(ikind)%acoef)
NULLIFY (lri(ikind)%v_int)
atomic_kind => atomic_kind_set(ikind)
CALL get_atomic_kind(atomic_kind=atomic_kind, natom=natom)
CALL get_qs_kind(qs_kind=qs_kind_set(ikind), nsgf=nsgf, basis_type="RI_AUX")
ALLOCATE (lri(ikind)%acoef(natom, nsgf))
lri(ikind)%acoef = 0._dp
ALLOCATE (lri(ikind)%v_int(natom, nsgf))
lri(ikind)%v_int = 0._dp
END DO
END SUBROUTINE make_lri_object
! **************************************************************************************************
!> \brief Read the external embedding potential, not to be optimized
!> \param qs_env ...
! **************************************************************************************************
SUBROUTINE given_embed_pot(qs_env)
TYPE(qs_environment_type), POINTER :: qs_env
LOGICAL :: open_shell_embed
TYPE(dft_control_type), POINTER :: dft_control
TYPE(pw_env_type), POINTER :: pw_env
TYPE(pw_pool_type), POINTER :: auxbas_pw_pool_subsys
TYPE(pw_r3d_rs_type), POINTER :: embed_pot, spin_embed_pot
TYPE(section_vals_type), POINTER :: input, qs_section
qs_env%given_embed_pot = .TRUE.
NULLIFY (input, dft_control, embed_pot, spin_embed_pot, embed_pot, spin_embed_pot, &
qs_section)
CALL get_qs_env(qs_env=qs_env, &
input=input, &
dft_control=dft_control, &
pw_env=pw_env)
qs_section => section_vals_get_subs_vals(input, "DFT%QS")
open_shell_embed = .FALSE.
IF (dft_control%nspins .EQ. 2) open_shell_embed = .TRUE.
! Prepare plane-waves pool
CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool_subsys)
! Create embedding potential
!CALL get_qs_env(qs_env=qs_env, &
! embed_pot=embed_pot)
ALLOCATE (embed_pot)
CALL auxbas_pw_pool_subsys%create_pw(embed_pot)
IF (open_shell_embed) THEN
! Create spin embedding potential
ALLOCATE (spin_embed_pot)
CALL auxbas_pw_pool_subsys%create_pw(spin_embed_pot)
END IF
! Read the cubes
CALL read_embed_pot_cube(embed_pot, spin_embed_pot, qs_section, open_shell_embed)
IF (.NOT. open_shell_embed) THEN
CALL set_qs_env(qs_env=qs_env, embed_pot=embed_pot)
ELSE
CALL set_qs_env(qs_env=qs_env, embed_pot=embed_pot, spin_embed_pot=spin_embed_pot)
END IF
END SUBROUTINE given_embed_pot
! **************************************************************************************************
!> \brief ...
!> \param qs_env ...
!> \param embed_pot ...
!> \param spin_embed_pot ...
!> \param section ...
!> \param opt_embed ...
! **************************************************************************************************
SUBROUTINE read_embed_pot(qs_env, embed_pot, spin_embed_pot, section, opt_embed)
TYPE(qs_environment_type), POINTER :: qs_env
TYPE(pw_r3d_rs_type), POINTER :: embed_pot, spin_embed_pot
TYPE(section_vals_type), POINTER :: section
TYPE(opt_embed_pot_type) :: opt_embed
! Read the potential as a vector in the auxiliary basis
IF (opt_embed%read_embed_pot) &
CALL read_embed_pot_vector(qs_env, embed_pot, spin_embed_pot, section, &
opt_embed%embed_pot_coef, opt_embed%open_shell_embed)
! Read the potential as a cube (two cubes for open shell)
IF (opt_embed%read_embed_pot_cube) &
CALL read_embed_pot_cube(embed_pot, spin_embed_pot, section, opt_embed%open_shell_embed)
END SUBROUTINE read_embed_pot
! **************************************************************************************************
!> \brief ...
!> \param embed_pot ...
!> \param spin_embed_pot ...
!> \param section ...
!> \param open_shell_embed ...
! **************************************************************************************************
SUBROUTINE read_embed_pot_cube(embed_pot, spin_embed_pot, section, open_shell_embed)
TYPE(pw_r3d_rs_type), INTENT(IN) :: embed_pot, spin_embed_pot
TYPE(section_vals_type), POINTER :: section
LOGICAL :: open_shell_embed
CHARACTER(LEN=default_path_length) :: filename
LOGICAL :: exist
REAL(KIND=dp) :: scaling_factor
exist = .FALSE.
CALL section_vals_val_get(section, "EMBED_CUBE_FILE_NAME", c_val=filename)
INQUIRE (FILE=filename, exist=exist)
IF (.NOT. exist) &
CPABORT("Embedding cube file not found. ")
scaling_factor = 1.0_dp
CALL cp_cube_to_pw(embed_pot, filename, scaling_factor)
! Spin-dependent part of the potential
IF (open_shell_embed) THEN
exist = .FALSE.
CALL section_vals_val_get(section, "EMBED_SPIN_CUBE_FILE_NAME", c_val=filename)
INQUIRE (FILE=filename, exist=exist)
IF (.NOT. exist) &
CPABORT("Embedding spin cube file not found. ")
scaling_factor = 1.0_dp
CALL cp_cube_to_pw(spin_embed_pot, filename, scaling_factor)
END IF
END SUBROUTINE read_embed_pot_cube
! **************************************************************************************************
!> \brief Read the embedding potential from the binary file
!> \param qs_env ...
!> \param embed_pot ...
!> \param spin_embed_pot ...
!> \param section ...
!> \param embed_pot_coef ...
!> \param open_shell_embed ...
! **************************************************************************************************
SUBROUTINE read_embed_pot_vector(qs_env, embed_pot, spin_embed_pot, section, embed_pot_coef, open_shell_embed)
TYPE(qs_environment_type), POINTER :: qs_env
TYPE(pw_r3d_rs_type), INTENT(INOUT) :: embed_pot
TYPE(pw_r3d_rs_type), INTENT(IN), POINTER :: spin_embed_pot
TYPE(section_vals_type), POINTER :: section
TYPE(cp_fm_type), INTENT(IN) :: embed_pot_coef
LOGICAL, INTENT(IN) :: open_shell_embed
CHARACTER(LEN=default_path_length) :: filename
INTEGER :: dimen_aux, dimen_restart_basis, &
dimen_var_aux, l_global, LLL, &
nrow_local, restart_unit
INTEGER, DIMENSION(:), POINTER :: row_indices
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: coef, coef_read
TYPE(cp_blacs_env_type), POINTER :: blacs_env
TYPE(cp_fm_struct_type), POINTER :: fm_struct
TYPE(cp_fm_type) :: my_embed_pot_coef
TYPE(mp_para_env_type), POINTER :: para_env
! Get the vector dimension
CALL find_aux_dimen(qs_env, dimen_aux)
IF (open_shell_embed) THEN
dimen_var_aux = dimen_aux*2
ELSE
dimen_var_aux = dimen_aux
END IF
! We need a temporary vector of coefficients
CALL get_qs_env(qs_env=qs_env, &
para_env=para_env)
NULLIFY (blacs_env)
NULLIFY (fm_struct)
CALL cp_blacs_env_create(blacs_env=blacs_env, para_env=para_env)
CALL cp_fm_struct_create(fm_struct, para_env=para_env, context=blacs_env, &
nrow_global=dimen_var_aux, ncol_global=1)
CALL cp_fm_create(my_embed_pot_coef, fm_struct, name="my_pot_coef")
CALL cp_fm_struct_release(fm_struct)
CALL cp_fm_set_all(my_embed_pot_coef, 0.0_dp)
! Read the coefficients vector
restart_unit = -1
! Allocate the attay to read the coefficients
ALLOCATE (coef(dimen_var_aux))
coef = 0.0_dp
IF (para_env%is_source()) THEN
! Get the restart file name
CALL embed_restart_file_name(filename, section)
CALL open_file(file_name=filename, &
file_action="READ", &
file_form="UNFORMATTED", &
file_status="OLD", &
unit_number=restart_unit)
READ (restart_unit) dimen_restart_basis
! Check the dimensions of the bases: the actual and the restart one
IF (.NOT. (dimen_restart_basis == dimen_aux)) &
CPABORT("Wrong dimension of the embedding basis in the restart file.")
ALLOCATE (coef_read(dimen_var_aux))
coef_read = 0.0_dp
READ (restart_unit) coef_read
coef(:) = coef_read(:)
DEALLOCATE (coef_read)
! Close restart file
CALL close_file(unit_number=restart_unit)
END IF
! Broadcast the coefficients on all processes
CALL para_env%bcast(coef)
! Copy to fm_type structure
! Information about full matrix gradient
CALL cp_fm_get_info(matrix=my_embed_pot_coef, &
nrow_local=nrow_local, &
row_indices=row_indices)
DO LLL = 1, nrow_local
l_global = row_indices(LLL)
my_embed_pot_coef%local_data(LLL, 1) = coef(l_global)
END DO
DEALLOCATE (coef)
! Copy to the my_embed_pot_coef to embed_pot_coef
CALL cp_fm_copy_general(my_embed_pot_coef, embed_pot_coef, para_env)
! Build the embedding potential on the grid
CALL update_embed_pot(embed_pot_coef, dimen_aux, embed_pot, spin_embed_pot, &
qs_env, .FALSE., open_shell_embed)
! Release my_embed_pot_coef
CALL cp_fm_release(my_embed_pot_coef)
! Release blacs environment
CALL cp_blacs_env_release(blacs_env)
END SUBROUTINE read_embed_pot_vector
! **************************************************************************************************
!> \brief Find the embedding restart file name
!> \param filename ...
!> \param section ...
! **************************************************************************************************
SUBROUTINE embed_restart_file_name(filename, section)
CHARACTER(LEN=default_path_length), INTENT(OUT) :: filename
TYPE(section_vals_type), POINTER :: section
LOGICAL :: exist
exist = .FALSE.
CALL section_vals_val_get(section, "EMBED_RESTART_FILE_NAME", c_val=filename)
INQUIRE (FILE=filename, exist=exist)
IF (.NOT. exist) &
CPABORT("Embedding restart file not found. ")
END SUBROUTINE embed_restart_file_name
! **************************************************************************************************
!> \brief Deallocate stuff for optimizing embedding potential
!> \param opt_embed ...
! **************************************************************************************************
SUBROUTINE release_opt_embed(opt_embed)
TYPE(opt_embed_pot_type) :: opt_embed
INTEGER :: i_dens, i_spin, ikind
IF (.NOT. opt_embed%grid_opt) THEN
CALL cp_fm_release(opt_embed%embed_pot_grad)
CALL cp_fm_release(opt_embed%embed_pot_coef)
CALL cp_fm_release(opt_embed%step)
CALL cp_fm_release(opt_embed%prev_step)
CALL cp_fm_release(opt_embed%embed_pot_hess)
CALL cp_fm_release(opt_embed%prev_embed_pot_grad)
CALL cp_fm_release(opt_embed%prev_embed_pot_coef)
CALL cp_fm_release(opt_embed%prev_embed_pot_hess)
CALL cp_fm_release(opt_embed%kinetic_mat)
DEALLOCATE (opt_embed%embed_pot_grad, opt_embed%embed_pot_coef, &
opt_embed%step, opt_embed%prev_step, opt_embed%embed_pot_hess, &
opt_embed%prev_embed_pot_grad, opt_embed%prev_embed_pot_coef, &
opt_embed%prev_embed_pot_hess, opt_embed%kinetic_mat)
DEALLOCATE (opt_embed%w_func)
DEALLOCATE (opt_embed%max_diff)
DEALLOCATE (opt_embed%int_diff)
DO ikind = 1, SIZE(opt_embed%lri)
DEALLOCATE (opt_embed%lri(ikind)%v_int)
DEALLOCATE (opt_embed%lri(ikind)%acoef)
END DO
DEALLOCATE (opt_embed%lri)
END IF
IF (ASSOCIATED(opt_embed%prev_subsys_dens)) THEN
DO i_dens = 1, SIZE(opt_embed%prev_subsys_dens)
CALL opt_embed%prev_subsys_dens(i_dens)%release()
END DO
DEALLOCATE (opt_embed%prev_subsys_dens)
END IF
DEALLOCATE (opt_embed%max_subsys_dens_diff)
DEALLOCATE (opt_embed%all_nspins)
IF (ASSOCIATED(opt_embed%const_pot)) THEN
CALL opt_embed%const_pot%release()
DEALLOCATE (opt_embed%const_pot)
END IF
IF (ASSOCIATED(opt_embed%pot_diff)) THEN
CALL opt_embed%pot_diff%release()
DEALLOCATE (opt_embed%pot_diff)
END IF
IF (ASSOCIATED(opt_embed%prev_embed_pot)) THEN
CALL opt_embed%prev_embed_pot%release()
DEALLOCATE (opt_embed%prev_embed_pot)
END IF
IF (ASSOCIATED(opt_embed%prev_spin_embed_pot)) THEN
CALL opt_embed%prev_spin_embed_pot%release()
DEALLOCATE (opt_embed%prev_spin_embed_pot)
END IF
IF (ASSOCIATED(opt_embed%v_w)) THEN
DO i_spin = 1, SIZE(opt_embed%v_w)
CALL opt_embed%v_w(i_spin)%release()
END DO
DEALLOCATE (opt_embed%v_w)
END IF
END SUBROUTINE release_opt_embed
! **************************************************************************************************
!> \brief Calculates subsystem Coulomb potential from the RESP charges of the total system
!> \param v_rspace ...
!> \param rhs ...
!> \param mapping_section ...
!> \param qs_env ...
!> \param nforce_eval ...
!> \param iforce_eval ...
!> \param eta ...
! **************************************************************************************************
SUBROUTINE Coulomb_guess(v_rspace, rhs, mapping_section, qs_env, nforce_eval, iforce_eval, eta)
TYPE(pw_r3d_rs_type) :: v_rspace
REAL(KIND=dp), DIMENSION(:), POINTER :: rhs
TYPE(section_vals_type), POINTER :: mapping_section
TYPE(qs_environment_type), POINTER :: qs_env
INTEGER :: nforce_eval, iforce_eval
REAL(KIND=dp) :: eta
INTEGER :: iparticle, jparticle, natom
INTEGER, DIMENSION(:), POINTER :: map_index
REAL(KIND=dp) :: dvol, normalize_factor
REAL(KIND=dp), DIMENSION(:), POINTER :: rhs_subsys
TYPE(particle_list_type), POINTER :: particles
TYPE(pw_c1d_gs_type) :: v_resp_gspace
TYPE(pw_env_type), POINTER :: pw_env
TYPE(pw_poisson_type), POINTER :: poisson_env
TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
TYPE(pw_r3d_rs_type) :: rho_resp, v_resp_rspace
TYPE(qs_subsys_type), POINTER :: subsys
! Get available particles
NULLIFY (subsys)
CALL get_qs_env(qs_env=qs_env, subsys=subsys, pw_env=pw_env)
CALL qs_subsys_get(subsys, particles=particles)
natom = particles%n_els
ALLOCATE (rhs_subsys(natom))
NULLIFY (map_index)
CALL get_subsys_map_index(mapping_section, natom, iforce_eval, nforce_eval, &
map_index, .TRUE.)
! Mapping particles from iforce_eval environment to the embed env
DO iparticle = 1, natom
jparticle = map_index(iparticle)
rhs_subsys(iparticle) = rhs(jparticle)
END DO
! Prepare plane waves
NULLIFY (auxbas_pw_pool)
CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool, &
poisson_env=poisson_env)
CALL auxbas_pw_pool%create_pw(v_resp_gspace)
CALL auxbas_pw_pool%create_pw(v_resp_rspace)
CALL auxbas_pw_pool%create_pw(rho_resp)
! Calculate charge density
CALL pw_zero(rho_resp)
CALL calculate_rho_resp_all(rho_resp, rhs_subsys, natom, eta, qs_env)
! Calculate potential
CALL pw_poisson_solve(poisson_env, rho_resp, &
vhartree=v_resp_rspace)
dvol = v_resp_rspace%pw_grid%dvol
CALL pw_scale(v_resp_rspace, dvol)
normalize_factor = SQRT((eta/pi)**3)
!normalize_factor = -2.0_dp
CALL pw_scale(v_resp_rspace, normalize_factor)
! Hard copy potential
CALL pw_copy(v_resp_rspace, v_rspace)
! Release plane waves
CALL v_resp_gspace%release()