pymatgen.io.cp2k.inputs module

This module defines the building blocks of a CP2K input file. The cp2k input structure is essentially a collection of “sections” which are similar to dictionary objects that activate modules of the cp2k executable, and then “keywords” which adjust variables inside of those modules. For example, FORCE_EVAL section will activate CP2K’s ability to calculate forces, and inside FORCE_EVAL, the Keyword “METHOD can be set to “QS” to set the method of force evaluation to be the quickstep (DFT) module.

A quick overview of the module:

– Section class defines the basis of Cp2k input and contains methods for manipulating these objects similarly to Dicts. – Keyword class defines the keywords used inside of Section objects that changes variables in Cp2k program – Cp2kInput class is special instantiation of Section that is used to represent the full cp2k calculation input. – The rest of the classes are children of Section intended to make initialization of common sections easier.

class BrokenSymmetry(l_alpha: int = - 1, n_alpha: int = 0, nel_alpha: int = - 1, l_beta: int = - 1, n_beta: int = 0, nel_beta: int = - 1)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Define the required atomic orbital occupation assigned in initialization of the density matrix, by adding or subtracting electrons from specific angular momentum channels. It works only with GUESS ATOMIC

Initialize the broken symmetry section

Parameters
  • l_alpha – Angular momentum quantum number of the orbitals whose occupation is changed

  • n_alpha – Principal quantum number of the orbitals whose occupation is changed. Default is the first not occupied

  • nel_alpha – Orbital occupation change per angular momentum quantum number. In unrestricted calculations applied to spin alpha

  • l_beta – Same as L_alpha for beta channel

  • n_beta – Same as N_alpha for beta channel

  • nel_beta – Same as NEL_alpha for beta channel

