Finite-element Poisson-Boltzmann calculations

Parameters for a finite-element polar solvation calculation.

class apbs.input_file.calculate.finite_element.FiniteElement(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:

a_priori_refinement

Strategy for refining the initial very coarse 8-tetrahedron initial finite element to a resolution suitable for the main solve-estimate-refine iteration.

Allowed values include:

  • geometric: geometry-based refinement at molecular surface and charges
  • uniform: uniform refinement of the mesh
Raises:
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
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:
domain_length

Length of rectangular prism domain.

Returns:

list of non-zero lengths for x, y, and z directions.

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_based_refinement

Specify error-based refinement strategy for simplices.

Can be assigned one of the following values:

  • global: this refines simplices until the global error is less than the amount specified by error_tolerance().
  • simplex: this refines simplices until the per-simplex error is less than the amount specified by error_tolerance().
  • fraction: this refines the specified fraction of simplices with the largest per-simplex error. The fraction is specified by error_tolerance().
Raises:
error_tolerance

Error tolerance for error-based refinement of simplices.

The meaning of this property changes based on the setting of error_based_refinement().

raises TypeError: if not a positive number.

from_dict(input_)[source]

Populate object from dictionary.

Parameters:input (dict) – dictionary with input data.
Raises:KeyError – if dictionary missing keys.
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
initial_mesh_resolution

Target spacing (in Å) of initial mesh.

Initial refinement will continue until this target spacing is met or until initial_mesh_vertices() is exceeded.

Raises:TypeError – if not set to a positive number
initial_mesh_vertices

Target number of vertices in initial mesh.

Initial refinement will continue until this number of vertices is exceeded or initial_mesh_resolution() is reached.

Raises:TypeError – if not set to a positive integer
ions

Descriptions of mobile ion species.

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

Maximum number of solve-estimate-refine iterations.

The solve-estimate-refine loop will continue until this number of iterations is reached or the mesh has more than maximum_vertices() vertices.

Raises:TypeError – if not sent to positive integer
maximum_vertices

Maximum number of vertices produced in solve-estimate-refine iterations.

The solve-estimate-refine loop will continue until the mesh has more than this number of vertices or maximum_refinement_iterations() has been reached.

Raises:TypeError – if not sent to positive integer
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
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]

Dump dictionary from object.

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 this object.

Assumes that all attributes have been set via setters.

Raises:ValueError – if object is invalid.
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