1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255
|
#ifndef KERNEL_FUNCTION_H
#define KERNEL_FUNCTION_H
#include <set>
#include "LammpsInterface.h"
#include "MatrixLibrary.h"
namespace ATC {
/**
* @class KernelFunctionMgr
* @brief Base class for managing kernels
*/
class KernelFunctionMgr {
public:
/** Static instance of this class */
static KernelFunctionMgr * instance();
class KernelFunction* function(char** arg, int nargs);
protected:
KernelFunctionMgr() {};
~KernelFunctionMgr();
private:
static KernelFunctionMgr * myInstance_;
std::set<KernelFunction*> pointerSet_;
};
/**
* @class KernelFunction
* @brief Base class for kernels for atom-continuum transfer
*/
class KernelFunction {
public:
// constructor
KernelFunction(int nparameters, double* parameters);
// destructor
virtual ~KernelFunction() {};
// function value
virtual double value(DENS_VEC& x_atom) const =0 ;
// function derivative
virtual void derivative(const DENS_VEC& x_atom, DENS_VEC& deriv) const =0 ;
// bond function value via quadrature
virtual double bond(DENS_VEC& xa, DENS_VEC&xb, double lam1, double lam2) const;
// localization-volume intercepts for bond calculation
virtual void bond_intercepts(DENS_VEC& xa,
DENS_VEC& xb, double &lam1, double &lam2) const;
virtual bool node_contributes(DENS_VEC node) const;
virtual bool in_support(DENS_VEC node) const;
double inv_vol(void) const { return invVol_; }
virtual double dimensionality_factor(void) const { return 1.; }
protected:
double invVol_; // normalization factor
double Rc_, invRc_; // cutoff radius
int nsd_ ; // number of dimensions
/** pointer to lammps interface class */
LammpsInterface * lammpsInterface_;
/** periodicity flags and lengths */
int periodicity[3];
double box_bounds[2][3];
double box_length[3];
};
/**
* @class KernelFunctionStep
* @brief Class for defining kernel function of a step with spherical support
*/
class KernelFunctionStep : public KernelFunction {
public:
// constructor
KernelFunctionStep(int nparameters, double* parameters);
// destructor
virtual ~KernelFunctionStep() {};
// function value
virtual double value(DENS_VEC& x_atom) const;
// function derivative
virtual void derivative(const DENS_VEC& x_atom, DENS_VEC& deriv) const;
// bond function value
virtual double bond(DENS_VEC& /* xa */, DENS_VEC& /* xb */, double lam1, double lam2) const
{ return lam2-lam1; }
};
/**
* @class KernelFunctionCell
* @brief Class for defining kernel function of a step with rectangular support
* suitable for a rectangular grid
*/
class KernelFunctionCell : public KernelFunction {
public:
// constructor
KernelFunctionCell(int nparameters, double* parameters);
// destructor
virtual ~KernelFunctionCell() {};
// function value
virtual double value(DENS_VEC& x_atom) const;
// function derivative
virtual void derivative(const DENS_VEC& x_atom, DENS_VEC& deriv) const;
// bond function value
virtual double bond(DENS_VEC& /* xa */, DENS_VEC& /* xb */, double lam1, double lam2) const
{return lam2 -lam1;}
// bond intercept values : origin is the node position
void bond_intercepts(DENS_VEC& xa, DENS_VEC& xb,
double &lam1, double &lam2) const;
bool node_contributes(DENS_VEC node) const;
bool in_support(DENS_VEC dx) const;
protected:
double hx, hy, hz;
DENS_VEC cellBounds_;
};
/**
* @class KernelFunctionCubicSphere
* @brief Class for defining kernel function of a cubic with spherical support
*/
class KernelFunctionCubicSphere : public KernelFunction {
public:
// constructor
KernelFunctionCubicSphere(int nparameters, double* parameters);
// destructor
virtual ~KernelFunctionCubicSphere() {};
// function value
virtual double value(DENS_VEC& x_atom) const;
// function derivative
virtual void derivative(const DENS_VEC& x_atom, DENS_VEC& deriv) const;
};
/**
* @class KernelFunctionQuarticSphere
* @brief Class for defining kernel function of a quartic with spherical support
*/
class KernelFunctionQuarticSphere : public KernelFunction {
public:
// constructor
KernelFunctionQuarticSphere(int nparameters, double* parameters);
// destructor
virtual ~KernelFunctionQuarticSphere() {};
// function value
virtual double value(DENS_VEC& x_atom) const;
// function derivative
virtual void derivative(const DENS_VEC& x_atom, DENS_VEC& deriv) const;
};
/**
* @class KernelFunctionCubicCyl
* @brief Class for defining kernel function of a cubic with cylindrical support
*/
class KernelFunctionCubicCyl : public KernelFunction {
public:
// constructor
KernelFunctionCubicCyl(int nparameters, double* parameters);
// destructor
virtual ~KernelFunctionCubicCyl() {};
// function value
virtual double value(DENS_VEC& x_atom) const;
// function derivative
virtual void derivative(const DENS_VEC& x_atom, DENS_VEC& deriv) const;
virtual double dimensionality_factor(void) const { return 0.5; }
};
/**
* @class KernelFunctionQuarticCyl
* @brief Class for defining kernel function of a quartic with cylindrical support
*/
class KernelFunctionQuarticCyl : public KernelFunction {
public:
// constructor
KernelFunctionQuarticCyl(int nparameters, double* parameters);
// destructor
virtual ~KernelFunctionQuarticCyl() {};
// function value
virtual double value(DENS_VEC& x_atom) const;
// function derivative
virtual void derivative(const DENS_VEC& x_atom, DENS_VEC& deriv) const;
virtual double dimensionality_factor(void) const { return 0.5; }
};
/**
* @class KernelFunctionCubicBar
* @brief Class for defining kernel function of a cubic with 1-dimensional (bar) support
*/
class KernelFunctionCubicBar : public KernelFunction {
public:
// constructor
KernelFunctionCubicBar(int nparameters, double* parameters);
// destructor
virtual ~KernelFunctionCubicBar() {};
// function value
virtual double value(DENS_VEC& x_atom) const;
// function derivative
virtual void derivative(const DENS_VEC& x_atom, DENS_VEC& deriv) const;
virtual double dimensionality_factor(void) const { return 0.25; }
};
/**
* @class KernelFunctionCubicBar
* @brief Class for defining kernel function of a cubic with 1-dimensional (bar) support
*/
class KernelFunctionLinearBar : public KernelFunction {
public:
// constructor
KernelFunctionLinearBar(int nparameters, double* parameters);
// destructor
virtual ~KernelFunctionLinearBar() {};
// function value
virtual double value(DENS_VEC& x_atom) const;
// function derivative
virtual void derivative(const DENS_VEC& x_atom, DENS_VEC& deriv) const;
virtual double dimensionality_factor(void) const { return 0.25; }
};
/**
* @class KernelFunctionQuarticBar
* @brief Class for defining kernel function of a quartic with 1-dimensional (bar) support
*/
class KernelFunctionQuarticBar : public KernelFunction {
public:
// constructor
KernelFunctionQuarticBar(int nparameters, double* parameters);
// destructor
virtual ~KernelFunctionQuarticBar() {};
// function value
virtual double value(DENS_VEC& x_atom) const;
// function derivative
virtual void derivative(const DENS_VEC& x_atom, DENS_VEC& deriv) const;
virtual double dimensionality_factor(void) const { return 0.25; }
};
};
#endif
|