class Cell(lattice: pymatgen.core.lattice.Lattice, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Defines the simulation cell (lattice)

Initialize the cell section.

Parameters

lattice – pymatgen lattice object

class Coord(structure: Union[pymatgen.core.structure.Structure, pymatgen.core.structure.Molecule], aliases: Optional[dict] = None, subsections: dict = None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Specifies the coordinates of the atoms using a pymatgen structure object.

Parameters
  • structure – Pymatgen structure object

  • alias (bool) – whether or not to identify the sites by Element + number so you can do things like assign unique magnetization do different elements.

class Cp2kInput(name: str = 'CP2K_INPUT', subsections: dict = None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Special instance of ‘Section’ class that is meant to represent the overall cp2k input. Distinguishes itself from Section by overriding get_string() to not print this section’s title and by implementing the file i/o

Initialize Cp2kInput by calling the super

static from_file(file: str)[source]

Initialize from a file

classmethod from_lines(lines: Union[List, tuple])[source]

Helper method to read lines of file

static from_string(s: str)[source]

Initialize from a string

get_string()[source]

Get string representation of the Cp2kInput

write_file(input_filename: str = 'cp2k.inp', output_dir: str = '.', make_dir_if_not_present: bool = True)[source]

Write input to a file.

class Davidson(new_prec_each: int = 20, preconditioner: str = 'FULL_SINGLE_INVERSE', **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Parameters for davidson diagonalization

Parameters
  • new_prec_each (int) – How often to recalculate the preconditioner

  • preconditioner (str) – Preconditioner to use

class Dft(basis_set_filenames='BASIS_MOLOPT', potential_filename='GTH_POTENTIALS', uks: bool = True, wfn_restart_file_name=None, subsections: dict = None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls the DFT parameters in Cp2k

Initialize the DFT section

Parameters
  • subsections – Any subsections to initialize with

  • basis_set_filename – Name of the file that contains the basis set information

  • potential_filename – Name of the file that contains the pseudopotential information

  • uks – Whether to run unrestricted Kohn Sham (spin polarized)

class DftPlusU(eps_u_ramping=1e-05, init_u_ramping_each_scf=False, l=- 1, u_minus_j=0, u_ramping=0)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls DFT+U for an atom kind

Initialize the DftPlusU section.

Parameters
  • eps_u_ramping – (float) SCF convergence threshold at which to start ramping the U value

  • init_u_ramping_each_scf – (bool) Whether or not to do u_ramping each scf cycle

  • l – (int) angular moment of the orbital to apply the +U correction

  • u_minus_j – (float) the effective U parameter, Ueff = U-J

  • u_ramping – (float) stepwise amount to increase during ramping until u_minus_j is reached

class Diagonalization(eps_adapt: float = 0, eps_iter: float = 1e-08, eps_jacobi: float = 0, jacobi_threshold: float = 1e-07, subsections: dict = None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls diagonalization settings (if using traditional diagonalization).

Initialize the diagronalization section

class E_Density_Cube(**kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls printing of the electron density cube file

Initialize the E_DENSITY_CUBE Section

class ForceEval(subsections: dict = None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls the calculation of energy and forces in Cp2k

Initialize the ForceEval section

class Global(project_name: str = 'CP2K', run_type: str = 'ENERGY_FORCE', **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls ‘global’ settings for cp2k execution such as RUN_TYPE and PROJECT_NAME

Initialize the global section

class Keyword(name: str, *values, description: str = None, units: str = None, verbose: bool = True, repeats: bool = False)[source]

Bases: monty.json.MSONable

Class representing a keyword argument in CP2K. Within CP2K Secitons, which activate features of the CP2K code, the keywords are arguments that control the functionality of that feature. For example, the section “FORCE_EVAL” activates the evaluation of forces/energies, but within “FORCE_EVAL” the keyword “METHOD” controls whether or not this will be done with, say, “Quickstep” (DFT) or “EIP” (empirical interatomic potential).

Initializes a keyword. These Keywords and the value passed to them are sometimes as simple as KEYWORD VALUE, but can also be more elaborate such as KEYWORD [UNITS] VALUE1 VALUE2, which is why this class exists: to handle many values and control easy printing to an input file.

Parameters
  • name (str) – The name of this keyword. Must match an acceptable keyword from CP2K

  • args – All non-keyword arguments after ‘name’ are interpreted as the values to set for this keyword. i.e: KEYWORD ARG1 ARG2 would provide two values to the keyword.

  • description (str) – The description for this keyword. This can make readability of input files easier for some. Default=None.

  • units (str) – The units for this keyword. If not specified, CP2K default units will be used. Consult manual for default units. Default=None.

  • repeats (bool) – Whether or not this keyword may be repeated. Default=False.

as_dict()[source]

Get a dictionary representation of the Keyword

classmethod from_dict(d)[source]

Initialise from dictonary

static from_string(s)[source]

Initialise from a string

get_string()[source]

String representation of Keyword

verbosity(v)[source]

Change the printing of this keyword’s description.

class KeywordList(*args, **kwds)[source]

Bases: monty.json.MSONable, collections.abc.Sequence, typing.Generic

Some keywords can be repeated, which makes accessing them via the normal dictionary methods a little unnatural. This class deals with this by defining a collection of same-named keywords that are accessed by one name.

Initializes a keyword. These Keywords and the value passed to them are sometimes as simple as KEYWORD VALUE, but can also be more elaborate such as KEYWORD [UNITS] VALUE1 VALUE2, which is why this class exists: to handle many values and control easy printing to an input file.

Parameters

keywords – A list of keywords. Must all have the same name (case-insensitive)

append(item)[source]

append the keyword list

extend(l)[source]

extend the keyword list

get_string(indent=0)[source]

String representation of Keyword

verbosity(verbosity)[source]

Silence all keywords in keyword list

class Kind(specie: str, alias: Optional[str] = None, magnetization: float = 0.0, subsections: dict = None, basis_set: str = 'GTH_BASIS', potential: str = 'GTH_POTENTIALS', ghost: bool = False, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Specifies the information for the different atom types being simulated.

Initialize a KIND section

Parameters
  • specie (Species or Element) – Object representing the atom.

  • alias (str) – Alias for the atom, can be used for specifying modifcations to certain atoms but not all, e.g. Mg_1 and Mg_2 to force difference oxidation states on the two atoms.

  • magnetization (float) – From the CP2K Manual: The magnetization used in the atomic initial guess. Adds magnetization/2 spin-alpha electrons and removes magnetization/2 spin-beta electrons.

  • basis_set (str) – Basis set for this atom, accessible from the basis set file specified

  • potential (str) – Pseudopotential for this atom, accessible from the potential file

  • kwargs – Additional kwargs to pass to Section()

class Kpoints(kpts: Union[Sequence, Sequence[Sequence[int]]], weights: Optional[Sequence] = None, eps_geo: float = 1e-06, full_grid: bool = False, parallel_group_size: int = - 1, scheme: str = 'MONKHORST-PACK', symmetry: bool = False, units: str = 'B_VECTOR', verbose: bool = False, wavefunctions: str = 'COMPLEX')[source]

Bases: pymatgen.io.cp2k.inputs.Section

Description of the k-points to use for the calculation.

Parameters
  • kpts (list, tuple) – a 2D array for the kpoints of the form [(1,1,1),]. If len(kpts) == 1. Then it is taken as subdivisions for automatic kpoint scheme. If it has more entries, it is taken as manual entries for kpoints.

  • weights (list, tuple) – a weight for each kpoint. Default is to weigh each by 1

  • eps_geo (float) – tolerance for symmetry. Default=1e-6

  • full_grid (bool) – use full (not reduced) kpoint grid. Default=False.

  • parallel_group_size (int) – from cp2k manual: Number of processors to be used for a single kpoint. This number must divide the total number of processes. The number of groups must divide the total number of kpoints. Value=-1 (smallest possible number of processes per group, satisfying the constraints). Value=0 (all processes). Value=n (exactly n processes). Default=-1.

  • scheme (str) – kpoint generation scheme. Default=’Monkhorst-Pack’

  • symmetry (bool) – Use symmetry to reduce the number of kpoints. Default=False.

  • units (str) – Units for the kpoint coordinates (reciprocal coordinates or cartesian). Default=’B_VECTOR’ (reciprocal)

  • verbose (bool) – verbose output for kpoints. Default=False

  • wavefunctions (str) – Whether to use complex or real valued wavefunctions (if available). Default=’complex’

classmethod from_kpoints(kpoints)[source]

Initialize the section from a Kpoints object (pymatgen.io.vasp.inputs).

Parameters

kpoints – A pymatgen kpoints object.

class LDOS(index: int = 1, alias: Optional[str] = None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls printing of the LDOS (List-Density of states). i.e. projects onto specific atoms.

Initialize the LDOS section

Parameters

index – Index of the atom to project onto

class MO_Cubes(write_cube: bool = False, nhomo: int = 1, nlumo: int = 1, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls printing of the molecular orbital eigenvalues

Initialize the MO_CUBES section

class Mgrid(cutoff: Union[int, float] = 1200, rel_cutoff: Union[int, float] = 80, ngrids: int = 5, progression_factor: int = 3, subsections: dict = None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls the multigrid for numerical integration

Initialize the MGRID section

Parameters
  • cutoff – Cutoff energy (in Rydbergs for historical reasons) defining how find of Gaussians will be used

  • rel_cutoff – The relative cutoff energy, which defines how to map the Gaussians onto the multigrid. If the the value is too low then, even if you have a high cutoff with sharp Gaussians, they will be mapped to the course part of the multigrid

  • ngrids – number of grids to use

  • progression_factor – divisor that decides how to map Gaussians the multigrid after the highest mapping is decided by rel_cutoff

class OrbitalTransformation(minimizer: str = 'CG', preconditioner: str = 'FULL_ALL', algorithm: str = 'STRICT', energy_gap: float = 0.01, linesearch: str = '2PNT', subsections: dict = None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Turns on the Orbital Transformation scheme for diagonalizing the Hamiltonian. Much faster and with guaranteed convergence compared to normal diagonalization, but requires the system to have a band gap.

NOTE: OT has poor convergence for metallic systems and cannot use SCF mixing or smearing. Therefore, you should not use it for metals or systems with ‘small’ band gaps. In that case, use normal diagonalization, which will be slower, but will converge properly.

Initialize the OT section

Parameters
  • minimizer (str) – The minimizer to use with the OT method. Default is conjugate gradient method, which is more robust, but more well-behaved systems should use DIIS, which can be as much as 50% faster.

  • preconditioner (str) – Preconditionar to use for OT, FULL_ALL tends to be most robust, but is not always most efficient. For difficult systems, FULL_SINGLE_INVERSE can be more robust, and is reasonably efficient with large systems. For huge, but well behaved, systems, where construction of the preconditioner can take a very long time, FULL_KINETIC can be a good choice.

  • energy_gap (float) – Guess for the band gap. For FULL_ALL, should be smaller than the actual band gap, so simply using 0.01 is a robust value. Choosing a larger value will help if you start with a bad initial guess though. For FULL_SINGLE_INVERSE, energy_gap is treated as a lower bound. Values lower than 0.05 in this case can lead to stability issues.

  • algorithm (str) – What algorithm to use for OT. ‘Strict’: Taylor or diagonalization based algorithm. IRAC: Orbital Transformation based Iterative Refinement of the Approximative Congruence transformation (OT/IR).

  • linesearch (str) – From the manual: 1D line search algorithm to be used with the OT minimizer, in increasing order of robustness and cost. MINIMIZER CG combined with LINESEARCH GOLD should always find an electronic minimum. Whereas the 2PNT minimizer is almost always OK, 3PNT might be needed for systems in which successive OT CG steps do not decrease the total energy.

class PBE(parameterization: str = 'ORIG', scale_c: Union[float, int] = 1, scale_x: Union[float, int] = 1)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Info about the PBE functional.

Parameters
  • parameterization (str) – ORIG: original PBE PBESOL: PBE for solids/surfaces REVPBE: revised PBE

  • scale_c (float) – scales the correlation part of the functional.

  • scale_x (float) – scales the exchange part of the functional.

class PDOS(nlumo: int = - 1, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls printing of projected density of states onto the different atom KINDS (elemental decomposed DOS).

Initialize the PDOS section

Parameters

nlumo – how many unoccupied orbitals to include (-1==ALL)

class QS(method: str = 'GPW', eps_default: float = 1e-07, extrapolation: str = 'PS', subsections: dict = None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls the quickstep settings (DFT driver)

Initialize the QS Section

Parameters
  • method – What dft methodology to use. Can be GPW (Gaussian Plane Waves) for DFT with pseudopotentials or GAPW (Gaussian Augmented Plane Waves) for all electron calculations

  • eps_default – The default level of convergence accuracy. NOTE: This is a global value for all the numerical value of all EPS_* values in QS module. It is not the same as EPS_SCF, which sets convergence accuracy of the SCF cycle alone.

  • extrapolation – Method use for extrapolation. If using gamma-point-only calculation, then one should use PS for relaxations and ASPC for MD. See the manual for other options.

  • subsections – Subsections to initialize with

class Scf(max_scf: int = 50, eps_scf: float = 1e-06, scf_guess: str = 'RESTART', subsections: dict = None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls the self consistent field loop

Initialize the Scf section

Parameters
  • max_scf – maximum number of SCF loops before terminating

  • eps_scf – convergence criteria for SCF loop

  • scf_guess – Initial guess for SCF loop (RESTART will switch to ATOMIC when no restart file is present)

class Section(name: str, subsections: dict = None, repeats: bool = False, description: Optional[str] = None, keywords: Dict = {}, section_parameters: Union[List, Tuple] = [], location: str = None, verbose: bool = True, alias: str = None, **kwargs)[source]

Bases: monty.json.MSONable

Basic input representation of input to Cp2k. Activates functionality inside of the Cp2k executable.

Basic object representing a CP2K Section. Sections activate different parts of the calculation. For example, FORCE_EVAL section will activate CP2K’s ability to calculate forces. Sections are described with the following:

Parameters
  • name (str) – The name of the section (must match name in CP2K)

  • subsections (dict) – A dictionary of subsections that are nested in this section. Format is {‘NAME’: Section(args, **kwargs). The name you chose for ‘NAME’ to index that subsection does not *have to be the same as the section’s true name, but we recommend matching them. You can specify a blank dictionary if there are no subsections, or if you want to insert the subsections later.

  • repeats (bool) – Whether or not this section can be repeated. Most sections cannot. Default=False.

  • description (str) – Description of this section for easier readability/more descriptiveness

  • keywords (list) – the keywords to be set for this section. Each element should be a Keyword object. This can be more cumbersome than simply using kwargs for building a class in a script, but is more convenient for the class instantiations of CP2K sections (see below).

  • section_parameters (list) – the section parameters for this section. Section parameters are specialized keywords that modify the behavior of the section overall. Most sections do not have section parameters, but some do. Unlike normal Keywords, these are specified as strings and not as Keyword objects.

  • location (str) – the path to the section in the form ‘SECTION/SUBSECTION1/SUBSECTION3’, example for QS module: ‘FORCE_EVAL/DFT/QS’. This location is used to automatically determine if a subsection requires a supersection to be activated.

  • verbose (str) – Controls how much is printed to Cp2k input files (Also see Keyword). If True, then a description of the section will be printed with it as a comment (if description is set). Default=True.

  • are interpreted as keyword (kwargs) –

  • pairs and added to the keywords array as Keyword objects (value) –

add(other)[source]

Add another keyword to the current section

by_path(path: str)[source]

Access a sub-section using a path. Used by the file parser.

Parameters

path (str) – Path to section of form ‘SUBSECTION1/SUBSECTION2/SUBSECTION_OF_INTEREST’

check(path: str)[source]

Check if section exists within the current using a path. Can be useful for cross-checking whether or not required dependencies have been satisfied, which CP2K does not enforce.

Parameters

path (str) – Path to section of form ‘SUBSECTION1/SUBSECTION2/SUBSECTION_OF_INTEREST’

get(d, default=None)[source]

Similar to get for dictionaries. This will attempt to retrieve the section or keyword matching d. Will not raise an error if d does not exist.

Parameters
  • d – the key to retrieve, if present

  • default – what to return if d is not found

get_string()[source]

Get string representation of Section

inc(d: dict)[source]

Mongo style dict modification. Include.

insert(d)[source]

Insert a new section as a subsection of the current one

required_keywords: tuple = ()[source]
required_sections: tuple = ()[source]
set(d: dict)[source]

Alias for update. Used by custodian.

silence()[source]

Recursively delete all print sections so that only defaults are printed out.

unset(d: dict)[source]

Dict based deletion. Used by custodian.

update(d: dict)[source]

Update the Section according to a dictionary argument. This is most useful for providing user-override settings to default parameters. As you pass a dictionary the class variables like “description”, “location”, or “repeats” are not included. Therefore, it is recommended that this be used to modify existing Section objects to a user’s needs, but not used for the creation of new Section child-classes.

Parameters

d (dict) –

A dictionary containing the update information. Should use nested dictionaries to specify the full path of the update. If a section or keyword does not exist, it will be created, but only with the values that are provided in “d”, not using default values from a Section object.

{‘SUBSECTION1’: {‘SUBSEC2’: {‘NEW_KEYWORD’, ‘NEW_VAL’},{‘NEW_SUBSEC’: {‘NEW_KWD’: ‘NEW_VAL’}}}

verbosity(verbosity)[source]

Change the section verbossity recursively by turning on/off the printing of descriptions. Turning off descriptions may reduce the appealing documentation of input files, but also helps de-clutter them.

class Smear(elec_temp: Union[int, float] = 300, method: str = 'FERMI_DIRAC', fixed_magnetic_moment: float = - 100.0, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Control electron smearing

Initialize the SMEAR section

class Subsys(subsections: dict = None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls the definition of the system to be simulated

Initialize the subsys section

class V_Hartree_Cube(keywords=None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Controls printing of the hartree potential as a cube file.

Initialize the V_HARTREE_CUBE section

class XC_FUNCTIONAL(functional: str, subsections: dict = None, **kwargs)[source]

Bases: pymatgen.io.cp2k.inputs.Section

Defines the XC functional to use

Initialize the XC_FUNCTIONAL class