-
Notifications
You must be signed in to change notification settings - Fork 1
/
input_cp2k_hfx.F
759 lines (658 loc) · 39.7 KB
/
input_cp2k_hfx.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
!--------------------------------------------------------------------------------------------------!
! 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 function that builds the hartree fock exchange section of the input
!> \par History
!> 09.2007 created
!> \author Manuel Guidon
! **************************************************************************************************
MODULE input_cp2k_hfx
USE bibliography, ONLY: Guidon2008,&
Guidon2009
USE cp_output_handling, ONLY: add_last_numeric,&
cp_print_key_section_create,&
high_print_level,&
medium_print_level
USE input_constants, ONLY: &
do_potential_coulomb, do_potential_gaussian, do_potential_id, do_potential_long, &
do_potential_mix_cl, do_potential_mix_cl_trunc, do_potential_mix_lg, do_potential_short, &
do_potential_truncated, ehrenfest, gaussian, hfx_ri_do_2c_cholesky, hfx_ri_do_2c_diag, &
hfx_ri_do_2c_iter
USE input_keyword_types, ONLY: keyword_create,&
keyword_release,&
keyword_type
USE input_section_types, ONLY: section_add_keyword,&
section_add_subsection,&
section_create,&
section_release,&
section_type
USE input_val_types, ONLY: real_t
USE kinds, ONLY: dp
USE string_utilities, ONLY: s2a
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .FALSE.
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_hfx'
INTEGER, PARAMETER, PUBLIC :: ri_mo = 1, ri_pmat = 2
PUBLIC :: create_hfx_section
CONTAINS
! **************************************************************************************************
!> \brief creates the input section for the hf part
!> \param section the section to create
!> \author Manuel Guidon
! **************************************************************************************************
SUBROUTINE create_hfx_section(section)
TYPE(section_type), POINTER :: section
TYPE(keyword_type), POINTER :: keyword
TYPE(section_type), POINTER :: print_key, subsection
CPASSERT(.NOT. ASSOCIATED(section))
CALL section_create(section, __LOCATION__, name="HF", &
description="Sets up the Hartree-Fock parameters if requested ", &
n_keywords=5, n_subsections=2, repeats=.TRUE., &
citations=(/Guidon2008, Guidon2009/))
NULLIFY (keyword, print_key, subsection)
CALL keyword_create(keyword, __LOCATION__, name="FRACTION", &
description="The fraction of Hartree-Fock to add to the total energy. "// &
"1.0 implies standard Hartree-Fock if used with XC_FUNCTIONAL NONE. "// &
"NOTE: In a mixed potential calculation this should be set to 1.0, otherwise "// &
"all parts are multiplied with this factor. ", &
usage="FRACTION 1.0", default_r_val=1.0_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="TREAT_LSD_IN_CORE", &
description="Determines how spin densities are taken into account. "// &
"If true, the beta spin density is included via a second in core call. "// &
"If false, alpha and beta spins are done in one shot ", &
usage="TREAT_LSD_IN_CORE TRUE", default_l_val=.FALSE.)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="PW_HFX", &
description="Compute the Hartree-Fock energy also in the plane wave basis. "// &
"The value is ignored, and intended for debugging only.", &
usage="PW_HFX FALSE", default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="PW_HFX_BLOCKSIZE", &
description="Improve the performance of pw_hfx at the cost of some additional memory "// &
"by storing the realspace representation of PW_HFX_BLOCKSIZE states.", &
usage="PW_HFX_BLOCKSIZE 20", default_i_val=20)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
NULLIFY (print_key)
CALL cp_print_key_section_create(print_key, __LOCATION__, "HF_INFO", &
description="Controls the printing basic info about hf method", &
print_level=medium_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
CALL section_add_subsection(section, print_key)
CALL section_release(print_key)
CALL create_hf_pbc_section(subsection)
CALL section_add_subsection(section, subsection)
CALL section_release(subsection)
CALL create_hf_screening_section(subsection)
CALL section_add_subsection(section, subsection)
CALL section_release(subsection)
CALL create_hf_potential_section(subsection)
CALL section_add_subsection(section, subsection)
CALL section_release(subsection)
CALL create_hf_load_balance_section(subsection)
CALL section_add_subsection(section, subsection)
CALL section_release(subsection)
CALL create_hf_memory_section(subsection)
CALL section_add_subsection(section, subsection)
CALL section_release(subsection)
CALL create_hf_ri_section(subsection)
CALL section_add_subsection(section, subsection)
CALL section_release(subsection)
END SUBROUTINE create_hfx_section
! **************************************************************************************************
!> \brief !****f* input_cp2k_dft/create_hf_load_balance_section [1.0] *
!>
!> creates the input section for the hf potential part
!> \param section the section to create
!> \author Manuel Guidon
! **************************************************************************************************
SUBROUTINE create_hf_load_balance_section(section)
TYPE(section_type), POINTER :: section
TYPE(keyword_type), POINTER :: keyword
TYPE(section_type), POINTER :: print_key
CPASSERT(.NOT. ASSOCIATED(section))
CALL section_create(section, __LOCATION__, name="LOAD_BALANCE", &
description="Parameters influencing the load balancing of the HF", &
n_keywords=1, n_subsections=0, repeats=.FALSE., &
citations=(/guidon2008/))
NULLIFY (keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="NBINS", &
description="Number of bins per process used to group atom quartets.", &
usage="NBINS 32", &
default_i_val=64)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="BLOCK_SIZE", &
description="Determines the blocking used for the atomic quartet loops. "// &
"A proper choice can speedup the calculation. The default (-1) is automatic.", &
usage="BLOCK_SIZE 4", &
default_i_val=-1)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
NULLIFY (keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="RANDOMIZE", &
description="This flag controls the randomization of the bin assignment to processes. "// &
"For highly ordered input structures with a bad load balance, setting "// &
"this flag to TRUE might improve.", &
usage="RANDOMIZE TRUE", &
default_l_val=.FALSE.)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
NULLIFY (print_key)
CALL cp_print_key_section_create(print_key, __LOCATION__, "PRINT", &
description="Controls the printing of info about load balance", &
print_level=medium_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
CALL section_add_subsection(section, print_key)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, &
name="LOAD_BALANCE_INFO", &
description="Activates the printing of load balance information ", &
default_l_val=.FALSE., &
lone_keyword_l_val=.TRUE.)
CALL section_add_keyword(print_key, keyword)
CALL keyword_release(keyword)
CALL section_release(print_key)
END SUBROUTINE create_hf_load_balance_section
! **************************************************************************************************
!> \brief !****f* input_cp2k_dft/create_hf_potential_section [1.0] *
!>
!> creates the input section for the hf potential part
!> \param section the section to create
!> \author Manuel Guidon
! **************************************************************************************************
SUBROUTINE create_hf_potential_section(section)
TYPE(section_type), POINTER :: section
TYPE(keyword_type), POINTER :: keyword
CPASSERT(.NOT. ASSOCIATED(section))
CALL section_create(section, __LOCATION__, name="INTERACTION_POTENTIAL", &
description="Sets up interaction potential if requested ", &
n_keywords=1, n_subsections=0, repeats=.FALSE., &
citations=(/guidon2008, guidon2009/))
NULLIFY (keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="POTENTIAL_TYPE", &
description="Which interaction potential should be used "// &
"(Coulomb, longrange or shortrange).", &
usage="POTENTIAL_TYPE SHORTRANGE", &
enum_c_vals=s2a("COULOMB", "SHORTRANGE", "LONGRANGE", "MIX_CL", "GAUSSIAN", &
"MIX_LG", "IDENTITY", "TRUNCATED", "MIX_CL_TRUNC"), &
enum_i_vals=(/do_potential_coulomb, do_potential_short, do_potential_long, &
do_potential_mix_cl, do_potential_gaussian, do_potential_mix_lg, &
do_potential_id, do_potential_truncated, do_potential_mix_cl_trunc/), &
enum_desc=s2a("Coulomb potential: 1/r", &
"Shortrange potential: erfc(omega*r)/r", &
"Longrange potential: erf(omega*r)/r", &
"Mix coulomb and longrange potential: 1/r + erf(omega*r)/r", &
"Damped Gaussian potential: exp(-omega^2*r^2)", &
"Mix Gaussian and longrange potential: erf(omega*r)/r + exp(-omega^2*r^2)", &
"Overlap", &
"Truncated coulomb potential: if (r < R_c) 1/r else 0", &
"Truncated Mix coulomb and longrange potential, assumes/requires that the erf has fully decayed at R_c"), &
default_i_val=do_potential_coulomb)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
NULLIFY (keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="OMEGA", &
description="Parameter for short/longrange interaction", &
usage="OMEGA 0.5", &
default_r_val=0.0_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="SCALE_COULOMB", &
description="Scales Hartree-Fock contribution arising from a coulomb potential. "// &
"Only valid when doing a mixed potential calculation", &
usage="SCALE_COULOMB 1.0", default_r_val=1.0_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="SCALE_LONGRANGE", &
description="Scales Hartree-Fock contribution arising from a longrange potential. "// &
"Only valid when doing a mixed potential calculation", &
usage="SCALE_LONGRANGE 1.0", default_r_val=1.0_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="SCALE_GAUSSIAN", &
description="Scales Hartree-Fock contribution arising from a gaussian potential. "// &
"Only valid when doing a mixed potential calculation", &
usage="SCALE_GAUSSIAN 1.0", default_r_val=1.0_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="CUTOFF_RADIUS", &
description="Determines cutoff radius (in Angstroms) for the truncated 1/r potential. "// &
"Only valid when doing truncated calculation", &
usage="CUTOFF_RADIUS 10.0", type_of_var=real_t, & ! default_r_val=10.0_dp,&
unit_str="angstrom")
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="T_C_G_DATA", &
description="Location of the file t_c_g.dat that contains the data for the "// &
"evaluation of the truncated gamma function ", &
usage="T_C_G_DATA /data/t_c_g.dat", &
default_c_val="t_c_g.dat")
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
END SUBROUTINE create_hf_potential_section
!****f* input_cp2k_dft/create_hf_screening_section [1.0] *
! **************************************************************************************************
!> \brief creates the input section for the hf screening part
!> \param section the section to create
!> \author Manuel Guidon
! **************************************************************************************************
SUBROUTINE create_hf_screening_section(section)
TYPE(section_type), POINTER :: section
TYPE(keyword_type), POINTER :: keyword
CPASSERT(.NOT. ASSOCIATED(section))
CALL section_create(section, __LOCATION__, name="SCREENING", &
description="Sets up screening parameters if requested ", &
n_keywords=1, n_subsections=0, repeats=.FALSE., &
citations=(/guidon2008, guidon2009/))
NULLIFY (keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="EPS_SCHWARZ", &
description="Screens the near field part of the electronic repulsion "// &
"integrals using the Schwarz inequality for the given "// &
"threshold.", &
usage="EPS_SCHWARZ 1.0E-6", &
default_r_val=1.0E-10_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
NULLIFY (keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="EPS_SCHWARZ_FORCES", &
description="Screens the near field part of the electronic repulsion "// &
"integrals using the Schwarz inequality for the given "// &
"threshold. This will be approximately the accuracy of the forces, "// &
"and should normally be similar to EPS_SCF. Default value is 100*EPS_SCHWARZ.", &
usage="EPS_SCHWARZ_FORCES 1.0E-5", &
default_r_val=1.0E-6_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
NULLIFY (keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="SCREEN_P_FORCES", &
description="Screens the electronic repulsion integrals for the forces "// &
"using the density matrix. Will be disabled for the "// &
"response part of forces in MP2/RPA/TDDFT. "// &
"This results in a significant speedup for large systems, "// &
"but might require a somewhat tigher EPS_SCHWARZ_FORCES.", &
usage="SCREEN_P_FORCES TRUE", &
default_l_val=.TRUE.)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
NULLIFY (keyword)
CALL keyword_create(keyword, __LOCATION__, name="SCREEN_ON_INITIAL_P", &
description="Screen on an initial density matrix. For the first MD step"// &
" this matrix must be provided by a Restart File.", &
usage="SCREEN_ON_INITIAL_P TRUE", default_l_val=.FALSE.)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
NULLIFY (keyword)
CALL keyword_create(keyword, __LOCATION__, name="P_SCREEN_CORRECTION_FACTOR", &
description="Recalculates integrals on the fly if the actual density matrix is"// &
" larger by a given factor than the initial one. If the factor is set"// &
" to 0.0_dp, this feature is disabled.", &
usage="P_SCREEN_CORRECTION_FACTOR 0.0_dp", default_r_val=0.0_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
END SUBROUTINE create_hf_screening_section
! **************************************************************************************************
!> \brief creates the input section for the hf-pbc part
!> \param section the section to create
!> \author Manuel Guidon
! **************************************************************************************************
SUBROUTINE create_hf_pbc_section(section)
TYPE(section_type), POINTER :: section
TYPE(keyword_type), POINTER :: keyword
CPASSERT(.NOT. ASSOCIATED(section))
CALL section_create(section, __LOCATION__, name="PERIODIC", &
description="Sets up periodic boundary condition parameters if requested ", &
n_keywords=1, n_subsections=0, repeats=.FALSE., &
citations=(/guidon2008, guidon2009/))
NULLIFY (keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="NUMBER_OF_SHELLS", &
description="Number of shells taken into account for periodicity. "// &
"By default, cp2k tries to automatically evaluate this number. "// &
"This algorithm might be to conservative, resulting in some overhead. "// &
"You can try to adjust this number in order to make a calculation cheaper. ", &
usage="NUMBER_OF_SHELLS 2", &
default_i_val=-1)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
END SUBROUTINE create_hf_pbc_section
! **************************************************************************************************
!> \brief creates the input section for the hf-memory part
!> \param section the section to create
!> \author Manuel Guidon
! **************************************************************************************************
SUBROUTINE create_hf_memory_section(section)
TYPE(section_type), POINTER :: section
TYPE(keyword_type), POINTER :: keyword
CPASSERT(.NOT. ASSOCIATED(section))
CALL section_create(section, __LOCATION__, name="MEMORY", &
description="Sets up memory parameters for the storage of the ERI's if requested ", &
n_keywords=1, n_subsections=0, repeats=.FALSE., &
citations=(/guidon2008/))
NULLIFY (keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="EPS_STORAGE_SCALING", &
variants=(/"EPS_STORAGE"/), &
description="Scaling factor to scale eps_schwarz. Storage threshold for compression "// &
"will be EPS_SCHWARZ*EPS_STORAGE_SCALING.", &
usage="EPS_STORAGE 1.0E-2", &
default_r_val=1.0E0_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="MAX_MEMORY", &
description="Defines the maximum amount of memory [MiB] to be consumed by the full HFX module. "// &
"All temporary buffers and helper arrays are subtracted from this number. "// &
"What remains will be used for storage of integrals. NOTE: This number "// &
"is assumed to represent the memory available to one MPI process. "// &
"When running a threaded version, cp2k automatically takes care of "// &
"distributing the memory among all the threads within a process.", &
usage="MAX_MEMORY 256", &
default_i_val=512)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="STORAGE_LOCATION", &
description="Loaction where ERI's are stored if MAX_DISK_SPACE /=0 "// &
"Expects a path to a directory. ", &
usage="STORAGE_LOCATION /data/scratch", &
default_c_val=".")
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="MAX_DISK_SPACE", &
description="Defines the maximum amount of disk space [MiB] used to store precomputed "// &
"compressed four-center integrals. If 0, nothing is stored to disk", &
usage="MAX_DISK_SPACE 256", &
default_i_val=0)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="TREAT_FORCES_IN_CORE", &
description="Determines whether the derivative ERI's should be stored to RAM or not. "// &
"Only meaningful when performing Ehrenfest MD. "// &
"Memory usage is defined via MAX_MEMORY, i.e. the memory is shared wit the energy ERI's.", &
usage="TREAT_FORCES_IN_CORE TRUE", default_l_val=.FALSE.)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
END SUBROUTINE create_hf_memory_section
! **************************************************************************************************
!> \brief ...
!> \param section ...
! **************************************************************************************************
SUBROUTINE create_hf_ri_section(section)
TYPE(section_type), POINTER :: section
TYPE(keyword_type), POINTER :: keyword
TYPE(section_type), POINTER :: print_key, subsection
NULLIFY (keyword, print_key, subsection)
CPASSERT(.NOT. ASSOCIATED(section))
CALL section_create(section, __LOCATION__, name="RI", &
description="Parameters for RI methods in HFX, including RI-HFXk with "// &
"k-point sampling. All keywords relevant to RI-HFXk have an "// &
"alias starting with KP_")
CALL keyword_create(keyword, __LOCATION__, name="_SECTION_PARAMETERS_", &
description="controls the activation of RI", &
usage="&RI T", &
default_l_val=.FALSE., &
lone_keyword_l_val=.TRUE.)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="EPS_FILTER", &
description="Filter threshold for DBT tensor contraction.", &
variants=(/"KP_EPS_FILTER"/), &
default_r_val=1.0E-09_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="EPS_FILTER_2C", &
description="Filter threshold for 2c integrals. Default should be kept.", &
default_r_val=1.0E-12_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, &
name="EPS_STORAGE_SCALING", &
description="Scaling factor to scale EPS_FILTER for storage of 3-center integrals. Storage threshold "// &
"will be EPS_FILTER*EPS_STORAGE_SCALING.", &
default_r_val=0.01_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="EPS_FILTER_MO", &
description="Filter threshold for contraction of 3-center integrals with MOs. "// &
"Default should be kept.", &
default_r_val=1.0E-12_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="OMEGA", &
description="The range parameter for the short range operator (in 1/a0). "// &
"Default is OMEGA from INTERACTION_POTENTIAL. ", &
variants=(/"KP_OMEGA"/), &
default_r_val=0.0_dp, &
repeats=.FALSE.)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="CUTOFF_RADIUS", &
description="The cutoff radius (in Angstroms) for the truncated Coulomb operator. "// &
"Default is CUTOFF_RADIUS from INTERACTION_POTENTIAL. ", &
variants=(/"KP_CUTOFF_RADIUS"/), &
default_r_val=0.0_dp, &
repeats=.FALSE., &
unit_str="angstrom")
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="KP_NGROUPS", &
description="The number of MPI subgroup that work in parallel during the SCF. "// &
"The default value is 1. Using N subgroups should speed up the "// &
"calculation by a factor ~N, at the cost of N times more memory usage.", &
variants=(/"NGROUPS"/), &
usage="KP_NGROUPS {int}", &
repeats=.FALSE., &
default_i_val=1)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="KP_USE_DELTA_P", &
description="This kweyword controls whether the KS matrix at each SCF cycle "// &
"is built by adding the contribution of the denisty difference (wrt to previous step) "// &
"to the KS matrix of the previous step. As the SCF converges, the density fluctuations "// &
"get smaller and sparsity increases, leading to faster SCF steps. Not always "// &
"numerically stable => turn off if SCF struggles to converge.", &
variants=s2a("USE_DELTA_P", "KP_USE_P_DIFF", "USE_P_DIFF"), &
usage="KP_USE_DELTA_P {logical}", &
repeats=.FALSE., &
default_l_val=.TRUE.)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="KP_STACK_SIZE", &
description="When doing contraction over periodic cells of the type: "// &
"T_mu^a,nu^b,P^c = (mu^a nu^b | Q^d) * (Q^d | P^c), with "// &
"a,b,c,d labeling cells, there are in principle Ncells "// &
"contractions taking place. Because a smaller number of "// &
"contractions involving larger tensors is more efficient, "// &
"the tensors can be stacked along the d direction. STCK_SIZE "// &
"controls the size of this stack. Larger stacks are more efficient, "// &
"but required more memory.", &
variants=(/"STACK_SIZE"/), &
usage="KP_STACK_SIZE {int}", &
repeats=.FALSE., &
default_i_val=32)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="KP_RI_BUMP_FACTOR", &
variants=s2a("RI_BUMP", "BUMP", "BUMP_FACTOR"), &
description="In KP-RI-HFX, the extended RI basis set has a bump radius. "// &
"All basis elements within that radius contribute with full weight. "// &
"All basis elements beyond that radius have decaying weight, from "// &
"1 at the bump radius, to zero at the RI extension radius. The "// &
"bump radius is calculated as a fraction of the RI extension radius: "// &
"bump radius = KP_RI_NUMP_FACTOR * RI extension radius", &
default_r_val=0.85_dp, &
repeats=.FALSE.)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="RI_METRIC", &
description="The type of RI operator. "// &
"Default is POTENTIAL_TYPE from INTERACTION_POTENTIAL. "// &
"The standard "// &
"Coulomb operator cannot be used in periodic systems.", &
usage="OPERATOR {string}", &
repeats=.FALSE., &
variants=(/"KP_RI_METRIC"/), &
default_i_val=0, &
enum_c_vals=s2a("HFX", "COULOMB", "IDENTITY", "TRUNCATED", "SHORTRANGE"), &
enum_desc=s2a("Same as HFX operator", &
"Standard Coulomb operator: 1/r", &
"Overlap", &
"Truncated Coulomb operator: 1/r if (r<R_c), 0 otherwise ", &
"Short range: erfc(omega*r)/r"), &
enum_i_vals=(/0, do_potential_coulomb, do_potential_id, do_potential_truncated, &
do_potential_short/))
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="2C_MATRIX_FUNCTIONS", &
description="Methods for matrix inverse and matrix square root.", &
default_i_val=hfx_ri_do_2c_cholesky, &
enum_c_vals=s2a("DIAG", "CHOLESKY", "ITER"), &
enum_desc=s2a("Diagonalization with eigenvalue quenching: stable", &
"Cholesky: not stable in case of ill-conditioned RI basis", &
"Iterative algorithms: linear scaling "// &
"Hotelling's method for inverse and Newton-Schulz iteration for matrix square root"), &
enum_i_vals=(/hfx_ri_do_2c_diag, hfx_ri_do_2c_cholesky, hfx_ri_do_2c_iter/))
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="EPS_EIGVAL", &
description="Throw away linear combinations of RI basis functions with a small eigenvalue, "// &
"this is applied only if 2C_MATRIX_FUNCTIONS DIAG", &
default_r_val=1.0e-7_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="CHECK_2C_MATRIX", &
description="Report accuracy for the inverse/sqrt of the 2-center integral matrix.", &
default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create( &
keyword, __LOCATION__, &
name="CALC_COND_NUM", &
variants=(/"CALC_CONDITION_NUMBER"/), &
description="Calculate the condition number of integral matrices.", &
usage="CALC_COND_NUM", &
default_l_val=.FALSE., &
lone_keyword_l_val=.TRUE.)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="SQRT_ORDER", &
description="Order of the iteration method for the calculation of "// &
"the sqrt of 2-center integral matrix.", &
default_i_val=3)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="EPS_LANCZOS", &
description="Threshold used for lanczos estimates.", &
default_r_val=1.0E-3_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="EPS_PGF_ORB", &
description="Sets precision of the integral tensors.", &
variants=(/"KP_EPS_PGF_ORB"/), &
default_r_val=1.0E-5_dp)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="MAX_ITER_LANCZOS", &
description="Maximum number of lanczos iterations.", &
usage="MAX_ITER_LANCZOS ", default_i_val=500)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="RI_FLAVOR", &
description="Flavor of RI: how to contract 3-center integrals", &
enum_c_vals=s2a("MO", "RHO"), &
enum_desc=s2a("with MO coefficients", "with density matrix"), &
enum_i_vals=(/ri_mo, ri_pmat/), &
default_i_val=ri_pmat)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="MIN_BLOCK_SIZE", &
description="Minimum tensor block size.", &
default_i_val=4)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="MAX_BLOCK_SIZE_MO", &
description="Maximum tensor block size for MOs.", &
default_i_val=64)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="MEMORY_CUT", &
description="Memory reduction factor. This keyword controls the batching of tensor "// &
"contractions into smaller, more manageable chunks. The details vary "// &
"depending on the RI_FLAVOR.", &
default_i_val=3)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL keyword_create(keyword, __LOCATION__, name="FLAVOR_SWITCH_MEMORY_CUT", &
description="Memory reduction factor to be applied upon RI_FLAVOR switching "// &
"from MO to RHO. The RHO flavor typically requires more memory, "// &
"and depending on the ressources available, a higher MEMORY_CUT.", &
default_i_val=3)
CALL section_add_keyword(section, keyword)
CALL keyword_release(keyword)
CALL section_create(subsection, __LOCATION__, name="PRINT", &
description="Section of possible print options in the RI-HFX code.", &
n_keywords=0, n_subsections=1, repeats=.FALSE.)
CALL cp_print_key_section_create(print_key, __LOCATION__, "RI_INFO", &
description="Controls the printing of DBCSR tensor log in RI HFX.", &
print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
CALL section_add_subsection(subsection, print_key)
CALL section_release(print_key)
CALL cp_print_key_section_create(print_key, __LOCATION__, "RI_DENSITY_COEFFS", &
description="Controls the printing of the projection of the elecontric "// &
"density on the RI_HFX basis. n(r) = sum_s c_s*phi_RI_s(r), "// &
"where c_s = sum_pqr P_pq (pq|r) (r|s)^-1 and | is the RI_METRIC", &
print_level=high_print_level, filename="RI_DENSITY_COEFFS", &
common_iter_levels=3)
CALL keyword_create(keyword, __LOCATION__, name="MULTIPLY_BY_RI_2C_INTEGRALS", &
variants=s2a("MULT_BY_RI", "MULT_BY_S", "MULT_BY_RI_INTS"), &
description="Whether the RI density coefficients to be printed should "// &
"be pre-multiplied by the RI_METRIC 2c-integrals: (r|s)*C_s.", &
default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
CALL section_add_keyword(print_key, keyword)
CALL keyword_release(keyword)
CALL section_add_subsection(subsection, print_key)
CALL section_release(print_key)
CALL cp_print_key_section_create(print_key, __LOCATION__, "RI_METRIC_2C_INTS", &
description="Controls the printing of RI 2-center integrals for the "// &
"HFX potential.", &
print_level=high_print_level, filename="RI_2C_INTS", &
common_iter_levels=3)
CALL section_add_subsection(subsection, print_key)
CALL section_release(print_key)
CALL section_add_subsection(section, subsection)
CALL section_release(subsection)
END SUBROUTINE
END MODULE input_cp2k_hfx