Finite-difference Poisson-Boltzmann calculations

Parameters for a finite-difference polar solvation calculation.

apbs.input_file.calculate.finite_difference.ERROR_TOLERANCE = 1e-06

Relative error tolerance for iterative solver.

class apbs.input_file.calculate.finite_difference.FiniteDifference(dict_=None, yaml=None, json=None)[source]

Bases: apbs.input_file.InputFile

Parameters for a finite-difference polar solvation Poisson-Boltzmann calculation.

Objects can be initialized with dictionary/JSON/YAML data with the following keys:

boundary_condition

Boundary condition for Poisson-Boltzmann equation.

This property can have one of the following values:

  • zero: Dirichlet condition where the potential at the boundary is set to zero. This condition is not commonly used and can result in large errors if used inappropriately.
  • single sphere: Dirichlet condition where the potential at the boundary is set to the values prescribed by a Debye-Hückel model for a single sphere with a point charge, dipole, and quadrupole. The sphere radius in this model is set to the radius of the biomolecule and the sphere charge, dipole, and quadrupole are set to the total moments of the protein. This condition works best when the boundary is sufficiently far (multiple Debye lengths) from the biomolecule.
  • multiple sphere: Dirichlet condition where the potential at the boundary is set to the values prescribed by a Debye-Hückel model for multiple, non-interacting spheres with a point charges. The radii of the non-interacting spheres are set to the atomic radii of and the sphere charges are set to the atomic charges. This condition works better than single sphere for closer boundaries but can be very slow for large biomolecules.
  • focus :c:var:`alias`: Dirichlet condition where the potential at the boundary is set to the values computed by a previous (usually lower-resolution) PB calculation with alias :c:var:`alias`. All of the boundary points should lie within the domain of the previous calculation for best accuracy; if any boundary points lie outside, their values are computed using the single sphere Debye-Hückel boundary condition (see above).
Raises:
  • ValueError – if set to an invalid boundary type
  • IndexError – if an insufficient number of words are present
calculate_energy

Indicate whether energy should be calculated.

Raises:TypeError – if not Boolean
calculate_forces

Indicate whether forces should be calculated.

Raises:TypeError – if not Boolean
calculation_parameters

Specify parameters specific to the calculation type.

The specific class is specific to the calculation type (see calculation_type()):

Note

The calculation_type() property must be set before setting this property (sorry…)

Raises:ValueError – if calculation parameter class doesn’t match calculation type or if
calculation_type
this uses multiple grids to generate high-resolution
solutions at a region of interest. See Focus for more information.
  • manual: perform a traditional non-focused calculation. See Manual for more information.
Raises:

ValueError – if invalid calculation type specified

Type:

Specify the type of finite difference calculation to perform

Type:
  • focus
charge_discretization

The method by which the biomolecular point charges (i.e., Dirac delta functions) by which charges are mapped to the grid used for the finite difference calculation.

As we are attempting to model delta functions, the support (domain) of these discretized charge distributions is always strongly dependent on the grid spacing.

The following types of discretization are supported:

  • linear: Traditional trilinear interpolation (linear splines). The charge is mapped onto the nearest-neighbor grid points. Resulting potentials are very sensitive to grid spacing, length, and position.
  • cubic: Cubic B-spline discretization. The charge is mapped onto the nearest- and next-nearest-neighbor grid points. Resulting potentials are somewhat less sensitive (than linear) to grid spacing, length, and position.
  • quintic: Quintic B-spline discretization. Similar to cubic, except the charge/multipole is additionally mapped to include next-next-nearest neighbors (125 grid points receive charge density).
Raises:
equation

Specifies which version of the Poisson-Boltzmann equation (PBE) to solve:

  • Most users should use one of these:
    • linearized pbe
    • nonlinear pbe
  • These versions are experimental and unstable:
    • linearized regularized pbe
    • nonlinear regularized pbe
Raises:
error_tolerance

Relative error tolerance for iterative solver.

If not specified, the default value is ERROR_TOLERANCE.

Raises:
from_dict(input_)[source]

Load object from dictionary.

Raises:
from_json(input_)

Parse JSON-format input string into this object.

Parameters:input (str) – JSON-format input string
from_yaml(input_)

Parse YAML-format input string into this object.

Parameters:input (str) – YAML-format input string
ions

Descriptions of mobile ion species.

Raises:TypeError – if not set to a Ions object
molecule

Specify which molecule to use for calculation.

Returns:alias for molecule read (see Data loading input file section (required))
Raises:TypeError – if not set to a string
noop

