CP2K

(contributed by Felipe Zapata, extended by Patrick Melix)

Settings

The cp2k input is rather complex one compared to other computational codes, but its input is structured as a set of nested block and sub-blocks that can be easily represented by the Settings class. Like the other interfaces, the CP2K input file is generated using the input branch of the job settings. For instance, a single point calculation for pentacene:

penta = Settings()

penta.input.force_eval.dft.basis_set_file_name = "/path/to/basis"
penta.input.force_eval.dft.potential_set_file_name = "/path/to/potential"
penta.input.force_eval.dft.mgrid.cutoff = 400
penta.input.force_eval.dft.mgrid.ngrids = 4
penta.input.force_eval.dft.poisson = Settings()
penta.input.force_eval.dft.localize._h = "T"


penta.input.force_eval.dft.qs.method = "GPW"
penta.input.force_eval.dft.scf.eps_scf = 1e-6
penta.input.force_eval.dft.scf.max_scf = 200
penta.input.force_eval.dft.xc.xc_functional._h = "PBE"

penta.input.force_eval.subsys.cell.A = '16.11886919 0.07814137 -0.697284243'
penta.input.force_eval.subsys.cell.B = '-0.215317662 4.389405268 1.408951791'
penta.input.force_eval.subsys.cell.C = '-0.216126961 1.732808365 9.74896108'
penta.input.force_eval.subsys.cell.periodic = 'XYZ'
penta.input.force_eval.subsys.kind.C.basis_set = "DZVP-MOLOPT-SR-GTH-Q4"
penta.input.force_eval.subsys.kind.C.potential = "GTH-PBE-Q4"
penta.input.force_eval.subsys.kind.H.basis_set = "DZVP-MOLOPT-SR-GTH-Q1"
penta.input.force_eval.subsys.kind.H.potential = "GTH-PBE-q1"
penta.input.force_eval.subsys.topology.coord_file_name = "./penta.xyz"
penta.input.force_eval.subsys.topology.coordinate = "xyz"

penta.input['global'].print_level = "low"
penta.input['global'].project  = "example"
penta.input['global'].run_type = "energy_force"

The input generated during the execution of the cp2k job is similar to:

&FORCE_EVAL
  &DFT
    BASIS_SET_FILE_NAME  /path/to/basis
    &MGRID
      CUTOFF  400
      NGRIDS  4
    &END
    &POISSON
    &END
    &LOCALIZE T
    &END
    POTENTIAL_FILE_NAME  /path/to/potential
    &QS
      METHOD  GPW
    &END
    &SCF
      EPS_SCF  1e-06
      MAX_SCF  200
    &END
    &XC
      &XC_FUNCTIONAL PBE
      &END
    &END
  &END
  &SUBSYS
    &CELL
      A  16.11886919 0.07814137 -0.697284243
      B  -0.215317662 4.389405268 1.408951791
      C  -0.216126961 1.732808365 9.748961085
      PERIODIC  XYZ
    &END
    &KIND  C
      BASIS_SET  DZVP-MOLOPT-SR-GTH-q4
      POTENTIAL  GTH-PBE-q4
    &END
    &KIND  H
      BASIS_SET  DZVP-MOLOPT-SR-GTH-q1
      POTENTIAL  GTH-PBE-q1
    &END
    &TOPOLOGY
      COORD_FILE_NAME  ./geometry.xyz
      COORDINATE  XYZ
    &END
  &END
&END

&GLOBAL
  PRINT_LEVEL  LOW
  PROJECT  example
  RUN_TYPE  ENERGY_FORCE
&END

PLAMS automatically creates the indented structure of the previous example together with the special character & at the beginning and end of each section, and finally the keyword END at the end of each section.

Notice that CP2K requires the explicit declaration of the basis set together with the charge and the name of the potential used for each one of the atoms. In the previous example the basis for the carbon is DZVP-MOLOPT-SR-GTH, while the potential is GTH-PBE and the charge q4.

The input parser also allows for header values like in this example &LOCALIZE T &END by using the <some_section>._h (for header) notation in the Settings instance (just like in ADFJob).

Inclusion of external files by using the @INCLUDE notation of Cp2k is supported. Also the @SET and @IF keys can be used, just replace the @ sign by at_ in the definition of your Settings. If you need some files to be copied to the actual execution directory, pass them to the constructor using the copy= option. See the API below.

Molecule parsing

Molecules can be parsed into the input automatically by passing them to the Cp2kJob constructor. Leaving the molecule unset is supported, make sure you then set your molecule by hand (e.g. using the TOPOLOGY section). Parsing of the molecule can be avoided by setting job.settings.ignore_molecule = True in the Settings of the job.

Also, the simulation cell can be specified using the x, y, z vectors like in the previous example. A cubic box can be easily specified by:

penta.input.force_eval.subsys.cell.ABC = "[angstrom] 50 50 50"

That results in a simulation cube of 50 cubic angstroms. For a more detailed description of the cp2k input see manual.

Loading jobs

Calculations done without PLAMS can be loaded using the load_external() functionality. The Cp2kResults class supports reading input files into Settings objects.

