Skip to content

Epsilon code input keywords (epsilon.inp)

Required keywords

Optional keywords

Keyword documentation

Parameters for full-frequency-dependent calculations

Three methods are available, depending on the value of frequency_dependence_method:


frequency_dependence_method [integer]

Full frequency dependence method for the polarizability, if frequency_dependence==2:

  • 0: Real-axis formalism, Adler-Wiser formula.
  • 1: Real-axis formalism, spectral method (PRB 74, 035101, (2006))
  • 2: Contour-deformation formalism, Adler-Wiser formula.

broadening [float]

Broadening parameter for each val->cond transition. It should correspond to the energy resolution due to k-point sampling, or a number as small as possible if you have a molecule. The default value is 0.1 for method 0 and 1, and 0.25 for method 2.


init_frequency [float]

Lower bound for the linear frequency grid.


frequency_low_cutoff [float]

Upper bound for the linear frequency grid. For methods 0 and 1, it is also the lower bound for the non-uniform frequency grid. Should be larger than the maximum transition, i.e., the energy difference between the highest conduction band and the lowest valence band. The default is 200 eV for method 0 and 1, and 10 eV for method 2. For method 2, you can increase frequency_low_cutoff if you wish to use the Sigma code and look into QP states deep in occupied manifold or high in the unoccupied manifold.


delta_frequency [float]

Frequency step for full-frequency integration for the linear grid Should be converged (the smaller, the better). For molecules, delta_frequency should be the same as broadening. Defaults to the value of broadening.


number_imaginary_freqs [integer]

Number of frequencies on the imaginary axis for method 2.


frequency_high_cutoff [float]

Upper limit of the non-uniform frequency grid. Defaults to 4*frequency_low_cutoff


delta_frequency_step [float]

Increase in the frequency step for the non-uniform frequency grid.


delta_sfrequency [float]

Frequency step for the linear grid for the spectral function method. Defaults to delta_frequency.


delta_sfrequency_step [float]

Increase in frequency step for the non-uniform grid for the spectral function method. Defaults to delta_frequency_step


sfrequency_low_cutoff [float]

Upper bound for the linear grid for the spectral function method and lower bound for the non-uniform frequency grid. Defaults to frequency_low_cutoff


sfrequency_high_cutoff [float]

Upper limit of the non-uniform frequency grid for the spectral function method. Defaults to frequency_low_cutoff


add_w0_freq

In case the frequency grid doesn't start at zero (non default value for init_frequency), it allows to add the zero frequency as the grid first value, second value being then the one specified by init_frequency.

Static subspace approximation

Input parameters controlling the full-frequency static subspace approximation method. The method speeds up full-frequency calculations by expanding the frequency-dependent part of the polarizability in the subspace basis formed by the lowest eigenvectors of the static polarizability (see Phys. Rev. B 99, 125128, 2019). The method is implemented for the case of full-frequency contour-deformation formalism, i.e., frequency_dependence==2 and frequency_dependence_method==2.


chi_eigenvalue_cutoff [float]

Activate the static subspace approximation and define the screening threshold for the eigenvalues of the static polarizability (alias eps_trunc_eigen).


nbasis_subspace [integer]

Define a maximum fixed number of eigenvectors to be used (recommended). Set for example to 25% of the number of G-vectors employed for the expansion of chi0.


write_subspace_epsinv

Write frequency-dependent epsilon matrices in eps0mat[.h5] and epsmat[.h5] files using the subspace basis instead of the full G-vector basis (recommended). This flag need to be specified to use the full-frequency static subspace approximation in sigma.


subspace_dont_keep_full_eps_omega0

Discharge the output of the full static epsilon matrix in the eps0mat[.h5] and epsmat[.h5] files.


subspace_use_elpa

If ELPA library is linked, specify to use ELPA for the diagonalization of the static polarizability.


dont_keep_fix_buffers

Using comm_nonblocking_cyclic force the algorithm to use buffers for communication with variable size, mainly used for debugging.


sub_collective_eigen_redistr

Replicate eigenvectors using collective operation in the basis transformation step, performs in general worse than point to point, mainly used for debugging.


