Main Page | Modules | Data Structures | File List | Data Fields | Globals | Related Pages

Vpmg class

A wrapper for Mike Holst's PMG multigrid code. More...


Files

file  vpmg.h
 Contains declarations for class Vpmg.


Data Structures

struct  sVpmg
 Contains public data members for Vpmg class/module. More...


Defines

#define VPMGMAXPART   2000

Typedefs

typedef sVpmg Vpmg
 Declaration of the Vpmg class as the Vpmg structure.


Functions

unsigned long int Vpmg_memChk (Vpmg *thee)
 Return the memory used by this structure (and its contents) in bytes.

VpmgVpmg_ctor (Vpmgp *parms, Vpbe *pbe, int focusFlag, Vpmg *pmgOLD, MGparm *mgparm, PBEparm_calcEnergy energyFlag)
 Constructor for the Vpmg class (allocates new memory).

int Vpmg_ctor2 (Vpmg *thee, Vpmgp *parms, Vpbe *pbe, int focusFlag, Vpmg *pmgOLD, MGparm *mgparm, PBEparm_calcEnergy energyFlag)
 FORTRAN stub constructor for the Vpmg class (uses previously-allocated memory).

void Vpmg_dtor (Vpmg **thee)
 Object destructor.

void Vpmg_dtor2 (Vpmg *thee)
 FORTRAN stub object destructor.

int Vpmg_fillco (Vpmg *thee, Vsurf_Meth surfMeth, double splineWin, Vchrg_Meth chargeMeth, int useDielXMap, Vgrid *dielXMap, int useDielYMap, Vgrid *dielYMap, int useDielZMap, Vgrid *dielZMap, int useKappaMap, Vgrid *kappaMap, int useChargeMap, Vgrid *chargeMap)
 Fill the coefficient arrays prior to solving the equation.

int Vpmg_solve (Vpmg *thee)
 Solve the PBE using PMG.

int Vpmg_solveLaplace (Vpmg *thee)
 Solve Poisson's equation with a homogeneous Laplacian operator using the solvent dielectric constant. This solution is performed by a sine wave decomposition.

double Vpmg_energy (Vpmg *thee, int extFlag)
 Get the total electrostatic energy.

double Vpmg_qfEnergy (Vpmg *thee, int extFlag)
 Get the "fixed charge" contribution to the electrostatic energy.

double Vpmg_qfAtomEnergy (Vpmg *thee, Vatom *atom)
 Get the per-atom "fixed charge" contribution to the electrostatic energy.

double Vpmg_qmEnergy (Vpmg *thee, int extFlag)
 Get the "mobile charge" contribution to the electrostatic energy.

double Vpmg_dielEnergy (Vpmg *thee, int extFlag)
 Get the "polarization" contribution to the electrostatic energy.

double Vpmg_npEnergy (Vpmg *thee, int extFlag)
 Get the "apolar" energy.

double Vpmg_dielGradNorm (Vpmg *thee)
 Get the integral of the gradient of the dielectric function.

int Vpmg_force (Vpmg *thee, double *force, int atomID, Vsurf_Meth srfm, Vchrg_Meth chgm)
 Calculate the total force on the specified atom in units of k_B T/AA.

int Vpmg_qfForce (Vpmg *thee, double *force, int atomID, Vchrg_Meth chgm)
 Calculate the "charge-field" force on the specified atom in units of k_B T/AA.

int Vpmg_dbnpForce (Vpmg *thee, double *dbForce, double *npForce, int atomID, Vsurf_Meth srfm)
 Calculate the dielectric boundary and apolar forces on the specified atom in units of k_B T/AA.

int Vpmg_ibForce (Vpmg *thee, double *force, int atomID, Vsurf_Meth srfm)
 Calculate the osmotic pressure on the specified atom in units of k_B T/AA.

void Vpmg_setPart (Vpmg *thee, double lowerCorner[3], double upperCorner[3], int bflags[6])
 Set partition information which restricts the calculation of observables to a (rectangular) subset of the problem domain.

