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

Vfetk class

FEtk master class (interface between FEtk and APBS). More...


Files

file  vfetk.h
 Contains declarations for class Vfetk.


Data Structures

struct  sVfetk
 Contains public data members for Vfetk class/module. More...

struct  sVfetk_LocalVar
 Vfetk LocalVar subclass Contains variables used when solving the PDE with FEtk. More...


Typedefs

typedef enum eVfetk_LsolvType Vfetk_LsolvType
 Declare FEMparm_LsolvType type.

typedef enum eVfetk_NsolvType Vfetk_NsolvType
 Declare FEMparm_NsolvType type.

typedef enum eVfetk_GuessType Vfetk_GuessType
 Declare FEMparm_GuessType type.

typedef enum eVfetk_PrecType Vfetk_PrecType
 Declare FEMparm_GuessType type.

typedef sVfetk_LocalVar Vfetk_LocalVar
 Declaration of the Vfetk_LocalVar subclass as the Vfetk_LocalVar structure.

typedef sVfetk Vfetk
 Declaration of the Vfetk class as the Vfetk structure.


Enumerations

enum  eVfetk_LsolvType {
  VLT_SLU = 0,
  VLT_MG = 1,
  VLT_CG = 2,
  VLT_BCG = 3
}
 Linear solver type. More...

enum  eVfetk_NsolvType {
  VNT_NEW = 0,
  VNT_INC = 1,
  VNT_ARC = 2
}
 Non-linear solver type. More...

enum  eVfetk_GuessType {
  VGT_ZERO = 0,
  VGT_DIRI = 1,
  VGT_PREV = 2
}
 Initial guess type. More...

enum  eVfetk_PrecType {
  VPT_IDEN = 0,
  VPT_DIAG = 1,
  VPT_MG = 2
}
 Preconditioner type. More...


Functions

Gem * Vfetk_getGem (Vfetk *thee)
 Get a pointer to the Gem (grid manager) object.

AM * Vfetk_getAM (Vfetk *thee)
 Get a pointer to the AM (algebra manager) object.

VpbeVfetk_getVpbe (Vfetk *thee)
 Get a pointer to the Vpbe (PBE manager) object.

VcsmVfetk_getVcsm (Vfetk *thee)
 Get a pointer to the Vcsm (charge-simplex map) object.

int Vfetk_getAtomColor (Vfetk *thee, int iatom)
 Get the partition information for a particular atom.

VfetkVfetk_ctor (Vpbe *pbe, Vhal_PBEType type)
 Constructor for Vfetk object.

int Vfetk_ctor2 (Vfetk *thee, Vpbe *pbe, Vhal_PBEType type)
 FORTRAN stub constructor for Vfetk object.

void Vfetk_dtor (Vfetk **thee)
 Object destructor.

void Vfetk_dtor2 (Vfetk *thee)
 FORTRAN stub object destructor.

double * Vfetk_getSolution (Vfetk *thee, int *length)
 Create an array containing the solution (electrostatic potential in units of $k_B T/e$) at the finest mesh level.

void Vfetk_setParameters (Vfetk *thee, PBEparm *pbeparm, FEMparm *feparm)
 Set the parameter objects.

double Vfetk_energy (Vfetk *thee, int color, int nonlin)
 Return the total electrostatic energy.

double Vfetk_dqmEnergy (Vfetk *thee, int color)
 Get the "mobile charge" and "polarization" contributions to the electrostatic energy.

double Vfetk_qfEnergy (Vfetk *thee, int color)
 Get the "fixed charge" contribution to the electrostatic energy.

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

void Vfetk_setAtomColors (Vfetk *thee)
 Transfer color (partition ID) information frmo a partitioned mesh to the atoms.

void Bmat_printHB (Bmat *thee, char *fname)
 Writes a Bmat to disk in Harwell-Boeing sparse matrix format.

int Vfetk_genCube (Vfetk *thee, double center[VAPBS_DIM], double length[VAPBS_DIM])
 Generates a simple cubic tetrahedral mesh.

PDE * Vfetk_PDE_ctor (Vfetk *fetk)
 Constructs the FEtk PDE object.

int Vfetk_PDE_ctor2 (PDE *thee, Vfetk *fetk)
 Intializes the FEtk PDE object.

void Vfetk_PDE_dtor (PDE **thee)
 Destroys FEtk PDE object.

