-
Notifications
You must be signed in to change notification settings - Fork 1
/
energy_corrections.F
3682 lines (3177 loc) · 164 KB
/
energy_corrections.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 Routines for an energy correction on top of a Kohn-Sham calculation
!> \par History
!> 03.2014 created
!> 09.2019 Moved from KG to Kohn-Sham
!> 08.2022 Add Density-Corrected DFT methods
!> 04.2023 Add hybrid functionals for DC-DFT
!> 10.2024 Add external energy method
!> \author JGH
! **************************************************************************************************
MODULE energy_corrections
USE admm_dm_methods, ONLY: admm_dm_calc_rho_aux
USE admm_methods, ONLY: admm_mo_calc_rho_aux
USE atomic_kind_types, ONLY: atomic_kind_type,&
get_atomic_kind
USE basis_set_types, ONLY: get_gto_basis_set,&
gto_basis_set_type
USE bibliography, ONLY: Belleflamme2023,&
cite_reference
USE cell_types, ONLY: cell_type,&
pbc
USE core_ae, ONLY: build_core_ae
USE core_ppl, ONLY: build_core_ppl
USE core_ppnl, ONLY: build_core_ppnl
USE cp_blacs_env, ONLY: cp_blacs_env_type
USE cp_control_types, ONLY: dft_control_type
USE cp_dbcsr_api, ONLY: &
dbcsr_add, dbcsr_copy, dbcsr_create, dbcsr_desymmetrize, dbcsr_distribution_type, &
dbcsr_dot, dbcsr_filter, dbcsr_get_info, dbcsr_init_p, dbcsr_multiply, dbcsr_p_type, &
dbcsr_release, dbcsr_scale, dbcsr_set, dbcsr_type, dbcsr_type_no_symmetry, &
dbcsr_type_symmetric
USE cp_dbcsr_cp2k_link, ONLY: cp_dbcsr_alloc_block_from_nbl
USE cp_dbcsr_operations, ONLY: copy_dbcsr_to_fm,&
copy_fm_to_dbcsr,&
cp_dbcsr_sm_fm_multiply,&
dbcsr_allocate_matrix_set,&
dbcsr_deallocate_matrix_set
USE cp_fm_basic_linalg, ONLY: cp_fm_column_scale,&
cp_fm_triangular_invert
USE cp_fm_cholesky, ONLY: cp_fm_cholesky_decompose
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_create,&
cp_fm_init_random,&
cp_fm_release,&
cp_fm_set_all,&
cp_fm_to_fm_triangular,&
cp_fm_type
USE cp_log_handling, ONLY: cp_get_default_logger,&
cp_logger_get_default_unit_nr,&
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_result_methods, ONLY: cp_results_erase,&
put_results
USE cp_result_types, ONLY: cp_result_type
USE cp_units, ONLY: cp_unit_from_cp2k
USE distribution_1d_types, ONLY: distribution_1d_type
USE distribution_2d_types, ONLY: distribution_2d_type
USE dm_ls_scf, ONLY: calculate_w_matrix_ls,&
post_scf_sparsities
USE dm_ls_scf_methods, ONLY: apply_matrix_preconditioner,&
density_matrix_sign,&
density_matrix_tc2,&
density_matrix_trs4,&
ls_scf_init_matrix_s
USE dm_ls_scf_qs, ONLY: matrix_ls_create,&
matrix_ls_to_qs,&
matrix_qs_to_ls
USE dm_ls_scf_types, ONLY: ls_scf_env_type
USE ec_efield_local, ONLY: ec_efield_integrals,&
ec_efield_local_operator
USE ec_env_types, ONLY: energy_correction_type
USE ec_external, ONLY: ec_ext_energy
USE ec_methods, ONLY: create_kernel,&
ec_mos_init
USE external_potential_types, ONLY: get_potential,&
gth_potential_type,&
sgp_potential_type
USE hfx_exx, ONLY: add_exx_to_rhs,&
calculate_exx
USE input_constants, ONLY: &
ec_diagonalization, ec_functional_dc, ec_functional_ext, ec_functional_harris, &
ec_matrix_sign, ec_matrix_tc2, ec_matrix_trs4, ec_ot_atomic, ec_ot_diag, ec_ot_gs, &
ot_precond_full_single_inverse, ot_precond_solver_default, vdw_pairpot_dftd3, &
vdw_pairpot_dftd3bj, xc_vdw_fun_pairpot
USE input_section_types, ONLY: section_get_ival,&
section_get_lval,&
section_vals_get,&
section_vals_get_subs_vals,&
section_vals_type,&
section_vals_val_get
USE kinds, ONLY: default_path_length,&
default_string_length,&
dp
USE mao_basis, ONLY: mao_generate_basis
USE mathlib, ONLY: det_3x3
USE message_passing, ONLY: mp_para_env_type
USE molecule_types, ONLY: molecule_type
USE moments_utils, ONLY: get_reference_point
USE particle_types, ONLY: particle_type
USE periodic_table, ONLY: ptable
USE physcon, ONLY: bohr,&
debye,&
pascal
USE preconditioner, ONLY: make_preconditioner
USE preconditioner_types, ONLY: destroy_preconditioner,&
init_preconditioner,&
preconditioner_type
USE pw_env_types, ONLY: pw_env_get,&
pw_env_type
USE pw_grid_types, ONLY: pw_grid_type
USE pw_methods, ONLY: pw_axpy,&
pw_copy,&
pw_integral_ab,&
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_p_type,&
pw_pool_type
USE pw_types, ONLY: pw_c1d_gs_type,&
pw_r3d_rs_type
USE qs_atomic_block, ONLY: calculate_atomic_block_dm
USE qs_collocate_density, ONLY: calculate_rho_elec
USE qs_core_energies, ONLY: calculate_ecore_overlap,&
calculate_ptrace
USE qs_density_matrices, ONLY: calculate_density_matrix,&
calculate_w_matrix
USE qs_dispersion_pairpot, ONLY: calculate_dispersion_pairpot
USE qs_dispersion_types, ONLY: qs_dispersion_type
USE qs_energy_types, ONLY: qs_energy_type
USE qs_environment_types, ONLY: get_qs_env,&
qs_environment_type
USE qs_force_types, ONLY: qs_force_type,&
total_qs_force,&
zero_qs_force
USE qs_integrate_potential, ONLY: integrate_v_core_rspace,&
integrate_v_rspace
USE qs_kind_types, ONLY: get_qs_kind,&
get_qs_kind_set,&
qs_kind_type
USE qs_kinetic, ONLY: build_kinetic_matrix
USE qs_ks_methods, ONLY: calc_rho_tot_gspace
USE qs_ks_reference, ONLY: ks_ref_potential
USE qs_ks_types, ONLY: qs_ks_env_type
USE qs_mo_methods, ONLY: calculate_subspace_eigenvalues,&
make_basis_sm
USE qs_mo_types, ONLY: deallocate_mo_set,&
get_mo_set,&
mo_set_type
USE qs_moments, ONLY: build_local_moment_matrix
USE qs_neighbor_list_types, ONLY: neighbor_list_set_p_type
USE qs_neighbor_lists, ONLY: atom2d_build,&
atom2d_cleanup,&
build_neighbor_lists,&
local_atoms_type,&
pair_radius_setup
USE qs_ot_eigensolver, ONLY: ot_eigensolver
USE qs_overlap, ONLY: build_overlap_matrix
USE qs_rho_types, ONLY: qs_rho_get,&
qs_rho_type
USE qs_vxc, ONLY: qs_vxc_create
USE response_solver, ONLY: response_calculation,&
response_force
USE rtp_admm_methods, ONLY: rtp_admm_calc_rho_aux
USE string_utilities, ONLY: uppercase
USE task_list_methods, ONLY: generate_qs_task_list
USE task_list_types, ONLY: allocate_task_list,&
deallocate_task_list
USE virial_methods, ONLY: one_third_sum_diag,&
write_stress_tensor,&
write_stress_tensor_components
USE virial_types, ONLY: symmetrize_virial,&
virial_type
USE voronoi_interface, ONLY: entry_voronoi_or_bqb
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
! Global parameters
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'energy_corrections'
PUBLIC :: energy_correction
CONTAINS
! **************************************************************************************************
!> \brief Energy Correction to a Kohn-Sham simulation
!> Available energy corrections: (1) Harris energy functional
!> (2) Density-corrected DFT
!> (3) Energy from external source
!>
!> \param qs_env ...
!> \param ec_init ...
!> \param calculate_forces ...
!> \par History
!> 03.2014 created
!> \author JGH
! **************************************************************************************************
SUBROUTINE energy_correction(qs_env, ec_init, calculate_forces)
TYPE(qs_environment_type), POINTER :: qs_env
LOGICAL, INTENT(IN), OPTIONAL :: ec_init, calculate_forces
CHARACTER(len=*), PARAMETER :: routineN = 'energy_correction'
INTEGER :: handle, unit_nr
LOGICAL :: my_calc_forces
TYPE(cp_logger_type), POINTER :: logger
TYPE(energy_correction_type), POINTER :: ec_env
TYPE(qs_energy_type), POINTER :: energy
TYPE(qs_force_type), DIMENSION(:), POINTER :: ks_force
CALL timeset(routineN, handle)
logger => cp_get_default_logger()
IF (logger%para_env%is_source()) THEN
unit_nr = cp_logger_get_default_unit_nr(logger, local=.TRUE.)
ELSE
unit_nr = -1
END IF
CALL cite_reference(Belleflamme2023)
NULLIFY (ec_env)
CALL get_qs_env(qs_env, ec_env=ec_env)
! Skip energy correction if ground-state is NOT converged
IF (.NOT. ec_env%do_skip) THEN
ec_env%should_update = .TRUE.
IF (PRESENT(ec_init)) ec_env%should_update = ec_init
my_calc_forces = .FALSE.
IF (PRESENT(calculate_forces)) my_calc_forces = calculate_forces
IF (ec_env%should_update) THEN
ec_env%old_etotal = 0.0_dp
ec_env%etotal = 0.0_dp
ec_env%eband = 0.0_dp
ec_env%ehartree = 0.0_dp
ec_env%ex = 0.0_dp
ec_env%exc = 0.0_dp
ec_env%vhxc = 0.0_dp
ec_env%edispersion = 0.0_dp
ec_env%exc_aux_fit = 0.0_dp
! Save total energy of reference calculation
CALL get_qs_env(qs_env, energy=energy)
ec_env%old_etotal = energy%total
END IF
IF (my_calc_forces) THEN
IF (unit_nr > 0) THEN
WRITE (unit_nr, '(T2,A,A,A,A,A)') "!", REPEAT("-", 25), &
" Energy Correction Forces ", REPEAT("-", 26), "!"
END IF
CALL get_qs_env(qs_env, force=ks_force)
CALL zero_qs_force(ks_force)
ELSE
IF (unit_nr > 0) THEN
WRITE (unit_nr, '(T2,A,A,A,A,A)') "!", REPEAT("-", 29), &
" Energy Correction ", REPEAT("-", 29), "!"
END IF
END IF
! Perform the energy correction
CALL energy_correction_low(qs_env, ec_env, my_calc_forces, unit_nr)
! Update total energy in qs environment and amount fo correction
IF (ec_env%should_update) THEN
energy%nonscf_correction = ec_env%etotal - ec_env%old_etotal
energy%total = ec_env%etotal
END IF
IF (.NOT. my_calc_forces .AND. unit_nr > 0) THEN
WRITE (unit_nr, '(T3,A,T56,F25.15)') "Energy Correction ", energy%nonscf_correction
END IF
IF (unit_nr > 0) THEN
WRITE (unit_nr, '(T2,A,A,A)') "!", REPEAT("-", 77), "!"
END IF
ELSE
! Ground-state energy calculation did not converge,
! do not calculate energy correction
IF (unit_nr > 0) THEN
WRITE (unit_nr, '(T2,A,A,A)') "!", REPEAT("-", 77), "!"
WRITE (unit_nr, '(T2,A,A,A,A,A)') "!", REPEAT("-", 26), &
" Skip Energy Correction ", REPEAT("-", 27), "!"
WRITE (unit_nr, '(T2,A,A,A)') "!", REPEAT("-", 77), "!"
END IF
END IF
CALL timestop(handle)
END SUBROUTINE energy_correction
! **************************************************************************************************
!> \brief Energy Correction to a Kohn-Sham simulation
!>
!> \param qs_env ...
!> \param ec_env ...
!> \param calculate_forces ...
!> \param unit_nr ...
!> \par History
!> 03.2014 created
!> \author JGH
! **************************************************************************************************
SUBROUTINE energy_correction_low(qs_env, ec_env, calculate_forces, unit_nr)
TYPE(qs_environment_type), POINTER :: qs_env
TYPE(energy_correction_type), POINTER :: ec_env
LOGICAL, INTENT(IN) :: calculate_forces
INTEGER, INTENT(IN) :: unit_nr
INTEGER :: ispin, nspins
LOGICAL :: debug_f
REAL(KIND=dp) :: exc
TYPE(dft_control_type), POINTER :: dft_control
TYPE(pw_env_type), POINTER :: pw_env
TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
IF (ec_env%should_update) THEN
CALL ec_build_neighborlist(qs_env, ec_env)
CALL ks_ref_potential(qs_env, &
ec_env%vh_rspace, &
ec_env%vxc_rspace, &
ec_env%vtau_rspace, &
ec_env%vadmm_rspace, &
ec_env%ehartree, exc)
SELECT CASE (ec_env%energy_functional)
CASE (ec_functional_harris)
CALL ec_build_core_hamiltonian(qs_env, ec_env)
CALL ec_build_ks_matrix(qs_env, ec_env)
IF (ec_env%mao) THEN
! MAO basis
IF (ASSOCIATED(ec_env%mao_coef)) CALL dbcsr_deallocate_matrix_set(ec_env%mao_coef)
NULLIFY (ec_env%mao_coef)
CALL mao_generate_basis(qs_env, ec_env%mao_coef, ref_basis_set="HARRIS", molecular=.TRUE., &
max_iter=ec_env%mao_max_iter, eps_grad=ec_env%mao_eps_grad, &
eps1_mao=ec_env%mao_eps1, iolevel=ec_env%mao_iolevel, unit_nr=unit_nr)
END IF
CALL ec_ks_solver(qs_env, ec_env)
CALL evaluate_ec_core_matrix_traces(qs_env, ec_env)
CASE (ec_functional_dc)
! Prepare Density-corrected DFT (DC-DFT) calculation
CALL ec_dc_energy(qs_env, ec_env, calculate_forces=.FALSE.)
! Rebuild KS matrix with DC-DFT XC functional evaluated in ground-state density.
! KS matrix might contain unwanted contributions
! Calculate Hartree and XC related energies here
CALL ec_build_ks_matrix(qs_env, ec_env)
CASE (ec_functional_ext)
CALL ec_ext_energy(qs_env, ec_env, calculate_forces=.FALSE.)
CASE DEFAULT
CPABORT("unknown energy correction")
END SELECT
! dispersion through pairpotentials
CALL ec_disp(qs_env, ec_env, calculate_forces=.FALSE.)
! Calculate total energy
CALL ec_energy(ec_env, unit_nr)
END IF
IF (calculate_forces) THEN
debug_f = ec_env%debug_forces .OR. ec_env%debug_stress
CALL ec_disp(qs_env, ec_env, calculate_forces=.TRUE.)
SELECT CASE (ec_env%energy_functional)
CASE (ec_functional_harris)
CALL ec_build_core_hamiltonian_force(qs_env, ec_env, &
ec_env%matrix_p, &
ec_env%matrix_s, &
ec_env%matrix_w)
CALL ec_build_ks_matrix_force(qs_env, ec_env)
CASE (ec_functional_dc)
! Prepare Density-corrected DFT (DC-DFT) calculation
! by getting ground-state matrices
CALL ec_dc_energy(qs_env, ec_env, calculate_forces=.TRUE.)
CALL ec_build_core_hamiltonian_force(qs_env, ec_env, &
ec_env%matrix_p, &
ec_env%matrix_s, &
ec_env%matrix_w)
CALL ec_dc_build_ks_matrix_force(qs_env, ec_env)
CASE (ec_functional_ext)
CALL ec_ext_energy(qs_env, ec_env, calculate_forces=.TRUE.)
CASE DEFAULT
CPABORT("unknown energy correction")
END SELECT
CALL response_calculation(qs_env, ec_env)
! Allocate response density on real space grid for use in properties
! Calculated in response_force
CALL get_qs_env(qs_env=qs_env, dft_control=dft_control, pw_env=pw_env)
nspins = dft_control%nspins
CPASSERT(ASSOCIATED(pw_env))
CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
ALLOCATE (ec_env%rhoz_r(nspins))
DO ispin = 1, nspins
CALL auxbas_pw_pool%create_pw(ec_env%rhoz_r(ispin))
END DO
CALL response_force(qs_env, &
vh_rspace=ec_env%vh_rspace, &
vxc_rspace=ec_env%vxc_rspace, &
vtau_rspace=ec_env%vtau_rspace, &
vadmm_rspace=ec_env%vadmm_rspace, &
matrix_hz=ec_env%matrix_hz, &
matrix_pz=ec_env%matrix_z, &
matrix_pz_admm=ec_env%z_admm, &
matrix_wz=ec_env%matrix_wz, &
rhopz_r=ec_env%rhoz_r, &
zehartree=ec_env%ehartree, &
zexc=ec_env%exc, &
zexc_aux_fit=ec_env%exc_aux_fit, &
p_env=ec_env%p_env, &
debug=debug_f)
CALL ec_properties(qs_env, ec_env)
! Deallocate Harris density and response density on grid
IF (ASSOCIATED(ec_env%rhoout_r)) THEN
DO ispin = 1, nspins
CALL auxbas_pw_pool%give_back_pw(ec_env%rhoout_r(ispin))
END DO
DEALLOCATE (ec_env%rhoout_r)
END IF
IF (ASSOCIATED(ec_env%rhoz_r)) THEN
DO ispin = 1, nspins
CALL auxbas_pw_pool%give_back_pw(ec_env%rhoz_r(ispin))
END DO
DEALLOCATE (ec_env%rhoz_r)
END IF
! Deallocate matrices
IF (ASSOCIATED(ec_env%matrix_ks)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_ks)
IF (ASSOCIATED(ec_env%matrix_h)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_h)
IF (ASSOCIATED(ec_env%matrix_s)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_s)
IF (ASSOCIATED(ec_env%matrix_t)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_t)
IF (ASSOCIATED(ec_env%matrix_p)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_p)
IF (ASSOCIATED(ec_env%matrix_w)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_w)
IF (ASSOCIATED(ec_env%matrix_hz)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_hz)
IF (ASSOCIATED(ec_env%matrix_wz)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_wz)
IF (ASSOCIATED(ec_env%matrix_z)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_z)
END IF
END SUBROUTINE energy_correction_low
! **************************************************************************************************
!> \brief Calculates the traces of the core matrices and the density matrix.
!> \param qs_env ...
!> \param ec_env ...
!> \author Ole Schuett
!> adapted for energy correction fbelle
! **************************************************************************************************
SUBROUTINE evaluate_ec_core_matrix_traces(qs_env, ec_env)
TYPE(qs_environment_type), POINTER :: qs_env
TYPE(energy_correction_type), POINTER :: ec_env
CHARACTER(LEN=*), PARAMETER :: routineN = 'evaluate_ec_core_matrix_traces'
INTEGER :: handle
TYPE(dft_control_type), POINTER :: dft_control
TYPE(qs_energy_type), POINTER :: energy
CALL timeset(routineN, handle)
NULLIFY (energy)
CALL get_qs_env(qs_env, dft_control=dft_control, energy=energy)
! Core hamiltonian energy
CALL calculate_ptrace(ec_env%matrix_h, ec_env%matrix_p, energy%core, dft_control%nspins)
! kinetic energy
CALL calculate_ptrace(ec_env%matrix_t, ec_env%matrix_p, energy%kinetic, dft_control%nspins)
CALL timestop(handle)
END SUBROUTINE evaluate_ec_core_matrix_traces
! **************************************************************************************************
!> \brief Prepare DC-DFT calculation by copying unaffected ground-state matrices (core Hamiltonian,
!> density matrix) into energy correction environment and rebuild the overlap matrix
!>
!> \param qs_env ...
!> \param ec_env ...
!> \param calculate_forces ...
!> \par History
!> 07.2022 created
!> \author fbelle
! **************************************************************************************************
SUBROUTINE ec_dc_energy(qs_env, ec_env, calculate_forces)
TYPE(qs_environment_type), POINTER :: qs_env
TYPE(energy_correction_type), POINTER :: ec_env
LOGICAL, INTENT(IN) :: calculate_forces
CHARACTER(LEN=*), PARAMETER :: routineN = 'ec_dc_energy'
CHARACTER(LEN=default_string_length) :: headline
INTEGER :: handle, ispin, nspins
TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_h, matrix_p, matrix_s, matrix_w
TYPE(dft_control_type), POINTER :: dft_control
TYPE(qs_ks_env_type), POINTER :: ks_env
TYPE(qs_rho_type), POINTER :: rho
CALL timeset(routineN, handle)
NULLIFY (dft_control, ks_env, matrix_h, matrix_p, matrix_s, matrix_w, rho)
CALL get_qs_env(qs_env=qs_env, &
dft_control=dft_control, &
ks_env=ks_env, &
matrix_h_kp=matrix_h, &
matrix_s_kp=matrix_s, &
matrix_w_kp=matrix_w, &
rho=rho)
CALL qs_rho_get(rho, rho_ao_kp=matrix_p)
nspins = dft_control%nspins
! For density-corrected DFT only the ground-state matrices are required
! Comply with ec_env environment for property calculations later
CALL build_overlap_matrix(ks_env, matrixkp_s=ec_env%matrix_s, &
matrix_name="OVERLAP MATRIX", &
basis_type_a="HARRIS", &
basis_type_b="HARRIS", &
sab_nl=ec_env%sab_orb)
! Core Hamiltonian matrix
IF (ASSOCIATED(ec_env%matrix_h)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_h)
CALL dbcsr_allocate_matrix_set(ec_env%matrix_h, 1, 1)
headline = "CORE HAMILTONIAN MATRIX"
ALLOCATE (ec_env%matrix_h(1, 1)%matrix)
CALL dbcsr_create(ec_env%matrix_h(1, 1)%matrix, name=TRIM(headline), &
template=matrix_s(1, 1)%matrix, matrix_type=dbcsr_type_symmetric)
CALL cp_dbcsr_alloc_block_from_nbl(ec_env%matrix_h(1, 1)%matrix, ec_env%sab_orb)
CALL dbcsr_copy(ec_env%matrix_h(1, 1)%matrix, matrix_h(1, 1)%matrix)
! Density matrix
IF (ASSOCIATED(ec_env%matrix_p)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_p)
CALL dbcsr_allocate_matrix_set(ec_env%matrix_p, nspins, 1)
headline = "DENSITY MATRIX"
DO ispin = 1, nspins
ALLOCATE (ec_env%matrix_p(ispin, 1)%matrix)
CALL dbcsr_create(ec_env%matrix_p(ispin, 1)%matrix, name=TRIM(headline), &
template=matrix_s(1, 1)%matrix, matrix_type=dbcsr_type_symmetric)
CALL cp_dbcsr_alloc_block_from_nbl(ec_env%matrix_p(ispin, 1)%matrix, ec_env%sab_orb)
CALL dbcsr_copy(ec_env%matrix_p(ispin, 1)%matrix, matrix_p(ispin, 1)%matrix)
END DO
IF (calculate_forces) THEN
! Energy-weighted density matrix
IF (ASSOCIATED(ec_env%matrix_w)) CALL dbcsr_deallocate_matrix_set(ec_env%matrix_w)
CALL dbcsr_allocate_matrix_set(ec_env%matrix_w, nspins, 1)
headline = "ENERGY-WEIGHTED DENSITY MATRIX"
DO ispin = 1, nspins
ALLOCATE (ec_env%matrix_w(ispin, 1)%matrix)
CALL dbcsr_create(ec_env%matrix_w(ispin, 1)%matrix, name=TRIM(headline), &
template=matrix_s(1, 1)%matrix, matrix_type=dbcsr_type_symmetric)
CALL cp_dbcsr_alloc_block_from_nbl(ec_env%matrix_w(ispin, 1)%matrix, ec_env%sab_orb)
CALL dbcsr_copy(ec_env%matrix_w(ispin, 1)%matrix, matrix_w(ispin, 1)%matrix)
END DO
END IF
! External field (nonperiodic case)
ec_env%efield_nuclear = 0.0_dp
ec_env%efield_elec = 0.0_dp
CALL ec_efield_local_operator(qs_env, ec_env, calculate_forces=.FALSE.)
CALL timestop(handle)
END SUBROUTINE ec_dc_energy
! **************************************************************************************************
!> \brief Kohn-Sham matrix contributions to force in DC-DFT
!> also calculate right-hand-side matrix B for response equations AX=B
!> \param qs_env ...
!> \param ec_env ...
!> \par History
!> 08.2022 adapted from qs_ks_build_kohn_sham_matrix
!> \author fbelle
! **************************************************************************************************
SUBROUTINE ec_dc_build_ks_matrix_force(qs_env, ec_env)
TYPE(qs_environment_type), POINTER :: qs_env
TYPE(energy_correction_type), POINTER :: ec_env
CHARACTER(LEN=*), PARAMETER :: routineN = 'ec_dc_build_ks_matrix_force'
CHARACTER(LEN=default_string_length) :: unit_string
INTEGER :: handle, i, iounit, ispin, natom, nspins
LOGICAL :: debug_forces, debug_stress, do_ec_hfx, &
use_virial
REAL(dp) :: dummy_real, dummy_real2(2), ehartree, &
eovrl, exc, fconv
REAL(dp), ALLOCATABLE, DIMENSION(:, :) :: ftot
REAL(dp), DIMENSION(3) :: fodeb, fodeb2
REAL(KIND=dp), DIMENSION(3, 3) :: h_stress, pv_loc, stdeb, sttot
TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
TYPE(cell_type), POINTER :: cell
TYPE(cp_logger_type), POINTER :: logger
TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_ks, matrix_s, scrm
TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_p
TYPE(dft_control_type), POINTER :: dft_control
TYPE(mp_para_env_type), POINTER :: para_env
TYPE(neighbor_list_set_p_type), DIMENSION(:), &
POINTER :: sab_orb
TYPE(pw_c1d_gs_type) :: rho_tot_gspace, v_hartree_gspace
TYPE(pw_env_type), POINTER :: pw_env
TYPE(pw_grid_type), POINTER :: pw_grid
TYPE(pw_poisson_type), POINTER :: poisson_env
TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
TYPE(pw_r3d_rs_type) :: v_hartree_rspace
TYPE(pw_r3d_rs_type), DIMENSION(:), POINTER :: rho_r, v_rspace, v_rspace_in, &
v_tau_rspace
TYPE(qs_force_type), DIMENSION(:), POINTER :: force
TYPE(qs_ks_env_type), POINTER :: ks_env
TYPE(qs_rho_type), POINTER :: rho
TYPE(section_vals_type), POINTER :: ec_hfx_sections
TYPE(virial_type), POINTER :: virial
CALL timeset(routineN, handle)
debug_forces = ec_env%debug_forces
debug_stress = ec_env%debug_stress
logger => cp_get_default_logger()
IF (logger%para_env%is_source()) THEN
iounit = cp_logger_get_default_unit_nr(logger, local=.TRUE.)
ELSE
iounit = -1
END IF
NULLIFY (atomic_kind_set, cell, dft_control, force, ks_env, matrix_ks, &
matrix_p, matrix_s, para_env, pw_env, rho, sab_orb, virial)
CALL get_qs_env(qs_env=qs_env, &
cell=cell, &
dft_control=dft_control, &
force=force, &
ks_env=ks_env, &
matrix_ks=matrix_ks, &
matrix_s=matrix_s, &
para_env=para_env, &
pw_env=pw_env, &
rho=rho, &
sab_orb=sab_orb, &
virial=virial)
CPASSERT(ASSOCIATED(pw_env))
nspins = dft_control%nspins
use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
fconv = 1.0E-9_dp*pascal/cell%deth
IF (debug_stress .AND. use_virial) THEN
sttot = virial%pv_virial
END IF
! Get density matrix of reference calculation
CALL qs_rho_get(rho, rho_ao_kp=matrix_p)
NULLIFY (auxbas_pw_pool, poisson_env)
! gets the tmp grids
CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool, &
poisson_env=poisson_env)
! Calculate the Hartree potential
CALL auxbas_pw_pool%create_pw(v_hartree_gspace)
CALL auxbas_pw_pool%create_pw(rho_tot_gspace)
CALL auxbas_pw_pool%create_pw(v_hartree_rspace)
! Get the total input density in g-space [ions + electrons]
CALL calc_rho_tot_gspace(rho_tot_gspace, qs_env, rho)
! v_H[n_in]
IF (use_virial) THEN
! Stress tensor - Volume and Green function contribution
h_stress(:, :) = 0.0_dp
CALL pw_poisson_solve(poisson_env, &
density=rho_tot_gspace, &
ehartree=ehartree, &
vhartree=v_hartree_gspace, &
h_stress=h_stress)
virial%pv_ehartree = virial%pv_ehartree + h_stress/REAL(para_env%num_pe, dp)
virial%pv_virial = virial%pv_virial + h_stress/REAL(para_env%num_pe, dp)
IF (debug_stress) THEN
stdeb = fconv*(h_stress/REAL(para_env%num_pe, dp))
CALL para_env%sum(stdeb)
IF (iounit > 0) WRITE (UNIT=iounit, FMT="(T2,A,T41,2(1X,ES19.11))") &
'STRESS| GREEN 1st V_H[n_in]*n_in ', one_third_sum_diag(stdeb), det_3x3(stdeb)
END IF
ELSE
CALL pw_poisson_solve(poisson_env, rho_tot_gspace, ehartree, &
v_hartree_gspace)
END IF
CALL pw_transfer(v_hartree_gspace, v_hartree_rspace)
CALL pw_scale(v_hartree_rspace, v_hartree_rspace%pw_grid%dvol)
! Save density on real space grid for use in properties
CALL qs_rho_get(rho, rho_r=rho_r)
ALLOCATE (ec_env%rhoout_r(nspins))
DO ispin = 1, nspins
CALL auxbas_pw_pool%create_pw(ec_env%rhoout_r(ispin))
CALL pw_copy(rho_r(ispin), ec_env%rhoout_r(ispin))
END DO
! Getting nuclear force contribution from the core charge density
! Vh(rho_c + rho_in)
IF (debug_forces) fodeb(1:3) = force(1)%rho_core(1:3, 1)
IF (debug_stress .AND. use_virial) stdeb = virial%pv_ehartree
CALL integrate_v_core_rspace(v_hartree_rspace, qs_env)
IF (debug_forces) THEN
fodeb(1:3) = force(1)%rho_core(1:3, 1) - fodeb(1:3)
CALL para_env%sum(fodeb)
IF (iounit > 0) WRITE (iounit, "(T3,A,T33,3F16.8)") "DEBUG:: Vtot*dncore", fodeb
END IF
IF (debug_stress .AND. use_virial) THEN
stdeb = fconv*(virial%pv_ehartree - stdeb)
CALL para_env%sum(stdeb)
IF (iounit > 0) WRITE (UNIT=iounit, FMT="(T2,A,T41,2(1X,ES19.11))") &
'STRESS| Vtot*dncore', one_third_sum_diag(stdeb), det_3x3(stdeb)
END IF
! v_XC[n_in]_DC
! v_rspace and v_tau_rspace are generated from the auxbas pool
NULLIFY (v_rspace, v_tau_rspace)
! only activate stress calculation if
IF (use_virial) virial%pv_calculate = .TRUE.
! Exchange-correlation potential
CALL qs_vxc_create(ks_env=ks_env, rho_struct=rho, xc_section=ec_env%xc_section, &
vxc_rho=v_rspace, vxc_tau=v_tau_rspace, exc=exc, just_energy=.FALSE.)
IF (.NOT. ASSOCIATED(v_rspace)) THEN
ALLOCATE (v_rspace(nspins))
DO ispin = 1, nspins
CALL auxbas_pw_pool%create_pw(v_rspace(ispin))
CALL pw_zero(v_rspace(ispin))
END DO
END IF
IF (use_virial) THEN
virial%pv_exc = virial%pv_exc - virial%pv_xc
virial%pv_virial = virial%pv_virial - virial%pv_xc
! virial%pv_xc will be zeroed in the xc routines
END IF
! initialize srcm matrix
NULLIFY (scrm)
CALL dbcsr_allocate_matrix_set(scrm, nspins)
DO ispin = 1, nspins
ALLOCATE (scrm(ispin)%matrix)
CALL dbcsr_create(scrm(ispin)%matrix, template=ec_env%matrix_ks(ispin, 1)%matrix)
CALL dbcsr_copy(scrm(ispin)%matrix, ec_env%matrix_ks(ispin, 1)%matrix)
CALL dbcsr_set(scrm(ispin)%matrix, 0.0_dp)
END DO
pw_grid => v_hartree_rspace%pw_grid
ALLOCATE (v_rspace_in(nspins))
DO ispin = 1, nspins
CALL v_rspace_in(ispin)%create(pw_grid)
END DO
! v_rspace_in = v_H[n_in] + v_xc[n_in] calculated in ks_ref_potential
DO ispin = 1, nspins
! v_xc[n_in]_GS
CALL pw_transfer(ec_env%vxc_rspace(ispin), v_rspace_in(ispin))
! add v_H[n_in]
CALL pw_axpy(ec_env%vh_rspace, v_rspace_in(ispin))
END DO
!------------------------------------------------
! If hybrid functional in DC-DFT
ec_hfx_sections => section_vals_get_subs_vals(qs_env%input, "DFT%ENERGY_CORRECTION%XC%HF")
CALL section_vals_get(ec_hfx_sections, explicit=do_ec_hfx)
IF (do_ec_hfx) THEN
IF (debug_forces) fodeb(1:3) = force(1)%fock_4c(1:3, 1)
IF (debug_forces) fodeb2(1:3) = force(1)%overlap_admm(1:3, 1)
! Calculate direct HFX forces here
! Virial contribution (fock_4c) done inside calculate_exx
dummy_real = 0.0_dp
CALL calculate_exx(qs_env=qs_env, &
unit_nr=iounit, &
hfx_sections=ec_hfx_sections, &
x_data=ec_env%x_data, &
do_gw=.FALSE., &
do_admm=ec_env%do_ec_admm, &
calc_forces=.TRUE., &
reuse_hfx=ec_env%reuse_hfx, &
do_im_time=.FALSE., &
E_ex_from_GW=dummy_real, &
E_admm_from_GW=dummy_real2, &
t3=dummy_real)
IF (debug_forces) THEN
fodeb(1:3) = force(1)%fock_4c(1:3, 1) - fodeb(1:3)
CALL para_env%sum(fodeb)
IF (iounit > 0) WRITE (iounit, "(T3,A,T33,3F16.8)") "DEBUG:: P*hfx_DC ", fodeb
fodeb2(1:3) = force(1)%overlap_admm(1:3, 1) - fodeb2(1:3)
CALL para_env%sum(fodeb2)
IF (iounit > 0) WRITE (iounit, "(T3,A,T33,3F16.8)") "DEBUG:: P*hfx_DC*S ", fodeb2
END IF
IF (debug_stress .AND. use_virial) THEN
stdeb = -1.0_dp*fconv*virial%pv_fock_4c
CALL para_env%sum(stdeb)
IF (iounit > 0) WRITE (UNIT=iounit, FMT="(T2,A,T41,2(1X,ES19.11))") &
'STRESS| P*hfx_DC ', one_third_sum_diag(stdeb), det_3x3(stdeb)
END IF
END IF
!------------------------------------------------
! Stress-tensor contribution derivative of integrand
! int v_Hxc[n^în]*n^out
IF (use_virial) THEN
pv_loc = virial%pv_virial
END IF
IF (debug_forces) fodeb(1:3) = force(1)%rho_elec(1:3, 1)
IF (debug_stress .AND. use_virial) stdeb = virial%pv_virial
DO ispin = 1, nspins
! Add v_H[n_in] + v_xc[n_in] = v_rspace
CALL pw_scale(v_rspace(ispin), v_rspace(ispin)%pw_grid%dvol)
CALL pw_axpy(v_hartree_rspace, v_rspace(ispin))
! integrate over potential <a|V|b>
CALL integrate_v_rspace(v_rspace=v_rspace(ispin), &
hmat=scrm(ispin), &
pmat=matrix_p(ispin, 1), &
qs_env=qs_env, &
calculate_forces=.TRUE., &
basis_type="HARRIS", &
task_list_external=ec_env%task_list)
END DO
IF (debug_forces) THEN
fodeb(1:3) = force(1)%rho_elec(1:3, 1) - fodeb(1:3)
CALL para_env%sum(fodeb)
IF (iounit > 0) WRITE (iounit, "(T3,A,T33,3F16.8)") "DEBUG:: Pout*dVhxc ", fodeb
END IF
IF (debug_stress .AND. use_virial) THEN
stdeb = fconv*(virial%pv_virial - stdeb)
CALL para_env%sum(stdeb)
IF (iounit > 0) WRITE (UNIT=iounit, FMT="(T2,A,T41,2(1X,ES19.11))") &
'STRESS| INT Pout*dVhxc ', one_third_sum_diag(stdeb), det_3x3(stdeb)
END IF
IF (ASSOCIATED(v_tau_rspace)) THEN
IF (debug_forces) fodeb(1:3) = force(1)%rho_elec(1:3, 1)
IF (debug_stress .AND. use_virial) stdeb = virial%pv_virial
DO ispin = 1, nspins
CALL pw_scale(v_tau_rspace(ispin), v_tau_rspace(ispin)%pw_grid%dvol)
! integrate over Tau-potential <nabla.a|V|nabla.b>
CALL integrate_v_rspace(v_rspace=v_tau_rspace(ispin), &
hmat=scrm(ispin), &
pmat=matrix_p(ispin, 1), &
qs_env=qs_env, &
calculate_forces=.TRUE., &
compute_tau=.TRUE., &
basis_type="HARRIS", &
task_list_external=ec_env%task_list)
END DO
IF (debug_forces) THEN
fodeb(1:3) = force(1)%rho_elec(1:3, 1) - fodeb(1:3)
CALL para_env%sum(fodeb)
IF (iounit > 0) WRITE (iounit, "(T3,A,T33,3F16.8)") "DEBUG:: Pout*dVhxc_tau ", fodeb
END IF
IF (debug_stress .AND. use_virial) THEN
stdeb = fconv*(virial%pv_virial - stdeb)
CALL para_env%sum(stdeb)
IF (iounit > 0) WRITE (UNIT=iounit, FMT="(T2,A,T41,2(1X,ES19.11))") &
'STRESS| INT Pout*dVhxc_tau ', one_third_sum_diag(stdeb), det_3x3(stdeb)
END IF
END IF
! Stress-tensor
IF (use_virial) THEN
virial%pv_ehartree = virial%pv_ehartree + (virial%pv_virial - pv_loc)
END IF
! delete scrm matrix
CALL dbcsr_deallocate_matrix_set(scrm)
!----------------------------------------------------
! Right-hand-side matrix B for linear response equations AX = B
!----------------------------------------------------
! RHS = int v_Hxc[n]_DC - v_Hxc[n]_GS dr + alpha_DC * E_X[n] - alpha_gs * E_X[n]
! = int v_Hxc[n]_DC - v_Hxc[n]_GS dr + alpha_DC / alpha_GS * E_X[n]_GS - E_X[n]_GS
!
! with v_Hxc[n] = v_H[n] + v_xc[n]
!
! Actually v_H[n_in] same for DC and GS, just there for convenience
! v_xc[n_in]_GS = 0 if GS is HF BUT =/0 if hybrid
! so, we keep this general form
NULLIFY (ec_env%matrix_hz)
CALL dbcsr_allocate_matrix_set(ec_env%matrix_hz, nspins)
DO ispin = 1, nspins
ALLOCATE (ec_env%matrix_hz(ispin)%matrix)
CALL dbcsr_create(ec_env%matrix_hz(ispin)%matrix, template=matrix_s(1)%matrix)
CALL dbcsr_copy(ec_env%matrix_hz(ispin)%matrix, matrix_s(1)%matrix)
CALL dbcsr_set(ec_env%matrix_hz(ispin)%matrix, 0.0_dp)
END DO
DO ispin = 1, nspins
! v_rspace = v_rspace - v_rspace_in
! = v_Hxc[n_in]_DC - v_Hxc[n_in]_GS
CALL pw_axpy(v_rspace_in(ispin), v_rspace(ispin), -1.0_dp)
END DO
DO ispin = 1, nspins
CALL integrate_v_rspace(v_rspace=v_rspace(ispin), &
hmat=ec_env%matrix_hz(ispin), &
pmat=matrix_p(ispin, 1), &
qs_env=qs_env, &
calculate_forces=.FALSE., &
basis_type="HARRIS", &
task_list_external=ec_env%task_list)
END DO
! Check if mGGA functionals are used
IF (dft_control%use_kinetic_energy_density) THEN
! If DC-DFT without mGGA functional, this needs to be allocated now.
IF (.NOT. ASSOCIATED(v_tau_rspace)) THEN
ALLOCATE (v_tau_rspace(nspins))
DO ispin = 1, nspins
CALL auxbas_pw_pool%create_pw(v_tau_rspace(ispin))
CALL pw_zero(v_tau_rspace(ispin))
END DO
END IF
DO ispin = 1, nspins
! v_tau_rspace = v_Hxc_tau[n_in]_DC - v_Hxc_tau[n_in]_GS
IF (ASSOCIATED(ec_env%vtau_rspace)) THEN
CALL pw_axpy(ec_env%vtau_rspace(ispin), v_tau_rspace(ispin), -1.0_dp)
END IF
! integrate over Tau-potential <nabla.a|V|nabla.b>
CALL integrate_v_rspace(v_rspace=v_tau_rspace(ispin), &
hmat=ec_env%matrix_hz(ispin), &
pmat=matrix_p(ispin, 1), &
qs_env=qs_env, &
calculate_forces=.FALSE., compute_tau=.TRUE., &
basis_type="HARRIS", &
task_list_external=ec_env%task_list)
END DO
END IF