Just do Cp2kJob.load_external(path) to get the settings from the file.

Molecule loading

A Molecule can be recreated from a Settings instance using the Cp2kSettings2Mol function.

API

class Cp2kJob[source]

A class representing a single computational job with CP2K.

In addition to the arguments of SingleJob, Cp2kJob takes a copy argument. copy can be a list or string, containing paths to files to be copied to the jobs directory. This might e.g. be a molecule, further input files etc.

_result_type

alias of scm.plams.interfaces.thirdparty.cp2k.Cp2kResults

__init__(copy=None, **kwargs)[source]

Initialize self. See help(type(self)) for accurate signature.

_get_ready()[source]

Copy files to execution dir if self.copy_files is set.

get_input()[source]

Transform all contents of input branch of Settings into string with blocks, subblocks, keys and values.

get_runscript()[source]

Run a parallel version of CP2K.

The exact CP2K executable, and whether or not one wants to use srun or mpirun can be specified under self.settings.executable. Currently supported executables are:

  • "sdbg": Serial single core testing and debugging

  • "sopt": Serial general single core usage

  • "ssmp": Parallel (only OpenMP), single node, multi core

  • "pdbg": Parallel (only MPI) multi-node testing and debugging

  • "popt": Parallel (only MPI) general usage, no threads

  • "psmp": parallel (MPI + OpenMP) general usage, threading might improve scalability and memory usage

For example:


>>> from scm.plams import Cp2kJob
>>> job = Cp2kJob(...)
>>> job.settings.executable = "cp2k.popt"
>>> job.settings.executable = "c2pk.ssmp"
>>> job.settings.executable = "mpirun -np 24 cp2k.psmp"
check()[source]

Look for the normal termination signal in Cp2k output.

class Cp2kResults[source]

A class for CP2K results.

recreate_settings()[source]

Recreate job for load_external().

If a keyword and a section with the same name appear, only the keyword is used. This happens e.g. when reading restart files where kind.symbol.potential is given as Potential ABCD and &POTENTIAL …. &END POTENTIAL.

Limited support of sections that have different formatting like KIND and COORD. Check the resulting Settings instance if the information you want is there. Be careful with reading restart files, since they are automatically generated and not every case is handled well here. You should get all the information but not sure if I know about all special cases of input.

get_runtime()[source]

Return runtime in seconds from output.

get_timings()[source]

Return the Timings from the End of the Output File

get_energy(index=- 1, unit='a.u.')[source]

Returns ‘Total energy:’ from the output file.

Set index to choose the n-th occurence of the total energy in the output, e.g. to choose an optimization step. Also supports slices. Defaults to the last occurence.

get_dispersion(index=- 1, unit='a.u.')[source]

Returns ‘Dispersion energy:’ from the output file.

Set index to choose the n-th occurence of the dispersion energy in the output, e.g. to choose an optimization step. Also supports slices. Defaults to the last occurence.

get_forces(file=None, index=- 1)[source]

Returns list of ndarrays with forces for each atom.

Set file to use other than the main output file.

Set index to choose the n-th occurence of the forces in the output, e.g. to choose an optimization step. Set to None to return all as a list. Defaults to the last occurence.

get_mulliken_charges(return_spin=False, index=- 1)[source]

Get Mulliken charges (and spin moments).

Set index to choose the n-th occurence of the Charges in the output, e.g. to choose an optimization step. Set to None to return all as a list. Defaults to the last occurence.

Returns list of charges. If return_spin is True returns tuple of charges and spins.

get_hirshfeld_charges(return_spin=False, index=- 1)[source]

Get Hirshfeld charges (and spin moments).

Set index to choose the n-th occurence of the Charges in the output, e.g. to choose an optimization step. Set to None to return all as a list. Defaults to the last occurence.

Returns list of charges. If return_spin is True returns tuple of charges and spins.

get_multigrid_info()[source]

Get Information on multigrids.

Usefull for converging cutoffs. Needs ‘Medium’ global print level.

Returns a dict with keys ‘counts’ and ‘cutoffs’.

get_md_infos(file=None, cache=False)[source]

Read the MD-info sections.

Returns a list with descriptors and a nested list containing the values for each timestep.

Set cache to save the results in self.md_infos to speed up analysis by avoiding I/O.

get_md_cell_volumes(file=None, unit='angstrom')[source]

Get cell Volumes using the get_md_infos() function.

get_md_pressure(file=None)[source]

Get pressures using the get_md_infos() function.

check_scf(file=None, return_n=False)[source]

Returns False if the string ‘SCF run NOT converged’ appears in file, otherwise True.

file defaults to self.job._filename('out').

Set return_n to recieve the number of occurences instead of a bool.

check_go(file=None)[source]

Returns False if the string ‘GEOMETRY OPTIMIZATION COMPLETED’ does not appear in file

file defaults to self.job._filename('out').

Cp2kSettings2Mol()[source]

Return a molecule from a Settings instance used for a Cp2kJob.

Loads coordinates from settings.input.force_eval.subsys.coord._h and cell information from settings.input.force_eval.subsys.cell.