void Vfetk_PDE_dtor2 (PDE *thee)
 FORTRAN stub: destroys FEtk PDE object.

void Vfetk_PDE_initAssemble (PDE *thee, int ip[], double rp[])
 Do once-per-assembly initialization.

void Vfetk_PDE_initElement (PDE *thee, int elementType, int chart, double tvx[][VAPBS_DIM], void *data)
 Do once-per-element initialization.

void Vfetk_PDE_initFace (PDE *thee, int faceType, int chart, double tnvec[])
 Do once-per-face initialization.

void Vfetk_PDE_initPoint (PDE *thee, int pointType, int chart, double txq[], double tU[], double tdU[][VAPBS_DIM])
 Do once-per-point initialization.

void Vfetk_PDE_Fu (PDE *thee, int key, double F[])
 Evaluate strong form of PBE. For interior points, this is:

\[ -\nabla \cdot \epsilon \nabla u + b(u) - f \]

where $b(u)$ is the (possibly nonlinear) mobile ion term and $f$ is the source charge distribution term (for PBE) or the induced surface charge distribution (for RPBE). For an interior-boundary (simplex face) point, this is:

\[ [\epsilon(x) \nabla u(x) \cdot n(x)]_{x=0^+} - [\epsilon(x) \nabla u(x) \cdot n(x)]_{x=0^-} \]

where $n(x)$ is the normal to the simplex face and the term represents the jump in dielectric displacement across the face. There is no outer-boundary contribution for this problem.

double Vfetk_PDE_Fu_v (PDE *thee, int key, double V[], double dV[][VAPBS_DIM])
 This is the weak form of the PBE; i.e. the strong form integrated with a test function to give:

\[ \int_\Omega \left[ \epsilon \nabla u \cdot \nabla v + b(u) v - f v \right] dx \]

where $b(u)$ denotes the mobile ion term.

double Vfetk_PDE_DFu_wv (PDE *thee, int key, double W[], double dW[][VAPBS_DIM], double V[], double dV[][VAPBS_DIM])
 This is the linearization of the weak form of the PBE; e.g., for use in a Newton iteration. This is the functional linearization of the strong form integrated with a test function to give:

