39#include "./base/base_uses.f90"
44 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'input_cp2k_transport'
60 cpassert(.NOT.
ASSOCIATED(section))
62 description=
"Specifies the parameters for transport, sets parameters for the OMEN code, "// &
63 "see also <https://nano-tcad.ee.ethz.ch>.", &
65 n_keywords=29, n_subsections=4, repeats=.false.)
70 keyword, __location__, name=
"TRANSPORT_METHOD", &
71 description=
"Preferred method for transport calculations.", &
72 usage=
"TRANSPORT_METHOD <method>", &
74 enum_c_vals=
s2a(
"TRANSPORT",
"LOCAL_SCF",
"TRANSMISSION"), &
75 enum_desc=
s2a(
"self-consistent CP2K and OMEN transport calculations", &
76 "CP2K valence Hamiltonian + OMEN self-consistent calculations on conduction electrons", &
77 "self-consistent transmission calculations without applied bias voltage"), &
83 keyword, __location__, name=
"QT_FORMALISM", &
84 description=
"Preferred quantum transport formalism to compute the current and density.", &
85 usage=
"QT_FORMALISM <method>", &
87 enum_c_vals=
s2a(
"NEGF",
"QTBM"), &
88 enum_desc=
s2a(
"The non-equilibrium Green's function formalism.", &
89 "The quantum transmitting boundary method / wave-function formalism."), &
95 description=
"The number of terms in the PEXSI's pole expansion method.", &
96 usage=
"NUM_POLE <integer>", default_i_val=64)
101 description=
"The number of k points for determination of the singularities.", &
102 usage=
"N_KPOINTS <integer>", default_i_val=64)
107 description=
"Max number of energy points per small interval.", &
108 usage=
"NUM_INTERVAL <integer>", default_i_val=10)
112 CALL keyword_create(keyword, __location__, name=
"TASKS_PER_ENERGY_POINT", &
113 description=
"Number of tasks per energy point. The value should be a divisor of the total "// &
114 "number of MPI ranks.", &
115 usage=
"TASKS_PER_ENERGY_POINT <integer>", default_i_val=1)
119 CALL keyword_create(keyword, __location__, name=
"TASKS_PER_POLE", &
120 description=
"Number of tasks per pole in the pole expansion method. The value should be a "// &
121 "divisor of the total number of MPI ranks.", &
122 usage=
"TASKS_PER_POLE <integer>", default_i_val=1)
126 CALL keyword_create(keyword, __location__, name=
"GPUS_PER_POINT", &
127 description=
"Number of GPUs per energy point for SplitSolve. Needs to be a power of two", &
128 usage=
"GPUS_PER_POINT <integer>", default_i_val=2)
132 CALL keyword_create(keyword, __location__, name=
"COLZERO_THRESHOLD", &
133 description=
"The smallest number that is not zero in the full diagonalization part.", &
134 usage=
"COLZERO_THRESHOLD <real>", default_r_val=1.0e-12_dp)
139 description=
"The smallest eigenvalue that is kept.", &
140 usage=
"EPS_LIMIT <real>", default_r_val=1.0e-4_dp)
145 description=
"The smallest eigenvalue that is kept on the complex contour.", &
146 usage=
"EPS_LIMIT_CC <real>", default_r_val=1.0e-6_dp)
151 description=
"The smallest imaginary part that a decaying eigenvalue may have not to be "// &
152 "considered as propagating.", &
153 usage=
"EPS_DECAY <real>", default_r_val=1.0e-4_dp)
157 CALL keyword_create(keyword, __location__, name=
"EPS_SINGULARITY_CURVATURES", &
158 description=
"Filter for degenerate bands in the bandstructure.", &
159 usage=
"EPS_SINGULARITY_CURVATURES <real>", default_r_val=1.0e-12_dp)
164 description=
"Accuracy to which the Fermi level should be determined.", &
165 usage=
"EPS_MU <real>", default_r_val=1.0e-6_dp)
169 CALL keyword_create(keyword, __location__, name=
"EPS_EIGVAL_DEGEN", &
170 description=
"Filter for degenerate bands in the injection vector.", &
171 usage=
"EPS_EIGVAL_DEGEN <real>", default_r_val=1.0e-6_dp)
176 description=
"Cutoff for the tail of the Fermi function.", &
177 usage=
"EPS_FERMI <real>", default_r_val=0.0_dp)
181 CALL keyword_create(keyword, __location__, name=
"ENERGY_INTERVAL", &
182 description=
"Distance between energy points in eV.", &
183 usage=
"ENERGY_INTERVAL <real>", default_r_val=1.0e-3_dp)
188 description=
"Smallest enery distance in energy vector.", &
189 usage=
"MIN_INTERVAL <real>", default_r_val=1.0e-4_dp)
194 description=
"Temperature.", &
195 usage=
"TEMPERATURE [K] 300.0", &
201 CALL keyword_create(keyword, __location__, name=
"CSR_SCREENING", &
202 description=
"Whether distance screening should be applied to improve sparsity of CSR matrices.", &
203 default_l_val=.true., lone_keyword_l_val=.true.)
208 keyword, __location__, name=
"LINEAR_SOLVER", &
209 description=
"Preferred solver for solving the linear system of equations.", &
210 usage=
"LINEAR_SOLVER <solver>", &
212 enum_c_vals=
s2a(
"SplitSolve",
"SuperLU",
"MUMPS",
"Full",
"Banded",
"PARDISO",
"UMFPACK"), &
219 keyword, __location__, name=
"MATRIX_INVERSION_METHOD", &
220 description=
"Preferred matrix inversion method.", &
221 usage=
"MATRIX_INVERSION_METHOD <solver>", &
223 enum_c_vals=
s2a(
"Full",
"PEXSI",
"PARDISO",
"RGF"), &
228 CALL keyword_create(keyword, __location__, name=
"INJECTION_METHOD", &
229 description=
"Method to solve the eigenvalue problem for the open boundary conditions.", &
230 usage=
"INJECTION_METHOD <method>", &
232 enum_c_vals=
s2a(
"EVP",
"BEYN"), &
233 enum_desc=
s2a(
"Full eigenvalue solver.", &
234 "Beyn eigenvalue solver."), &
240 keyword, __location__, name=
"CUTOUT", &
241 description=
"The number of atoms at the beginning and the end of the structure where the density should "// &
243 usage=
"CUTOUT <integer> <integer>", &
244 n_var=2, default_i_vals=(/0, 0/))
248 CALL keyword_create(keyword, __location__, name=
"REAL_AXIS_INTEGRATION_METHOD", &
249 description=
"Integration method for the real axis.", &
250 usage=
"REAL_AXIS_INTEGRATION_METHOD <method>", &
252 enum_c_vals=
s2a(
"Gauss_Chebyshev",
"Trapezoidal_rule",
"Read"), &
253 enum_desc=
s2a(
"Gauss-Chebyshev integration between singularity points.", &
254 "Trapezoidal rule on the total range.", &
255 "Read integration points from a file (named E.dat)."), &
261 description=
"Number of integration points for the sigma solver on the complex contour.", &
262 usage=
"N_POINTS_INV <integer>", default_i_val=64)
266 CALL keyword_create(keyword, __location__, name=
"OBC_EQUILIBRIUM", &
267 description=
"Compute the equilibrium density with open boundary conditions.", &
268 default_l_val=.false., lone_keyword_l_val=.true.)
272 CALL keyword_create(keyword, __location__, name=
"CONTACT_FILLING", &
273 description=
"Determination of the contact Fermi levels. Note that this keyword "// &
274 "only works when the TRANSPORT_METHOD is specified as TRANSPORT.", &
276 enum_c_vals=
s2a(
"BAND_STRUCTURE",
"DOS"), &
277 enum_desc=
s2a(
"Determine the Fermi levels from the band structure.", &
278 "Determine the Fermi levels from the density of states."), &
283 CALL keyword_create(keyword, __location__, name=
"DENSITY_MIXING", &
284 description=
"Mixing parameter for a density mixing in OMEN.", &
285 usage=
"DENSITY_MIXING <real>", default_r_val=1.0_dp)
291 CALL create_contact_section(subsection)
295 CALL create_beyn_section(subsection)
299 CALL create_pexsi_section(subsection)
303 CALL create_transport_print_section(subsection)
313 SUBROUTINE create_contact_section(section)
318 cpassert(.NOT.
ASSOCIATED(section))
320 description=
"Parameters for defining device contacts.", &
321 n_keywords=5, n_subsections=0, repeats=.true.)
326 description=
"The number of neighboring unit cells that one unit cell interacts with.", &
327 usage=
"BANDWIDTH <integer>", default_i_val=0)
332 description=
"Index of the first atom in the contact unit cell. Set to 0 to define the contact "// &
333 "unit cell as the first/last N_ATOMS of the structure (after cutout)", &
334 usage=
"START <integer>", default_i_val=0)
339 description=
"Number of atoms in the contact unit cell.", &
340 usage=
"N_ATOMS <integer>", default_i_val=0)
344 CALL keyword_create(keyword, __location__, name=
"INJECTION_SIGN", &
345 description=
"Contact unit cell interacts with unit cells to the right (positive) or "// &
346 "to the left (negative).", &
347 usage=
"INJECTION_SIGN <integer>", &
349 enum_c_vals=
s2a(
"POSITIVE",
"NEGATIVE"), &
350 enum_desc=
s2a(
"When the contact unit cell is at the upper left corner of the Hamiltonian.", &
351 "When the contact unit cell is at the lower right corner of the Hamiltonian."), &
356 CALL keyword_create(keyword, __location__, name=
"INJECTING_CONTACT", &
357 description=
"whether or not the contact can inject electrons.", &
358 default_l_val=.true., lone_keyword_l_val=.true.)
362 END SUBROUTINE create_contact_section
368 SUBROUTINE create_beyn_section(section)
373 cpassert(.NOT.
ASSOCIATED(section))
375 description=
"Parameters for the Beyn eigenvalue solver.", &
376 n_keywords=6, n_subsections=0, repeats=.false.)
381 description=
"Number of random vectors as a fraction of the size of the unit cell.", &
382 usage=
"N_RAND <real>", default_r_val=1.0_dp)
387 description=
"Number of random vectors as a fraction of the size of the unit cell "// &
388 "for the complex contour.", &
389 usage=
"N_RAND_CC <real>", default_r_val=1.0_dp)
394 description=
"Cutoff for the singular values in the Beyn solver.", &
395 usage=
"SVD_CUTOFF <real>", default_r_val=1.0_dp)
399 CALL keyword_create(keyword, __location__, name=
"N_POINTS_BEYN", &
400 description=
"Number of integration points per circle in the Beyn solver.", &
401 usage=
"N_POINTS_BEYN <integer>", default_i_val=32)
406 description=
"Set to .TRUE. if only one circle instead of two should be used in the Beyn solver.", &
407 default_l_val=.false., lone_keyword_l_val=.true.)
411 CALL keyword_create(keyword, __location__, name=
"TASKS_PER_INTEGRATION_POINT", &
412 description=
"Number of tasks per integration point.", &
413 usage=
"TASKS_PER_INTEGRATION_POINT <integer>", default_i_val=0)
417 END SUBROUTINE create_beyn_section
423 SUBROUTINE create_pexsi_section(section)
428 cpassert(.NOT.
ASSOCIATED(section))
430 description=
"Parameters for the PEXSI solver to be used within OMEN.", &
431 n_keywords=4, n_subsections=0, repeats=.false., &
432 deprecation_notice=
"Support for the PEXSI library is slated for removal.")
437 description=
"Ordering strategy for factorization and selected inversion.", &
438 enum_c_vals=
s2a(
"PARALLEL",
"SEQUENTIAL",
"MULTIPLE_MINIMUM_DEGREE"), &
439 enum_desc=
s2a(
"Parallel ordering using ParMETIS/PT-SCOTCH (PARMETIS option in SuperLU_DIST)", &
440 "Sequential ordering using METIS (METIS_AT_PLUS_A option in SuperLU_DIST)", &
441 "Multiple minimum degree ordering (MMD_AT_PLUS_A option in SuperLU_DIST)"), &
442 enum_i_vals=(/0, 1, 2/), &
443 usage=
"ORDERING <integer>", default_i_val=1)
448 description=
"row permutation strategy for factorization and selected inversion.", &
449 enum_c_vals=
s2a(
"NOROWPERM",
"LARGEDIAG"), &
450 enum_desc=
s2a(
"No row permutation (NOROWPERM option in SuperLU_DIST)", &
451 "Make diagonal entry larger than off diagonal (LargeDiag option in SuperLU_DIST)"), &
452 enum_i_vals=(/0, 1/), &
453 usage=
"ROW_ORDERING <integer>", default_i_val=0)
458 description=
"The level of output information.", &
459 enum_c_vals=
s2a(
"SILENT",
"BASIC",
"DETAILED"), &
460 enum_i_vals=(/0, 1, 2/), &
461 usage=
"VERBOSITY <integer>", default_i_val=0)
466 description=
"Number of processors for PARMETIS/PT-SCOTCH. Only used if ORDERING is set to PARALLEL. "// &
467 "If 0, the number of processors for PARMETIS/PT-SCOTCH will be set equal to the number of "// &
468 "MPI ranks per pole. Note: if more than one processor is used, a segmentation fault may occur in the "// &
469 "symbolic factorization phase.", &
470 usage=
"NP_SYMB_FACT <integer>", default_i_val=1)
474 END SUBROUTINE create_pexsi_section
480 SUBROUTINE create_transport_print_section(section)
486 cpassert(.NOT.
ASSOCIATED(section))
488 description=
"Section of possible print options for transport calculations.", &
491 NULLIFY (keyword, print_key)
493 description=
"Controls the printing of current into cube files.", &
497 description=
"The stride (X,Y,Z) used to write the cube file "// &
498 "(larger values result in smaller cube files). You can provide 3 numbers (for X,Y,Z) or"// &
499 " 1 number valid for all components.", &
500 usage=
"STRIDE 2 2 2", n_var=-1, default_i_vals=(/2, 2, 2/), type_of_var=
integer_t)
504 description=
"append the cube files when they already exist", &
505 default_l_val=.false., lone_keyword_l_val=.true.)
512 END SUBROUTINE create_transport_print_section
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public bruck2014
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer, parameter, public high_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
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
Defines the basic variable types.
integer, parameter, public dp
Utilities for string manipulations.