Determine whether solver is run.

If set to True, then skip running the solver but still calculate coefficient maps, etc.

The default value of this property is False.

Raises:TypeError – if not set to bool
solute_dielectric

Solute dielectric.

The dielectric value of a solute is often chosen using the following rules of thumb:

  • 1: only used to compare calculation results with non-polarizable molecular simulation
  • 2-4: “molecular” dielectric value; used when conformational degrees of freedom are modeled explicitly
  • 4-8: used to mimic sidechain libration and other small-scale conformational degrees of freedom
  • 8-12: used to model larger-scale sidechain rearrangement
  • 20-40: used to model larger-scale macromolecular conformational changes and/or water penetration into interior of molecule

Note

What does the continuum dielectric value of a non-continuum molecule mean? Hard to say – this approximation can be very difficult to interpret and can significant affect your results.

Returns:

a floating point number greater than or equal to one

Raises:
solvent_dielectric

Solvent dielectric.

78.5 is a good choice for water at 298 K.

Returns:

a floating point number greater than or equal to one

Raises:
solvent_radius

Radius of the solvent molecules.

This parameter is used to define various solvent-related surfaces and volumes (see surface_method()). This value is usually set to 1.4 Å for a water-like molecular surface and set to 0 Å for a van der Waals surface.

Raises:ValueError – if value is not a non-negative number
surface_method

Method used to defined solute-solvent interface.

One of the following values:

  • molecular surface: The dielectric coefficient is defined based on a molecular surface definition. The problem domain is divided into two spaces. The “free volume” space is defined by the union of solvent-sized spheres (see solvent_radius()) which do not overlap with the solute atoms. This free volume is assigned bulk solvent dielectric values. The complement of this space is assigned solute dielectric values. When the solvent radius is set to zero, this method corresponds to a van der Waals surface definition. The ion-accessibility coefficient is defined by an “inflated” van der Waals model. Specifically, the radius of each biomolecular atom is increased by the radius of the ion species (as specified with the ion() property). The problem domain is then divided into two spaces. The space inside the union of these inflated atomic spheres is assigned an ion-accessibility value of 0; the complement space is assigned the bulk ion accessibility value. See Connolly ML, J Appl Crystallography 16 548-558, 1983 (10.1107/S0021889883010985).
  • smoothed molecular surface: The dielectric and ion-accessibility coefficients are defined as for the molecular surface (see above). However, they are then “smoothed” by a 9-point harmonic averaging to somewhat reduce sensitivity to the grid setup. See Bruccoleri et al. J Comput Chem 18 268-276, 1997 (10.1007/s00214-007-0397-0).
  • cubic spline: The dielectric and ion-accessibility coefficients are defined by a cubic-spline surface as described by Im et al, Comp Phys Commun 111 (1-3) 59-75, 1998 (10.1016/S0010-4655(98)00016-2). The width of the dielectric interface is controlled by the spline_window() property. These spline-based surface definitions are very stable with respect to grid parameters and therefore ideal for calculating forces. However, they require substantial reparameterization of the force field; interested users should consult Nina et al, Biophys Chem 78 (1-2) 89-96, 1999 (10.1016/S0301-4622(98)00236-1). Additionally, these surfaces can generate unphysical results with non-zero ionic strengths.
  • septic spline: The dielectric and ion-accessibility coefficients are defined by a 7th order polynomial. This surface definition has characteristics similar to the cubic spline, but provides higher order continuity necessary for stable force calculations with atomic multipole force fields (up to quadrupole).
Raises:
surface_spline_window

Window for spline-based surface definitions (not needed otherwise).

This is the distance (in Å) over which the spline transitions from the solvent dielectric value to the solute dielectric value. A typical value is 0.3 Å.

Returns:positive number
Raises:TypeError – if not a positive number
temperature

Temperature for the calculation in Kelvin.

Raises:ValueError – if not a positive number (no violations of the 3rd Law!)
to_dict() → dict[source]

Produce dictionary representation of self.

to_json() → str

Produce JSON representation of self.

to_yaml() → str

Produce YAML representation of self.

use_maps

Information for (optionally) using maps read into APBS.

Returns:list of UseMap objects
Raises:TypeError – if not a list of UseMap objects
validate()[source]

Validate the object.

Raises:ValueError – if object is not valid
write_atom_potentials

Write out the electrostatic potential at each atom location.

Write out text file with potential at center of atom in units of \(k_b \, T \, e_c^{-1}\).

Note

