(git:e7e05ae)
input_cp2k_mp2.F
Go to the documentation of this file.
1 !--------------------------------------------------------------------------------------------------!
2 ! CP2K: A general program to perform molecular dynamics simulations !
3 ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 ! !
5 ! SPDX-License-Identifier: GPL-2.0-or-later !
6 !--------------------------------------------------------------------------------------------------!
7 
8 ! **************************************************************************************************
9 !> \brief input section for MP2
10 !> \par History
11 !> 05.2011 created
12 !> \author MDB
13 ! **************************************************************************************************
15  USE bibliography, ONLY: &
26  USE cp_units, ONLY: cp_unit_to_cp2k
27  USE input_constants, ONLY: &
45  keyword_type
50  section_type
51  USE input_val_types, ONLY: integer_t, &
52  real_t
53  USE kinds, ONLY: dp
54  USE string_utilities, ONLY: newline, &
55  s2a
56 #include "./base/base_uses.f90"
57 
58  IMPLICIT NONE
59  PRIVATE
60 
61  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_mp2'
62 
63  PUBLIC :: create_mp2_section
64 
65 CONTAINS
66 
67 ! **************************************************************************************************
68 !> \brief creates the input section for the mp2 part
69 !> \param section the section to create
70 !> \author MDB
71 ! **************************************************************************************************
72  SUBROUTINE create_mp2_section(section)
73  TYPE(section_type), POINTER :: section
74 
75  TYPE(keyword_type), POINTER :: keyword
76  TYPE(section_type), POINTER :: print_key, subsection
77 
78  cpassert(.NOT. ASSOCIATED(section))
79  CALL section_create(section, __location__, name="WF_CORRELATION", &
80  description="Sets up the wavefunction-based correlation methods as MP2, "// &
81  "RI-MP2, RI-SOS-MP2, RI-RPA and GW (inside RI-RPA). ", &
82  n_keywords=4, n_subsections=7, repeats=.true., &
85 
86  NULLIFY (keyword, subsection)
87 
88  CALL keyword_create( &
89  keyword, __location__, &
90  name="MEMORY", &
91  description="Maximum allowed total memory usage during MP2 methods [MiB].", &
92  usage="MEMORY 1500 ", &
93  default_r_val=1.024e+3_dp)
94  CALL section_add_keyword(section, keyword)
95  CALL keyword_release(keyword)
96 
97  CALL keyword_create( &
98  keyword, __location__, &
99  name="E_GAP", &
100  description="Gap energy for integration grids in Hartree. Defaults to -1.0 (automatic determination). "// &
101  "Recommended to set if several RPA or SOS-MP2 gradient calculations are requested or to be restarted. "// &
102  "In this way, differences of integration grids across different runs are removed as CP2K "// &
103  "does not include derivatives thereof.", &
104  usage="E_GAP 0.5", &
105  default_r_val=-1.0_dp)
106  CALL section_add_keyword(section, keyword)
107  CALL keyword_release(keyword)
108 
109  CALL keyword_create( &
110  keyword, __location__, &
111  name="E_RANGE", &
112  description="Energy range (ratio of largest and smallest) energy difference "// &
113  "of unoccupied and occupied orbitals for integration grids. Defaults to 0.0 (automatic determination). "// &
114  "Recommended to set if several RPA or SOS-MP2 gradient calculations are requested or to be restarted. "// &
115  "In this way, differences of integration grids across different runs are removed as CP2K "// &
116  "does not include derivatives thereof.", &
117  usage="E_RANGE 10.0", &
118  default_r_val=-1.0_dp)
119  CALL section_add_keyword(section, keyword)
120  CALL keyword_release(keyword)
121 
122  CALL keyword_create( &
123  keyword, __location__, &
124  name="SCALE_S", &
125  description="Scaling factor of the singlet energy component (opposite spin, OS) of the "// &
126  "MP2, RI-MP2 and SOS-MP2 correlation energy. ", &
127  usage="SCALE_S 1.0", &
128  default_r_val=1.0_dp)
129  CALL section_add_keyword(section, keyword)
130  CALL keyword_release(keyword)
131 
132  CALL keyword_create( &
133  keyword, __location__, &
134  name="SCALE_T", &
135  description="Scaling factor of the triplet energy component (same spin, SS) of the MP2 "// &
136  "and RI-MP2 correlation energy.", &
137  usage="SCALE_T 1.0", &
138  default_r_val=1.0_dp)
139  CALL section_add_keyword(section, keyword)
140  CALL keyword_release(keyword)
141 
142  CALL keyword_create( &
143  keyword, __location__, &
144  name="GROUP_SIZE", &
145  variants=(/"NUMBER_PROC"/), &
146  description="Group size used in the computation of GPW and MME integrals and the MP2 correlation energy. "// &
147  "The group size must be a divisor of the total number of MPI ranks. "// &
148  "A smaller group size (for example the number of MPI ranks per node) "// &
149  "accelerates the computation of integrals but a too large group size increases communication costs. "// &
150  "A too small group size may lead to out of memory.", &
151  usage="GROUP_SIZE 2", &
152  default_i_val=1)
153  CALL section_add_keyword(section, keyword)
154  CALL keyword_release(keyword)
155 
156  NULLIFY (subsection)
157  CALL create_mp2_details_section(subsection)
158  CALL section_add_subsection(section, subsection)
159  CALL section_release(subsection)
160 
161  CALL create_ri_mp2(subsection)
162  CALL section_add_subsection(section, subsection)
163  CALL section_release(subsection)
164 
165  CALL create_ri_rpa(subsection)
166  CALL section_add_subsection(section, subsection)
167  CALL section_release(subsection)
168 
169  CALL create_ri_laplace(subsection)
170  CALL section_add_subsection(section, subsection)
171  CALL section_release(subsection)
172 
173  ! here we generate an imag. time subsection to use with RPA or Laplace-SOS-MP2
174  CALL create_low_scaling(subsection)
175  CALL section_add_subsection(section, subsection)
176  CALL section_release(subsection)
177 
178  CALL create_ri_section(subsection)
179  CALL section_add_subsection(section, subsection)
180  CALL section_release(subsection)
181 
182  CALL create_integrals_section(subsection)
183  CALL section_add_subsection(section, subsection)
184  CALL section_release(subsection)
185 
186  CALL create_canonical_gradients(subsection)
187  CALL section_add_subsection(section, subsection)
188  CALL section_release(subsection)
189 
190  NULLIFY (print_key)
191  CALL cp_print_key_section_create(print_key, __location__, "PRINT", &
192  description="Controls the printing basic info about WFC methods", &
193  print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
194  CALL section_add_subsection(section, print_key)
195  CALL section_release(print_key)
196 
197  END SUBROUTINE create_mp2_section
198 
199 ! **************************************************************************************************
200 !> \brief ...
201 !> \param section ...
202 ! **************************************************************************************************
203  SUBROUTINE create_mp2_details_section(section)
204  TYPE(section_type), POINTER :: section
205 
206  TYPE(keyword_type), POINTER :: keyword
207 
208  cpassert(.NOT. ASSOCIATED(section))
209  CALL section_create(section, __location__, name="MP2", &
210  description="Parameters influencing MP2 (non-RI).", &
211  n_keywords=3, n_subsections=0, repeats=.false.)
212 
213  NULLIFY (keyword)
214  CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
215  description="Activates MP2 calculations.", &
216  usage="&MP2 .TRUE.", &
217  default_l_val=.false., lone_keyword_l_val=.true.)
218  CALL section_add_keyword(section, keyword)
219  CALL keyword_release(keyword)
220 
221  CALL keyword_create( &
222  keyword, __location__, &
223  name="METHOD", &
224  citations=(/delben2012, delben2013/), &
225  description="Method that is used to compute the MP2 energy.", &
226  usage="METHOD MP2_GPW", &
227  enum_c_vals=s2a("NONE", "DIRECT_CANONICAL", "MP2_GPW"), &
229  enum_desc=s2a("Skip MP2 calculation.", &
230  "Use the direct mp2 canonical approach.", &
231  "Use the GPW approach to MP2."), &
232  default_i_val=mp2_method_direct)
233  CALL section_add_keyword(section, keyword)
234  CALL keyword_release(keyword)
235 
236  CALL keyword_create( &
237  keyword, __location__, &
238  name="BIG_SEND", &
239  description="Influencing the direct canonical MP2 method: Send big "// &
240  "messages between processes (useful for >48 processors).", &
241  usage="BIG_SEND", &
242  default_l_val=.true., &
243  lone_keyword_l_val=.true.)
244  CALL section_add_keyword(section, keyword)
245  CALL keyword_release(keyword)
246 
247  END SUBROUTINE create_mp2_details_section
248 
249 ! **************************************************************************************************
250 !> \brief ...
251 !> \param section ...
252 ! **************************************************************************************************
253  SUBROUTINE create_ri_mp2(section)
254  TYPE(section_type), POINTER :: section
255 
256  TYPE(keyword_type), POINTER :: keyword
257 
258  cpassert(.NOT. ASSOCIATED(section))
259  CALL section_create(section, __location__, name="RI_MP2", &
260  description="Parameters influencing the RI-MP2 method. RI-MP2 supports gradients.", &
261  n_keywords=3, n_subsections=1, repeats=.false., &
262  citations=(/delben2013/))
263 
264  NULLIFY (keyword)
265  CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
266  description="Putting the &RI_MP2 section activates RI-MP2 calculation.", &
267  usage="&RI_MP2 .TRUE.", &
268  default_l_val=.false., lone_keyword_l_val=.true.)
269  CALL section_add_keyword(section, keyword)
270  CALL keyword_release(keyword)
271 
272  CALL keyword_create(keyword, __location__, name="BLOCK_SIZE", &
273  variants=(/"MESSAGE_SIZE"/), &
274  description="Determines the blocking used for communication in RI-MP2. Larger BLOCK_SIZE "// &
275  "reduces communication but requires more memory. The default (-1) is automatic.", &
276  usage="BLOCK_SIZE 2", &
277  default_i_val=-1)
278  CALL section_add_keyword(section, keyword)
279  CALL keyword_release(keyword)
280 
281  CALL keyword_create(keyword, __location__, name="NUMBER_INTEGRATION_GROUPS", &
282  description="Sets the number of integration groups of the communication scheme in RI-MP2. "// &
283  "Integrals will be replicated such that each integration group has all integrals available. "// &
284  "Must be a divisor of the number of subgroups (see GROUP_SIZE keyword in the WF_CORRELATION "// &
285  "section. Smaller groups reduce the communication costs but increase the memory developments. "// &
286  "If the provided value is non-positive or not a divisor of the number of subgroups, "// &
287  "the number of integration groups is determined automatically (default).", &
288  usage="NUMBER_INTEGRATION_GROUPS 2", &
289  default_i_val=-1)
290  CALL section_add_keyword(section, keyword)
291  CALL keyword_release(keyword)
292 
293  CALL keyword_create( &
294  keyword, __location__, &
295  name="PRINT_DGEMM_INFO", &
296  description="Print details about all DGEMM calls.", &
297  lone_keyword_l_val=.true., &
298  default_l_val=.false.)
299  CALL section_add_keyword(section, keyword)
300  CALL keyword_release(keyword)
301 
302  END SUBROUTINE create_ri_mp2
303 
304 ! **************************************************************************************************
305 !> \brief ...
306 !> \param section ...
307 ! **************************************************************************************************
308  SUBROUTINE create_opt_ri_basis(section)
309  TYPE(section_type), POINTER :: section
310 
311  TYPE(keyword_type), POINTER :: keyword
312 
313  cpassert(.NOT. ASSOCIATED(section))
314  CALL section_create(section, __location__, name="OPT_RI_BASIS", &
315  description="Parameters influencing the optimization of the RI MP2 basis. "// &
316  "Only exponents of non-contracted auxiliary basis can be optimized. "// &
317  "An initial RI auxiliary basis has to be specified.", &
318  n_keywords=6, n_subsections=0, repeats=.false., &
319  citations=(/delben2013/))
320  NULLIFY (keyword)
321  CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
322  description="Putting the &OPT_RI_BASIS section activates optimization of RI basis.", &
323  usage="&OPT_RI_BASIS .TRUE.", &
324  default_l_val=.false., lone_keyword_l_val=.true.)
325  CALL section_add_keyword(section, keyword)
326  CALL keyword_release(keyword)
327 
328  CALL keyword_create(keyword, __location__, name="DELTA_I_REL", &
329  variants=(/"DI_REL"/), &
330  description="Target accuracy in the relative deviation of the amplitudes calculated with "// &
331  "and without RI approximation, (more details in Chem.Phys.Lett.294(1998)143).", &
332  usage="DELTA_I_REL 1.0E-6_dp", &
333  default_r_val=1.0e-6_dp)
334  CALL section_add_keyword(section, keyword)
335  CALL keyword_release(keyword)
336 
337  CALL keyword_create(keyword, __location__, name="DELTA_RI", &
338  variants=(/"DRI"/), &
339  description="Target accuracy in the absolute difference between the RI-MP2 "// &
340  "and the exact MP2 energy, DRI=ABS(E_MP2-E_RI-MP2).", &
341  usage="DELTA_RI 1.0E-6_dp", &
342  default_r_val=5.0e-6_dp)
343  CALL section_add_keyword(section, keyword)
344  CALL keyword_release(keyword)
345 
346  CALL keyword_create(keyword, __location__, name="EPS_DERIV", &
347  variants=(/"EPS_NUM_DERIV"/), &
348  description="The derivatives of the MP2 energy with respect to the "// &
349  "exponents of the basis are calculated numerically. "// &
350  "The change in the exponent a_i employed for the numerical evaluation "// &
351  "is defined as h_i=EPS_DERIV*a_i.", &
352  usage="EPS_DERIV 1.0E-3_dp", &
353  default_r_val=1.0e-3_dp)
354  CALL section_add_keyword(section, keyword)
355  CALL keyword_release(keyword)
356 
357  CALL keyword_create(keyword, __location__, name="MAX_ITER", &
358  variants=(/"MAX_NUM_ITER"/), &
359  description="Specifies the maximum number of steps in the RI basis optimization.", &
360  usage="MAX_ITER 100", &
361  default_i_val=50)
362  CALL section_add_keyword(section, keyword)
363  CALL keyword_release(keyword)
364 
365  CALL keyword_create(keyword, __location__, name="NUM_FUNC", &
366  description="Specifies the number of function, for each angular momentum (s, p, d ...), "// &
367  "employed in the automatically generated initial guess. "// &
368  "This will be effective only if RI_AUX_BASIS_SET in the KIND section is not specified.", &
369  usage="NUM_FUNC {number of s func.} {number of p func.} ...", &
370  n_var=-1, default_i_vals=(/-1/), type_of_var=integer_t)
371  CALL section_add_keyword(section, keyword)
372  CALL keyword_release(keyword)
373 
374  CALL keyword_create(keyword, __location__, name="BASIS_SIZE", &
375  description="Specifies the size of the auxiliary basis set automatically "// &
376  "generated as initial guess. This will be effective only if RI_AUX_BASIS_SET "// &
377  "in the KIND section and NUM_FUNC are not specified.", &
378  usage="BASIS_SIZE (MEDIUM|LARGE|VERY_LARGE)", &
379  enum_c_vals=s2a("MEDIUM", "LARGE", "VERY_LARGE"), &
380  enum_i_vals=(/0, 1, 2/), &
381  default_i_val=0)
382  CALL section_add_keyword(section, keyword)
383  CALL keyword_release(keyword)
384 
385  END SUBROUTINE create_opt_ri_basis
386 
387 ! **************************************************************************************************
388 !> \brief ...
389 !> \param section ...
390 ! **************************************************************************************************
391  SUBROUTINE create_ri_laplace(section)
392  TYPE(section_type), POINTER :: section
393 
394  TYPE(keyword_type), POINTER :: keyword
395 
396  cpassert(.NOT. ASSOCIATED(section))
397  CALL section_create(section, __location__, name="RI_SOS_MP2", &
398  description="Parameters influencing the RI-SOS-MP2-Laplace method", &
399  n_keywords=3, n_subsections=1, repeats=.false., &
400  citations=(/delben2013/))
401 
402  NULLIFY (keyword)
403  CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
404  description="Putting the &RI_SOS_MP2 section activates RI-SOS-MP2 calculation.", &
405  usage="&RI_SOS_MP2 .TRUE.", &
406  default_l_val=.false., lone_keyword_l_val=.true.)
407  CALL section_add_keyword(section, keyword)
408  CALL keyword_release(keyword)
409 
410  CALL keyword_create( &
411  keyword, __location__, name="QUADRATURE_POINTS", &
412  variants=(/"LAPLACE_NUM_QUAD_POINTS"/), &
413  description="Number of quadrature points for the numerical integration in the RI-SOS-MP2-Laplace method.", &
414  usage="QUADRATURE_POINTS 6", &
415  default_i_val=5)
416  CALL section_add_keyword(section, keyword)
417  CALL keyword_release(keyword)
418 
419  CALL keyword_create( &
420  keyword, __location__, name="NUM_INTEG_GROUPS", &
421  description="Number of groups for the integration in the Laplace method. Each groups processes "// &
422  "the same amount of quadrature points. It must be a divisor of the number of quadrature points and "// &
423  "NUM_INTEG_GROUPS*GROUP_SIZE must be a divisor of the total number of processes. The default (-1) is automatic.", &
424  usage="SIZE_INTEG_GROUP 2", &
425  default_i_val=-1)
426  CALL section_add_keyword(section, keyword)
427  CALL keyword_release(keyword)
428 
429  END SUBROUTINE create_ri_laplace
430 
431 ! **************************************************************************************************
432 !> \brief ...
433 !> \param section ...
434 ! **************************************************************************************************
435  SUBROUTINE create_canonical_gradients(section)
436  TYPE(section_type), POINTER :: section
437 
438  TYPE(keyword_type), POINTER :: keyword
439  TYPE(section_type), POINTER :: subsection
440 
441  cpassert(.NOT. ASSOCIATED(section))
442  CALL section_create(section, __location__, name="CANONICAL_GRADIENTS", &
443  description="Parameters influencing gradient calculations of canonical RI methods. "// &
444  "Ignored if the IM_TIME section is set.", &
445  n_keywords=3, n_subsections=1, repeats=.false.)
446 
447  NULLIFY (subsection, keyword)
448  CALL create_cphf(subsection)
449  CALL section_add_subsection(section, subsection)
450  CALL section_release(subsection)
451 
452  CALL keyword_create(keyword, __location__, name="EPS_CANONICAL", &
453  description="Threshold under which a given ij or ab pair is considered to be degenerate and "// &
454  "its contribution to the density matrix is calculated directly. "// &
455  "Ignored in case of energy-only calculation.", &
456  usage="EPS_CANONICAL 1.0E-8", type_of_var=real_t, &
457  default_r_val=1.0e-7_dp)
458  CALL section_add_keyword(section, keyword)
459  CALL keyword_release(keyword)
460 
461  CALL keyword_create( &
462  keyword, __location__, &
463  name="FREE_HFX_BUFFER", &
464  description="Free the buffer containing the 4 center integrals used in the Hartree-Fock exchange calculation. "// &
465  "This will be effective only for gradients calculations, since for the energy only "// &
466  "case, the buffers are released by default. (Right now debugging only).", &
467  usage="FREE_HFX_BUFFER", &
468  default_l_val=.true., &
469  lone_keyword_l_val=.true.)
470  CALL section_add_keyword(section, keyword)
471  CALL keyword_release(keyword)
472 
473  CALL keyword_create( &
474  keyword, __location__, &
475  name="USE_OLD_GRADIENT_CODE", &
476  description="Use the original RI-MP2 gradient code.", &
477  usage="USE_OLD_GRADIENT_CODE T", &
478  lone_keyword_l_val=.true., &
479  default_l_val=.true.)
480  CALL section_add_keyword(section, keyword)
481  CALL keyword_release(keyword)
482 
483  CALL keyword_create( &
484  keyword, __location__, &
485  name="DOT_PRODUCT_BLKSIZE", &
486  description="Dot products for the calculation of the RPA/SOS-MP2 density matrices "// &
487  "are calculated in batches of the size given by this keyword. Larger block sizes "// &
488  "improve the performance but reduce the numerical accuracy. Recommended block sizes are multiples of the number of "// &
489  "doubles per cache line (usually 8). Ignored with MP2 gradients. Set it to -1 to prevent blocking.", &
490  default_i_val=-1)
491  CALL section_add_keyword(section, keyword)
492  CALL keyword_release(keyword)
493 
494  CALL keyword_create( &
495  keyword, __location__, &
496  name="MAX_PARALLEL_COMM", &
497  description="Sets the maximum number of parallel communication steps of the non-blocking communication scheme. "// &
498  "The number of channels is determined from the available memory. If set to a value smaller than one, "// &
499  "CP2K will use all memory for communication. A value of one enforces the blocking communication scheme "// &
500  "increasing the communication costs.", &
501  default_i_val=2)
502  CALL section_add_keyword(section, keyword)
503  CALL keyword_release(keyword)
504 
505  END SUBROUTINE create_canonical_gradients
506 
507 ! **************************************************************************************************
508 !> \brief ...
509 !> \param section ...
510 ! **************************************************************************************************
511  SUBROUTINE create_ri_rpa(section)
512  TYPE(section_type), POINTER :: section
513 
514  TYPE(keyword_type), POINTER :: keyword
515  TYPE(section_type), POINTER :: subsection
516 
517  cpassert(.NOT. ASSOCIATED(section))
518  CALL section_create(section, __location__, name="RI_RPA", &
519  description="Parameters influencing RI-RPA and GW.", &
520  n_keywords=8, n_subsections=4, repeats=.false., &
521  citations=(/delben2013, delben2015/))
522 
523  NULLIFY (keyword, subsection)
524  CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
525  description="Putting the &RI_RPA section activates RI-RPA calculation.", &
526  usage="&RI_RPA .TRUE.", &
527  default_l_val=.false., lone_keyword_l_val=.true.)
528  CALL section_add_keyword(section, keyword)
529  CALL keyword_release(keyword)
530 
531  CALL keyword_create(keyword, __location__, name="QUADRATURE_POINTS", &
532  variants=(/"RPA_NUM_QUAD_POINTS"/), &
533  description="Number of quadrature points for the numerical integration in the RI-RPA method.", &
534  usage="QUADRATURE_POINTS 60", &
535  default_i_val=40)
536  CALL section_add_keyword(section, keyword)
537  CALL keyword_release(keyword)
538 
539  CALL keyword_create(keyword, __location__, name="NUM_INTEG_GROUPS", &
540  description="Number of groups for the integration in the Laplace method. Each groups processes "// &
541  "the same amount of quadrature points. It must be a divisor of the number of quadrature points and "// &
542  "NUM_INTEG_GROUPS*GROUP_SIZE must be a divisor of the total number of processes. "// &
543  "The default (-1) is automatic.", &
544  usage="SIZE_INTEG_GROUP 2", &
545  default_i_val=-1)
546  CALL section_add_keyword(section, keyword)
547  CALL keyword_release(keyword)
548 
549  CALL keyword_create(keyword, __location__, &
550  name="MM_STYLE", &
551  description="Matrix multiplication style for the Q matrix.", &
552  usage="MM_STYLE GEMM", &
553  enum_c_vals=s2a("GEMM", "SYRK"), &
554  enum_i_vals=(/wfc_mm_style_gemm, wfc_mm_style_syrk/), &
555  enum_desc=s2a("Use pdgemm: more flops, maybe faster.", &
556  "Use pdysrk: fewer flops, maybe slower."), &
557  default_i_val=wfc_mm_style_gemm)
558  CALL section_add_keyword(section, keyword)
559  CALL keyword_release(keyword)
560 
561  CALL keyword_create( &
562  keyword, __location__, &
563  name="MINIMAX_QUADRATURE", &
564  variants=(/"MINIMAX"/), &
565  description="Use the Minimax quadrature scheme for performing the numerical integration. "// &
566  "Maximum number of quadrature point limited to 20.", &
567  usage="MINIMAX_QUADRATURE", &
568  default_l_val=.false., &
569  lone_keyword_l_val=.true.)
570  CALL section_add_keyword(section, keyword)
571  CALL keyword_release(keyword)
572 
573  CALL keyword_create( &
574  keyword, __location__, &
575  name="RI_AXK", &
576  variants=(/"AXK"/), &
577  description="Decide whether to perform an RPA-AXK calculation.", &
578  usage="RI_AXK", &
579  default_l_val=.false., &
580  lone_keyword_l_val=.true.)
581  CALL section_add_keyword(section, keyword)
582  CALL keyword_release(keyword)
583 
584  CALL keyword_create( &
585  keyword, __location__, &
586  name="RSE", &
587  variants=(/"SE"/), &
588  description="Decide whether to add singles correction.", &
589  usage="RI_AXK", &
590  default_l_val=.false., &
591  lone_keyword_l_val=.true.)
592  CALL section_add_keyword(section, keyword)
593  CALL keyword_release(keyword)
594 
595  CALL keyword_create( &
596  keyword, __location__, &
597  name="ADMM", &
598  description="Decide whether to perform ADMM in the exact exchange calc. for RPA and/or GW. "// &
599  "The ADMM XC correction is governed by the AUXILIARY_DENSITY_MATRIX_METHOD section in &DFT. "// &
600  "In most cases, the Hartree-Fock exchange is not too expensive and there is no need for ADMM, "// &
601  "ADMM can however provide significant speedup and memory savings in case of diffuse basis sets. "// &
602  "If it is a GW bandgap calculations, RI_SIGMA_X can also be used. ", &
603  usage="ADMM", &
604  default_l_val=.false., &
605  lone_keyword_l_val=.true.)
606  CALL section_add_keyword(section, keyword)
607  CALL keyword_release(keyword)
608 
609  CALL keyword_create( &
610  keyword, __location__, &
611  name="SCALE_RPA", &
612  description="Scales RPA energy contributions (RPA, AXK).", &
613  usage="SCALE_RPA 1.0", &
614  default_r_val=1.0_dp)
615  CALL section_add_keyword(section, keyword)
616  CALL keyword_release(keyword)
617 
618  CALL keyword_create( &
619  keyword, __location__, &
620  name="PRINT_DGEMM_INFO", &
621  description="Print details about all DGEMM calls.", &
622  lone_keyword_l_val=.true., &
623  default_l_val=.false.)
624  CALL section_add_keyword(section, keyword)
625  CALL keyword_release(keyword)
626 
627  ! here we generate a hfx subsection to use in the case EXX has to be computed after RPA
628  CALL create_hfx_section(subsection)
629  CALL section_add_subsection(section, subsection)
630  CALL section_release(subsection)
631 
632  ! here we generate a G0W0 subsection to use if G0W0 is desired
633  CALL create_ri_g0w0(subsection)
634  CALL section_add_subsection(section, subsection)
635  CALL section_release(subsection)
636 
637  ! here we generate a RI_AXK subsection
638  CALL create_ri_axk(subsection)
639  CALL section_add_subsection(section, subsection)
640  CALL section_release(subsection)
641 
642  END SUBROUTINE create_ri_rpa
643 
644 ! **************************************************************************************************
645 !> \brief ...
646 !> \param section ...
647 ! **************************************************************************************************
648  SUBROUTINE create_ri_axk(section)
649  TYPE(section_type), POINTER :: section
650 
651  !TYPE(section_type), POINTER :: subsection
652 
653  cpassert(.NOT. ASSOCIATED(section))
654  CALL section_create(section, __location__, name="RI_AXK", &
655  description="Parameters influencing the RI-RPA-AXK method", &
656  n_keywords=0, n_subsections=0, repeats=.false., &
657  citations=(/bates2013/))
658 
659  !NULLIFY (keyword, subsection)
660 
661  END SUBROUTINE create_ri_axk
662 
663 ! **************************************************************************************************
664 !> \brief ...
665 !> \param section ...
666 ! **************************************************************************************************
667  SUBROUTINE create_ri_g0w0(section)
668  TYPE(section_type), POINTER :: section
669 
670  TYPE(keyword_type), POINTER :: keyword
671  TYPE(section_type), POINTER :: subsection
672 
673  cpassert(.NOT. ASSOCIATED(section))
674  CALL section_create(section, __location__, name="GW", &
675  description="Parameters influencing the RI-G0W0 method", &
676  n_keywords=24, n_subsections=1, repeats=.false.)
677 
678  NULLIFY (keyword, subsection)
679 
680  CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
681  description="Activates GW calculations.", &
682  usage="&GW .TRUE.", &
683  default_l_val=.false., lone_keyword_l_val=.true.)
684  CALL section_add_keyword(section, keyword)
685  CALL keyword_release(keyword)
686 
687  CALL keyword_create(keyword, __location__, name="CORR_MOS_OCC", &
688  variants=(/"CORR_OCC"/), &
689  description="Number of occupied MOs whose energies are corrected by RI-G0W0. "// &
690  "Counting beginning from HOMO, e.g. 3 corrected occ. MOs correspond "// &
691  "to correction of HOMO, HOMO-1 and HOMO-2. Numerical effort and "// &
692  "storage of RI-G0W0 increase linearly with this number. In case you "// &
693  "want to correct all occ. MOs, insert either a negative number or "// &
694  "a number larger than the number of occ. MOs.", &
695  usage="CORR_OCC 3", &
696  default_i_val=10)
697  CALL section_add_keyword(section, keyword)
698  CALL keyword_release(keyword)
699 
700  CALL keyword_create(keyword, __location__, name="CORR_MOS_VIRT", &
701  variants=(/"CORR_VIRT"/), &
702  description="Number of virtual MOs whose energies are corrected by RI-G0W0. "// &
703  "Counting beginning from LUMO, e.g. 3 corrected occ. MOs correspond "// &
704  "to correction of LUMO, LUMO+1 and LUMO+2. Numerical effort and "// &
705  "storage of RI-G0W0 increase linearly with this number. In case you "// &
706  "want to correct all virt. MOs, insert either a negative number or "// &
707  "a number larger than the number of virt. MOs.", &
708  usage="CORR_VIRT 3", &
709  default_i_val=10)
710  CALL section_add_keyword(section, keyword)
711  CALL keyword_release(keyword)
712 
713  CALL keyword_create(keyword, __location__, name="NUMB_POLES", &
714  description="Number of poles for the fitting. Usually, two poles are sufficient. ", &
715  usage="NUMB_POLES 2", &
716  default_i_val=2)
717  CALL section_add_keyword(section, keyword)
718  CALL keyword_release(keyword)
719 
720  CALL keyword_create(keyword, __location__, name="OMEGA_MAX_FIT", &
721  description="Determines fitting range for the self-energy on the imaginary axis: "// &
722  "[0, OMEGA_MAX_FIT] for virt orbitals, [-OMEGA_MAX_FIT,0] for occ orbitals. "// &
723  "Unit: Hartree. Default: 0.734996 H = 20 eV. ", &
724  usage="OMEGA_MAX_FIT 0.5", &
725  default_r_val=0.734996_dp)
726  CALL section_add_keyword(section, keyword)
727  CALL keyword_release(keyword)
728 
729  CALL keyword_create(keyword, __location__, name="CROSSING_SEARCH", &
730  description="Determines, how the self_energy is evaluated on the real axis.", &
731  usage="CROSSING_SEARCH Z_SHOT", &
732  enum_c_vals=s2a("Z_SHOT", "NEWTON", "BISECTION"), &
733  enum_i_vals=(/ri_rpa_g0w0_crossing_z_shot, &
735  enum_desc=s2a("Calculate the derivative of Sigma and out of it Z. Then extrapolate using Z.", &
736  "Make a Newton-Raphson fix point iteration.", &
737  "Make a bisection fix point iteration."), &
738  default_i_val=ri_rpa_g0w0_crossing_newton)
739  CALL section_add_keyword(section, keyword)
740  CALL keyword_release(keyword)
741 
742  CALL keyword_create(keyword, __location__, name="FERMI_LEVEL_OFFSET", &
743  description="Fermi level for occ. orbitals: e_HOMO + FERMI_LEVEL_OFFSET; "// &
744  "Fermi level for virt. orbitals: e_LUMO - FERMI_LEVEL_OFFSET. "// &
745  "In case e_homo + FERMI_LEVEL_OFFSET < e_lumo - FERMI_LEVEL_OFFSET, "// &
746  "we set Fermi level = (e_HOMO+e_LUMO)/2. For cubic-scaling GW, the Fermi level "// &
747  "is always equal to (e_HOMO+e_LUMO)/2 regardless of FERMI_LEVEL_OFFSET.", &
748  usage="FERMI_LEVEL_OFFSET 1.0E-2", &
749  default_r_val=2.0e-2_dp)
750  CALL section_add_keyword(section, keyword)
751  CALL keyword_release(keyword)
752 
753  CALL keyword_create(keyword, __location__, name="EV_GW_ITER", &
754  description="Maximum number of iterations for eigenvalue self-consistency cycle. The "// &
755  "computational effort of GW scales linearly with this number. ", &
756  usage="EV_GW_ITER 3", &
757  default_i_val=1)
758  CALL section_add_keyword(section, keyword)
759  CALL keyword_release(keyword)
760 
761  CALL keyword_create(keyword, __location__, name="SC_GW0_ITER", &
762  description="Maximum number of iterations for GW0 self-consistency cycle. The "// &
763  "computational effort of GW is not much affected by the number of scGW0 cycles. ", &
764  usage="SC_GW0_ITER 3", &
765  default_i_val=1)
766  CALL section_add_keyword(section, keyword)
767  CALL keyword_release(keyword)
768 
769  CALL keyword_create(keyword, __location__, name="EPS_ITER", &
770  description="Target accuracy for the eigenvalue self-consistency. "// &
771  "If the G0W0 HOMO-LUMO gap differs by less than the "// &
772  "target accuracy during the iteration, the eigenvalue "// &
773  "self-consistency cycle stops. Unit: Hartree.", &
774  usage="EPS_EV_SC_ITER 0.00005", &
775  default_r_val=cp_unit_to_cp2k(value=0.00136_dp, unit_str="eV"), &
776  unit_str="eV")
777 
778  CALL section_add_keyword(section, keyword)
779  CALL keyword_release(keyword)
780 
781  CALL keyword_create(keyword, __location__, name="PRINT_EXX", &
782  description="Print exchange self-energy minus exchange correlation potential for Gamma-only "// &
783  "calculation (PRINT). For a GW calculation with k-points we use this output as "// &
784  "exchange self-energy (READ). This is a temporary solution because the hybrid MPI/OMP "// &
785  "parallelization in the HFX by Manuel Guidon conflicts with the parallelization in "// &
786  "low-scaling GW k-points which is most efficient with maximum number of MPI tasks and "// &
787  "minimum number of OMP threads. For HFX by M. Guidon, the density matrix is "// &
788  "fully replicated on every MPI rank which necessitates a high number of OMP threads per MPI "// &
789  "rank for large systems to prevent out of memory. "// &
790  "Such a high number of OMP threads would slow down the GW calculation "// &
791  "severely. Therefore, it was decided to temporarily divide the GW k-point calculation in a "// &
792  "Gamma-only HF calculation with high number of OMP threads to prevent out of memory and "// &
793  "a GW k-point calculation with 1 OMP thread per MPI rank reading the previousHF output.", &
794  usage="PRINT_EXX TRUE", &
795  enum_c_vals=s2a("TRUE", "FALSE", "READ", "SKIP_FOR_REGTEST"), &
797  enum_desc=s2a("Please, put TRUE for Gamma only calculation to get the exchange self-energy. "// &
798  "If 'SIGMA_X' and the corresponding values for the exchange-energy are written, "// &
799  "the writing has been successful", &
800  "FALSE is needed if you want to do nothing here.", &
801  "Please, put READ for the k-point GW calculation to read the exact exchange. "// &
802  "You have to provide an output file including the exact exchange. This file "// &
803  "has to be named 'exx.dat'.", &
804  "SKIP_FOR_REGTEST is only used for the GW k-point regtest where no exchange "// &
805  "self-energy is computed."), &
806  default_i_val=gw_no_print_exx)
807  CALL section_add_keyword(section, keyword)
808  CALL keyword_release(keyword)
809 
810  CALL keyword_create(keyword, __location__, name="PRINT_SELF_ENERGY", &
811  description="If true, print the self-energy for all levels for real energy "// &
812  "together with the straight line to see the quasiparticle energy as intersection. "// &
813  "In addition, prints the self-energy for imaginary frequencies together with the Pade fit.", &
814  usage="SELF_ENERGY", &
815  default_l_val=.false., &
816  lone_keyword_l_val=.true.)
817  CALL section_add_keyword(section, keyword)
818  CALL keyword_release(keyword)
819 
820  CALL keyword_create(keyword, __location__, name="RI_SIGMA_X", &
821  description="If true, the exchange self-energy is calculated approximatively with RI. "// &
822  "If false, the Hartree-Fock implementation in CP2K is used.", &
823  usage="RI_SIGMA_X", &
824  default_l_val=.true., &
825  lone_keyword_l_val=.true.)
826  CALL section_add_keyword(section, keyword)
827  CALL keyword_release(keyword)
828 
829  CALL keyword_create(keyword, __location__, name="IC_CORR_LIST", &
830  description="List of image charge correction from a previous calculation to be applied in G0W0 "// &
831  "or evGW. Keyword is active, if the first entry is positive (since IC corrections are positive "// &
832  "occupied MOs. The start corresponds to the first corrected GW level.", &
833  usage="IC_CORR_LIST <REAL> ... <REAL>", &
834  default_r_vals=(/-1.0_dp/), &
835  type_of_var=real_t, n_var=-1, unit_str="eV")
836  CALL section_add_keyword(section, keyword)
837  CALL keyword_release(keyword)
838 
839  CALL keyword_create(keyword, __location__, name="IC_CORR_LIST_BETA", &
840  description="IC_CORR_LIST for beta spins in case of open shell calculation.", &
841  usage="IC_CORR_LIST_BETA <REAL> ... <REAL>", &
842  default_r_vals=(/-1.0_dp/), &
843  type_of_var=real_t, n_var=-1, unit_str="eV")
844  CALL section_add_keyword(section, keyword)
845  CALL keyword_release(keyword)
846 
847  CALL keyword_create(keyword, __location__, name="PERIODIC_CORRECTION", &
848  description="If true, the periodic correction scheme is used employing k-points. "// &
849  "Method is not recommended to use, use instead PERIODIC_LOW_SCALING which much "// &
850  "more accurate than the periodic correction.", &
851  usage="PERIODIC_CORRECTION", &
852  default_l_val=.false., &
853  lone_keyword_l_val=.true.)
854  CALL section_add_keyword(section, keyword)
855  CALL keyword_release(keyword)
856 
857  CALL keyword_create(keyword, __location__, name="IMAGE_CHARGE_MODEL", &
858  variants=(/"IC"/), &
859  description="If true, an image charge model is applied to mimic the renormalization of "// &
860  "electronic levels of a molecule at a metallic surface. For this calculation, the molecule "// &
861  "has to be reflected on the desired xy image plane. The coordinates of the reflected molecule "// &
862  "have to be added to the coord file as ghost atoms. For the ghost atoms, identical basis sets "// &
863  "the normal atoms have to be used.", &
864  usage="IC TRUE", &
865  default_l_val=.false., &
866  lone_keyword_l_val=.true.)
867  CALL section_add_keyword(section, keyword)
868  CALL keyword_release(keyword)
869 
870  CALL keyword_create(keyword, __location__, name="ANALYTIC_CONTINUATION", &
871  description="Defines which type of analytic continuation for the self energy is used", &
872  usage="ANALYTIC_CONTINUATION", &
873  enum_c_vals=s2a("TWO_POLE", "PADE"), &
874  enum_i_vals=(/gw_two_pole_model, gw_pade_approx/), &
875  enum_desc=s2a("Use 'two-pole' model.", &
876  "Use Pade approximation."), &
877  default_i_val=gw_pade_approx)
878  CALL section_add_keyword(section, keyword)
879  CALL keyword_release(keyword)
880 
881  CALL keyword_create(keyword, __location__, name="NPARAM_PADE", &
882  description="Number of parameters for the Pade approximation "// &
883  "when using the latter for the analytic continuation of the "// &
884  "self energy. 16 parameters (corresponding to 8 poles) are "// &
885  "are recommended.", &
886  usage="NPARAM_PADE 16", &
887  default_i_val=16)
888  CALL section_add_keyword(section, keyword)
889  CALL keyword_release(keyword)
890 
891  CALL keyword_create(keyword, __location__, name="GAMMA_ONLY_SIGMA", &
892  variants=(/"GAMMA"/), &
893  description="If true, the correlation self-energy is only computed at the Gamma point. "// &
894  "The Gamma point itself is obtained by averaging over all kpoints of the DFT mesh.", &
895  usage="GAMMA TRUE", &
896  default_l_val=.false., &
897  lone_keyword_l_val=.true.)
898  CALL section_add_keyword(section, keyword)
899  CALL keyword_release(keyword)
900 
901  CALL keyword_create(keyword, __location__, name="UPDATE_XC_ENERGY", &
902  description="If true, the Hartree-Fock and RPA total energy are printed and the total energy "// &
903  "is corrected using exact exchange and the RPA correlation energy.", &
904  usage="UPDATE_XC_ENERGY", &
905  default_l_val=.false., &
906  lone_keyword_l_val=.true.)
907  CALL section_add_keyword(section, keyword)
908  CALL keyword_release(keyword)
909 
910  CALL keyword_create(keyword, __location__, name="KPOINTS_SELF_ENERGY", &
911  description="Specify number of k-points for the k-point grid of the self-energy. Internally, a "// &
912  "Monkhorst-Pack grid is used. A dense k-point grid may be necessary to compute an accurate density "// &
913  "of state from GW. Large self-energy k-meshes do not cost much more computation time.", &
914  usage="KPOINTS nx ny nz", repeats=.true., &
915  n_var=3, type_of_var=integer_t, default_i_vals=(/0, 0, 0/))
916  CALL section_add_keyword(section, keyword)
917  CALL keyword_release(keyword)
918 
919  CALL keyword_create(keyword, __location__, name="REGULARIZATION_MINIMAX", &
920  description="Tikhonov regularization for computing weights of the Fourier transform "// &
921  "from imaginary time to imaginary frequency and vice versa. Needed for large minimax "// &
922  "grids with 20 or more points and a small range.", &
923  usage="REGULARIZATION_MINIMAX 1.0E-6", &
924  default_r_val=0.0_dp)
925  CALL section_add_keyword(section, keyword)
926  CALL keyword_release(keyword)
927 
928  CALL keyword_create(keyword, __location__, name="SOC", &
929  description="Calculate the spin-orbit splitting of the eigenvalues/band structure "// &
930  "using the spin-orbit part of the GTH pseudos parametrized in Hartwigsen, Goedecker, "// &
931  "Hutter, Phys. Rev. B 58, 3641 (1998), Eq. 19, "// &
932  "parameters in Table I.", &
933  usage="SOC", &
934  enum_c_vals=s2a("NONE", "LDA", "PBE"), &
935  enum_i_vals=(/soc_none, soc_lda, soc_pbe/), &
936  enum_desc=s2a("No SOC.", &
937  "Use parameters from LDA (PADE) pseudopotential.", &
938  "Use parameters from PBE pseudopotential."), &
939  default_i_val=soc_none)
940  CALL section_add_keyword(section, keyword)
941  CALL keyword_release(keyword)
942 
943  CALL keyword_create(keyword, __location__, name="SOC_ENERGY_WINDOW", &
944  description="For perturbative SOC calculation, only "// &
945  "take frontier levels in an energy window "// &
946  "[E_HOMO - SOC_ENERGY_WINDOW/2 , E_LUMO + SOC_ENERGY_WINDOW/2 "// &
947  "into account for the diagonalization of H^GW,SOC.", &
948  usage="SOC_ENERGY_WINDOW 20.0_eV", &
949  default_r_val=cp_unit_to_cp2k(value=50.0_dp, unit_str="eV"), &
950  unit_str="eV")
951  CALL section_add_keyword(section, keyword)
952  CALL keyword_release(keyword)
953 
954  ! here we generate a subsection for the periodic GW correction
955  CALL create_periodic_gw_correction_section(subsection)
956  CALL section_add_subsection(section, subsection)
957  CALL section_release(subsection)
958 
959  ! here we generate a subsection for Bethe-Salpeter
960  CALL create_bse_section(subsection)
961  CALL section_add_subsection(section, subsection)
962  CALL section_release(subsection)
963 
964  ! here we generate a subsection for image charge calculations
965  CALL create_ic_section(subsection)
966  CALL section_add_subsection(section, subsection)
967  CALL section_release(subsection)
968 
969  ! here we generate a subsection for calculating the GW band structures
970  CALL create_kpoint_set_section(subsection)
971  CALL section_add_subsection(section, subsection)
972  CALL section_release(subsection)
973 
974  ! here we generate a subsection for additional printing
975  CALL create_print_section(subsection)
976  CALL section_add_subsection(section, subsection)
977  CALL section_release(subsection)
978 
979  END SUBROUTINE create_ri_g0w0
980 
981 ! **************************************************************************************************
982 !> \brief ...
983 !> \param section ...
984 ! **************************************************************************************************
985  SUBROUTINE create_print_section(section)
986  TYPE(section_type), POINTER :: section
987 
988  TYPE(keyword_type), POINTER :: keyword
989  TYPE(section_type), POINTER :: gw_dos_section, print_key
990 
991  cpassert(.NOT. ASSOCIATED(section))
992  NULLIFY (print_key, keyword)
993  NULLIFY (gw_dos_section, keyword)
994  CALL section_create(section, __location__, name="PRINT", &
995  description="Section of possible print options specific for the GW code.", &
996  n_keywords=0, n_subsections=2, repeats=.false.)
997 
998  CALL cp_print_key_section_create(print_key, __location__, "LOCAL_BANDGAP", &
999  description="Prints a local bandgap E_gap(r), derived from the local density of "// &
1000  "states rho(r,E). Details and formulae in the SI of the periodic GW paper (2023).", &
1001  print_level=high_print_level, add_last=add_last_numeric, &
1002  filename="LOCAL_BANDGAP", &
1003  common_iter_levels=3)
1004 
1005  CALL keyword_create(keyword, __location__, name="ENERGY_WINDOW", &
1006  description="Energy window in the LDOS for searching the gap.", &
1007  usage="ENERGY_WINDOW 6.0", &
1008  default_r_val=cp_unit_to_cp2k(value=6.0_dp, unit_str="eV"), &
1009  unit_str="eV")
1010  CALL section_add_keyword(print_key, keyword)
1011  CALL keyword_release(keyword)
1012 
1013  CALL keyword_create(keyword, __location__, name="ENERGY_SPACING", &
1014  description="Energy spacing of the LDOS for searching the gap.", &
1015  usage="ENERGY_SPACING 0.03", &
1016  default_r_val=cp_unit_to_cp2k(value=0.03_dp, unit_str="eV"), &
1017  unit_str="eV")
1018  CALL section_add_keyword(print_key, keyword)
1019  CALL keyword_release(keyword)
1020 
1021  CALL keyword_create(keyword, __location__, name="LDOS_THRESHOLD_GAP", &
1022  description="Relative LDOS threshold that determines the local bandgap.", &
1023  usage="LDOS_THRESHOLD_GAP 0.1", &
1024  default_r_val=0.1_dp)
1025  CALL section_add_keyword(print_key, keyword)
1026  CALL keyword_release(keyword)
1027 
1028  CALL keyword_create(keyword, __location__, name="STRIDE", &
1029  description="The stride (X,Y,Z) used to write the cube file "// &
1030  "(larger values result in smaller cube files). You can provide 3 numbers (for X,Y,Z) or"// &
1031  " 1 number valid for all components.", &
1032  usage="STRIDE 2 2 2", n_var=-1, default_i_vals=(/2, 2, 2/), type_of_var=integer_t)
1033  CALL section_add_keyword(print_key, keyword)
1034  CALL keyword_release(keyword)
1035 
1036  CALL section_add_subsection(section, print_key)
1037  CALL section_release(print_key)
1038 
1039  CALL section_create(gw_dos_section, __location__, name="GW_DOS", &
1040  description="Section for printing the spectral function.", &
1041  n_keywords=6, n_subsections=0, repeats=.false.)
1042 
1043  CALL keyword_create(keyword, __location__, name="LOWER_BOUND", &
1044  description="Lower bound for GW-DOS in eV.", &
1045  usage="GW_DOS_LOWER_BOUND -20.0", &
1046  default_r_val=cp_unit_to_cp2k(value=-20.0_dp, unit_str="eV"), &
1047  unit_str="eV")
1048  CALL section_add_keyword(gw_dos_section, keyword)
1049  CALL keyword_release(keyword)
1050 
1051  CALL keyword_create(keyword, __location__, name="UPPER_BOUND", &
1052  description="Upper bound for GW-DOS in eV.", &
1053  usage="GW_DOS_UPPER_BOUND 5.0", &
1054  default_r_val=cp_unit_to_cp2k(value=0.0_dp, unit_str="eV"), &
1055  unit_str="eV")
1056  CALL section_add_keyword(gw_dos_section, keyword)
1057  CALL keyword_release(keyword)
1058 
1059  CALL keyword_create(keyword, __location__, name="STEP", &
1060  description="Difference of two consecutive energy levels for GW-DOS.", &
1061  usage="GW_DOS_PRECISION 0.1", &
1062  default_r_val=cp_unit_to_cp2k(value=0.0_dp, unit_str="eV"), &
1063  unit_str="eV")
1064  CALL section_add_keyword(gw_dos_section, keyword)
1065  CALL keyword_release(keyword)
1066 
1067  CALL keyword_create(keyword, __location__, name="MIN_LEVEL_SPECTRAL", &
1068  description="Lowest energy level to print the self energy to files.", &
1069  usage="MIN_LEVEL_SELF_ENERGY 3", &
1070  default_i_val=1)
1071  CALL section_add_keyword(gw_dos_section, keyword)
1072  CALL keyword_release(keyword)
1073 
1074  CALL keyword_create(keyword, __location__, name="MAX_LEVEL_SPECTRAL", &
1075  description="Highest energy level to print the self energy to files.", &
1076  usage="MAX_LEVEL_SELF_ENERGY 6", &
1077  default_i_val=0)
1078  CALL section_add_keyword(gw_dos_section, keyword)
1079  CALL keyword_release(keyword)
1080 
1081  CALL keyword_create(keyword, __location__, name="MIN_LEVEL_SELF_ENERGY", &
1082  description="Lowest energy level to print the self energy to files.", &
1083  usage="MIN_LEVEL_SELF_ENERGY 3", &
1084  default_i_val=1)
1085  CALL section_add_keyword(gw_dos_section, keyword)
1086  CALL keyword_release(keyword)
1087 
1088  CALL keyword_create(keyword, __location__, name="MAX_LEVEL_SELF_ENERGY", &
1089  description="Highest energy level to print the self energy to files.", &
1090  usage="MAX_LEVEL_SELF_ENERGY 6", &
1091  default_i_val=0)
1092  CALL section_add_keyword(gw_dos_section, keyword)
1093  CALL keyword_release(keyword)
1094 
1095  CALL keyword_create(keyword, __location__, name="BROADENING", &
1096  description="Broadening parameter for spectral function.", &
1097  usage="BROADENING 0.001", &
1098  default_r_val=cp_unit_to_cp2k(value=0.0_dp, unit_str="eV"), &
1099  unit_str="eV")
1100  CALL section_add_keyword(gw_dos_section, keyword)
1101  CALL keyword_release(keyword)
1102 
1103  CALL section_add_subsection(section, gw_dos_section)
1104  CALL section_release(gw_dos_section)
1105 
1106  END SUBROUTINE
1107 
1108 ! **************************************************************************************************
1109 !> \brief ...
1110 !> \param section ...
1111 ! **************************************************************************************************
1112  SUBROUTINE create_periodic_gw_correction_section(section)
1113  TYPE(section_type), POINTER :: section
1114 
1115  TYPE(keyword_type), POINTER :: keyword
1116 
1117  cpassert(.NOT. ASSOCIATED(section))
1118  CALL section_create(section, __location__, name="PERIODIC_CORRECTION", &
1119  description="Parameters influencing correction for periodic GW. Old method, "// &
1120  "not recommended to use", &
1121  n_keywords=12, n_subsections=1, repeats=.false.)
1122 
1123  NULLIFY (keyword)
1124 
1125  CALL keyword_create(keyword, __location__, name="KPOINTS", &
1126  description="Specify number of k-points for a single k-point grid. Internally, a "// &
1127  "Monkhorst-Pack grid is used. Typically, even numbers are chosen such that the Gamma "// &
1128  "point is excluded from the k-point mesh.", &
1129  usage="KPOINTS nx ny nz", repeats=.true., &
1130  n_var=3, type_of_var=integer_t, default_i_vals=(/16, 16, 16/))
1131  CALL section_add_keyword(section, keyword)
1132  CALL keyword_release(keyword)
1133 
1134  CALL keyword_create(keyword, __location__, name="NUM_KP_GRIDS", &
1135  description="Number of k-point grids around the Gamma point with different resolution. "// &
1136  "E.g. for KPOINTS 4 4 4 and NUM_KP_GRIDS 3, there will be a 3x3x3 Monkhorst-Pack (MP) k-point "// &
1137  "grid for the whole Brillouin zone (excluding Gamma), another 3x3x3 MP grid with smaller "// &
1138  "spacing around Gamma (again excluding Gamma) and a very fine 4x4x4 MP grid around Gamma.", &
1139  usage="NUM_KP_GRIDS 5", &
1140  default_i_val=1)
1141  CALL section_add_keyword(section, keyword)
1142  CALL keyword_release(keyword)
1143 
1144  CALL keyword_create(keyword, __location__, name="EPS_KPOINT", &
1145  description="If the absolute value of a k-point is below EPS_KPOINT, this kpoint is "// &
1146  "neglected since the Gamma point is not included in the periodic correction.", &
1147  usage="EPS_KPOINT 1.0E-4", &
1148  default_r_val=1.0e-05_dp)
1149  CALL section_add_keyword(section, keyword)
1150  CALL keyword_release(keyword)
1151 
1152  CALL keyword_create(keyword, __location__, name="MO_COEFF_GAMMA", &
1153  description="If true, only the MO coefficients at the Gamma point are used for the periodic "// &
1154  "correction. Otherwise, the MO coeffs are computed at every k-point which is much more "// &
1155  "expensive. It should be okay to use the Gamma MO coefficients.", &
1156  usage="MO_COEFF_GAMMA", &
1157  default_l_val=.true., &
1158  lone_keyword_l_val=.true.)
1159  CALL section_add_keyword(section, keyword)
1160  CALL keyword_release(keyword)
1161 
1162  CALL keyword_create(keyword, __location__, name="AVERAGE_DEGENERATE_LEVELS", &
1163  variants=(/"ADL"/), &
1164  description="If true, the correlation self-energy of degenerate levels is averaged.", &
1165  usage="AVERAGE_DEGENERATE_LEVELS", &
1166  default_l_val=.true., &
1167  lone_keyword_l_val=.true.)
1168  CALL section_add_keyword(section, keyword)
1169  CALL keyword_release(keyword)
1170 
1171  CALL keyword_create(keyword, __location__, name="EPS_EIGENVAL", &
1172  description="Threshold for considering levels as degenerate. Unit: Hartree.", &
1173  usage="EPS_EIGENVAL 1.0E-5", &
1174  default_r_val=2.0e-04_dp)
1175  CALL section_add_keyword(section, keyword)
1176  CALL keyword_release(keyword)
1177 
1178  CALL keyword_create(keyword, __location__, name="EXTRAPOLATE_KPOINTS", &
1179  variants=(/"EXTRAPOLATE"/), &
1180  description="If true, extrapolates the k-point mesh. Only working if k-point mesh numbers are "// &
1181  "divisible by 4, e.g. 8x8x8 or 12x12x12 is recommended.", &
1182  usage="EXTRAPOLATE_KPOINTS FALSE", &
1183  default_l_val=.true., &
1184  lone_keyword_l_val=.true.)
1185  CALL section_add_keyword(section, keyword)
1186  CALL keyword_release(keyword)
1187 
1188  CALL keyword_create(keyword, __location__, name="DO_AUX_BAS_GW", &
1189  description="If true, use a different basis for the periodic correction. This can be necessary "// &
1190  "in case a diffused basis is used for GW to converge the HOMO-LUMO gap. In this case, "// &
1191  "numerical problems may occur due to diffuse functions in the basis. This keyword only works if "// &
1192  "AUX_GW <basis set> is specified in the kind section for every atom kind.", &
1193  usage="AUX_BAS_GW TRUE", &
1194  default_l_val=.false., &
1195  lone_keyword_l_val=.true.)
1196  CALL section_add_keyword(section, keyword)
1197  CALL keyword_release(keyword)
1198 
1199  CALL keyword_create(keyword, __location__, name="FRACTION_AUX_MOS", &
1200  description="Fraction how many MOs are used in the auxiliary basis.", &
1201  usage="FRACTION_AUX_MOS 0.6", &
1202  default_r_val=0.5_dp)
1203  CALL section_add_keyword(section, keyword)
1204  CALL keyword_release(keyword)
1205 
1206  CALL keyword_create(keyword, __location__, name="NUM_OMEGA_POINTS", &
1207  description="Number of Clenshaw-Curtis integration points for the periodic correction in cubic- "// &
1208  "scaling GW. This variable is a dummy variable for canonical N^4 GW calculations.", &
1209  usage="NUM_OMEGA_POINTS 200", &
1210  default_i_val=300)
1211  CALL section_add_keyword(section, keyword)
1212  CALL keyword_release(keyword)
1213 
1214  END SUBROUTINE
1215 
1216 ! **************************************************************************************************
1217 !> \brief ...
1218 !> \param section ...
1219 ! **************************************************************************************************
1220  SUBROUTINE create_bse_section(section)
1221  TYPE(section_type), POINTER :: section
1222 
1223  TYPE(keyword_type), POINTER :: keyword
1224  TYPE(section_type), POINTER :: subsection
1225 
1226  cpassert(.NOT. ASSOCIATED(section))
1227  CALL section_create(section, __location__, name="BSE", &
1228  description="Parameters for a calculation solving the Bethe-Salpeter equation "// &
1229  "(BSE) for electronic excitations. The full BSE "// &
1230  "$\left( \begin{array}{cc}A & B\\B & A\end{array} \right)$ "// &
1231  "$\left( \begin{array}{cc}\mathbf{X}^{(n)}\\\mathbf{Y}^{(n)}\end{array} \right) = "// &
1232  "\Omega^{(n)}\left(\begin{array}{cc}1&0\\0&-1\end{array}\right)$ "// &
1233  "$\left(\begin{array}{cc}\mathbf{X}^{(n)}\\\mathbf{Y}^{(n)}\end{array}\right)$ "// &
1234  "enables of electronic excitation energies $\Omega^{(n)}$. The BSE can be solved by diagonalizing "// &
1235  "the full ABBA-matrix or by setting B=0, i.e. within the Tamm-Dancoff approximation (TDA). "// &
1236  "Preliminary reference: Eq. (35) in PRB 92, 045209 (2015); http://dx.doi.org/10.1103/PhysRevB.92.045209", &
1237  n_keywords=8, n_subsections=2, repeats=.false.)
1238 
1239  NULLIFY (keyword)
1240 
1241  CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
1242  description="Activates BSE calculations.", &
1243  usage="&BSE .TRUE.", &
1244  default_l_val=.false., lone_keyword_l_val=.true.)
1245  CALL section_add_keyword(section, keyword)
1246  CALL keyword_release(keyword)
1247 
1248  CALL keyword_create(keyword, __location__, name="SPIN_CONFIG", &
1249  description="Choose between calculation of singlet or triplet excitation (cf. given Reference above).", &
1250  usage="SPIN_CONFIG TRIPLET", &
1251  enum_c_vals=s2a("SINGLET", "TRIPLET"), &
1252  enum_i_vals=(/bse_singlet, bse_triplet/), &
1253  enum_desc=s2a("Computes singlet excitations.", &
1254  "Computes triplet excitations."), &
1255  default_i_val=bse_singlet)
1256  CALL section_add_keyword(section, keyword)
1257  CALL keyword_release(keyword)
1258 
1259  CALL keyword_create(keyword, __location__, name="BSE_DIAG_METHOD", &
1260  description="Method for BSE calculations. "// &
1261  "Choose between full or iterative diagonalization.", &
1262  usage="&BSE_DIAG_METHOD FULLDIAG", &
1263  enum_c_vals=s2a("FULLDIAG", "ITERDIAG"), &
1264  enum_i_vals=(/bse_fulldiag, bse_iterdiag/), &
1265  enum_desc=s2a("Fully diagonalizes the BSE matrices within the chosen level of approximation.", &
1266  "Iterative diagonalization has not been implemented yet."), &
1267  default_i_val=bse_fulldiag)
1268  CALL section_add_keyword(section, keyword)
1269  CALL keyword_release(keyword)
1270 
1271  CALL keyword_create(keyword, __location__, name="BSE_APPROX", &
1272  description="Level of approximation applied to BSE calculations. "// &
1273  "Choose between Tamm Dancoff approximation (TDA) and diagonalization of the full ABBA-matrix.", &
1274  usage="&BSE_APPROX TDA", &
1275  enum_c_vals=s2a("TDA", "ABBA", "BOTH"), &
1276  enum_i_vals=(/bse_tda, bse_abba, bse_both/), &
1277  enum_desc=s2a("The TDA is applied, i.e. B=0.", &
1278  "The ABBA-matrix is diagonalized, i.e. the TDA is not applied.", &
1279  "The BSE is solved within the TDA (B=0) as well as for the full ABBA-matrix."), &
1280  default_i_val=bse_tda)
1281  CALL section_add_keyword(section, keyword)
1282  CALL keyword_release(keyword)
1283 
1284  CALL keyword_create(keyword, __location__, name="ENERGY_CUTOFF_OCC", &
1285  description="Removing all orbitals with indices i,j from A_ia,jb and B_ia,jb with energy difference "// &
1286  "to HOMO level larger than the given energy cutoff, i.e. "// &
1287  "$\epsilon_{\mathrm{HOMO}}-\epsilon_a>\mathtt{ENERGY\_CUTOFF\_OCC}$. "// &
1288  "Can be used to accelerate runtime and reduce memory consumption.", &
1289  usage="ENERGY_CUTOFF_OCC 10.0", unit_str="eV", &
1290  type_of_var=real_t, default_r_val=-1.0_dp)
1291  CALL section_add_keyword(section, keyword)
1292  CALL keyword_release(keyword)
1293 
1294  CALL keyword_create(keyword, __location__, name="ENERGY_CUTOFF_VIRT", &
1295  description="Removing all orbitals with indices a,b from A_ia,jb and B_ia,jb with energy difference "// &
1296  "to LUMO level larger than the given energy cutoff, i.e. "// &
1297  "$\epsilon_i-\epsilon_{\mathrm{LUMO}}>\mathtt{ENERGY\_CUTOFF\_VIRT}$. "// &
1298  "Can be used to accelerate runtime and reduce memory consumption.", &
1299  usage="ENERGY_CUTOFF_VIRT 10.0", unit_str="eV", &
1300  type_of_var=real_t, default_r_val=-1.0_dp)
1301  CALL section_add_keyword(section, keyword)
1302  CALL keyword_release(keyword)
1303 
1304  CALL keyword_create(keyword, __location__, name="BSE_DEBUG_PRINT", &
1305  description="Activates debug print statements in the BSE calculation.", &
1306  usage="&BSE_DEBUG_PRINT .TRUE.", &
1307  default_l_val=.false., lone_keyword_l_val=.true.)
1308  CALL section_add_keyword(section, keyword)
1309  CALL keyword_release(keyword)
1310 
1311  CALL keyword_create(keyword, __location__, name="NUM_PRINT_EXC", &
1312  description="Number of printed excitation energies. Does not affect computation time.", &
1313  usage="NUM_PRINT_EXC 10", &
1314  default_i_val=10)
1315  CALL section_add_keyword(section, keyword)
1316  CALL keyword_release(keyword)
1317 
1318  CALL keyword_create(keyword, __location__, name="EPS_X", &
1319  description="Threshold for printing contributions of singleparticle "// &
1320  "transitions, i.e. elements of the eigenvector $X_{ia}^{(n)}$.", &
1321  usage="EPS_X 0.1", &
1322  type_of_var=real_t, default_r_val=0.1_dp)
1323  CALL section_add_keyword(section, keyword)
1324  CALL keyword_release(keyword)
1325 
1326  NULLIFY (subsection)
1327  CALL create_bse_iterat_section(subsection)
1328  CALL section_add_subsection(section, subsection)
1329  CALL section_release(subsection)
1330 
1331  END SUBROUTINE
1332 
1333 ! **************************************************************************************************
1334 !> \brief ...
1335 !> \param section ...
1336 ! **************************************************************************************************
1337  SUBROUTINE create_bse_iterat_section(section)
1338  TYPE(section_type), POINTER :: section
1339 
1340  TYPE(keyword_type), POINTER :: keyword
1341 
1342  cpassert(.NOT. ASSOCIATED(section))
1343  CALL section_create(section, __location__, name="BSE_ITERAT", &
1344  description="Parameters influencing the iterative Bethe-Salpeter calculation. "// &
1345  "The iterative solver has not been fully implemented yet.", &
1346  n_keywords=9, n_subsections=0, repeats=.false.)
1347 
1348  NULLIFY (keyword)
1349 
1350  CALL keyword_create(keyword, __location__, name="DAVIDSON_ABORT_COND", &
1351  description="Desired abortion condition for Davidson solver", &
1352  usage="DAVIDSON_ABORT_COND OR", &
1353  enum_c_vals=s2a("EN", "RES", "OR"), &
1355  enum_desc=s2a("Uses energy threshold for successfully exiting solver.", &
1356  "Uses residual threshold for successfully exiting solver.", &
1357  "Uses either energy or residual threshold for successfully exiting solver."), &
1358  default_i_val=bse_iter_en_cond)
1359  CALL section_add_keyword(section, keyword)
1360  CALL keyword_release(keyword)
1361 
1362  CALL keyword_create(keyword, __location__, name="NUM_EXC_EN", &
1363  description="Number of lowest excitation energies to be computed.", &
1364  usage="NUM_EXC_EN 3", &
1365  default_i_val=3)
1366  CALL section_add_keyword(section, keyword)
1367  CALL keyword_release(keyword)
1368 
1369  CALL keyword_create(keyword, __location__, name="NUM_ADD_START_Z_SPACE", &
1370  description="Determines the initial dimension of the subspace as "// &
1371  "dim = (NUM_EXC_EN+NUM_ADD_START_Z_SPACE)", &
1372  usage="NUM_ADD_START_Z_SPACE 1", &
1373  default_i_val=0)
1374  CALL section_add_keyword(section, keyword)
1375  CALL keyword_release(keyword)
1376 
1377  CALL keyword_create(keyword, __location__, name="FAC_MAX_Z_SPACE", &
1378  description="Factor to determine maximum dimension of the Davidson subspace. "// &
1379  "dimension = (NUM_EXC_EN+NUM_ADD_START_Z_SPACE)*FAC_MAX_Z_SPACE", &
1380  usage="FAC_MAX_Z_SPACE 5", &
1381  default_i_val=5)
1382  CALL section_add_keyword(section, keyword)
1383  CALL keyword_release(keyword)
1384 
1385  CALL keyword_create(keyword, __location__, name="NUM_NEW_T", &
1386  description="Number of new t vectors added. "// &
1387  "Must be smaller/equals (NUM_EXC_EN+NUM_ADD_START_Z_SPACE)", &
1388  usage="NUM_NEW_T 4", &
1389  default_i_val=1)
1390  CALL section_add_keyword(section, keyword)
1391  CALL keyword_release(keyword)
1392 
1393  CALL keyword_create(keyword, __location__, name="EPS_RES", &
1394  description="Threshold for stopping the iteration for computing the transition energies. "// &
1395  "If the residuals inside the Davidson space change by less than EPS_RES (in eV), the iteration "// &
1396  "stops.", &
1397  usage="EPS_RES 0.001", unit_str="eV", &
1398  type_of_var=real_t, default_r_val=0.001_dp)
1399  CALL section_add_keyword(section, keyword)
1400  CALL keyword_release(keyword)
1401 
1402  CALL keyword_create(keyword, __location__, name="EPS_EXC_EN", &
1403  description="Threshold for stopping the iteration for computing the transition energies. "// &
1404  "If the desired excitation energies change by less than EPS_EXC_EN (in eV), the iteration "// &
1405  "stops.", &
1406  usage="EPS_EXC_EN 0.001", unit_str="eV", &
1407  type_of_var=real_t, default_r_val=0.001_dp)
1408  CALL section_add_keyword(section, keyword)
1409  CALL keyword_release(keyword)
1410 
1411  CALL keyword_create(keyword, __location__, name="NUM_DAVIDSON_ITER", &
1412  description="Maximum number of iterations for determining the transition energies.", &
1413  usage="MAX_ITER 100", &
1414  default_i_val=100)
1415  CALL section_add_keyword(section, keyword)
1416  CALL keyword_release(keyword)
1417 
1418  CALL keyword_create(keyword, __location__, name="Z_SPACE_ENERGY_CUTOFF", &
1419  description="Cutoff (in eV) for maximal energy difference entering the A matrix. "// &
1420  "Per default and for negative values, there is no cutoff applied.", &
1421  usage="Z_SPACE_ENERGY_CUTOFF 60", unit_str="eV", &
1422  type_of_var=real_t, default_r_val=-1.0_dp)
1423  CALL section_add_keyword(section, keyword)
1424  CALL keyword_release(keyword)
1425  END SUBROUTINE
1426 
1427 ! **************************************************************************************************
1428 !> \brief ...
1429 !> \param section ...
1430 ! **************************************************************************************************
1431  SUBROUTINE create_ic_section(section)
1432  TYPE(section_type), POINTER :: section
1433 
1434  TYPE(keyword_type), POINTER :: keyword
1435 
1436  cpassert(.NOT. ASSOCIATED(section))
1437  CALL section_create(section, __location__, name="IC", &
1438  description="Parameters influencing the image charge correction. "// &
1439  "The image plane is always an xy plane, so adjust the molecule according "// &
1440  "to that. ", &
1441  n_keywords=3, n_subsections=1, repeats=.false.)
1442 
1443  NULLIFY (keyword)
1444 
1445  CALL keyword_create(keyword, __location__, name="PRINT_IC_LIST", &
1446  description="If true, the image charge correction values are printed in a list, "// &
1447  "such that it can be used as input for a subsequent evGW calculation.", &
1448  usage="PRINT_IC_VALUES", &
1449  default_l_val=.false., &
1450  lone_keyword_l_val=.true.)
1451  CALL section_add_keyword(section, keyword)
1452  CALL keyword_release(keyword)
1453 
1454  CALL keyword_create(keyword, __location__, name="EPS_DIST", &
1455  description="Threshold where molecule and image molecule have to coincide. ", &
1456  usage="EPS_DIST 0.1", unit_str="angstrom", &
1457  type_of_var=real_t, default_r_val=3.0e-02_dp, repeats=.false.)
1458  CALL section_add_keyword(section, keyword)
1459  CALL keyword_release(keyword)
1460 
1461  END SUBROUTINE
1462 
1463 ! **************************************************************************************************
1464 !> \brief ...
1465 !> \param section ...
1466 ! **************************************************************************************************
1467  SUBROUTINE create_low_scaling(section)
1468  TYPE(section_type), POINTER :: section
1469 
1470  TYPE(keyword_type), POINTER :: keyword
1471  TYPE(section_type), POINTER :: subsection
1472 
1473  cpassert(.NOT. ASSOCIATED(section))
1474  CALL section_create( &
1475  section, __location__, name="LOW_SCALING", &
1476  description="Cubic scaling RI-RPA, GW and Laplace-SOS-MP2 method using the imaginary time formalism. "// &
1477  "EPS_GRID in WFC_GPW section controls accuracy / req. memory for 3-center integrals. "// &
1478  "SORT_BASIS EXP should be specified in DFT section.", &
1479  n_keywords=12, n_subsections=2, repeats=.false.)
1480 
1481  NULLIFY (keyword)
1482  CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
1483  description="Activates cubic-scaling RPA, GW and Laplace-SOS-MP2 calculations.", &
1484  usage="&LOW_SCALING .TRUE.", &
1485  default_l_val=.false., lone_keyword_l_val=.true.)
1486  CALL section_add_keyword(section, keyword)
1487  CALL keyword_release(keyword)
1488 
1489  CALL keyword_create(keyword, __location__, name="MEMORY_CUT", &
1490  description="Reduces memory for sparse tensor contractions by this factor. "// &
1491  "A high value leads to some loss of performance. "// &
1492  "This memory reduction factor applies to storage of the tensors 'M occ' / 'M virt' "// &
1493  "but does not reduce storage of '3c ints'.", &
1494  usage="MEMORY_CUT 16", &
1495  default_i_val=5)
1496  CALL section_add_keyword(section, keyword)
1497  CALL keyword_release(keyword)
1498 
1499  CALL keyword_create(keyword, __location__, name="MEMORY_INFO", &
1500  description="Decide whether to print memory info on the sparse matrices.", &
1501  usage="MEMORY_INFO", &
1502  default_l_val=.false., &
1503  lone_keyword_l_val=.true.)
1504  CALL section_add_keyword(section, keyword)
1505  CALL keyword_release(keyword)
1506 
1507  CALL keyword_create( &
1508  keyword, __location__, name="EPS_FILTER", &
1509  description="Determines a threshold for the DBCSR based multiply. "// &
1510  "Normally, this EPS_FILTER determines accuracy and timing of low-scaling RPA and GW calculations.", &
1511  usage="EPS_FILTER 1.0E-10 ", type_of_var=real_t, &
1512  default_r_val=1.0e-9_dp)
1513  CALL section_add_keyword(section, keyword)
1514  CALL keyword_release(keyword)
1515 
1516  CALL keyword_create( &
1517  keyword, __location__, name="EPS_FILTER_FACTOR", &
1518  description="Multiply EPS_FILTER with this factor to determine filter epsilon "// &
1519  "for DBCSR based multiply P(it)=(Mocc(it))^T*Mvirt(it) "// &
1520  "Default should be kept.", &
1521  type_of_var=real_t, &
1522  default_r_val=10.0_dp)
1523  CALL section_add_keyword(section, keyword)
1524  CALL keyword_release(keyword)
1525 
1526  CALL keyword_create( &
1527  keyword, __location__, &
1528  name="EPS_STORAGE_SCALING", &
1529  variants=(/"EPS_STORAGE"/), &
1530  description="Scaling factor to scale EPS_FILTER. Storage threshold for compression "// &
1531  "will be EPS_FILTER*EPS_STORAGE_SCALING.", &
1532  default_r_val=1.0e-3_dp)
1533  CALL section_add_keyword(section, keyword)
1534  CALL keyword_release(keyword)
1535 
1536  CALL keyword_create( &
1537  keyword, __location__, &
1538  name="DO_KPOINTS", &
1539  description="Besides in DFT, this keyword has to be switched on if one wants to do kpoints in. "// &
1540  "cubic RPA.", &
1541  usage="DO_KPOINTS", &
1542  default_l_val=.false., &
1543  lone_keyword_l_val=.true.)
1544  CALL section_add_keyword(section, keyword)
1545  CALL keyword_release(keyword)
1546 
1547  CALL keyword_create( &
1548  keyword, __location__, name="KPOINTS", &
1549  description="Keyword activates periodic, low-scaling GW calculations (&LOW_SCALING section also needed). "// &
1550  "For periodic calculations, kpoints are used for the density response, the "// &
1551  "Coulomb interaction and the screened Coulomb interaction. For 2d periodic systems, e.g. xz "// &
1552  "periodicity, please also specify KPOINTS, e.g. N_x 1 N_z.", &
1553  usage="KPOINTS N_x N_y N_z", &
1554  n_var=3, type_of_var=integer_t, default_i_vals=(/0, 0, 0/))
1555  CALL section_add_keyword(section, keyword)
1556  CALL keyword_release(keyword)
1557 
1558  CALL keyword_create( &
1559  keyword, __location__, &
1560  name="KPOINT_WEIGHTS_W", &
1561  description="For kpoints in low-scaling GW, a Monkhorst-Pack mesh is used. The screened Coulomb "// &
1562  "interaction W(k) needs special care near the Gamma point (e.g. in 3d, W(k) diverges at the "// &
1563  "Gamma point with W(k) ~ k^alpha). KPOINT_WEIGHTS_W decides how the weights of the "// &
1564  "Monkhorst-Pack mesh are chosen to compute W(R) = int_BZ W(k) exp(ikR) dk (BZ=Brllouin zone). ", &
1565  usage="KPOINT_WEIGHTS_W AUTO", &
1566  enum_c_vals=s2a("TAILORED", "AUTO", "UNIFORM"), &
1568  enum_desc=s2a("Choose k-point integration weights such that the function f(k)=k^alpha is "// &
1569  "exactly integrated. alpha is specified using EXPONENT_TAILORED_WEIGHTS.", &
1570  "As 'TAILORED', but alpha is chosen automatically according to dimensionality "// &
1571  "(3D: alpha = -2 for 3D, 2D: alpha = -1 for exchange self-energy, uniform "// &
1572  "weights for correlation self-energy).", &
1573  "Choose the same weight for every k-point (original Monkhorst-Pack method)."), &
1574  default_i_val=kp_weights_w_uniform)
1575  CALL section_add_keyword(section, keyword)
1576  CALL keyword_release(keyword)
1577 
1578  CALL keyword_create( &
1579  keyword, __location__, &
1580  name="EXPONENT_TAILORED_WEIGHTS", &
1581  description="Gives the exponent of exactly integrated function in case 'KPOINT_WEIGHTS_W "// &
1582  "TAILORED' is chosen.", &
1583  usage="EXPONENT_TAILORED_WEIGHTS -2", &
1584  default_r_val=-2.0_dp)
1585  CALL section_add_keyword(section, keyword)
1586  CALL keyword_release(keyword)
1587 
1588  CALL keyword_create( &
1589  keyword, __location__, &
1590  name="REGULARIZATION_RI", &
1591  description="Parameter to reduce the expansion coefficients in RI for periodic GW. Larger parameter "// &
1592  "means smaller expansion coefficients that leads to a more stable calculation at the price "// &
1593  "of a slightly worse RI approximation. In case the parameter 0.0 is chosen, ordinary RI is used.", &
1594  usage="REGULARIZATION_RI 1.0E-4", &
1595  default_r_val=0.0_dp)
1596  CALL section_add_keyword(section, keyword)
1597  CALL keyword_release(keyword)
1598 
1599  CALL keyword_create( &
1600  keyword, __location__, &
1601  name="EPS_EIGVAL_S", &
1602  description="Parameter to reduce the expansion coefficients in RI for periodic GW. Removes all "// &
1603  "eigenvectors and eigenvalues of S_PQ(k) that are smaller than EPS_EIGVAL_S. ", &
1604  usage="EPS_EIGVAL_S 1.0E-3", &
1605  default_r_val=0.0_dp)
1606  CALL section_add_keyword(section, keyword)
1607  CALL keyword_release(keyword)
1608 
1609  CALL keyword_create( &
1610  keyword, __location__, &
1611  name="EPS_EIGVAL_S_GAMMA", &
1612  description="Parameter to reduce the expansion coefficients in RI for periodic GW. Removes all "// &
1613  "eigenvectors and eigenvalues of M_PQ(k=0) that are smaller than EPS_EIGVAL_S. ", &
1614  usage="EPS_EIGVAL_S 1.0E-3", &
1615  default_r_val=0.0_dp)
1616  CALL section_add_keyword(section, keyword)
1617  CALL keyword_release(keyword)
1618 
1619  CALL keyword_create( &
1620  keyword, __location__, &
1621  name="MAKE_CHI_POS_DEFINITE", &
1622  description="If true, makes eigenvalue decomposition of chi(iw,k) and removes negative "// &
1623  "eigenvalues. May increase computational cost significantly. Only recommended to try in case "// &
1624  "Cholesky decomposition of epsilon(iw,k) fails.", &
1625  usage="MAKE_CHI_POS_DEFINITE", &
1626  default_l_val=.true., &
1627  lone_keyword_l_val=.true.)
1628  CALL section_add_keyword(section, keyword)
1629  CALL keyword_release(keyword)
1630 
1631  CALL keyword_create( &
1632  keyword, __location__, &
1633  name="MAKE_OVERLAP_MAT_AO_POS_DEFINITE", &
1634  description="If true, makes eigenvalue decomposition of S_mu,nu(k) and removes negative "// &
1635  "eigenvalues. Slightly increases computational cost. Only recommended to try in case "// &
1636  "Cholesky decomposition of S_mu,nu(k) fails (error message: Cholesky decompose failed: "// &
1637  "matrix is not positive definite or ill-conditioned; when calling create_kp_and_calc_kp_orbitals).", &
1638  usage="MAKE_OVERLAP_MAT_AO_POS_DEFINITE", &
1639  default_l_val=.false., &
1640  lone_keyword_l_val=.true.)
1641  CALL section_add_keyword(section, keyword)
1642  CALL keyword_release(keyword)
1643 
1644  CALL keyword_create( &
1645  keyword, __location__, &
1646  name="DO_EXTRAPOLATE_KPOINTS", &
1647  description="If true, use a larger k-mesh to extrapolate the k-point integration of W. "// &
1648  "For example, in 2D, when using KPOINTS 4 4 1, an additional 6x6x1 mesh will be used to "// &
1649  "extrapolate the k-point integration of W with N_k^-0.5, where Nk is the number of k-points.", &
1650  usage="DO_EXTRAPOLATE_KPOINTS FALSE", &
1651  default_l_val=.true., &
1652  lone_keyword_l_val=.true.)
1653  CALL section_add_keyword(section, keyword)
1654  CALL keyword_release(keyword)
1655 
1656  CALL keyword_create( &
1657  keyword, __location__, &
1658  name="TRUNC_COULOMB_RI_X", &
1659  description="If true, use the truncated Coulomb operator for the exchange-self-energy in "// &
1660  "periodic GW.", &
1661  usage="TRUNC_COULOMB_RI_X", &
1662  default_l_val=.true., &
1663  lone_keyword_l_val=.true.)
1664  CALL section_add_keyword(section, keyword)
1665  CALL keyword_release(keyword)
1666 
1667  CALL keyword_create( &
1668  keyword, __location__, &
1669  name="REL_CUTOFF_TRUNC_COULOMB_RI_X", &
1670  description="Only active in case TRUNC_COULOMB_RI_X = True. Normally, relative cutoff = 0.5 is "// &
1671  "good choice; still needs to be evaluated for RI schemes. ", &
1672  usage="REL_CUTOFF_TRUNC_COULOMB_RI_X 0.3", &
1673  default_r_val=0.5_dp)
1674  CALL section_add_keyword(section, keyword)
1675  CALL keyword_release(keyword)
1676 
1677  CALL keyword_create( &
1678  keyword, __location__, &
1679  name="KEEP_QUADRATURE", &
1680  variants=s2a("KEEP_WEIGHTS", "KEEP_QUAD", "KEEP_WEIGHT"), &
1681  description="Keep the Laplace quadrature defined at the first energy evaluations throughout "// &
1682  "the run. Allows to have consistent force evaluations.", &
1683  usage="KEEP_QUADRATURE", &
1684  default_l_val=.true., &
1685  lone_keyword_l_val=.true.)
1686  CALL section_add_keyword(section, keyword)
1687  CALL keyword_release(keyword)
1688 
1689  CALL keyword_create( &
1690  keyword, __location__, &
1691  name="K_MESH_G_FACTOR", &
1692  description="The k-mesh for the Green's function can be chosen to be larger than the k-mesh for "// &
1693  "W (without much higher computational cost). The factor given here multiplies the mesh for W to obtain "// &
1694  "the k-mesh for G. Example: factor 4, k-mesh for W: 4x4x1 -> k-mesh for G: 16x16x1 (z-dir. is "// &
1695  "non-periodic).", &
1696  default_i_val=1)
1697  CALL section_add_keyword(section, keyword)
1698  CALL keyword_release(keyword)
1699 
1700  CALL keyword_create( &
1701  keyword, __location__, &
1702  name="MIN_BLOCK_SIZE", &
1703  description="Minimum tensor block size. Adjusting this value may have minor effect on "// &
1704  "performance but default should be good enough.", &
1705  default_i_val=5)
1706  CALL section_add_keyword(section, keyword)
1707  CALL keyword_release(keyword)
1708 
1709  CALL keyword_create( &
1710  keyword, __location__, &
1711  name="MIN_BLOCK_SIZE_MO", &
1712  description="Tensor block size for MOs. Only relevant for GW calculations. "// &
1713  "The memory consumption of GW scales as O(MIN_BLOCK_SIZE_MO). It is recommended to "// &
1714  "set this parameter to a smaller number if GW runs out of memory. "// &
1715  "Otherwise the default should not be changed.", &
1716  default_i_val=64)
1717  CALL section_add_keyword(section, keyword)
1718  CALL keyword_release(keyword)
1719 
1720  NULLIFY (subsection)
1721  CALL create_low_scaling_cphf(subsection)
1722  CALL section_add_subsection(section, subsection)
1723  CALL section_release(subsection)
1724 
1725  END SUBROUTINE
1726 
1727 ! **************************************************************************************************
1728 !> \brief ...
1729 !> \param section ...
1730 ! **************************************************************************************************
1731  SUBROUTINE create_wfc_gpw(section)
1732  TYPE(section_type), POINTER :: section
1733 
1734  TYPE(keyword_type), POINTER :: keyword
1735 
1736  cpassert(.NOT. ASSOCIATED(section))
1737  CALL section_create(section, __location__, name="WFC_GPW", &
1738  description="Parameters for the GPW approach in Wavefunction-based Correlation methods", &
1739  n_keywords=5, n_subsections=0, repeats=.false.)
1740 
1741  NULLIFY (keyword)
1742  CALL keyword_create(keyword, __location__, name="EPS_GRID", &
1743  description="Determines a threshold for the GPW based integration", &
1744  usage="EPS_GRID 1.0E-9 ", type_of_var=real_t, &
1745  default_r_val=1.0e-8_dp)
1746  CALL section_add_keyword(section, keyword)
1747  CALL keyword_release(keyword)
1748 
1749  CALL keyword_create( &
1750  keyword, __location__, name="EPS_FILTER", &
1751  description="Determines a threshold for the DBCSR based multiply (usually 10 times smaller than EPS_GRID). "// &
1752  "Normally, this EPS_FILTER determines accuracy and timing of cubic-scaling RPA calculation.", &
1753  usage="EPS_FILTER 1.0E-10 ", type_of_var=real_t, &
1754  default_r_val=1.0e-9_dp)
1755  CALL section_add_keyword(section, keyword)
1756  CALL keyword_release(keyword)
1757 
1758  CALL keyword_create(keyword, __location__, name="CUTOFF", &
1759  description="The cutoff of the finest grid level in the MP2 gpw integration.", &
1760  usage="CUTOFF 300", type_of_var=real_t, &
1761  default_r_val=300.0_dp)
1762  CALL section_add_keyword(section, keyword)
1763  CALL keyword_release(keyword)
1764 
1765  CALL keyword_create(keyword, __location__, name="REL_CUTOFF", &
1766  variants=(/"RELATIVE_CUTOFF"/), &
1767  description="Determines the grid at which a Gaussian is mapped.", &
1768  usage="REL_CUTOFF 50", type_of_var=real_t, &
1769  default_r_val=50.0_dp)
1770  CALL section_add_keyword(section, keyword)
1771  CALL keyword_release(keyword)
1772 
1773  CALL keyword_create(keyword, __location__, name="PRINT_LEVEL", &
1774  variants=(/"IOLEVEL"/), &
1775  description="How much output is written by the individual groups.", &
1776  usage="PRINT_LEVEL HIGH", &
1777  default_i_val=silent_print_level, enum_c_vals= &
1778  s2a("SILENT", "LOW", "MEDIUM", "HIGH", "DEBUG"), &
1779  enum_desc=s2a("Almost no output", &
1780  "Little output", "Quite some output", "Lots of output", &
1781  "Everything is written out, useful for debugging purposes only"), &
1784  CALL section_add_keyword(section, keyword)
1785  CALL keyword_release(keyword)
1786 
1787  CALL keyword_create( &
1788  keyword, __location__, name="EPS_PGF_ORB_S", &
1789  description="Screening for overlap matrix in RI. Usually, it is best to choose this parameter "// &
1790  "to be very small since the inversion of overlap matrix might be ill-conditioned.", &
1791  usage="EPS_PGF_ORB_S 1.0E-10 ", type_of_var=real_t, &
1792  default_r_val=1.0e-10_dp)
1793  CALL section_add_keyword(section, keyword)
1794  CALL keyword_release(keyword)
1795 
1796  END SUBROUTINE create_wfc_gpw
1797 
1798 ! **************************************************************************************************
1799 !> \brief ...
1800 !> \param section ...
1801 ! **************************************************************************************************
1802  SUBROUTINE create_cphf(section)
1803  TYPE(section_type), POINTER :: section
1804 
1805  TYPE(keyword_type), POINTER :: keyword
1806 
1807  cpassert(.NOT. ASSOCIATED(section))
1808  CALL section_create( &
1809  section, __location__, name="CPHF", &
1810  description="Parameters influencing the solution of the Z-vector equations in MP2 gradients calculations.", &
1811  n_keywords=2, n_subsections=0, repeats=.false., &
1812  citations=(/delben2013/))
1813 
1814  NULLIFY (keyword)
1815 
1816  CALL keyword_create(keyword, __location__, name="MAX_ITER", &
1817  variants=(/"MAX_NUM_ITER"/), &
1818  description="Maximum number of iterations allowed for the solution of the Z-vector equations.", &
1819  usage="MAX_ITER 50", &
1820  default_i_val=30)
1821  CALL section_add_keyword(section, keyword)
1822  CALL keyword_release(keyword)
1823 
1824  CALL keyword_create(keyword, __location__, name="RESTART_EVERY", &
1825  description="Restart iteration every given number of steps.", &
1826  usage="RESTART_EVERY 5", &
1827  default_i_val=5)
1828  CALL section_add_keyword(section, keyword)
1829  CALL keyword_release(keyword)
1830 
1831  CALL keyword_create(keyword, __location__, name="SOLVER_METHOD", &
1832  description="Chose solver of the z-vector equations.", &
1833  usage="SOLVER_METHOD POPLE", enum_c_vals= &
1834  s2a("POPLE", "CG", "RICHARDSON", "SD"), &
1835  enum_desc=s2a("Pople's method (Default).", &
1836  "Conjugated gradient method (equivalent to Pople).", &
1837  "Richardson iteration", &
1838  "Steepest Descent iteration"), &
1840  default_i_val=z_solver_pople)
1841  CALL section_add_keyword(section, keyword)
1842  CALL keyword_release(keyword)
1843 
1844  CALL keyword_create(keyword, __location__, name="EPS_CONV", &
1845  description="Convergence threshold for the solution of the Z-vector equations. "// &
1846  "The Z-vector equations have the form of a linear system of equations Ax=b, "// &
1847  "convergence is achieved when |Ax-b|<=EPS_CONV.", &
1848  usage="EPS_CONV 1.0E-6", type_of_var=real_t, &
1849  default_r_val=1.0e-4_dp)
1850  CALL section_add_keyword(section, keyword)
1851  CALL keyword_release(keyword)
1852 
1853  CALL keyword_create(keyword, __location__, name="SCALE_STEP_SIZE", &
1854  description="Scaling factor of each step.", &
1855  usage="SCALE_STEP_SIZE 1.0", &
1856  default_r_val=1.0_dp)
1857  CALL section_add_keyword(section, keyword)
1858  CALL keyword_release(keyword)
1859 
1860  CALL keyword_create(keyword, __location__, name="ENFORCE_DECREASE", &
1861  description="Restarts if residual does not decrease.", &
1862  usage="ENFORCE_DECREASE T", &
1863  lone_keyword_l_val=.true., &
1864  default_l_val=.false.)
1865  CALL section_add_keyword(section, keyword)
1866  CALL keyword_release(keyword)
1867 
1868  CALL keyword_create(keyword, __location__, name="DO_POLAK_RIBIERE", &
1869  description="Use a Polak-Ribiere update of the search vector in CG instead of the Fletcher "// &
1870  "Reeves update. Improves the convergence with modified step sizes. "// &
1871  "Ignored with other methods than CG.", &
1872  usage="ENFORCE_DECREASE T", &
1873  lone_keyword_l_val=.true., &
1874  default_l_val=.false.)
1875  CALL section_add_keyword(section, keyword)
1876  CALL keyword_release(keyword)
1877 
1878  CALL keyword_create(keyword, __location__, name="RECALC_RESIDUAL", &
1879  description="Recalculates residual in every step.", &
1880  usage="RECALC_RESIDUAL T", &
1881  lone_keyword_l_val=.true., &
1882  default_l_val=.false.)
1883  CALL section_add_keyword(section, keyword)
1884  CALL keyword_release(keyword)
1885 
1886  END SUBROUTINE create_cphf
1887 
1888 ! **************************************************************************************************
1889 !> \brief ...
1890 !> \param section ...
1891 ! **************************************************************************************************
1892  SUBROUTINE create_low_scaling_cphf(section)
1893  TYPE(section_type), POINTER :: section
1894 
1895  TYPE(keyword_type), POINTER :: keyword
1896 
1897  NULLIFY (keyword)
1898 
1899  cpassert(.NOT. ASSOCIATED(section))
1900  CALL section_create(section, __location__, name="CPHF", &
1901  description="Parameters influencing the solution of the Z-vector equations "// &
1902  "in low-scaling Laplace-SOS-MP2 gradients calculations.", &
1903  n_keywords=5, n_subsections=0, repeats=.false.)
1904 
1905  CALL keyword_create(keyword, __location__, name="EPS_CONV", &
1906  description="Target accuracy for Z-vector euation solution.", &
1907  usage="EPS 1.e-6", default_r_val=1.e-6_dp)
1908  CALL section_add_keyword(section, keyword)
1909  CALL keyword_release(keyword)
1910 
1911  CALL keyword_create(keyword, __location__, name="MAX_ITER", &
1912  description="Maximum number of conjugate gradient iteration to be performed for one optimization.", &
1913  usage="MAX_ITER 200", default_i_val=50)
1914  CALL section_add_keyword(section, keyword)
1915  CALL keyword_release(keyword)
1916 
1917  CALL keyword_create( &
1918  keyword, __location__, name="PRECONDITIONER", &
1919  description="Type of preconditioner to be used with all minimization schemes. "// &
1920  "They differ in effectiveness, cost of construction, cost of application. "// &
1921  "Properly preconditioned minimization can be orders of magnitude faster than doing nothing.", &
1922  usage="PRECONDITIONER FULL_ALL", &
1923  default_i_val=ot_precond_full_all, &
1924  enum_c_vals=s2a("FULL_ALL", "FULL_SINGLE_INVERSE", "FULL_SINGLE", "FULL_KINETIC", "FULL_S_INVERSE", &
1925  "NONE"), &
1926  enum_desc=s2a("Most effective state selective preconditioner based on diagonalization, "// &
1927  "requires the ENERGY_GAP parameter to be an underestimate of the HOMO-LUMO gap. "// &
1928  "This preconditioner is recommended for almost all systems, except very large systems where "// &
1929  "make_preconditioner would dominate the total computational cost.", &
1930  "Based on H-eS cholesky inversion, similar to FULL_SINGLE in preconditioning efficiency "// &
1931  "but cheaper to construct, "// &
1932  "might be somewhat less robust. Recommended for large systems.", &
1933  "Based on H-eS diagonalisation, not as good as FULL_ALL, but somewhat cheaper to apply. ", &
1934  "Cholesky inversion of S and T, fast construction, robust, and relatively good, "// &
1935  "use for very large systems.", &
1936  "Cholesky inversion of S, not as good as FULL_KINETIC, yet equally expensive.", &
1937  "skip preconditioning"), &
1940  CALL section_add_keyword(section, keyword)
1941  CALL keyword_release(keyword)
1942 
1943  CALL keyword_create(keyword, __location__, name="ENERGY_GAP", &
1944  description="Energy gap estimate [a.u.] for preconditioning", &
1945  usage="ENERGY_GAP 0.1", &
1946  default_r_val=0.2_dp)
1947  CALL section_add_keyword(section, keyword)
1948  CALL keyword_release(keyword)
1949 
1950  END SUBROUTINE create_low_scaling_cphf
1951 
1952 ! **************************************************************************************************
1953 !> \brief ...
1954 !> \param section ...
1955 ! **************************************************************************************************
1956  SUBROUTINE create_mp2_potential(section)
1957  TYPE(section_type), POINTER :: section
1958 
1959  TYPE(keyword_type), POINTER :: keyword
1960 
1961  cpassert(.NOT. ASSOCIATED(section))
1962  CALL section_create(section, __location__, name="INTERACTION_POTENTIAL", &
1963  description="Parameters the interaction potential in computing the biel integrals", &
1964  n_keywords=4, n_subsections=0, repeats=.false.)
1965 
1966  NULLIFY (keyword)
1967  CALL keyword_create( &
1968  keyword, __location__, &
1969  name="POTENTIAL_TYPE", &
1970  description="Which interaction potential should be used "// &
1971  "(Coulomb, TShPSC operator).", &
1972  usage="POTENTIAL_TYPE TSHPSC", &
1973  enum_c_vals=s2a("COULOMB", "TShPSC", "LONGRANGE", "SHORTRANGE", "TRUNCATED", "MIX_CL", "IDENTITY"), &
1974  enum_i_vals=(/do_potential_coulomb, &
1980  do_potential_id/), &
1981  enum_desc=s2a("Coulomb potential: 1/r", &
1982  "| Range | TShPSC |"//newline// &
1983  "| ----- | ------ |"//newline// &
1984  "| $ x \leq R_c $ | $ 1/x - s/R_c $ |"//newline// &
1985  "| $ R_c < x \leq nR_c $ | "// &
1986  "$ (1 - s)/R_c - (x - R_c)/R_c^2 + (x - R_c)^2/R_c^3 - "// &
1987  "(2n^2 - 7n + 9 - 4s)(x - R_c)^3/(R_c^4(n^2 - 2n + 1)(n - 1)) + "// &
1988  "(6-3s - 4n + n^2)(x - R_c)^4/(R_c^5(n^4 - 4n^3 + 6n^2 - 4n + 1)) $ "// &
1989  "(4th order polynomial) | "//newline// &
1990  "| $ x > nR_c $ | $ 0 $ | "//newline, &
1991  "Longrange Coulomb potential: $ \operatorname{erf}(wr)/r $", &
1992  "Shortrange Coulomb potential: $ \operatorname{erfc}(wr)/r $", &
1993  "Truncated Coulomb potential", &
1994  "Mixed Coulomb/Longrange Coulomb potential", &
1995  "Delta potential"), &
1996  default_i_val=do_potential_coulomb)
1997  CALL section_add_keyword(section, keyword)
1998  CALL keyword_release(keyword)
1999 
2000  CALL keyword_create(keyword, __location__, name="TRUNCATION_RADIUS", &
2001  variants=(/"CUTOFF_RADIUS"/), &
2002  description="Determines truncation radius for the truncated potentials. "// &
2003  "Only valid when doing truncated calculations", &
2004  usage="TRUNCATION_RADIUS 10.0", type_of_var=real_t, &
2005  default_r_val=10.0_dp, &
2006  unit_str="angstrom")
2007  CALL section_add_keyword(section, keyword)
2008  CALL keyword_release(keyword)
2009 
2010  CALL keyword_create( &
2011  keyword, __location__, &
2012  name="POTENTIAL_DATA", &
2013  variants=s2a("TShPSC_DATA", "T_C_G_DATA"), &
2014  description="Location of the file TShPSC.dat or t_c_g.dat that contains the data for the "// &
2015  "evaluation of the evaluation of the truncated potentials", &
2016  usage="TShPSC_DATA t_sh_p_s_c.dat", &
2017  default_c_val="t_sh_p_s_c.dat")
2018  CALL section_add_keyword(section, keyword)
2019  CALL keyword_release(keyword)
2020 
2021  CALL keyword_create( &
2022  keyword, __location__, &
2023  name="OMEGA", &
2024  description="Range separation parameter for the longrange or shortrange potential. "// &
2025  "Only valid when longrange or shortrange potential is requested.", &
2026  usage="OMEGA 0.5", type_of_var=real_t, &
2027  default_r_val=0.5_dp)
2028  CALL section_add_keyword(section, keyword)
2029  CALL keyword_release(keyword)
2030 
2031  CALL keyword_create( &
2032  keyword, __location__, &
2033  name="SCALE_COULOMB", &
2034  description="Scaling factor of (truncated) Coulomb potential in mixed (truncated) Coulomb/Longrange potential. "// &
2035  "Only valid when mixed potential is requested.", &
2036  usage="OMEGA 0.5", type_of_var=real_t, &
2037  default_r_val=1.0_dp)
2038  CALL section_add_keyword(section, keyword)
2039  CALL keyword_release(keyword)
2040 
2041  CALL keyword_create( &
2042  keyword, __location__, &
2043  name="SCALE_LONGRANGE", &
2044  description="Scaling factor of longrange Coulomb potential in mixed (truncated) Coulomb/Longrange potential. "// &
2045  "Only valid when mixed potential is requested.", &
2046  usage="OMEGA 0.5", type_of_var=real_t, &
2047  default_r_val=1.0_dp)
2048  CALL section_add_keyword(section, keyword)
2049  CALL keyword_release(keyword)
2050 
2051  END SUBROUTINE create_mp2_potential
2052 
2053 ! **************************************************************************************************
2054 !> \brief ...
2055 !> \param section ...
2056 ! **************************************************************************************************
2057  SUBROUTINE create_ri_section(section)
2058  TYPE(section_type), POINTER :: section
2059 
2060  TYPE(keyword_type), POINTER :: keyword
2061  TYPE(section_type), POINTER :: subsection
2062 
2063  cpassert(.NOT. ASSOCIATED(section))
2064  CALL section_create(section, __location__, name="RI", &
2065  description="Parameters influencing resolution of the identity (RI) that is "// &
2066  "used in RI-MP2, RI-RPA, RI-SOS-MP2 and GW (inside RI-RPA).", &
2067  n_keywords=6, n_subsections=2, repeats=.false.)
2068 
2069  NULLIFY (subsection)
2070  CALL create_ri_metric_section(subsection)
2071  CALL section_add_subsection(section, subsection)
2072  CALL section_release(subsection)
2073 
2074  CALL create_opt_ri_basis(subsection)
2075  CALL section_add_subsection(section, subsection)
2076  CALL section_release(subsection)
2077 
2078  NULLIFY (keyword)
2079  CALL keyword_create( &
2080  keyword, __location__, &
2081  name="ROW_BLOCK", &
2082  variants=(/"ROW_BLOCK_SIZE"/), &
2083  description="Size of the row block used in the SCALAPACK block cyclic data distribution. "// &
2084  "Default is (ROW_BLOCK=-1) is automatic. "// &
2085  "A proper choice can speedup the parallel matrix multiplication in the case of RI-RPA and RI-SOS-MP2-Laplace.", &
2086  usage="ROW_BLOCK 512", &
2087  default_i_val=-1)
2088  CALL section_add_keyword(section, keyword)
2089  CALL keyword_release(keyword)
2090 
2091  CALL keyword_create( &
2092  keyword, __location__, &
2093  name="COL_BLOCK", &
2094  variants=(/"COL_BLOCK_SIZE"/), &
2095  description="Size of the column block used in the SCALAPACK block cyclic data distribution. "// &
2096  "Default is (COL_BLOCK=-1) is automatic. "// &
2097  "A proper choice can speedup the parallel matrix multiplication in the case of RI-RPA and RI-SOS-MP2-Laplace.", &
2098  usage="COL_BLOCK 512", &
2099  default_i_val=-1)
2100  CALL section_add_keyword(section, keyword)
2101  CALL keyword_release(keyword)
2102 
2103  CALL keyword_create( &
2104  keyword, __location__, &
2105  name="CALC_COND_NUM", &
2106  variants=(/"CALC_CONDITION_NUMBER"/), &
2107  description="Calculate the condition number of the (P|Q) matrix for the RI methods.", &
2108  usage="CALC_COND_NUM", &
2109  default_l_val=.false., &
2110  lone_keyword_l_val=.true.)
2111  CALL section_add_keyword(section, keyword)
2112  CALL keyword_release(keyword)
2113 
2114  CALL keyword_create(keyword, __location__, name="DO_SVD", &
2115  description="Wether to perform a singular value decomposition instead of the Cholesky decomposition "// &
2116  "of the potential operator in the RI basis. Computationally expensive but numerically more stable. "// &
2117  "It reduces the computational costs of some subsequent steps. Recommended when a longrange Coulomb "// &
2118  "potential is employed.", &
2119  usage="DO_SVD .TRUE.", &
2120  default_l_val=.false., lone_keyword_l_val=.true.)
2121  CALL section_add_keyword(section, keyword)
2122  CALL keyword_release(keyword)
2123 
2124  CALL keyword_create(keyword, __location__, name="ERI_BLKSIZE", &
2125  description="block sizes for tensors (only used if ERI_METHOD=MME). First value "// &
2126  "is the block size for ORB basis, second value is the block size for RI_AUX basis.", &
2127  usage="ERI_BLKSIZE", &
2128  n_var=2, &
2129  default_i_vals=(/4, 16/))
2130  CALL section_add_keyword(section, keyword)
2131  CALL keyword_release(keyword)
2132 
2133  END SUBROUTINE create_ri_section
2134 
2135 ! **************************************************************************************************
2136 !> \brief ...
2137 !> \param section ...
2138 ! **************************************************************************************************
2139  SUBROUTINE create_integrals_section(section)
2140  TYPE(section_type), POINTER :: section
2141 
2142  TYPE(keyword_type), POINTER :: keyword
2143  TYPE(section_type), POINTER :: subsection
2144 
2145  cpassert(.NOT. ASSOCIATED(section))
2146  CALL section_create(section, __location__, name="INTEGRALS", &
2147  description="Parameters controlling how to compute integrals that are needed "// &
2148  "in MP2, RI-MP2, RI-RPA, RI-SOS-MP2 and GW (inside RI-RPA).", &
2149  n_keywords=2, n_subsections=3, repeats=.false.)
2150 
2151  NULLIFY (subsection)
2152  CALL create_eri_mme_section(subsection)
2153  CALL section_add_subsection(section, subsection)
2154  CALL section_release(subsection)
2155 
2156  CALL create_wfc_gpw(subsection)
2157  CALL section_add_subsection(section, subsection)
2158  CALL section_release(subsection)
2159 
2160  CALL create_mp2_potential(subsection)
2161  CALL section_add_subsection(section, subsection)
2162  CALL section_release(subsection)
2163 
2164  NULLIFY (keyword)
2165  CALL keyword_create(keyword, __location__, name="ERI_METHOD", &
2166  description="Method for calculating periodic electron repulsion integrals "// &
2167  "(MME method is faster but experimental, forces not yet implemented). "// &
2168  "Obara-Saika (OS) for the Coulomb operator can only be used for non-periodic calculations.", &
2169  usage="ERI_METHOD MME", &
2170  enum_c_vals=s2a("DEFAULT", "GPW", "MME", "OS"), &
2171  enum_i_vals=(/eri_default, do_eri_gpw, do_eri_mme, do_eri_os/), &
2172  enum_desc=s2a("Use default ERI method (for periodic systems: GPW, for molecules: OS, "// &
2173  "for MP2 and RI-MP2: GPW in any case).", &
2174  "Uses Gaussian Plane Wave method [DelBen2013].", &
2175  "Uses MiniMax-Ewald method (experimental, ERI_MME subsection, only for fully periodic "// &
2176  "systems with orthorhombic cells).", &
2177  "Use analytical Obara-Saika method."), &
2178  default_i_val=eri_default)
2179  CALL section_add_keyword(section, keyword)
2180  CALL keyword_release(keyword)
2181 
2182  CALL keyword_create(keyword, __location__, name="SIZE_LATTICE_SUM", &
2183  description="Size of sum range L. ", &
2184  usage="SIZE_LATTICE_SUM 10", &
2185  default_i_val=5)
2186  CALL section_add_keyword(section, keyword)
2187  CALL keyword_release(keyword)
2188 
2189  END SUBROUTINE create_integrals_section
2190 
2191 ! **************************************************************************************************
2192 !> \brief ...
2193 !> \param section ...
2194 ! **************************************************************************************************
2195  SUBROUTINE create_ri_metric_section(section)
2196  TYPE(section_type), POINTER :: section
2197 
2198  TYPE(keyword_type), POINTER :: keyword
2199 
2200  cpassert(.NOT. ASSOCIATED(section))
2201  CALL section_create(section, __location__, name="RI_METRIC", &
2202  description="Sets up RI metric", &
2203  repeats=.false.)
2204 
2205  NULLIFY (keyword)
2206  CALL keyword_create( &
2207  keyword, __location__, &
2208  name="POTENTIAL_TYPE", &
2209  description="Decides which operator/metric is used for resolution of the identity (RI).", &
2210  usage="POTENTIAL_TYPE DEFAULT", &
2211  enum_c_vals=s2a("DEFAULT", "COULOMB", "IDENTITY", "LONGRANGE", "SHORTRANGE", "TRUNCATED"), &
2214  enum_desc=s2a("Use Coulomb metric for RI-MP2 and normal-scaling RI-SOS-MP2, RI-RPA and GW. "// &
2215  "Use Overlap metric for low-scaling RI-SOS-MP2, RI-RPA and GW for periodic systems. "// &
2216  "Use truncated Coulomb metric for low-scaling RI-SOS-MP2, RI-RPA and GW for non-periodic systems.", &
2217  "Coulomb metric: 1/r. Recommended for RI-MP2,", &
2218  "Overlap metric: delta(r).", &
2219  "Longrange metric: erf(omega*r)/r. Not recommended with DO_SVD .TRUE.", &
2220  "Shortrange metric: erfc(omega*r)/r", &
2221  "Truncated Coulomb metric: if (r &lt; R_c) 1/r else 0. More "// &
2222  "accurate than IDENTITY for non-periodic systems. Recommended for low-scaling methods."), &
2223  default_i_val=ri_default)
2224  CALL section_add_keyword(section, keyword)
2225  CALL keyword_release(keyword)
2226 
2227  NULLIFY (keyword)
2228  CALL keyword_create( &
2229  keyword, __location__, &
2230  name="OMEGA", &
2231  description=operator"The range parameter for the short/long range (in 1/a0).", &
2232  usage="OMEGA 0.5", &
2233  default_r_val=0.0_dp)
2234  CALL section_add_keyword(section, keyword)
2235  CALL keyword_release(keyword)
2236 
2237  CALL keyword_create(keyword, __location__, name="CUTOFF_RADIUS", &
2238  description="The cutoff radius (in Angstrom) for the truncated Coulomb operator.", &
2239  usage="CUTOFF_RADIUS 3.0", default_r_val=cp_unit_to_cp2k(value=3.0_dp, unit_str="angstrom"), &
2240  type_of_var=real_t, unit_str="angstrom")
2241  CALL section_add_keyword(section, keyword)
2242  CALL keyword_release(keyword)
2243 
2244  CALL keyword_create( &
2245  keyword, __location__, &
2246  name="T_C_G_DATA", &
2247  description="Location of the file t_c_g.dat that contains the data for the "// &
2248  "evaluation of the truncated gamma function ", &
2249  default_c_val="t_c_g.dat")
2250  CALL section_add_keyword(section, keyword)
2251  CALL keyword_release(keyword)
2252 
2253  CALL keyword_create(keyword, __location__, name="EPS_RANGE", &
2254  description="The threshold to determine the effective range of the short range "// &
2255  "RI metric: erfc(omega*eff_range)/eff_range = EPS_RANGE", &
2256  default_r_val=1.0e-08_dp, &
2257  repeats=.false.)
2258  CALL section_add_keyword(section, keyword)
2259  CALL keyword_release(keyword)
2260 
2261  END SUBROUTINE create_ri_metric_section
2262 
2263 END MODULE input_cp2k_mp2
collects all references to literature in CP2K as new algorithms / method are included from literature...
Definition: bibliography.F:28
integer, save, public wilhelm2016b
Definition: bibliography.F:43
integer, save, public delben2015
Definition: bibliography.F:43
integer, save, public wilhelm2016a
Definition: bibliography.F:43
integer, save, public wilhelm2018
Definition: bibliography.F:43
integer, save, public rybkin2016
Definition: bibliography.F:43
integer, save, public delben2013
Definition: bibliography.F:43
integer, save, public delben2012
Definition: bibliography.F:43
integer, save, public delben2015b
Definition: bibliography.F:43
integer, save, public wilhelm2017
Definition: bibliography.F:43
integer, save, public bates2013
Definition: bibliography.F:43
Interface to Minimax-Ewald method for periodic ERI's to be used in CP2K.
subroutine, public create_eri_mme_section(section, default_n_minimax)
Create main input section.
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer, parameter, public debug_print_level
integer, parameter, public low_print_level
integer, parameter, public medium_print_level
integer, parameter, public high_print_level
integer, parameter, public add_last_numeric
integer, parameter, public silent_print_level
subroutine, public cp_print_key_section_create(print_key_section, location, name, description, print_level, each_iter_names, each_iter_values, add_last, filename, common_iter_levels, citations, unit_str)
creates a print_key section
unit conversion facility
Definition: cp_units.F:30
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
Definition: cp_units.F:1150
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public wfc_mm_style_syrk
integer, parameter, public soc_lda
integer, parameter, public bse_iterdiag
integer, parameter, public gw_pade_approx
integer, parameter, public z_solver_pople
integer, parameter, public bse_iter_both_cond
integer, parameter, public do_eri_os
integer, parameter, public gw_skip_for_regtest
integer, parameter, public mp2_method_direct
integer, parameter, public kp_weights_w_auto
integer, parameter, public kp_weights_w_uniform
integer, parameter, public bse_iter_res_cond
integer, parameter, public do_eri_mme
integer, parameter, public ri_rpa_g0w0_crossing_bisection
integer, parameter, public wfc_mm_style_gemm
integer, parameter, public ot_precond_full_kinetic
integer, parameter, public gw_print_exx
integer, parameter, public do_potential_mix_cl
integer, parameter, public bse_singlet
integer, parameter, public z_solver_cg
integer, parameter, public bse_fulldiag
integer, parameter, public bse_triplet
integer, parameter, public ot_precond_full_single
integer, parameter, public bse_tda
integer, parameter, public bse_both
integer, parameter, public do_potential_truncated
integer, parameter, public bse_iter_en_cond
integer, parameter, public z_solver_sd
integer, parameter, public ot_precond_none
integer, parameter, public mp2_method_gpw
integer, parameter, public ot_precond_full_single_inverse
integer, parameter, public do_potential_id
integer, parameter, public mp2_method_none
integer, parameter, public soc_pbe
integer, parameter, public gw_read_exx
integer, parameter, public eri_default
integer, parameter, public ri_default
integer, parameter, public bse_abba
integer, parameter, public do_potential_coulomb
integer, parameter, public gaussian
integer, parameter, public ri_rpa_g0w0_crossing_z_shot
integer, parameter, public do_potential_short
integer, parameter, public z_solver_richardson
integer, parameter, public gw_no_print_exx
integer, parameter, public kp_weights_w_tailored
integer, parameter, public ot_precond_s_inverse
integer, parameter, public do_potential_long
integer, parameter, public gw_two_pole_model
integer, parameter, public do_eri_gpw
integer, parameter, public do_potential_tshpsc
integer, parameter, public ri_rpa_g0w0_crossing_newton
integer, parameter, public numerical
integer, parameter, public ot_precond_full_all
integer, parameter, public soc_none
function that builds the hartree fock exchange section of the input
subroutine, public create_hfx_section(section)
creates the input section for the hf part
function that build the kpoints section of the input
subroutine, public create_kpoint_set_section(section, section_name)
...
input section for MP2
subroutine, public create_mp2_section(section)
creates the input section for the mp2 part
represents keywords in an input
subroutine, public keyword_release(keyword)
releases the given keyword (see doc/ReferenceCounting.html)
subroutine, public keyword_create(keyword, location, name, description, usage, type_of_var, n_var, repeats, variants, default_val, default_l_val, default_r_val, default_lc_val, default_c_val, default_i_val, default_l_vals, default_r_vals, default_c_vals, default_i_vals, lone_keyword_val, lone_keyword_l_val, lone_keyword_r_val, lone_keyword_c_val, lone_keyword_i_val, lone_keyword_l_vals, lone_keyword_r_vals, lone_keyword_c_vals, lone_keyword_i_vals, enum_c_vals, enum_i_vals, enum, enum_strict, enum_desc, unit_str, citations, deprecation_notice, removed)
creates a keyword object
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_create(section, location, name, description, n_keywords, n_subsections, repeats, citations)
creates a list of keywords
subroutine, public section_add_keyword(section, keyword)
adds a keyword to the given section
subroutine, public section_add_subsection(section, subsection)
adds a subsection to the given section
recursive subroutine, public section_release(section)
releases the given keyword list (see doc/ReferenceCounting.html)
a wrapper for basic fortran types.
integer, parameter, public real_t
integer, parameter, public integer_t
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
Utilities for string manipulations.
character(len=1), parameter, public newline