void Vpmg_unsetPart (Vpmg *thee)
 Remove partition restrictions.

int Vpmg_fillArray (Vpmg *thee, double *vec, Vdata_Type type, double parm, Vhal_PBEType pbetype)
 Fill the specified array with accessibility values.

void Vpmg_printColComp (Vpmg *thee, char path[72], char title[72], char mxtype[3], int flag)
 Print out a column-compressed sparse matrix in Harwell-Boeing format.


Detailed Description

A wrapper for Mike Holst's PMG multigrid code.

Note:
Many of the routines and macros are borrowed from the main.c driver (written by Mike Holst) provided with the PMG code.

Define Documentation

#define VPMGMAXPART   2000
 

The maximum number of partitions the mesh can be divided into


Function Documentation

Vpmg* Vpmg_ctor Vpmgp parms,
Vpbe pbe,
int  focusFlag,
Vpmg pmgOLD,
MGparm mgparm,
PBEparm_calcEnergy  energyFlag
 

Constructor for the Vpmg class (allocates new memory).

Author:
Nathan Baker
Returns:
Pointer to newly allocated Vpmg object
Parameters:
pbe  PMG parameter object
focusFlag  PBE-specific variables
pmgOLD  1 for focusing, 0 otherwise
mgparm  Old Vpmg object to use for boundary conditions
energyFlag  MGparm parameter object for boundary conditions What types of energies to calculate

int Vpmg_ctor2 Vpmg thee,
Vpmgp parms,
Vpbe pbe,
int  focusFlag,
Vpmg pmgOLD,
MGparm mgparm,
PBEparm_calcEnergy  energyFlag
 

FORTRAN stub constructor for the Vpmg class (uses previously-allocated memory).

Author:
Nathan Baker
Returns:
1 if successful, 0 otherwise
Parameters:
parms  Memory location for object
pbe  PMG parameter object
focusFlag  PBE-specific variables
pmgOLD  1 for focusing, 0 otherwise
mgparm  Old Vpmg object to use for boundary conditions (can be VNULL if focusFlag = 0)
energyFlag  MGparm parameter object for boundary conditions (can be VNULL if focusFlag = 0) What types of energies to calculate (ignored if focusFlag = 0)

int Vpmg_dbnpForce Vpmg thee,
double *  dbForce,
double *  npForce,
int  atomID,
Vsurf_Meth  srfm
 

Calculate the dielectric boundary and apolar forces on the specified atom in units of k_B T/AA.

Author:
Nathan Baker
Note:
  • Using the force evaluation methods of Im et al (Roux group), Comput Phys Commun, 111, 59--75 (1998). However, this gives the whole (self-interactions included) force -- reaction field forces will have to be calculated at higher level.
  • No contributions are made from higher levels of focusing.
Returns:
1 if successful, 0 otherwise
Parameters:
dbForce  Vpmg object
npForce  3*sizeof(double) space to hold the dielectric boundary force in units of k_B T/AA
atomID  3*sizeof(double) space to hold the apolar boundary force in units of k_B T/AA
srfm  Valist ID of desired atom Surface discretization method

double Vpmg_dielEnergy Vpmg thee,
int  extFlag
 

Get the "polarization" contribution to the electrostatic energy.

Using the solution at the finest mesh level, get the electrostatic energy due to the interaction of the mobile charges with the potential:

\[ G = \frac{1}{2} \int \epsilon (\nabla u)^2 dx \]

where epsilon is the dielectric parameter and u(x) is the dimensionless electrostatic potential. The energy is scaled to units of k_b T.

Author:
Nathan Baker
Note:
The value of this observable may be modified by setting restrictions on the subdomain over which it is calculated. Such limits can be set via Vpmg_setPart and are generally useful for parallel runs.
Returns:
The polarization electrostatic energy in units of k_B T.
Parameters:
extFlag  Vpmg object If this was a focused calculation, include (1 -- for serial calculations) or ignore (0 -- for parallel calculations) energy contributions from outside the focusing domain

double Vpmg_dielGradNorm Vpmg thee  ) 
 