These numbers are meaningless by themselves due to the presence of “self-energy” terms that are sensitive to grid spacing and position. These numbers should be evaluated with respect to a reference calculation: the potentials from that reference calculation should be subtracted from the target system. For example, one calculation might include a molecule with a heterogeneous dielectric coefficient and the reference system might be exactly the same system setup but with a homeogeneous dielectric coefficient. If the results from the reference calculation are substracted from the first calculation, then the result will be a physically meaningful reaction field potential. However, the results from the first and reference calculations are meaningless by themselves.

Returns:path to text file for writing atom potential values.
Raises:TypeError – if not set to string
write_maps

Write out maps related to computed properties.

See WriteMap for more information.

Raises:TypeError – if set to wrong type
class apbs.input_file.calculate.finite_difference.Focus(dict_=None, yaml=None, json=None)[source]

Bases: apbs.input_file.InputFile

Parameters specific to a focused-grid finite-difference polar solvation Poisson-Boltzmann calculation.

Focusing provides automatically configured finite difference Poisson-Boltzmann calculations.

This multigrid calculation automatically sets up and performs a string of single-point PBE calculations to “focus” on a region of interest (binding site, etc.) in a system. It is basically an automated way to set parameters in Manual, allowing for (hopefully) easier use. Most users should use this Focus rather than Manual.

Focusing is a method for solving the Poisson-Boltzmann equation in a finite difference setting. Some of the earliest references to this method are from Gilson and Honig (DOI:10.1038/330084a0). The method starts by solving the equation on a coarse grid (i.e., few grid points) with large dimensions (i.e., grid lengths). The solution on this coarse grid is then used to set the Dirichlet boundary condition values for a smaller problem domain – and therefore a finer grid – surrounding the region of interest. The finer grid spacing in the smaller problem domain often provides greater accuracy in the solution.

Note

During focusing calculations, you may encounter the message WARNING! Unusually large potential values detected on the focusing boundary! for some highly charged systems based on location of the focusing boundary. First, you should determine if you received any other warning or error messages as part of this calculation, particularly those referring to exceeded number of iterations or error tolerance. Next, you should check if the calculation converged to a reasonable answer. In particular, you should check sensitivity to the grid spacing by making small changes to the grid lengths and see if the changes in energies are correspondingly small. If so, then this warning can be safely ignored.

Objects can be initialized with dictionary/JSON/YAML data with the following keys:

coarse_grid_center

The center of the coarse grid in a focusing calculation.

Raises:TypeError – if not GridCenter object
coarse_grid_dimensions

Dimensions of the coarse grid in a focusing finite-difference Poisson-Boltzmann calculation.

This is the starting mesh, so it should be large enough to completely enclose the biomolecule and ensure that the chosen boundary condition is appropriate.

Raises:TypeError – if the value is not GridDimensions.
fine_grid_center

The center of the fine grid in a focusing calculation.

Raises:TypeError – if not GridCenter object
fine_grid_dimensions

Dimensions of the coarse grid in a focusing finite-difference Poisson-Boltzmann calculation.

This should enclose the region of interest in the biomolecule.

Raises:TypeError – if the value is not GridDimensions.
from_dict(input_)[source]

Load object from dictionary.

Raises:KeyError – if missing items.
from_json(input_)

Parse JSON-format input string into this object.

Parameters:input (str) – JSON-format input string
from_yaml(input_)

Parse YAML-format input string into this object.

Parameters:input (str) – YAML-format input string
parallel

Indicate whether a parallel calculation should be performed.

If set, the parallel_parameters() property must also be set. Set the parallel_parameters() property before setting this value to True (sorry).

Raises:TypeError – if set to something other than bool.
parallel_parameters

Provide parameters for a parallel focusing calculation.

This property is optional and only used if parallel() is True.

Raises:TypeError – if not ParallelFocus class
to_dict() → dict[source]

Produce dictionary representation of self.

to_json() → str

Produce JSON representation of self.

to_yaml() → str

Produce YAML representation of self.

validate()[source]

Validate the object.

Raises:ValueError – if object is not valid
class apbs.input_file.calculate.finite_difference.GridCenter(dict_=None, yaml=None, json=None)[source]

Bases: apbs.input_file.InputFile

Parameters for specifying the center of a finite difference grid.

Objects can be initialized with dictionary/JSON/YAML data with one of the following keys:

  • molecule: use the center of the specified molecule. See molecule().
  • position: use a specific position (coordinates). See position().
from_dict(input_)[source]

Parse dictionary-format input into this object.

Parameters:input (dict) – input dictionary
Raises:KeyError – when input is missing
from_json(input_)