cluster_size_subspace_diagonalization [integer]

Allows to control the cluster size when performing the diagonalization with ScaLAPACK. In case there are too many eigenvalues clustered together and the diagonalization fails, try increasing this parameter (default is 250).

GPUs Parameters

Input parameters controlling the execution when running with GPU support. For the default algorithm (keys ending with _algo) the hierarchy depends on which programming models are compiled in (OpenACC, OpenMP-target and/or CPU (threaded or not)): if all possible implementations have been compiled then the hierarchy is OPENACC_ALGO > OMP_TARGET_ALGO > CPU_ALGO.


mtxel_algo

Controlling which algorithm to use for the GPU offload of the mtxel kernel calculating transition matrix elements. Possible values are CPU_ALGO, OPENACC_ALGO and OMP_TARGET_ALGO (See section header for default values).


chi_summation_algo

Controlling which algorithm to use for the GPU offload of the chi_summation kernel calculating the polarizability. Possible values are CPU_ALGO, OPENACC_ALGO and OMP_TARGET_ALGO (See section header for default values).


n_ffts_per_batch [integer]

Control how many FFTs to perform simultaneously in the batched mtxel kernel (Default 20). If set to -1 the band by band mtxel kernel will be used, see accel_mtxel_band_block_size.


accel_mtxel_band_block_size [integer]

Control the inner block size in the band by band mtxel kernel, the larger the better, but uses more memory (Default 20).


accel_chi_summation_no_full_offload

Avoid offloading the whole gme matrix elements to GPU to accelerate the preparation step in chi_summation kernel. Generally it slows down considerably the execution but uses less memory and it's useful in case the prep execution fails. To better control the memory usage also have a look to max_mem_nv_block_algo.

Scissors operator

Scissors operator (linear fit of the quasiparticle energy corrections) for the bands in WFN and WFNq. For valence-band energies:

  • ev_cor = ev_in + evs + evdel * (ev_in - ev0)

For conduction-band energies:

  • ec_cor = ec_in + ecs + ecdel * (ec_in - ec0)

Defaults is zero for all entries, i.e., no scissors corrections. evs, ev0, ecs, ec0 are in eV. If you have eqp.dat and eqp_q.dat files this information is ignored in favor of the eigenvalues in eqp.dat and eqp_q.dat. One can specify all parameters for scissors operator in a single line with cvfit evs ev0 evdel ecs ec0 ecdel


evs [float]

ev0 [float]

evdel [float]

ecs [float]

ec0 [float]

ecdel [float]

cvfit [array of integers]

Truncation schemes for the Coulomb potential

Since BerkerleyGW is a plane-wave-based code, one must truncate the Coulomb potential to avoid spurious interactions between repeated supercells when dealing with systems with reduced dimensionality. Make sure you understand how to setup your mean-field calculation so that the supercell is large enough to perform a truncation of the Coulomb potential.


cell_box_truncation

Truncate the Coulomb potential based on the Wigner-Seitz cell. This is the recommended truncation for 0D systems.


cell_wire_truncation

Truncation scheme for 1D systems, such as carbon nanotubes. The z direction is assumed to be periodic, and x and y confined.


cell_slab_truncation

Truncation scheme for 2D systems, such as graphene or monolayer MoS2. The z direction is assumed to be confined, and x and y periodic.


spherical_truncation

Truncate the Coulomb potential based on an analytical scheme. This is ok for quasi-spherical systems, such as CH4 molecule or C60, but the cell_box_truncation is the most general and recommended scheme. When using spherical truncation, you must also specify the radius for the truncation in spherical_truncation.


coulomb_truncation_radius [float]

This specifies the radius of for spherical truncation, in Bohr, so that the Coulomb potential v(r) is zero for r larger than these values. This flag is to be used together with spherical_truncation.

Misc. parameters

epsilon_cutoff [float]

