APBS 3.0.0
|
Contains declarations for class Vfetk. More...
#include "apbscfg.h"
#include "maloc/maloc.h"
#include "mc/mc.h"
#include "generic/vhal.h"
#include "generic/vatom.h"
#include "generic/vpbe.h"
#include "generic/vunit.h"
#include "generic/vgreen.h"
#include "generic/vcap.h"
#include "generic/pbeparm.h"
#include "generic/femparm.h"
#include "fem/vcsm.h"
Go to the source code of this file.
Data Structures | |
struct | sVfetk |
Contains public data members for Vfetk class/module. More... | |
struct | sVfetk_LocalVar |
Vfetk LocalVar subclass. More... | |
Typedefs | |
typedef enum eVfetk_LsolvType | Vfetk_LsolvType |
Declare FEMparm_LsolvType type. | |
typedef enum eVfetk_MeshLoad | Vfetk_MeshLoad |
Declare FEMparm_GuessType 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 struct sVfetk | Vfetk |
Declaration of the Vfetk class as the Vfetk structure. | |
typedef struct sVfetk_LocalVar | Vfetk_LocalVar |
Declaration of the Vfetk_LocalVar subclass as the Vfetk_LocalVar structure. | |
Enumerations | |
enum | eVfetk_LsolvType { VLT_SLU =0 , VLT_MG =1 , VLT_CG =2 , VLT_BCG =3 } |
Linear solver type. More... | |
enum | eVfetk_MeshLoad { VML_DIRICUBE , VML_NEUMCUBE , VML_EXTERNAL } |
Mesh loading operation. 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 | |
VEXTERNC Gem * | Vfetk_getGem (Vfetk *thee) |
Get a pointer to the Gem (grid manager) object. | |
VEXTERNC AM * | Vfetk_getAM (Vfetk *thee) |
Get a pointer to the AM (algebra manager) object. | |
VEXTERNC Vpbe * | Vfetk_getVpbe (Vfetk *thee) |
Get a pointer to the Vpbe (PBE manager) object. | |
VEXTERNC Vcsm * | Vfetk_getVcsm (Vfetk *thee) |
Get a pointer to the Vcsm (charge-simplex map) object. | |
VEXTERNC int | Vfetk_getAtomColor (Vfetk *thee, int iatom) |
Get the partition information for a particular atom. | |
VEXTERNC Vfetk * | Vfetk_ctor (Vpbe *pbe, Vhal_PBEType type) |
Constructor for Vfetk object. | |
VEXTERNC int | Vfetk_ctor2 (Vfetk *thee, Vpbe *pbe, Vhal_PBEType type) |
FORTRAN stub constructor for Vfetk object. | |
VEXTERNC void | Vfetk_dtor (Vfetk **thee) |
Object destructor. | |
VEXTERNC void | Vfetk_dtor2 (Vfetk *thee) |
FORTRAN stub object destructor. | |
VEXTERNC double * | Vfetk_getSolution (Vfetk *thee, int *length) |
Create an array containing the solution (electrostatic potential in units of ![]() | |
VEXTERNC void | Vfetk_setParameters (Vfetk *thee, PBEparm *pbeparm, FEMparm *feparm) |
Set the parameter objects. | |
VPUBLIC double | Vfetk_energy (Vfetk *thee, int color, int nonlin) |
Return the total electrostatic energy. | |
VEXTERNC double | Vfetk_dqmEnergy (Vfetk *thee, int color) |
Get the "mobile charge" and "polarization" contributions to the electrostatic energy. | |
VEXTERNC double | Vfetk_qfEnergy (Vfetk *thee, int color) |
Get the "fixed charge" contribution to the electrostatic energy. | |
VEXTERNC unsigned long int | Vfetk_memChk (Vfetk *thee) |
Return the memory used by this structure (and its contents) in bytes. | |
VEXTERNC void | Vfetk_setAtomColors (Vfetk *thee) |
Transfer color (partition ID) information frmo a partitioned mesh to the atoms. | |
VEXTERNC void | Bmat_printHB (Bmat *thee, char *fname) |
Writes a Bmat to disk in Harwell-Boeing sparse matrix format. | |
VEXTERNC Vrc_Codes | Vfetk_genCube (Vfetk *thee, double center[3], double length[3], Vfetk_MeshLoad meshType) |
Construct a rectangular mesh (in the current Vfetk object) | |
VEXTERNC Vrc_Codes | Vfetk_loadMesh (Vfetk *thee, double center[3], double length[3], Vfetk_MeshLoad meshType, Vio *sock) |
Loads a mesh into the Vfetk (and associated) object(s). | |
VEXTERNC PDE * | Vfetk_PDE_ctor (Vfetk *fetk) |
Constructs the FEtk PDE object. | |
VEXTERNC int | Vfetk_PDE_ctor2 (PDE *thee, Vfetk *fetk) |
Intializes the FEtk PDE object. | |
VEXTERNC void | Vfetk_PDE_dtor (PDE **thee) |
Destroys FEtk PDE object. | |
VEXTERNC void | Vfetk_PDE_dtor2 (PDE *thee) |
FORTRAN stub: destroys FEtk PDE object. | |
VEXTERNC void | Vfetk_PDE_initAssemble (PDE *thee, int ip[], double rp[]) |
Do once-per-assembly initialization. | |
VEXTERNC void | Vfetk_PDE_initElement (PDE *thee, int elementType, int chart, double tvx[][VAPBS_DIM], void *data) |
Do once-per-element initialization. | |
VEXTERNC void | Vfetk_PDE_initFace (PDE *thee, int faceType, int chart, double tnvec[]) |
Do once-per-face initialization. | |
VEXTERNC void | Vfetk_PDE_initPoint (PDE *thee, int pointType, int chart, double txq[], double tU[], double tdU[][VAPBS_DIM]) |
Do once-per-point initialization. | |
VEXTERNC void | Vfetk_PDE_Fu (PDE *thee, int key, double F[]) |
Evaluate strong form of PBE. For interior points, this is: | |
VEXTERNC 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: | |
VEXTERNC 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: | |
VEXTERNC 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. | |
VEXTERNC void | Vfetk_PDE_u_D (PDE *thee, int type, int chart, double txq[], double F[]) |
Evaluate the Dirichlet boundary condition at the given point. | |
VEXTERNC 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. | |
VEXTERNC void | Vfetk_PDE_bisectEdge (int dim, int dimII, int edgeType, int chart[], double vx[][VAPBS_DIM]) |
Define the way manifold edges are bisected. | |
VEXTERNC 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. | |
VEXTERNC 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. | |
VEXTERNC 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. | |
VEXTERNC double | Vfetk_PDE_Ju (PDE *thee, int key) |
Energy functional. This returns the energy (less delta function terms) in the form: | |
VEXTERNC 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) | |
VEXTERNC 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. | |
VEXTERNC 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. | |
VEXTERNC void | Vfetk_readMesh (Vfetk *thee, int skey, Vio *sock) |
Read in mesh and initialize associated internal structures. | |
VEXTERNC void | Vfetk_dumpLocalVar () |
Debugging routine to print out local variables used by PDE object. | |
VEXTERNC int | Vfetk_fillArray (Vfetk *thee, Bvec *vec, Vdata_Type type) |
Fill an array with the specified data. | |
VEXTERNC 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. | |
VEXTERNC Vrc_Codes | Vfetk_loadGem (Vfetk *thee, Gem *gm) |
Load a Gem geometry manager object into Vfetk. | |
Contains declarations for class Vfetk.
* * APBS -- Adaptive Poisson-Boltzmann Solver * * Nathan A. Baker (nathan.baker@pnnl.gov) * Pacific Northwest National Laboratory * * Additional contributing authors listed in the code documentation. * * Copyright (c) 2010-2020 Battelle Memorial Institute. Developed at the * Pacific Northwest National Laboratory, operated by Battelle Memorial * Institute, Pacific Northwest Division for the U.S. Department of Energy. * * Portions Copyright (c) 2002-2010, Washington University in St. Louis. * Portions Copyright (c) 2002-2010, Nathan A. Baker. * Portions Copyright (c) 1999-2002, The Regents of the University of * California. * Portions Copyright (c) 1995, Michael Holst. * All rights reserved. * * Redistribution and use in source and binary forms, with or without * modification, are permitted provided that the following conditions are met: * * Redistributions of source code must retain the above copyright notice, this * list of conditions and the following disclaimer. * * Redistributions in binary form must reproduce the above copyright notice, * this list of conditions and the following disclaimer in the documentation * and/or other materials provided with the distribution. * * Neither the name of the developer nor the names of its contributors may be * used to endorse or promote products derived from this software without * specific prior written permission. * * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF * THE POSSIBILITY OF SUCH DAMAGE. * *
Definition in file vfetk.h.