Parse JSON-format input string into this object.

Parameters:input (str) – JSON-format input string
from_yaml(input_)

Parse YAML-format input string into this object.

Parameters:input (str) – YAML-format input string
molecule

Alias for molecule used to center the finite difference grid. See Data loading input file section (required) for more information on reading molecules into APBS.

Returns:alias for molecule at center
Raises:TypeError – if alias is not string
position

Coordinates for the grid center.

Returns:

3-element list of float positions in Å units

Raises:
to_dict() → dict[source]

Produce dictionary representation of self.

to_json() → str

Produce JSON representation of self.

to_yaml() → str

Produce YAML representation of self.

validate()[source]

Validate the object.

Raises:ValueError – if object is not valid
class apbs.input_file.calculate.finite_difference.GridDimensions(dict_=None, yaml=None, json=None)[source]

Bases: apbs.input_file.InputFile

Parameters for controlling the number of grid points and the grid spacing.

For a given direction \(i \in [x, y, z]\), the grid spacing \(h_i\) is related to the number of grid points \(n_i\) and the grid length \(l_i\) by

\[l_i = h_i (n_i - 1)\]

Therefore, only two of the following three properties are needed. When all three are provided, the spacings() will be inferred from the counts() and lengths().

  • counts: the number of grid points in each direction; see counts() for more information.
  • spacings: the spacing of grid points in each direction; see spacings() for more information.
  • lengths: the span of the grid in each dimension; see lengths() for more information.

The following properties are read-only and cannot be set:

  • levels: the number of multigrid levels in the calculation; see levels() for more information.