\[ \int_\Omega \left[ \epsilon \nabla w \cdot \nabla v + b'(u) w v - f v \right] dx \]

where $b'(u)$ denotes the functional derivation of the mobile ion term.

void Vfetk_PDE_delta (PDE *thee, int type, int chart, double txq[], void *user, double F[])
 Evaluate a (discretized) delta function source term at the given point.

void Vfetk_PDE_u_D (PDE *thee, int type, int chart, double txq[], double F[])
 Evaluate the Dirichlet boundary condition at the given point.

void Vfetk_PDE_u_T (PDE *thee, int type, int chart, double txq[], double F[])
 Evaluate the "true solution" at the given point for comparison with the numerical solution.

void Vfetk_PDE_bisectEdge (int dim, int dimII, int edgeType, int chart[], double vx[][VAPBS_DIM])
 Define the way manifold edges are bisected.

void Vfetk_PDE_mapBoundary (int dim, int dimII, int vertexType, int chart, double vx[VAPBS_DIM])
 Map a boundary point to some pre-defined shape.

int Vfetk_PDE_markSimplex (int dim, int dimII, int simplexType, int faceType[VAPBS_NVS], int vertexType[VAPBS_NVS], int chart[], double vx[][VAPBS_DIM], void *simplex)
 User-defined error estimator -- in our case, a geometry-based refinement method; forcing simplex refinement at the dielectric boundary and (for non-regularized PBE) the charges.

void Vfetk_PDE_oneChart (int dim, int dimII, int objType, int chart[], double vx[][VAPBS_DIM], int dimV)
 Unify the chart for different coordinate systems -- a no-op for us.

double Vfetk_PDE_Ju (PDE *thee, int key)
 Energy functional. This returns the energy (less delta function terms) in the form:

\[ c^{-1}/2 \int (\epsilon (\nabla u)^2 + \kappa^2 (cosh u - 1)) dx \]

for a 1:1 electrolyte where $c$ is the output from Vpbe_getZmagic.

void Vfetk_externalUpdateFunction (SS **simps, int num)
 External hook to simplex subdivision routines in Gem. Called each time a simplex is subdivided (we use it to update the charge-simplex map).

int Vfetk_PDE_simplexBasisInit (int key, int dim, int comp, int *ndof, int dof[])
 Initialize the bases for the trial or the test space, for a particular component of the system, at all quadrature points on the master simplex element.

void Vfetk_PDE_simplexBasisForm (int key, int dim, int comp, int pdkey, double xq[], double basis[])
 Evaluate the bases for the trial or test space, for a particular component of the system, at all quadrature points on the master simplex element.

void Vfetk_readMesh (Vfetk *thee, int skey, Vio *sock)
 Read in mesh and initialize associated internal structures.

void Vfetk_dumpLocalVar ()
 Debugging routine to print out local variables used by PDE object.

int Vfetk_fillArray (Vfetk *thee, Bvec *vec, Vdata_Type type)
 Fill an array with the specified data.

int Vfetk_write (Vfetk *thee, const char *iodev, const char *iofmt, const char *thost, const char *fname, Bvec *vec, Vdata_Format format)
 Write out data.


Detailed Description

FEtk master class (interface between FEtk and APBS).


Enumeration Type Documentation

enum eVfetk_GuessType
 

Initial guess type.

Note:
Do not change these values; they correspond to settings in FEtk
Enumeration values:
VGT_ZERO  Zero initial guess
VGT_DIRI  Dirichlet boundary condition initial guess
VGT_PREV  Previous level initial guess

enum eVfetk_LsolvType
 

Linear solver type.

Note:
Do not change these values; they correspond to settings in FEtk
Enumeration values:
VLT_SLU  SuperLU direct solve
VLT_MG  Multigrid
VLT_CG  Conjugate gradient
VLT_BCG  BiCGStab

enum eVfetk_NsolvType
 

Non-linear solver type.

Note:
Do not change these values; they correspond to settings in FEtk
Enumeration values:
VNT_NEW  Newton solver
VNT_INC  Incremental
VNT_ARC  Psuedo-arclength

enum eVfetk_PrecType
 

Preconditioner type.

Note:
Do not change these values; they correspond to settings in FEtk
Enumeration values:
VPT_IDEN  Identity matrix
VPT_DIAG  Diagonal scaling
VPT_MG  Multigrid


Function Documentation

void Bmat_printHB Bmat *  thee,
char *  fname
 

Writes a Bmat to disk in Harwell-Boeing sparse matrix format.

Author:
Stephen Bond
Note:
This is a friend function of Bmat
Bug:
Hardwired to only handle the single block symmetric case.
Parameters:
fname  The matrix to write Filename for output

Vfetk* Vfetk_ctor Vpbe pbe,
Vhal_PBEType  type
 

Constructor for Vfetk object.

Author:
Nathan Baker
Returns:
Pointer to newly allocated Vfetk object
Note:
This sets up the Gem, AM, and Aprx FEtk objects but does not create a mesh. The easiest way to create a mesh is to then call Vfetk_genCube
Parameters:
type  Vpbe (PBE manager object) Version of PBE to solve

int Vfetk_ctor2 Vfetk thee,
Vpbe pbe,
Vhal_PBEType  type
 

FORTRAN stub constructor for Vfetk object.

Author:
Nathan Baker
Returns:
1 if successful, 0 otherwise
Note:
This sets up the Gem, AM, and Aprx FEtk objects but does not create a mesh. The easiest way to create a mesh is to then call Vfetk_genCube
Parameters:
pbe  Vfetk object memory
type  PBE manager object Version of PBE to solve

double Vfetk_dqmEnergy Vfetk thee,
int  color
 

Get the "mobile charge" and "polarization" contributions 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 and polarization of the dielectric medium:

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

for the NPBE and

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

for the LPBE. Here $i$ denotes the counterion species, $I_s$ is the bulk ionic strength, $\overline{\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$, $\epsilon$ is the dielectric function, and $u(x)$ is the dimensionless electrostatic potential. The energy is scaled to units of $k_b T$.

Author:
Nathan Baker
Parameters:
thee Vfetk object
color Partition restriction for energy evaluation, only used if non-negative
Returns:
The "mobile charge" and "polarization" contributions to the electrostatic energy in units of $k_B T$.
Parameters:
color  The Vfetk object Partition restriction for energy calculation; if non-negative, energy calculation is restricted to the specified partition (indexed by simplex and atom colors

void Vfetk_dtor Vfetk **  thee  ) 
 

Object destructor.

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

void Vfetk_dtor2 Vfetk thee  ) 
 

FORTRAN stub object destructor.

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

void Vfetk_dumpLocalVar  ) 
 

Debugging routine to print out local variables used by PDE object.

Author:
Nathan Baker
Bug:
This function is not thread-safe

double Vfetk_energy Vfetk thee,
int  color,
int  nonlin
 

Return the total electrostatic energy.

Using the solution at the finest mesh level, get the electrostatic energy using the free energy functional for the Poisson-Boltzmann equation without removing any self-interaction terms (i.e., removing the reference state of isolated charges present in an infinite dielectric continuum with the same relative permittivity as the interior of the protein) and return the result in units of $k_B T$. The argument color allows the user to control the partition on which this energy is calculated; if (color == -1) no restrictions are used. The solution is obtained from the finest level of the passed AM object, but atomic data from the Vfetk object is used to calculate the energy.

Author:
Nathan Baker
Returns:
Total electrostatic energy in units of $k_B T$.
Parameters:
color  THe Vfetk object
nonlin  Partition restriction for energy calculation; if non-negative, energy calculation is restricted to the specified partition (indexed by simplex and atom colors If 1, the NPBE energy functional is used; otherwise, the LPBE energy functional is used.

void Vfetk_externalUpdateFunction SS **  simps,
int  num
 

External hook to simplex subdivision routines in Gem. Called each time a simplex is subdivided (we use it to update the charge-simplex map).

Author:
Nathan Baker
Bug:
This function is not thread-safe.
Parameters:
num  List of parent (simps[0]) and children (remainder) simplices Number of simplices in list

int Vfetk_fillArray Vfetk thee,
Bvec *  vec,
Vdata_Type  type
 

Fill an array with the specified data.

Author:
Nathan Baker
Note:
This function is thread-safe
Bug:
Several values of type are not implemented
Returns:
1 if successful, 0 otherwise
Parameters:
vec  The Vfetk object with the data
type  The vector to hold the data THe type of data to write

int Vfetk_genCube Vfetk thee,
double  center[VAPBS_DIM],
double  length[VAPBS_DIM]
 

Generates a simple cubic tetrahedral mesh.

Author:
Nathan Baker (based on Mike Holst's Gem_makeCube code)
Returns:
1 if successful, 0 otherwise
Parameters:
center  The Vfetk object
length  Desired mesh center (in Å) Desired mesh length (in Å)

AM* Vfetk_getAM Vfetk thee  ) 
 

Get a pointer to the AM (algebra manager) object.

Author:
Nathan Baker
Returns:
Pointer to the AM (algebra manager) object
Parameters:
thee  The Vfetk object

int Vfetk_getAtomColor Vfetk thee,
int  iatom
 

Get the partition information for a particular atom.

Author:
Nathan Baker
Note:
Friend function of Vatom
Parameters:
thee Vfetk object
iatom Valist atom ID
Returns:
Partition ID
Parameters:
iatom  The Vfetk object Valist atom index

Gem* Vfetk_getGem Vfetk thee  ) 
 

Get a pointer to the Gem (grid manager) object.

Author:
Nathan Baker
Returns:
Pointer to the Gem (grid manager) object
Parameters:
thee  Vfetk object

double* Vfetk_getSolution Vfetk thee,
int *  length
 

Create an array containing the solution (electrostatic potential in units of $k_B T/e$) at the finest mesh level.

Author:
Nathan Baker and Michael Holst
Note:
The user is responsible for destroying the newly created array
Returns:
Newly created array of length "length" (see above); the user is responsible for destruction
Parameters:
length  Vfetk object with solution Ste to length of the newly created solution array

Vcsm* Vfetk_getVcsm Vfetk thee  ) 
 