Energy cutoff for the dielectric matrix, in Ry. The dielectric matrix \varepsilon_{GG'} will contain all G-vectors with kinetic energy |q+G|^2 up to this cutoff.

number_bands [integer]

Total number of bands (valence+conduction) to sum over. Defaults to the number of bands in the WFN file minus 1.

use_wfn_hdf5

Read input wavefunction files in HDF5 format: WFN.h5 instead of WFN and, if required, WFNq.h5 instead of WFNq.

frequency_dependence [integer]

This flags specifies the frequency dependence of the inverse dielectric matrix:

  • Set to 0 to compute the static inverse dielectric matrix (default).
  • Set to 2 to compute the full frequency dependent inverse dielectric matrix.
  • Set to 3 to compute the two frequencies needed for Godby-Needs GPP model.

plasma_freq [float]

Plasma frequency (eV) needed for the contour-deformation method (i.e., frequency_dependence==2). The exact value is unimportant, especially if you have enough imaginary frequency points. We recommend you keep this value fixed at 2 Ry.

imaginary_frequency [float]

For frequency_dependence==3, the value of the purely imaginary frequency, in eV:

comm_nonblocking_cyclic

Within gcomm_matrix, employs a non-blocking cyclic communication scheme overlapping computation and communication in the evaluation of the polarizability (drastically reduce the time spent in communication for large runs, but require more memory).

full_chi_conv_log [integer]

Logging convergence of the head & tail of polarizability matrix with respect to conduction bands:

  • Set to -1 for no convergence test
  • Set to 0 for the 5 column format including the extrapolated values (default).
  • Set to 1 for the 2 column format, real part only.
  • Set to 2 for the 2 column format, real and imaginary parts.

begin qpoints ... end

qx qy qz 1/scale_factor is_q0

scale_factor is for specifying values such as 1/3 is_q0 = 0 for regular, non-zero q-vectors (read val WFNs from WFN) is_q0 = 1 for a small q-vector in semiconductors (read val WFNs from WFNq) is_q0 = 2 for a small q-vector in metals (read val WFNs from WFN) if present the small q-vector should be first in the list You can generate this list with kgrid.x: just set the shifts to zero and use same grid numbers as for WFN. Then replace the zero vector with q0. In the case of RPA energy calculations (do_rpa) an additional value at the end of each q-point line has to be specified with the weight of the q-point for the BZ integration, do_rpa flag needs to be specified before the begin qpoints block.

fermi_level [float]

Specify the Fermi level (in eV), if you want implicit doping Note that value refers to energies after scissor shift or eqp corrections. See also fermi_level_absolute and fermi_level_relative to control the meaning of the Fermi level.

The Fermi level in keyword fermi_level can be treated as an absolute value or relative to that found from the mean field (default). Using fermi_level_absolute will force the code to recompute the Fermi level regardless of existing occupations/ifmax array.


fermi_level_absolute

fermi_level_relative

dont_use_hdf5

Read from traditional simple binary format for epsmat/eps0mat instead of HDF5 file format. Relevant only if code is compiled with HDF5 support.

verbosity [integer]

Verbosity level, options are:

  • 1: default
  • 2: medium - info about k-points, symmetries, and eqp corrections.
  • 3: high - full dump of the reduced and unfolded k-points.
  • 4: log - log of various function calls. Use to debug code.
  • 5: debug - extra debug statements. Use to debug code.
  • 6: max - only use if instructed to, severe performance downgrade. Note that verbosity levels are cumulative. Most users will want to stick with level 1 and, at most, level 3. Only use level 4+ if debugging the code.

eqp_corrections

Set this to use eigenvalues in eqp.dat and eqp_q.dat If not set, these files will be ignored.

write_vcoul

Write the bare Coulomb potential v(q+G) to file

Matrix Element Communication Method (Chi Sum Comm). Default is gcomm_matrix which is good if nk, nc, nv > nmtx, nfreq. If nk, nc, nv < nfreq, nmtx (nk, nv < nfreq since nc \sim nmtx), use gcomm_elements. Only gcomm_elements is supported with the spectral method.


gcomm_matrix

gcomm_elements

number_valence_pools [integer]

Number of pools for distribution of valence bands The default is chosen to minimize memory in calculation

By default, the code computes the polarizability matrix, constructs the dielectric matrix, inverts it and writes the result to file epsmat. Use keyword skip_epsilon to compute the polarizability matrix and write it to file chimat. Use keyword skip_chi to read the polarizability matrix from file chimat, construct the dielectric matrix, invert it and write the result to file epsmat.


skip_epsilon

skip_chi

stop_after_qpt [integer]

Allows to stop the epsilon calculation after the specified number of q-points are calculated (default -1, meaning all q-points calculated).

read_chi_add

Read in screening files (chimat and chi0mat) for system outside of the electronic system in your current calculation to add to the electronic polarizability you calculate. For example, you could include the fluid polarizability from the JDFTx computational package, or the electronic polarizabilities for the substrate for your electronic system that you calculated previously. Works only with hdf5 and you have to have the following files in your directory: chi0mat_add.h5 and chimat_add.h5 (only the latter if you are working on a metallic system).

chi_add_fact

If reading in additional screening files, you can scale the read-in polarizability by this factor. Default is 1.0.

read_chi_maxmem_per_block [float]

When using skip_chi allows to control the max memory (in Mb) to use when reading in block of columns of chimat in h5 format.

nv_block [integer]

Input parameters to activate the valence bands blocked algorithm. Particularly useful in case of memory issues, the entire epsilon computation is broken into blocks of valence bands, and block size can be given in input such that memory limit is not exceeded in the allocation of the matrix elements array (pol%gme). You can either specify the block size (nv_block) or the max memory (in Gb) available per MPI task (max_mem_nv_block_algo).


max_mem_nv_block_algo [float]

nfreq_group [integer]

(Full Frequency only) Calculates several frequencies in parallel. No "new" processors are used here, the chi summation is simply done in another order way to decrease communication. This also allows the inversion of multiple dielectric matrices simultaneously via ScaLAPACK, circumventing ScaLAPACK's scaling problems. Can be very efficient for system with lots of G-vectors and when you have many frequencies. In general gives speedup. However, in order to calculate N frequencies in parallel, the memory to store pol%gme is currently multiplied by N as well.

partial_occ_degeneracy_tolerance [float]

Allows to define the threshold for determine the degeneracy between partially occupied valence and conduction eigenvalues in the computation of transition matrix elements (default 1.0E-6).

do_rpa

Enable the calculation of the RPA energy. IMPORTANT NOTES: (1) Works only in combination with Contour-deformation formalism with or without subspace approximation. (2) The flag has to be specified before the begin qpoints block, additionally an extra value at the end of each q-point line has to be specified with the weight of the q-point for the BZ integration. (3) The imaginary frequency grid will be different than the one used in the standard GW Contour-deformation calculation so if the eps0mat[.h5] and epsmat[.h5] files will be used for a subsequent sigma calculation the calculated eqp will be different.

'unfolded BZ' is from the kpoints in the WFN file 'full BZ' is generated from the kgrid parameters in the WFN file See comments in Common/checkbz.f90 for more details


fullbz_replace

Replace unfolded BZ with full BZ


fullbz_write

Write unfolded BZ and full BZ to files

degeneracy_check_override

The requested number of bands cannot break degenerate subspace Use the following keyword to suppress this check Note that you must still provide one more band in wavefunction file in order to assess degeneracy Automatically enabled if pseudobands are detected.

no_min_fftgrid

Instead of using the RHO FFT box to perform convolutions, we automatically determine (and use) the smallest box that is compatible with your epsilon cutoff. This also reduces the amount of memory needed for the FFTs. Although this optimization is safe, you can disable it by uncommenting the following line:

restart

Use this flag if you would like to restart your Epsilon calculation instead of starting it from scratch. Note that we can only reuse q-points that were fully calculated. This flag is ignored unless you are running the code with HDF5.

qgrid [array of integers]

Q-grid for the epsmat file. Defaults to the WFN k-grid.

subsample

Use this option to generate an eps0mat file suitable for subsampled Sigma calculations. The only thing this flag does is to allow a calculation with more than one q\rightarrow0 points without complaining.

dont_check_norms

Whether we want to check WFN norms. Automatically enabled if pseudobands are detected.

occ_broadening [float]

Apply the first order Methfessel-Paxton smearing scheme to the band occupations. Note that this keyword must be used with fermi_level_absolute.