The number of levels \(p\) is a constant for the entire domain (i.e., it doesn’t vary by grid direction. The number of levels is related to the number of grid points \(n_i\) in the direction \(i\) by the formula

\[n_i = c \, 2^{p} + 1\]

where \(c\) is a non-zero integer. The value of \(p\) is chosen as the smallest that satisfies the formula above for all \(n_i\). Calculations become more efficient for larger \(p\) so APBS will adjust the numbers of grid points \(n_i\) downwards to ensure that \(p\) is equal to or greater than MIN_LEVEL, likely resulting in lower grid resolution (and accuracy).

The most common values for \(n_i\) are 65, 97, 129, and 161 (they can be different in each direction); these are all compatible with \(p = 4\).

adjust_counts() → bool[source]

Adjust numbers of grid points to acheive multigrid level of at least MIN_LEVEL.

This routine alters self._counts and (re)sets self._level. The new counts will always be less than or equal to the old counts.

Parameters:targets (list) – list of target numbers of grid points.
Returns:True if counts were changes from targets
counts

Target number of grid points in a finite-difference Poisson-Boltzmann calculation.

Note

the actual number of grid points might be adjusted to achieve a better multigrid level. The most common values for grid dimensions are 65, 97, 129, and 161 (they can be different in each direction); these are all compatible with a levels() value of 4.

Returns:

3-vector of integers greater than 32 indicating number of grid points in the x-, y-, and z-directions. The numbers may be different in each direction.

Raises:
  • IndexError – if the length of the value is not 3.
  • TypeError – if the counts are not integers greater than 2**``MIN_LEVEL``.
  • ValueError – if not enough information is present to calculate counts
static find_count(target, min_level=4, max_level=None) → tuple[source]

Find the number of grid points closest to the given target that satisfies

\[n_i = \arg \max_p c \, 2^p + 1\]

where \(n_i\) is the number of grid points (always less than or equal to target), \(c\) is a positive integer constant, and \(p\) is the number of levels, greater than or equal to MIN_LEVEL.

Note

assumes target is greater than or equal to 2 ** MIN_LEVEL

Parameters:
  • target (int) – target number of grid points
  • max_level (int) – maximum level for multigrid
  • min_level (int) – minimum level for multigrid
Returns:

(\(n_i\), \(c\), \(p\))

from_dict(input_)[source]

Initialize object from dictionary.

from_json(input_)

Parse JSON-format input string into this object.

Parameters:input (str) – JSON-format input string
from_yaml(input_)

Parse YAML-format input string into this object.

Parameters:input (str) – YAML-format input string
lengths

Length of the grid in a finite-difference Poisson-Boltzmann calculation.

Returns:

3-vector of positive numbers indicating grid lengths in the x-, y-, and z-directions in Å. The lengths may be different in each direction.

Raises:
  • IndexError – if the length of the value is not 3.
  • TypeError – if the lengths are not positive numbers.
  • ValueError – if not enough information is available to calculate lengths
levels

Number of multigrid levels, an integer greater than 2.

Raises:
spacings

Spacings of the grid in a finite-difference Poisson-Boltzmann calculation.

Returns:

3-vector of positive numbers indicating grid spacings in the x-, y-, and z-directions in Å. The spacings may be different in each direction.

Raises:
  • IndexError – if the length of the value is not 3.
  • TypeError – if the spacings are not positive numbers.
  • ValueError – if not enough information is available to calculate spacings
to_dict()[source]

Produce dictionary representation of self.

to_json() → str

Produce JSON representation of self.

to_yaml() → str

Produce YAML representation of self.

validate()[source]

Validate the object.

Raises:ValueError – if object is not valid
apbs.input_file.calculate.finite_difference.MIN_LEVEL = 4

Mininum multigrid level in calculations.

class apbs.input_file.calculate.finite_difference.Manual(dict_=None, yaml=None, json=None)[source]

Bases: apbs.input_file.InputFile

Parameters specific to a manual finite-difference polar solvation Poisson-Boltzmann calculation.

Manually-configured finite difference multigrid Poisson-Boltzmann calculations.

This is a standard single-point multigrid PBE calculation without focusing or additional refinement. This FiniteDifference.calculation_type() offers the most control of parameters to the user. Several of these calculations can be strung together to perform focusing calculations by judicious choice of the FiniteDifference.boundary_condition() property; however, the setup of the focusing is not automated as it is in the Focus FiniteDifference.calculation_type(). Therefore, this command should primarily be used by more experienced users.

Objects can be initialized with dictionary/JSON/YAML data with the following keys:

from_dict(input_)[source]

Load object from dictionary.

Raises:KeyError – if missing entries
from_json(input_)

Parse JSON-format input string into this object.

Parameters:input (str) – JSON-format input string
from_yaml(input_)

Parse YAML-format input string into this object.

Parameters:input (str) – YAML-format input string
grid_center

The center of the grid in a finite difference Poisson-Boltzmann calculation.

Raises:TypeError – if not GridCenter object
grid_dimensions

Dimensions of the grid in a focusing finite-difference Poisson-Boltzmann calculation.

Raises:TypeError – if the value is not GridDimensions.
to_dict() → dict[source]

Produce dictionary representation of self.

to_json() → str

Produce JSON representation of self.

to_yaml() → str

Produce YAML representation of self.

validate()[source]

Validate the object.

Raises:ValueError – if object is not valid
class apbs.input_file.calculate.finite_difference.ParallelFocus(dict_=None, yaml=None, json=None)[source]

Bases: apbs.input_file.InputFile

Parameters specific to a parallel focusing finite-difference polar solvation Poisson-Boltzmann calculation.

Objects can be initialized with dictionary/JSON/YAML data with the following keys:

asynchronous_rank

Integer rank of processor in processor array.

This should be a positive integer.

Raises:TypeError – if not a positive integer
from_dict(input_)[source]

Load object from dictionary.

Raises:KeyError – for missing items
from_json(input_)

Parse JSON-format input string into this object.

Parameters:input (str) – JSON-format input string
from_yaml(input_)

Parse YAML-format input string into this object.

Parameters:input (str) – YAML-format input string
overlap_fraction

Fraction of size of parallel focusing domains that overlap.

This is a positive floating point number less than 1.

Raises:
  • TypeError – if the value is not a positive number
  • ValueError – if the value is greater than 1
processor_array

Distribution of processors over problem domain.

Returns:

a 3-vector of positive integers: [npx npy npz] the integer number of processors to be used in the x-, y- and z-directions of the system. The product npx × npy × npz should be less than or equal to the total number of processors with which APBS was invoked (usually via mpirun). If more processors are provided at invocation than actually used during the run, the extra processors are not used in the calculation. The processors are tiled across the domain in a Cartesian fashion with a specified amount of overlap (see overlap_fraction()) between each processor to ensure continuity of the solution. Each processor’s subdomain will contain the number of grid points specified by the dime keyword. For broad spatial support of the splines, every charge included in partition needs to be at least 1 grid space (first-order spline), 2 grid spaces (third-order spline), or 3 grid spaces (fifth-order spline) away from the partition boundary.

Raises:
  • TypeError – if not a list or list of positive integers.
  • IndexError – if vector doesn’t have length 3.
to_dict() → dict[source]

Produce dictionary representation of self.

to_json() → str

Produce JSON representation of self.

to_yaml() → str

Produce YAML representation of self.

validate()[source]

Validate the object.

Raises:ValueError – if object is not valid