Get a pointer to the Vcsm (charge-simplex map) object.

Author:
Nathan Baker
Returns:
Pointer to the Vcsm (charge-simplex map) object
Parameters:
thee  The Vfetk object

Vpbe* Vfetk_getVpbe Vfetk thee  ) 
 

Get a pointer to the Vpbe (PBE manager) object.

Author:
Nathan Baker
Returns:
Pointer to the Vpbe (PBE manager) object
Parameters:
thee  The Vfetk object

unsigned long int Vfetk_memChk Vfetk 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  THe Vfetk object

void Vfetk_PDE_bisectEdge int  dim,
int  dimII,
int  edgeType,
int  chart[],
double  vx[][VAPBS_DIM]
 

Define the way manifold edges are bisected.

Author:
Nathan Baker and Mike Holst
Note:
This function is thread-safe.
Parameters:
dimII  Intrinsic dimension of manifold
edgeType  Embedding dimension of manifold
chart  Type of edge being refined
vx  Chart for edge vertices, used here as accessibility bitfields Edge vertex coordindates

PDE* Vfetk_PDE_ctor Vfetk fetk  ) 
 

Constructs the FEtk PDE object.

Author:
Nathan Baker
Returns:
Newly-allocated PDE object
Bug:
Not thread-safe
Parameters:
fetk  The Vfetk object