Get the integral of the gradient of the dielectric function.

Using the dielectric map at the finest mesh level, calculate the integral of the norm of the dielectric function gradient routines of Im et al (see Vpmg_dbnpForce for reference):

\[ \int \| \nabla \epsilon \| dx \]

where epsilon is the dielectric parameter. The integral is returned in units of A^2.

Author:
Nathan Baker restrictions on the subdomain over which it is calculated. Such limits can be set via Vpmg_setPart and are generally useful for parallel runs.
Returns:
The integral in units of A^2.
Parameters:
thee  Vpmg object

void Vpmg_dtor Vpmg **  thee  ) 
 

Object destructor.

Author:
Nathan Baker
Parameters:
thee  Pointer to memory location of object to be destroyed

void Vpmg_dtor2 Vpmg thee  ) 
 

FORTRAN stub object destructor.

Author:
Nathan Baker
Parameters:
thee  Pointer to object to be destroyed

double Vpmg_energy Vpmg thee,
int  extFlag
 

Get the total electrostatic energy.

Author:
Nathan Baker
Note:
The value of this observable may be modified by setting restrictions on the subdomain over which it is calculated. Such limits can be set via Vpmg_setPart and are generally useful for parallel runs.
Returns:
The electrostatic energy in units of k_B T.
Parameters:
extFlag  Vpmg object If this was a focused calculation, include (1 -- for serial calculations) or ignore (0 -- for parallel calculations) energy contributions from outside the focusing domain

int Vpmg_fillArray Vpmg thee,
double *  vec,
Vdata_Type  type,
double  parm,
Vhal_PBEType  pbetype
 

Fill the specified array with accessibility values.

Author:
Nathan Baker
Returns:
1 if successful, 0 otherwise
Parameters:
vec  Vpmg object
type  A nx*ny*nz*sizeof(double) array to contain the values to be written
parm  What to write
pbetype  Parameter for data type definition (if needed) Parameter for PBE type (if needed)

int Vpmg_fillco Vpmg thee,
Vsurf_Meth  surfMeth,
double  splineWin,
Vchrg_Meth  chargeMeth,
int  useDielXMap,
Vgrid dielXMap,
int  useDielYMap,
Vgrid dielYMap,
int  useDielZMap,
Vgrid dielZMap,
int  useKappaMap,
Vgrid kappaMap,
int  useChargeMap,
Vgrid chargeMap
 

Fill the coefficient arrays prior to solving the equation.

Author:
Nathan Baker
Returns:
1 if successful, 0 otherwise
Parameters:
surfMeth  Vpmg object
splineWin  Surface discretization method
chargeMeth  Spline window (in A) for surfMeth = VSM_SPLINE
useDielXMap  Charge discretization method
dielXMap  Boolean to use dielectric map argument
useDielYMap  External dielectric map
dielYMap  Boolean to use dielectric map argument
useDielZMap  External dielectric map
dielZMap  Boolean to use dielectric map argument
useKappaMap  External dielectric map
kappaMap  Boolean to use kappa map argument
useChargeMap  External kappa map
chargeMap  Boolean to use charge map argument External charge map

int Vpmg_force Vpmg thee,
double *  force,
int  atomID,
Vsurf_Meth  srfm,
Vchrg_Meth  chgm
 

Calculate the total force on the specified atom in units of k_B T/AA.

Author:
Nathan Baker
Note:
  • Using the force evaluation methods of Im et al (Roux group), Comput Phys Commun, 111, 59--75 (1998). However, this gives the whole (self-interactions included) force -- reaction field forces will have to be calculated at higher level.
  • No contributions are made from higher levels of focusing.
Returns:
1 if successful, 0 otherwise
Parameters:
force  Vpmg object
atomID  3*sizeof(double) space to hold the force in units of k_B T/AA
srfm  Valist ID of desired atom
chgm  Surface discretization method Charge discretization method

int Vpmg_ibForce Vpmg thee,
double *  force,
int  atomID,
Vsurf_Meth  srfm
 

