#include "apbscfg.h"
#include "apbs/vfetk.h"
Include dependency graph for vfetk.c:
Defines | |
#define | VRINGMAX 1000 |
Maximum number of simplices in a simplex ring. | |
#define | VATOMMAX 1000000 |
Maximum number of atoms associated with a vertex. | |
Functions | |
VPUBLIC Gem * | Vfetk_getGem (Vfetk *thee) |
Get a pointer to the Gem (grid manager) object. | |
VPUBLIC AM * | Vfetk_getAM (Vfetk *thee) |
Get a pointer to the AM (algebra manager) object. | |
VPUBLIC Vpbe * | Vfetk_getVpbe (Vfetk *thee) |
Get a pointer to the Vpbe (PBE manager) object. | |
VPUBLIC Vcsm * | Vfetk_getVcsm (Vfetk *thee) |
Get a pointer to the Vcsm (charge-simplex map) object. | |
VPUBLIC int | Vfetk_getAtomColor (Vfetk *thee, int iatom) |
Get the partition information for a particular atom. | |
VPUBLIC Vfetk * | Vfetk_ctor (Vpbe *pbe, Vhal_PBEType type) |
Constructor for Vfetk object. | |
VPUBLIC int | Vfetk_ctor2 (Vfetk *thee, Vpbe *pbe, Vhal_PBEType type) |
FORTRAN stub constructor for Vfetk object. | |
VPUBLIC void | Vfetk_setParameters (Vfetk *thee, PBEparm *pbeparm, FEMparm *feparm) |
Set the parameter objects. | |
VPUBLIC void | Vfetk_dtor (Vfetk **thee) |
Object destructor. | |
VPUBLIC void | Vfetk_dtor2 (Vfetk *thee) |
FORTRAN stub object destructor. | |
VPUBLIC double * | Vfetk_getSolution (Vfetk *thee, int *length) |
Create an array containing the solution (electrostatic potential in units of ![]() | |
VPUBLIC double | Vfetk_energy (Vfetk *thee, int color, int nonlin) |
Return the total electrostatic energy. | |
VPUBLIC double | Vfetk_qfEnergy (Vfetk *thee, int color) |
Get the "fixed charge" contribution to the electrostatic energy. | |
VPUBLIC double | Vfetk_dqmEnergy (Vfetk *thee, int color) |
Get the "mobile charge" and "polarization" contributions to the electrostatic energy. | |
VPUBLIC void | Vfetk_setAtomColors (Vfetk *thee) |
Transfer color (partition ID) information frmo a partitioned mesh to the atoms. | |
VPUBLIC unsigned long int | Vfetk_memChk (Vfetk *thee) |
Return the memory used by this structure (and its contents) in bytes. | |
VPUBLIC void | Vfetk_readMesh (Vfetk *thee, int skey, Vio *sock) |
Read in mesh and initialize associated internal structures. | |
VPUBLIC void | Bmat_printHB (Bmat *thee, char *fname) |
Writes a Bmat to disk in Harwell-Boeing sparse matrix format. | |
VPUBLIC PDE * | Vfetk_PDE_ctor (Vfetk *fetk) |
Constructs the FEtk PDE object. | |
VPUBLIC int | Vfetk_PDE_ctor2 (PDE *thee, Vfetk *fetk) |
Intializes the FEtk PDE object. | |
VPUBLIC void | Vfetk_PDE_dtor (PDE **thee) |
Destroys FEtk PDE object. | |
VPUBLIC void | Vfetk_PDE_dtor2 (PDE *thee) |
FORTRAN stub: destroys FEtk PDE object. | |
VPUBLIC void | Vfetk_PDE_initAssemble (PDE *thee, int ip[], double rp[]) |
Do once-per-assembly initialization. | |
VPUBLIC void | Vfetk_PDE_initFace (PDE *thee, int faceType, int chart, double tnvec[]) |
Do once-per-face initialization. | |
VPUBLIC void | Vfetk_PDE_Fu (PDE *thee, int key, double F[]) |
Evaluate strong form of PBE. For interior points, this is:
where
where | |
VPUBLIC 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:
| |
VPUBLIC 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. | |
VPUBLIC void | Vfetk_PDE_u_D (PDE *thee, int type, int chart, double txq[], double F[]) |
Evaluate the Dirichlet boundary condition at the given point. | |
VPUBLIC 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. | |
VPUBLIC double | Vfetk_PDE_Ju (PDE *thee, int key) |
Energy functional. This returns the energy (less delta function terms) in the form:
for a 1:1 electrolyte where | |
VPUBLIC 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). | |
VPUBLIC 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. | |
VPUBLIC 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. | |
VPUBLIC void | Vfetk_dumpLocalVar () |
Debugging routine to print out local variables used by PDE object. | |
VPUBLIC int | Vfetk_fillArray (Vfetk *thee, Bvec *vec, Vdata_Type type) |
Fill an array with the specified data. | |
VPUBLIC 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. |
* * APBS -- Adaptive Poisson-Boltzmann Solver * * Nathan A. Baker (baker@biochem.wustl.edu) * Dept. of Biochemistry and Molecular Biophysics * Center for Computational Biology * Washington University in St. Louis * * Additional contributing authors listed in the code documentation. * * Copyright (c) 2002-2007. Washington University in St. Louis. * All Rights Reserved. * Portions Copyright (c) 1999-2002. The Regents of the University of * California. * Portions Copyright (c) 1995. Michael Holst. * * This file is part of APBS. * * APBS is free software; you can redistribute it and/or modify * it under the terms of the GNU General Public License as published by * the Free Software Foundation; either version 2 of the License, or * (at your option) any later version. * * APBS is distributed in the hope that it will be useful, * but WITHOUT ANY WARRANTY; without even the implied warranty of * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * GNU General Public License for more details. * * You should have received a copy of the GNU General Public License * along with APBS; if not, write to the Free Software * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA * * Linking APBS statically or dynamically with other modules is making a * combined work based on APBS. Thus, the terms and conditions of the GNU * General Public License cover the whole combination. * * SPECIAL GPL EXCEPTION * In addition, as a special exception, the copyright holders of APBS * give you permission to combine the APBS program with free software * programs and libraries that are released under the GNU LGPL or with * code included in releases of ISIM, Ion Simulator Interface, PMV, PyMOL * SMOL, VMD, and Vision. Such combined software may be linked with APBS and * redistributed together in original or modified form as mere aggregation * without requirement that the entire work be under the scope of the GNU * General Public License. This special exception permission is also extended * to any software listed in the SPECIAL GPL EXCEPTION clauses by the PMG, * FEtk, MC, or MALOC libraries. * * Note that people who make modified versions of APBS are not obligated * to grant this special exception for their modified versions; it is * their choice whether to do so. The GNU General Public License gives * permission to release a modified version without this exception; this * exception also makes it possible to release a modified version which * carries forward this exception. * *