int Vfetk_PDE_ctor2 PDE *  thee,
Vfetk fetk
 

Intializes the FEtk PDE object.

Author:
Nathan Baker (with code by Mike Holst)
Returns:
1 if successful, 0 otherwise
Bug:
Not thread-safe
Parameters:
fetk  The newly-allocated PDE object The parent Vfetk object

void Vfetk_PDE_delta PDE *  thee,
int  type,
int  chart,
double  txq[],
void *  user,
double  F[]
 

Evaluate a (discretized) delta function source term at the given point.

Author:
Nathan Baker
Bug:
This function is not thread-safe
Parameters:
type  PDE object
chart  Vertex type
txq  Chart for point coordinates
user  Point coordinates
F  Vertex object pointer Set to delta function value

double Vfetk_PDE_DFu_wv PDE *  thee,
int  key,
double  W[],
double  dW[][VAPBS_DIM],
double  V[],
double  dV[][VAPBS_DIM]
 

This is the linearization of the weak form of the PBE; e.g., for use in a Newton iteration. This is the functional linearization of the strong form integrated with a test function to give:

\[ \int_\Omega \left[ \epsilon \nabla w \cdot \nabla v + b'(u) w v - f v \right] dx \]

where $b'(u)$ denotes the functional derivation of the mobile ion term.

Author:
Nathan Baker and Mike Holst
Returns:
Integrand value
Bug:
This function is not thread-safe
Parameters:
key  The PDE object
W  Integrand to evaluate (0 = interior weak form, 1 = boundary weak form)
dW  Trial function value at current point
V  Trial function gradient at current point
dV  Test function value at current point Test function gradient

void Vfetk_PDE_dtor PDE **  thee  ) 
 

Destroys FEtk PDE object.

Author:
Nathan Baker
Note:
Thread-safe
Parameters:
thee  Pointer to PDE object memory

void Vfetk_PDE_dtor2 PDE *  thee  ) 
 

FORTRAN stub: destroys FEtk PDE object.

Author:
Nathan Baker
Note:
Thread-safe
Parameters:
thee  PDE object memory

void Vfetk_PDE_Fu PDE *  thee,
int  key,
double  F[]
 

Evaluate strong form of PBE. For interior points, this is:

\[ -\nabla \cdot \epsilon \nabla u + b(u) - f \]

where $b(u)$ is the (possibly nonlinear) mobile ion term and $f$ is the source charge distribution term (for PBE) or the induced surface charge distribution (for RPBE). For an interior-boundary (simplex face) point, this is:

\[ [\epsilon(x) \nabla u(x) \cdot n(x)]_{x=0^+} - [\epsilon(x) \nabla u(x) \cdot n(x)]_{x=0^-} \]

where $n(x)$ is the normal to the simplex face and the term represents the jump in dielectric displacement across the face. There is no outer-boundary contribution for this problem.

Author:
Nathan Baker
Bug:
This function is not thread-safe

This function is not implemented (sets error to zero)