Calculate the osmotic pressure on the specified atom in units of k_B T/AA.

Author:
Nathan Baker
Note:
  • Using the force evaluation methods of Im et al (Roux group), Comput Phys Commun, 111, 59--75 (1998). However, this gives the whole (self-interactions included) force -- reaction field forces will have to be calculated at higher level.
  • No contributions are made from higher levels of focusing.
Returns:
1 if successful, 0 otherwise
Parameters:
force  Vpmg object
atomID  3*sizeof(double) space to hold the boundary force in units of k_B T/AA
srfm  Valist ID of desired atom Surface discretization method

unsigned long int Vpmg_memChk Vpmg thee  ) 
 

Return the memory used by this structure (and its contents) in bytes.

Author:
Nathan Baker
Returns:
The memory used by this structure and its contents in bytes
Parameters:
thee  Object for memory check

double Vpmg_npEnergy Vpmg thee,
int  extFlag
 

Get the "apolar" energy.

Using the dielectric map at the finest mesh level, calculate the surface area in a manner consistent with the force evaluation routines of Im et al (see Vpmg_dbnpForce and Vpmg_dielGradNorm):

\[ A = \frac{1}{\epsilon_s-\epsilon_p} \int \| \nabla \epsilon \| dx \]

where epsilon is the dielectric parameter, epsilon_s is the dielectric constant for the solvent and epsilon_p is the dielectric constant for the protein. The apolar energy is then,

\[G_{np} = \gamma S \]

where gamma is the apolar coefficient set in Vpbe (see Vpbe_ctor). The energy is returned in units of k_b T.

Author:
Nathan Baker
Note:
I personally feel that this routine should not find its way into the main APBS driver. In this case, the apolar energy is calculated in a manner consistent with the force evaluation, but it is not the only possible apolar energy definition... The value of this observable may be modified by setting restrictions on the subdomain over which it is calculated. Such limits can be set via Vpmg_setPart and are generally useful for parallel runs.
Returns:
The apolar energy in units of k_B T.
Parameters:
extFlag  Vpmg object If this was a focused calculation, include (1 -- for serial calculations) or ignore (0 -- for parallel calculations) energy contributions from outside the focusing domain

void Vpmg_printColComp Vpmg thee,
char  path[72],
char  title[72],
char  mxtype[3],
int  flag
 

Print out a column-compressed sparse matrix in Harwell-Boeing format.

Author:
Nathan Baker
Bug:
Can this path variable be replaced with a Vio socket?
Parameters:
path  Vpmg object
title  The file to which the matrix is to be written
mxtype  The title of the matrix
flag  The type of REAL-valued matrix, a 3-character string of the form "R_A" where the '_' can be one of:
  • S: symmetric matrix
  • U: unsymmetric matrix
  • H: Hermitian matrix
  • Z: skew-symmetric matrix
  • R: rectangular matrix The operator to compress:
  • 0: Poisson operator
  • 1: Linearization of the full Poisson-Boltzmann operator around the current solution

double Vpmg_qfAtomEnergy Vpmg thee,
Vatom atom
 

Get the per-atom "fixed charge" contribution to the electrostatic energy.

Using the solution at the finest mesh level, get the electrostatic energy due to the interaction of the fixed charges with the potential:

\[ G = q u(r), \]

where q$ is the charge and r is the location of the atom of interest. The result is returned in units of k_B T. Clearly, no self-interaction terms are removed. A factor a 1/2 has to be included to convert this to a real energy.

Author:
Nathan Baker
Note:
The value of this observable may be modified by setting restrictions on the subdomain over which it is calculated. Such limits can be set via Vpmg_setPart and are generally useful for parallel runs.
Returns:
The fixed charge electrostatic energy in units of k_B T.
Parameters:
atom  The Vpmg object The atom for energy calculations

double Vpmg_qfEnergy Vpmg thee,
int  extFlag
 

Get the "fixed charge" contribution to the electrostatic energy.

Using the solution at the finest mesh level, get the electrostatic energy due to the interaction of the fixed charges with the potential:

