(git:374b731)
Loading...
Searching...
No Matches
input_cp2k_transport.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 NEGF based quantum transport calculations
10!> (integration with the quantum transport code OMEN)
11!>
12!> \par History
13!> 07.2013 created [Hossein Bani-Hashemian]
14!> \author Hossein Bani-Hashemian
15! **************************************************************************************************
17 USE bibliography, ONLY: bruck2014
20 USE cp_units, ONLY: cp_unit_to_cp2k
21 USE input_constants, ONLY: &
36 USE input_val_types, ONLY: integer_t
37 USE kinds, ONLY: dp
38 USE string_utilities, ONLY: s2a
39#include "./base/base_uses.f90"
40
41 IMPLICIT NONE
42 PRIVATE
43
44 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_transport'
45
47
48CONTAINS
49
50! **************************************************************************************************
51!> \brief creates the TRANSPORT section
52!> \param[inout] section the section to be created
53! **************************************************************************************************
54 SUBROUTINE create_transport_section(section)
55 TYPE(section_type), POINTER :: section
56
57 TYPE(keyword_type), POINTER :: keyword
58 TYPE(section_type), POINTER :: subsection
59
60 cpassert(.NOT. ASSOCIATED(section))
61 CALL section_create(section, __location__, name="TRANSPORT", &
62 description="Specifies the parameters for transport, sets parameters for the OMEN code, "// &
63 "see also <https://nano-tcad.ee.ethz.ch>.", &
64 citations=(/bruck2014/), &
65 n_keywords=29, n_subsections=4, repeats=.false.)
66
67 NULLIFY (keyword)
68
69 CALL keyword_create( &
70 keyword, __location__, name="TRANSPORT_METHOD", &
71 description="Preferred method for transport calculations.", &
72 usage="TRANSPORT_METHOD <method>", &
73 default_i_val=transport_negf, &
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"), &
79 CALL section_add_keyword(section, keyword)
80 CALL keyword_release(keyword)
81
82 CALL keyword_create( &
83 keyword, __location__, name="QT_FORMALISM", &
84 description="Preferred quantum transport formalism to compute the current and density.", &
85 usage="QT_FORMALISM <method>", &
86 default_i_val=rho_qtbm, &
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."), &
90 enum_i_vals=(/rho_negf, rho_qtbm/))
91 CALL section_add_keyword(section, keyword)
92 CALL keyword_release(keyword)
93
94 CALL keyword_create(keyword, __location__, name="NUM_POLE", &
95 description="The number of terms in the PEXSI's pole expansion method.", &
96 usage="NUM_POLE <integer>", default_i_val=64)
97 CALL section_add_keyword(section, keyword)
98 CALL keyword_release(keyword)
99
100 CALL keyword_create(keyword, __location__, name="N_KPOINTS", &
101 description="The number of k points for determination of the singularities.", &
102 usage="N_KPOINTS <integer>", default_i_val=64)
103 CALL section_add_keyword(section, keyword)
104 CALL keyword_release(keyword)
105
106 CALL keyword_create(keyword, __location__, name="NUM_INTERVAL", &
107 description="Max number of energy points per small interval.", &
108 usage="NUM_INTERVAL <integer>", default_i_val=10)
109 CALL section_add_keyword(section, keyword)
110 CALL keyword_release(keyword)
111
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)
116 CALL section_add_keyword(section, keyword)
117 CALL keyword_release(keyword)
118
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)
123 CALL section_add_keyword(section, keyword)
124 CALL keyword_release(keyword)
125
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)
129 CALL section_add_keyword(section, keyword)
130 CALL keyword_release(keyword)
131
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)
135 CALL section_add_keyword(section, keyword)
136 CALL keyword_release(keyword)
137
138 CALL keyword_create(keyword, __location__, name="EPS_LIMIT", &
139 description="The smallest eigenvalue that is kept.", &
140 usage="EPS_LIMIT <real>", default_r_val=1.0e-4_dp)
141 CALL section_add_keyword(section, keyword)
142 CALL keyword_release(keyword)
143
144 CALL keyword_create(keyword, __location__, name="EPS_LIMIT_CC", &
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)
147 CALL section_add_keyword(section, keyword)
148 CALL keyword_release(keyword)
149
150 CALL keyword_create(keyword, __location__, name="EPS_DECAY", &
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)
154 CALL section_add_keyword(section, keyword)
155 CALL keyword_release(keyword)
156
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)
160 CALL section_add_keyword(section, keyword)
161 CALL keyword_release(keyword)
162
163 CALL keyword_create(keyword, __location__, name="EPS_MU", &
164 description="Accuracy to which the Fermi level should be determined.", &
165 usage="EPS_MU <real>", default_r_val=1.0e-6_dp)
166 CALL section_add_keyword(section, keyword)
167 CALL keyword_release(keyword)
168
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)
172 CALL section_add_keyword(section, keyword)
173 CALL keyword_release(keyword)
174
175 CALL keyword_create(keyword, __location__, name="EPS_FERMI", &
176 description="Cutoff for the tail of the Fermi function.", &
177 usage="EPS_FERMI <real>", default_r_val=0.0_dp)
178 CALL section_add_keyword(section, keyword)
179 CALL keyword_release(keyword)
180
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)
184 CALL section_add_keyword(section, keyword)
185 CALL keyword_release(keyword)
186
187 CALL keyword_create(keyword, __location__, name="MIN_INTERVAL", &
188 description="Smallest enery distance in energy vector.", &
189 usage="MIN_INTERVAL <real>", default_r_val=1.0e-4_dp)
190 CALL section_add_keyword(section, keyword)
191 CALL keyword_release(keyword)
192
193 CALL keyword_create(keyword, __location__, name="TEMPERATURE", &
194 description="Temperature.", &
195 usage="TEMPERATURE [K] 300.0", &
196 default_r_val=cp_unit_to_cp2k(value=300.0_dp, unit_str="K"), &
197 unit_str="K")
198 CALL section_add_keyword(section, keyword)
199 CALL keyword_release(keyword)
200
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.)
204 CALL section_add_keyword(section, keyword)
205 CALL keyword_release(keyword)
206
207 CALL keyword_create( &
208 keyword, __location__, name="LINEAR_SOLVER", &
209 description="Preferred solver for solving the linear system of equations.", &
210 usage="LINEAR_SOLVER <solver>", &
211 default_i_val=linsolver_full, &
212 enum_c_vals=s2a("SplitSolve", "SuperLU", "MUMPS", "Full", "Banded", "PARDISO", "UMFPACK"), &
215 CALL section_add_keyword(section, keyword)
216 CALL keyword_release(keyword)
217
218 CALL keyword_create( &
219 keyword, __location__, name="MATRIX_INVERSION_METHOD", &
220 description="Preferred matrix inversion method.", &
221 usage="MATRIX_INVERSION_METHOD <solver>", &
222 default_i_val=matrixinv_full, &
223 enum_c_vals=s2a("Full", "PEXSI", "PARDISO", "RGF"), &
225 CALL section_add_keyword(section, keyword)
226 CALL keyword_release(keyword)
227
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>", &
231 default_i_val=injmethod_beyn, &
232 enum_c_vals=s2a("EVP", "BEYN"), &
233 enum_desc=s2a("Full eigenvalue solver.", &
234 "Beyn eigenvalue solver."), &
235 enum_i_vals=(/injmethod_evp, injmethod_beyn/))
236 CALL section_add_keyword(section, keyword)
237 CALL keyword_release(keyword)
238
239 CALL keyword_create( &
240 keyword, __location__, name="CUTOUT", &
241 description="The number of atoms at the beginning and the end of the structure where the density should "// &
242 "not be changed.", &
243 usage="CUTOUT <integer> <integer>", &
244 n_var=2, default_i_vals=(/0, 0/))
245 CALL section_add_keyword(section, keyword)
246 CALL keyword_release(keyword)
247
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>", &
251 default_i_val=rlaxisint_gausschebyshev, &
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)."), &
257 CALL section_add_keyword(section, keyword)
258 CALL keyword_release(keyword)
259
260 CALL keyword_create(keyword, __location__, name="N_POINTS_INV", &
261 description="Number of integration points for the sigma solver on the complex contour.", &
262 usage="N_POINTS_INV <integer>", default_i_val=64)
263 CALL section_add_keyword(section, keyword)
264 CALL keyword_release(keyword)
265
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.)
269 CALL section_add_keyword(section, keyword)
270 CALL keyword_release(keyword)
271
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.", &
275 default_i_val=neutlead_bs, &
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."), &
279 enum_i_vals=(/neutlead_bs, neutlead_dos/))
280 CALL section_add_keyword(section, keyword)
281 CALL keyword_release(keyword)
282
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)
286 CALL section_add_keyword(section, keyword)
287 CALL keyword_release(keyword)
288
289 NULLIFY (subsection)
290
291 CALL create_contact_section(subsection)
292 CALL section_add_subsection(section, subsection)
293 CALL section_release(subsection)
294
295 CALL create_beyn_section(subsection)
296 CALL section_add_subsection(section, subsection)
297 CALL section_release(subsection)
298
299 CALL create_pexsi_section(subsection)
300 CALL section_add_subsection(section, subsection)
301 CALL section_release(subsection)
302
303 CALL create_transport_print_section(subsection)
304 CALL section_add_subsection(section, subsection)
305 CALL section_release(subsection)
306
307 END SUBROUTINE create_transport_section
308
309! **************************************************************************************************
310!> \brief Creates the section for creating contacts.
311!> \param[inout] section the section to be created
312! **************************************************************************************************
313 SUBROUTINE create_contact_section(section)
314 TYPE(section_type), POINTER :: section
315
316 TYPE(keyword_type), POINTER :: keyword
317
318 cpassert(.NOT. ASSOCIATED(section))
319 CALL section_create(section, __location__, name="CONTACT", &
320 description="Parameters for defining device contacts.", &
321 n_keywords=5, n_subsections=0, repeats=.true.)
322
323 NULLIFY (keyword)
324
325 CALL keyword_create(keyword, __location__, name="BANDWIDTH", &
326 description="The number of neighboring unit cells that one unit cell interacts with.", &
327 usage="BANDWIDTH <integer>", default_i_val=0)
328 CALL section_add_keyword(section, keyword)
329 CALL keyword_release(keyword)
330
331 CALL keyword_create(keyword, __location__, name="START", &
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)
335 CALL section_add_keyword(section, keyword)
336 CALL keyword_release(keyword)
337
338 CALL keyword_create(keyword, __location__, name="N_ATOMS", &
339 description="Number of atoms in the contact unit cell.", &
340 usage="N_ATOMS <integer>", default_i_val=0)
341 CALL section_add_keyword(section, keyword)
342 CALL keyword_release(keyword)
343
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>", &
348 default_i_val=injsign_positive, &
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."), &
352 enum_i_vals=(/injsign_positive, injsign_negative/))
353 CALL section_add_keyword(section, keyword)
354 CALL keyword_release(keyword)
355
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.)
359 CALL section_add_keyword(section, keyword)
360 CALL keyword_release(keyword)
361
362 END SUBROUTINE create_contact_section
363
364! **************************************************************************************************
365!> \brief Creates the section for the Beyn eigenvalue solver.
366!> \param[inout] section the section to be created
367! **************************************************************************************************
368 SUBROUTINE create_beyn_section(section)
369 TYPE(section_type), POINTER :: section
370
371 TYPE(keyword_type), POINTER :: keyword
372
373 cpassert(.NOT. ASSOCIATED(section))
374 CALL section_create(section, __location__, name="BEYN", &
375 description="Parameters for the Beyn eigenvalue solver.", &
376 n_keywords=6, n_subsections=0, repeats=.false.)
377
378 NULLIFY (keyword)
379
380 CALL keyword_create(keyword, __location__, name="N_RAND", &
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)
383 CALL section_add_keyword(section, keyword)
384 CALL keyword_release(keyword)
385
386 CALL keyword_create(keyword, __location__, name="N_RAND_CC", &
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)
390 CALL section_add_keyword(section, keyword)
391 CALL keyword_release(keyword)
392
393 CALL keyword_create(keyword, __location__, name="SVD_CUTOFF", &
394 description="Cutoff for the singular values in the Beyn solver.", &
395 usage="SVD_CUTOFF <real>", default_r_val=1.0_dp)
396 CALL section_add_keyword(section, keyword)
397 CALL keyword_release(keyword)
398
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)
402 CALL section_add_keyword(section, keyword)
403 CALL keyword_release(keyword)
404
405 CALL keyword_create(keyword, __location__, name="ONE_CIRCLE", &
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.)
408 CALL section_add_keyword(section, keyword)
409 CALL keyword_release(keyword)
410
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)
414 CALL section_add_keyword(section, keyword)
415 CALL keyword_release(keyword)
416
417 END SUBROUTINE create_beyn_section
418
419! **************************************************************************************************
420!> \brief Creates the section for the PEXSI solver.
421!> \param[inout] section the section to be created
422! **************************************************************************************************
423 SUBROUTINE create_pexsi_section(section)
424 TYPE(section_type), POINTER :: section
425
426 TYPE(keyword_type), POINTER :: keyword
427
428 cpassert(.NOT. ASSOCIATED(section))
429 CALL section_create(section, __location__, name="PEXSI", &
430 description="Parameters for the PEXSI solver to be used within OMEN.", &
431 n_keywords=4, n_subsections=0, repeats=.false.)
432
433 NULLIFY (keyword)
434
435 CALL keyword_create(keyword, __location__, name="ORDERING", &
436 description="Ordering strategy for factorization and selected inversion.", &
437 enum_c_vals=s2a("PARALLEL", "SEQUENTIAL", "MULTIPLE_MINIMUM_DEGREE"), &
438 enum_desc=s2a("Parallel ordering using ParMETIS/PT-SCOTCH (PARMETIS option in SuperLU_DIST)", &
439 "Sequential ordering using METIS (METIS_AT_PLUS_A option in SuperLU_DIST)", &
440 "Multiple minimum degree ordering (MMD_AT_PLUS_A option in SuperLU_DIST)"), &
441 enum_i_vals=(/0, 1, 2/), &
442 usage="ORDERING <integer>", default_i_val=1)
443 CALL section_add_keyword(section, keyword)
444 CALL keyword_release(keyword)
445
446 CALL keyword_create(keyword, __location__, name="ROW_ORDERING", &
447 description="row permutation strategy for factorization and selected inversion.", &
448 enum_c_vals=s2a("NOROWPERM", "LARGEDIAG"), &
449 enum_desc=s2a("No row permutation (NOROWPERM option in SuperLU_DIST)", &
450 "Make diagonal entry larger than off diagonal (LargeDiag option in SuperLU_DIST)"), &
451 enum_i_vals=(/0, 1/), &
452 usage="ROW_ORDERING <integer>", default_i_val=0)
453 CALL section_add_keyword(section, keyword)
454 CALL keyword_release(keyword)
455
456 CALL keyword_create(keyword, __location__, name="VERBOSITY", &
457 description="The level of output information.", &
458 enum_c_vals=s2a("SILENT", "BASIC", "DETAILED"), &
459 enum_i_vals=(/0, 1, 2/), &
460 usage="VERBOSITY <integer>", default_i_val=0)
461 CALL section_add_keyword(section, keyword)
462 CALL keyword_release(keyword)
463
464 CALL keyword_create(keyword, __location__, name="NP_SYMB_FACT", &
465 description="Number of processors for PARMETIS/PT-SCOTCH. Only used if ORDERING is set to PARALLEL. "// &
466 "If 0, the number of processors for PARMETIS/PT-SCOTCH will be set equal to the number of "// &
467 "MPI ranks per pole. Note: if more than one processor is used, a segmentation fault may occur in the "// &
468 "symbolic factorization phase.", &
469 usage="NP_SYMB_FACT <integer>", default_i_val=1)
470 CALL section_add_keyword(section, keyword)
471 CALL keyword_release(keyword)
472
473 END SUBROUTINE create_pexsi_section
474
475!**************************************************************************************************
476!> \brief Creates print section for transport calculations.
477!> \param[inout] section the section to be created
478! **************************************************************************************************
479 SUBROUTINE create_transport_print_section(section)
480 TYPE(section_type), POINTER :: section
481
482 TYPE(keyword_type), POINTER :: keyword
483 TYPE(section_type), POINTER :: print_key
484
485 cpassert(.NOT. ASSOCIATED(section))
486 CALL section_create(section, __location__, name="PRINT", &
487 description="Section of possible print options for transport calculations.", &
488 repeats=.false.)
489
490 NULLIFY (keyword, print_key)
491 CALL cp_print_key_section_create(print_key, __location__, "CURRENT", &
492 description="Controls the printing of current into cube files.", &
493 print_level=high_print_level, filename="current")
494
495 CALL keyword_create(keyword, __location__, name="stride", &
496 description="The stride (X,Y,Z) used to write the cube file "// &
497 "(larger values result in smaller cube files). You can provide 3 numbers (for X,Y,Z) or"// &
498 " 1 number valid for all components.", &
499 usage="STRIDE 2 2 2", n_var=-1, default_i_vals=(/2, 2, 2/), type_of_var=integer_t)
500 CALL section_add_keyword(print_key, keyword)
501 CALL keyword_release(keyword)
502 CALL keyword_create(keyword, __location__, name="APPEND", &
503 description="append the cube files when they already exist", &
504 default_l_val=.false., lone_keyword_l_val=.true.)
505 CALL section_add_keyword(print_key, keyword)
506 CALL keyword_release(keyword)
507
508 CALL section_add_subsection(section, print_key)
509 CALL section_release(print_key)
510
511 END SUBROUTINE create_transport_print_section
512
513END MODULE input_cp2k_transport
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
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 transport_negf
integer, parameter, public injmethod_beyn
integer, parameter, public matrixinv_full
integer, parameter, public linsolver_full
integer, parameter, public linsolver_umfpack
integer, parameter, public injmethod_evp
integer, parameter, public linsolver_superlu
integer, parameter, public neutlead_bs
integer, parameter, public matrixinv_rgf
integer, parameter, public linsolver_mumps
integer, parameter, public rho_qtbm
integer, parameter, public transport_localscf
integer, parameter, public matrixinv_pardiso
integer, parameter, public injsign_positive
integer, parameter, public neutlead_dos
integer, parameter, public injsign_negative
integer, parameter, public rlaxisint_gausschebyshev
integer, parameter, public matrixinv_pexsi
integer, parameter, public linsolver_splitsolve
integer, parameter, public linsolver_pardiso
integer, parameter, public rho_negf
integer, parameter, public linsolver_banded
integer, parameter, public rlaxisint_readfromfile
integer, parameter, public transport_transmission
integer, parameter, public rlaxisint_trapezoidal
input section for NEGF based quantum transport calculations (integration with the quantum transport c...
subroutine, public create_transport_section(section)
creates the TRANSPORT section
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 integer_t
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Utilities for string manipulations.
represent a keyword in the input
represent a section of the input file