Parameters:
key  The PDE object
F  Type of point (0 = interior, 1 = boundary, 2 = interior boundary Set to value of residual

double Vfetk_PDE_Fu_v PDE *  thee,
int  key,
double  V[],
double  dV[][VAPBS_DIM]
 

This is the weak form of the PBE; i.e. the strong form integrated with a test function to give:

\[ \int_\Omega \left[ \epsilon \nabla u \cdot \nabla v + b(u) v - f v \right] dx \]

where $b(u)$ denotes the mobile ion term.

Author:
Nathan Baker and Mike Holst
Returns:
Integrand value
Bug:
This function is not thread-safe
Parameters:
key  The PDE object
V  Integrand to evaluate (0 = interior weak form, 1 = boundary weak form
dV  Test function at current point Test function derivative at current point

void Vfetk_PDE_initAssemble PDE *  thee,
int  ip[],
double  rp[]
 

Do once-per-assembly initialization.

Author:
Nathan Baker and Mike Holst
Note:
Thread-safe
Parameters:
ip  PDE object
rp  Integer parameter array (not used) Double parameter array (not used)

void Vfetk_PDE_initElement PDE *  thee,
int  elementType,
int  chart,
double  tvx[][VAPBS_DIM],
void *  data
 

Do once-per-element initialization.

Author:
Nathan Baker and Mike Holst
Todo:
Jump term is not implemented
Bug:
This function is not thread-safe
Parameters:
elementType  PDE object
chart  Material type (not used)
tvx  Chart in which the vertex coordinates are provided, used here as a bitfield to store molecular accessibility
data  Vertex coordinates Simplex pointer (hack)

void Vfetk_PDE_initFace PDE *  thee,
int  faceType,
int  chart,
double  tnvec[]
 

Do once-per-face initialization.

Author:
Nathan Baker and Mike Holst
Bug:
This function is not thread-safe
Parameters:
faceType  THe PDE object
chart  Simplex face type (interior or various boundary types)
tnvec  Chart in which the vertex coordinates are provided, used here as a bitfield for molecular accessibility Coordinates of outward normal vector for face

void Vfetk_PDE_initPoint PDE *  thee,
int  pointType,
int  chart,
double  txq[],
double  tU[],
double  tdU[][VAPBS_DIM]
 

Do once-per-point initialization.

Author:
Nathan Baker
Bug:
This function is not thread-safe

This function uses pre-defined boudnary definitions for the molecular surface.

Parameters:
pointType  The PDE object
chart  The type of point -- interior or various faces
txq  The chart in which the point coordinates are provided, used here as bitfield for molecular accessibility
tU  Point coordinates
tdU  Solution value at point Solution derivative at point

double Vfetk_PDE_Ju PDE *  thee,
int  key
 

Energy functional. This returns the energy (less delta function terms) in the form:

\[ c^{-1}/2 \int (\epsilon (\nabla u)^2 + \kappa^2 (cosh u - 1)) dx \]

for a 1:1 electrolyte where $c$ is the output from Vpbe_getZmagic.

Author:
Nathan Baker
Returns:
Energy value (in kT)
Bug:
This function is not thread-safe.
Parameters:
key  The PDE object What to evluate: interior (0) or boundary (1)?

void Vfetk_PDE_mapBoundary int  dim,
int  dimII,
int  vertexType,
int  chart,
double  vx[VAPBS_DIM]
 

Map a boundary point to some pre-defined shape.

Author:
Nathan Baker and Mike Holst
Note:
This function is thread-safe and is a no-op
Parameters:
dimII  Intrinsic dimension of manifold
vertexType  Embedding dimension of manifold
chart  Type of vertex
vx  Chart for vertex coordinates Vertex coordinates

int Vfetk_PDE_markSimplex int  dim,
int  dimII,
int  simplexType,
int  faceType[VAPBS_NVS],
int  vertexType[VAPBS_NVS],
int  chart[],
double  vx[][VAPBS_DIM],
void *  simplex
 

User-defined error estimator -- in our case, a geometry-based refinement method; forcing simplex refinement at the dielectric boundary and (for non-regularized PBE) the charges.

Author:
Nathan Baker
Returns:
1 if mark simplex for refinement, 0 otherwise
Bug:
This function is not thread-safe
Parameters:
dimII  Intrinsic manifold dimension
simplexType  Embedding manifold dimension
faceType  Type of simplex being refined
vertexType  Types of faces in simplex
chart  Types of vertices in simplex
vx  Charts for vertex coordinates
simplex  Vertex coordinates Simplex pointer

void Vfetk_PDE_oneChart int  dim,
int  dimII,
int  objType,
int  chart[],
double  vx[][VAPBS_DIM],
int  dimV
 

Unify the chart for different coordinate systems -- a no-op for us.

Author:
Nathan Baker
Note:
Thread-safe; a no-op
Parameters:
dimII  Intrinsic manifold dimension
objType  Embedding manifold dimension
chart  ???
vx  Charts of vertices' coordinates
dimV  Vertices' coordinates Number of vertices

void Vfetk_PDE_simplexBasisForm int  key,
int  dim,
int  comp,
int  pdkey,
double  xq[],
double  basis[]
 

Evaluate the bases for the trial or test space, for a particular component of the system, at all quadrature points on the master simplex element.

Author:
Mike Holst
Parameters:
dim  Basis type to evaluate (0 = trial, 1 = test, 2 = trialB, 3 = testB)
comp  Spatial dimension Which component of elliptic system to produce basis for
xq  Basis partial differential equation evaluation key:
  • 0 = evaluate basis(x,y,z)
  • 1 = evaluate basis_x(x,y,z)
  • 2 = evaluate basis_y(x,y,z)
  • 3 = evaluate basis_z(x,y,z)
  • 4 = evaluate basis_xx(x,y,z)
  • 5 = evaluate basis_yy(x,y,z)
  • 6 = evaluate basis_zz(x,y,z)
  • 7 = etc...
basis  Set to quad pt coordinate Set to all basis functions evaluated at all quadrature pts

int Vfetk_PDE_simplexBasisInit int  key,
int  dim,
int  comp,
int *  ndof,
int  dof[]
 

Initialize the bases for the trial or the test space, for a particular component of the system, at all quadrature points on the master simplex element.

Author:
Mike Holst
Note:
 *   The basis ordering is important.  For a fixed quadrature
 *   point iq, you must follow the following ordering in p[iq][],
 *   based on how you specify the degrees of freedom in dof[]:
 *   
 *   <v_0 vDF_0>,      <v_1 vDF_0>,      ..., <v_{nv} vDF_0>
 *   <v_0 vDF_1>,      <v_1 vDF_1>,      ..., <v_{nv} vDF_1>
 *                           ...
 *   <v_0 vDF_{nvDF}>, <v_0 vDF_{nvDF}>, ..., <v_{nv} vDF_{nvDF}>
 *   
 *   <e_0 eDF_0>,      <e_1 eDF_0>,      ..., <e_{ne} eDF_0>
 *   <e_0 eDF_1>,      <e_1 eDF_1>,      ..., <e_{ne} eDF_1>
 *                           ...
 *   <e_0 eDF_{neDF}>, <e_1 eDF_{neDF}>, ..., <e_{ne} eDF_{neDF}>
 *
 *   <f_0 fDF_0>,      <f_1 fDF_0>,      ..., <f_{nf} fDF_0>
 *   <f_0 fDF_1>,      <f_1 fDF_1>,      ..., <f_{nf} fDF_1>
 *                           ...
 *   <f_0 fDF_{nfDF}>, <f_1 fDF_{nfDF}>, ..., <f_{nf} fDF_{nfDF}>
 *
 *   <s_0 sDF_0>,      <s_1 sDF_0>,      ..., <s_{ns} sDF_0>
 *   <s_0 sDF_1>,      <s_1 sDF_1>,      ..., <s_{ns} sDF_1>
 *                           ...
 *   <s_0 sDF_{nsDF}>, <s_1 sDF_{nsDF}>, ..., <s_{ns} sDF_{nsDF}>
 *
 *   For example, linear elements in R^3, with one degree of freedom at each *
 *   vertex, would use the following ordering: 
 *
 *     <v_0 vDF_0>, <v_1 vDF_0>, <v_2 vDF_0>, <v_3 vDF_0> 
 *
 *   Quadratic elements in R^2, with one degree of freedom at each vertex and
 *   edge, would use the following ordering: 
 * 
 *     <v_0 vDF_0>, <v_1 vDF_0>, <v_2 vDF_0> 
 *     <e_0 eDF_0>, <e_1 eDF_0>, <e_2 eDF_0> 
 *
 *   You can use different trial and test spaces for each component of the
 *   elliptic system, thereby allowing for the use of Petrov-Galerkin methods.
 *   You MUST then tag the bilinear form symmetry entries as nonsymmetric in
 *   your PDE constructor to reflect that DF(u)(w,v) will be different from
 *   DF(u)(v,w), even if your form acts symmetrically when the same basis is
 *   used for w and v.
 * 
 *   You can also use different trial spaces for each component of the elliptic
 *   system, and different test spaces for each component of the elliptic
 *   system.  This allows you to e.g.  use a basis which is vertex-based for 
 *   one component, and a basis which is edge-based for another.  This is
 *   useful in fluid mechanics, eletromagnetics, or simply to play around with
 *   different elements.  
 *   
 *   This function is called by MC to build new master elements whenever it
 *   reads in a new mesh.  Therefore, this function does not have to be all
 *   that fast, and e.g.  could involve symbolic computation.
 *   
Parameters:
dim  Basis type to evaluate (0 = trial, 1 = test, 2 = trialB, 3 = testB)
comp  Spatial dimension
ndof  Which component of elliptic system to produce basis for?
dof  Set to the number of degrees of freedom Set to degree of freedom per v/e/f/s

void Vfetk_PDE_u_D PDE *  thee,
int  type,
int  chart,
double  txq[],
double  F[]
 

Evaluate the Dirichlet boundary condition at the given point.

Author:
Nathan Baker

Bug:
This function is hard-coded to call only multiple-sphere Debye-Hü functions.

This function is not thread-safe.

Parameters:
type  PDE object
chart  Vertex boundary type
txq  Chart for point coordinates
F  Point coordinates Set to boundary values

void Vfetk_PDE_u_T PDE *  thee,
int  type,
int  chart,
double  txq[],
double  F[]
 

Evaluate the "true solution" at the given point for comparison with the numerical solution.

Author:
Nathan Baker
Note:
This function only returns zero.
Bug:
This function is not thread-safe.
Parameters:
type  PDE object
chart  Point type
txq  Chart for point coordinates
F  Point coordinates Set to value at point

double Vfetk_qfEnergy Vfetk thee,
int  color
 

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
Parameters:
thee Vfetk object
color Partition restriction for energy evaluation, only used if non-negative
Returns:
The fixed charge electrostatic energy in units of $k_B T$.
Parameters:
color  The Vfetk object Partition restriction for energy evaluation, only used if non-negative

void Vfetk_readMesh Vfetk thee,
int  skey,
Vio *  sock
 

Read in mesh and initialize associated internal structures.

Author:
Nathan Baker
Note:
See also:
Vfetk_genCube
Parameters:
skey  THe Vfetk object
sock  The sock format key (0 = MCSF simplex format) Socket object ready for reading

void Vfetk_setAtomColors Vfetk thee  ) 
 

Transfer color (partition ID) information frmo a partitioned mesh to the atoms.

Transfer color information from partitioned mesh to the atoms. In the case that a charge is shared between two partitions, the partition color of the first simplex is selected. Due to the arbitrary nature of this selection, THIS METHOD SHOULD ONLY BE USED IMMEDIATELY AFTER PARTITIONING!!!

Warning:
This function should only be used immediately after mesh partitioning
Author:
Nathan Baker
Note:
This is a friend function of Vcsm
Parameters:
thee  THe Vfetk object

void Vfetk_setParameters Vfetk thee,
PBEparm pbeparm,
FEMparm feparm
 

Set the parameter objects.

Author:
Nathan Baker
Parameters:
pbeparm  The Vfetk object
feparm  Parameters for solution of the PBE FEM-speecific solution parameters

int Vfetk_write Vfetk thee,
const char *  iodev,
const char *  iofmt,
const char *  thost,
const char *  fname,
Bvec *  vec,
Vdata_Format  format
 

Write out data.

Author:
Nathan Baker
Parameters:
thee Vfetk object
vec FEtk Bvec vector to use
format Format for data
iodev Output device type (FILE/BUFF/UNIX/INET)
iofmt Output device format (ASCII/XDR)
thost Output hostname (for sockets)
fname Output FILE/BUFF/UNIX/INET name
Note:
This function is thread-safe
Bug:
Some values of format are not implemented
Returns:
1 if successful, 0 otherwise
Parameters:
iodev  The Vfetk object
iofmt  Output device type (FILE = file, BUFF = buffer, UNIX = unix pipe, INET = network socket)
thost  Output device format (ASCII = ascii/plaintext, XDR = xdr)
fname  Output hostname for sockets
vec  Output filename for other
format  Data vector Data format


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