\[ G = \sum_i q_i u(r_i) \]

and return the result in units of k_B T. Clearly, no self-interaction terms are removed. A factor a 1/2 has to be included to convert this to a real energy.

Author:
Nathan Baker
Note:
The value of this observable may be modified by setting restrictions on the subdomain over which it is calculated. Such limits can be set via Vpmg_setPart and are generally useful for parallel runs.
Returns:
The fixed charge electrostatic energy in units of k_B T.
Parameters:
extFlag  Vpmg object If this was a focused calculation, include (1 -- for serial calculations) or ignore (0 -- for parallel calculations) energy contributions from outside the focusing domain

int Vpmg_qfForce Vpmg thee,
double *  force,
int  atomID,
Vchrg_Meth  chgm
 

Calculate the "charge-field" force on the specified atom in units of k_B T/AA.

Author:
Nathan Baker
Note:
  • Using the force evaluation methods of Im et al (Roux group), Comput Phys Commun, 111, 59--75 (1998). However, this gives the whole (self-interactions included) force -- reaction field forces will have to be calculated at higher level.
  • No contributions are made from higher levels of focusing.
Returns:
1 if sucessful, 0 otherwise
Parameters:
force  Vpmg object
atomID  3*sizeof(double) space to hold the force in units of k_B T/A
chgm  Valist ID of desired atom Charge discretization method

double Vpmg_qmEnergy Vpmg thee,
int  extFlag
 

Get the "mobile charge" contribution to the electrostatic energy.

Using the solution at the finest mesh level, get the electrostatic energy due to the interaction of the mobile charges with the potential:

\[ G = \frac{1}{4 I_s} \sum_i c_i q_i^2 \int \kappa^2(x) e^{-q_i u(x)} dx \]

for the NPBE and

\[ G = \frac{1}{2} \int \overline{\kappa}^2(x) u^2(x) dx \]

for the LPBE. Here i denotes the counterion species, I_s is the bulk ionic strength, kappa^2(x) is the modified Debye-Huckel parameter, c_i is the concentration of species i, q_i is the charge of species i, and u(x) is the dimensionless electrostatic potential. The energy is scaled to units of k_b T.

Author:
Nathan Baker
Note:
The value of this observable may be modified by setting restrictions on the subdomain over which it is calculated. Such limits can be set via Vpmg_setPart and are generally useful for parallel runs.
Returns:
The mobile charge electrostatic energy in units of k_B T.
Parameters:
extFlag  Vpmg object If this was a focused calculation, include (1 -- for serial calculations) or ignore (0 -- for parallel calculations) energy contributions from outside the focusing domain

void Vpmg_setPart Vpmg thee,
double  lowerCorner[3],
double  upperCorner[3],
int  bflags[6]
 

Set partition information which restricts the calculation of observables to a (rectangular) subset of the problem domain.

Author:
Nathan Baker
Parameters:
lowerCorner  Vpmg object
upperCorner  Partition lower corner
bflags  Partition upper corner Booleans indicating whether a particular processor is on the boundary with another partition. 0 if the face is not bounded (next to) another partition, and 1 otherwise.

int Vpmg_solve Vpmg thee  ) 
 

Solve the PBE using PMG.

Author:
Nathan Baker
Returns:
1 if successful, 0 otherwise
Parameters:
thee  Vpmg object

int Vpmg_solveLaplace Vpmg thee  ) 
 

Solve Poisson's equation with a homogeneous Laplacian operator using the solvent dielectric constant. This solution is performed by a sine wave decomposition.

Author:
Nathan Baker
Returns:
1 if successful, 0 otherwise
Note:
This function is really only for testing purposes as the PMG multigrid solver can solve the homogeneous system much more quickly. Perhaps we should implement an FFT version at some point...
Parameters:
thee  Vpmg object

void Vpmg_unsetPart Vpmg thee  ) 
 

Remove partition restrictions.

Author:
Nathan Baker
Parameters:
thee  Vpmg object


Generated on Tue Dec 6 10:05:59 2005 for APBS by doxygen 1.3.5