JDFTx  1.5.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Data structures

Files

file  ManagedMemory.h
 Base class and operators for managed-memory objects.
 
file  matrix.h
 
file  matrix3.h
 3x3 matrices with CPU and GPU operators
 
file  RadialFunction.h
 
file  scalar.h
 Complex numbers with CPU and GPU operators.
 
file  ScalarField.h
 Real and complex scalar fields in real and reciprocal space.
 
file  ScalarFieldArray.h
 Variable length arrays of ScalarField and ScalarFieldTildeArray, and their operators.
 
file  tensor3.h
 Symmetric traceless tensor with CPU and GPU operators.
 
file  vector3.h
 3-vector with CPU and GPU operators
 
file  VectorField.h
 Fixed-length multiplet of ScalarField and ScalarFieldTilde, and specialization to vector fields (3-component)
 
file  ColumnBundle.h
 

Classes

class  ManagedMemoryBase
 Base class for managed-memory objects (that could potentially live on GPUs as well) with unspecified data type. More...
 
class  ManagedMemory< T >
 Base class for managed memory of a specified data type. More...
 
struct  ManagedArray< T >
 
struct  IndexArray
 Managed array of integers (indices) More...
 
struct  IndexVecArray
 Managed array of integer vectors. More...
 
class  diagMatrix
 Real diagonal matrix. More...
 
class  matrix
 General complex matrix. More...
 
struct  matrixScaledTransOp
 Matrix with pending scale, submatrix and/or transpose operations. More...
 
class  tiledBlockMatrix
 A block matrix formed by repeating (tiling) a dense matrix along the diagonal. More...
 
class  matrix3< scalar >
 3x3 matrix More...
 
struct  SpaceGroupOp
 Space group operation r -> rot * r + a in real-space lattice coordinates. More...
 
struct  RadialFunctionG
 G-space radial function stored on a uniform grid (of |G|) More...
 
struct  RadialFunctionR
 A function on a non-uniform real-space radial grid. More...
 
struct  array< T, N >
 
struct  complex
 Complex number (need to define our own because we need operators for gpu code as well) More...
 
struct  FieldData< T >
 ManagedMemory wrapper with gridInfo and pending scale factor for ScalarField* classes. More...
 
struct  ScalarFieldData
 Real space real scalar field data Do not use this data structure directly or from a simple pointer ScalarFieldData*; work only with ScalarField's. The public functions of ScalarFieldData can be accessed with -> from the ScalarField. More...
 
struct  ScalarFieldTildeData
 Reciprocal space real scalar field data Do not use this data structure directly or from a simple pointer ScalarFieldTildeData*; work only with ScalarFieldTilde's. The public functions of ScalarFieldTildeData can be accessed with -> from the ScalarFieldTilde. More...
 
struct  complexScalarFieldData
 Real space complex scalar field data Do not use this data structure directly or from a simple pointer complexScalarFieldData*; work only with complexScalarField's. The public functions of complexScalarFieldData can be accessed with -> from the complexScalarField. More...
 
struct  complexScalarFieldTildeData
 Reciprocal space complex scalar field data Do not use this data structure directly or from a simple pointer complexScalarFieldTildeData*; work only with complexScalarFieldTilde's. The public functions of complexScalarFieldTildeData can be accessed with -> from the complexScalarFieldTilde. More...
 
struct  RealKernel
 Special class for storing real reciprocal-space kernels encountered ever so often for convolutions. More...
 
struct  scaled< T >
 Template to avoid (delay) scaling operations on linear objects. More...
 
class  tensor3< scalar >
 Symmetric traceless rank-2 tensor in 3D. More...
 
class  vector3< scalar >
 Generic 3-vector. More...
 
struct  ScalarFieldMultiplet< T, N >
 Generic multiplet object with overloaded arithmetic. More...
 
class  ColumnBundle
 Wavefunction data structure. More...
 
struct  ColumnBundleMatrixProduct
 ColumnBundle with a pending matrix multiply (on the right side) More...
 

Macros

#define MUL_MAT_VEC(retType)
 
#define MUL_VEC_MAT(retType)
 
#define MUL_MAT_MAT(retType)
 
#define __hostanddev__   inline
 
#define TptrCollection   std::vector<std::shared_ptr<T> >
 shorthand for templates below (undef'd at end of file)
 
#define LOOP5(code)   { for(int k=0; k<5; k++) { code } }
 
#define LOOP3(code)   { for(int k=0; k<3; k++) { code } }
 
#define Tptr   std::shared_ptr<T>
 shorthand for writing the template operators (undef'd at end of header)
 
#define TptrMul   ScalarFieldMultiplet<T,N>
 shorthand for the template operators/functions (undef'd at end of file)
 
#define RptrMul   ScalarFieldMultiplet<ScalarFieldData,N>
 shorthand for real-space-only template operators/functions (undef'd at end of file)
 
#define GptrMul   ScalarFieldMultiplet<ScalarFieldTildeData,N>
 shorthand for reciprocal-space-only template operators/functions (undef'd at end of file)
 

Typedefs

typedef std::shared_ptr
< ScalarFieldData
ScalarField
 A smart reference-counting pointer to ScalarFieldData.
 
typedef std::shared_ptr
< ScalarFieldTildeData
ScalarFieldTilde
 A smart reference-counting pointer to ScalarFieldTildeData.
 
typedef std::shared_ptr
< complexScalarFieldData
complexScalarField
 A smart reference-counting pointer to complexScalarFieldData.
 
typedef std::shared_ptr
< complexScalarFieldTildeData
complexScalarFieldTilde
 A smart reference-counting pointer to complexScalarFieldTildeData.
 
typedef std::vector< ScalarFieldScalarFieldArray
 dynamic size collection of real space scalar fields
 
typedef std::vector
< ScalarFieldTilde
ScalarFieldTildeArray
 dynamic size collection of reciprocal space scalar fields
 
typedef ScalarFieldMultiplet
< ScalarFieldData, 3 > 
VectorField
 Real space vector field.
 
typedef ScalarFieldMultiplet
< ScalarFieldTildeData, 3 > 
VectorFieldTilde
 Reciprocal space vector field.
 
typedef ScalarFieldMultiplet
< ScalarFieldData, 5 > 
TensorField
 Symmetric traceless tensor: real space field.
 
typedef ScalarFieldMultiplet
< ScalarFieldTildeData, 5 > 
TensorFieldTilde
 Symmetric traceless tensor: reciprocal space field.
 

Functions

template<typename T >
void memcpy (ManagedMemory< T > &, const ManagedMemory< T > &)
 copy entire object over
 
void scale (double alpha, ManagedMemory< double > &y)
 
void scale (double alpha, ManagedMemory< complex > &y)
 scale y *= alpha
 
void scale (complex alpha, ManagedMemory< complex > &y)
 scale y *= alpha
 
void scale (const ManagedMemory< complex > &x, ManagedMemory< complex > &y)
 scale y *= alpha
 
void axpy (complex alpha, const ManagedMemory< complex > &x, ManagedMemory< complex > &y)
 scale y *= x elementwise More...
 
double nrm2 (const ManagedMemory< complex > &)
 2-norm, pretending it is a vector
 
complex dotc (const ManagedMemory< complex > &a, const ManagedMemory< complex > &b)
 return a^H b
 
matrix conj (const scaled< matrix > &A)
 return element-wise complex conjugate of A
 
matrixScaledTransOp dagger (const matrixScaledTransOp &A)
 return hermitian adjoint of A
 
matrixScaledTransOp transpose (const matrixScaledTransOp &A)
 return transpose of A
 
matrix dagger_symmetrize (const matrixScaledTransOp &A)
 
matrix transpose_symmetrize (const matrixScaledTransOp &A)
 return adjoint symmetric part: (A + Adag)/2
 
double nrm2 (const matrixScaledTransOp &A)
 return transpose symmetric part: (A + AT)/2 More...
 
matrixoperator*= (matrix &m, double s)
 
scaled< matrixoperator* (double s, const matrix &m)
 
scaled< matrixoperator* (const matrix &m, double s)
 
scaled< matrixoperator- (const matrix &m)
 
matrixoperator*= (matrix &m, complex s)
 
matrix operator* (complex s, const matrix &m)
 
matrix operator* (const matrix &m, complex s)
 
diagMatrix operator* (double s, const diagMatrix &m)
 
diagMatrix operator* (const diagMatrix &m, double s)
 
diagMatrix operator- (const diagMatrix &m)
 
matrix operator* (const matrixScaledTransOp &, const matrixScaledTransOp &)
 
matrix operator* (const diagMatrix &, const matrix &)
 
matrix operator* (const matrix &, const diagMatrix &)
 
matrix operator* (const std::vector< complex > &d, const matrix &m)
 
matrix operator* (const matrix &m, const std::vector< complex > &d)
 
diagMatrix operator* (const diagMatrix &, const diagMatrix &)
 
matrixoperator+= (matrix &m1, const matrix &m2)
 
matrixoperator-= (matrix &m1, const matrix &m2)
 
matrix operator+ (const matrix &m1, const matrix &m2)
 
matrix operator- (const matrix &m1, const matrix &m2)
 
void axpy (double alpha, const diagMatrix &x, matrix &y)
 
matrixoperator+= (matrix &m, const diagMatrix &d)
 
matrixoperator-= (matrix &m, const diagMatrix &d)
 
matrix operator+ (const matrix &m, const diagMatrix &d)
 
matrix operator+ (const diagMatrix &d, const matrix &m)
 
matrix operator- (const matrix &m, const diagMatrix &d)
 
matrix operator- (const diagMatrix &d, const matrix &m)
 
void axpy (double alpha, const diagMatrix &x, diagMatrix &y)
 
diagMatrixoperator+= (diagMatrix &d1, const diagMatrix &d2)
 
diagMatrixoperator-= (diagMatrix &d1, const diagMatrix &d2)
 
diagMatrix operator+ (const diagMatrix &d1, const diagMatrix &d2)
 
diagMatrix operator- (const diagMatrix &d1, const diagMatrix &d2)
 
diagMatrix clone (const diagMatrix &x)
 
matrix clone (const matrix &x)
 
double dot (const diagMatrix &x, const diagMatrix &y)
 
double dot (const matrix &x, const matrix &y)
 
void randomize (diagMatrix &x)
 
void randomize (matrix &x)
 
matrix inv (const matrix &A)
 Compute inverse of an arbitrary matrix A (via LU decomposition) More...
 
diagMatrix inv (const diagMatrix &A)
 inverse of diagonal matrix
 
matrix invApply (const matrix &A, const matrix &b)
 return inv(A) * b (A must be hermitian, positive-definite)
 
matrix LU (const matrix &A)
 Compute the LU decomposition of the matrix.
 
complex det (const matrix &A)
 
double det (const diagMatrix &A)
 Compute the determinant of an diagonal matrix A.
 
matrix pow (const matrix &A, double exponent, matrix *Aevecs=0, diagMatrix *Aeigs=0, bool *isSingular=0)
 
matrix invsqrt (const matrix &A, matrix *Aevecs=0, diagMatrix *Aeigs=0, bool *isSingular=0)
 
matrix cis (const matrix &A, matrix *Aevecs=0, diagMatrix *Aeigs=0)
 Compute cis(A) = exp(iota A) and optionally the eigensystem of A (if non-null)
 
matrix cisInv (const matrix &A, matrix *Bevecs=0, diagMatrix *Beigs=0)
 
matrix sqrt_grad (const matrix &grad_sqrtA, const matrix &Aevecs, const diagMatrix &Aeigs)
 Return gradient w.r.t A given gradient w.r.t sqrt(A) and A's eigensystem.
 
matrix cis_grad (const matrix &grad_cisA, const matrix &Aevecs, const diagMatrix &Aeigs)
 Return gradient w.r.t A given gradient w.r.t cis(A) and A's eigensystem.
 
complex trace (const matrix &m)
 trace of matrix
 
double trace (const diagMatrix &m)
 trace of diagonal matrix
 
double nrm2 (const diagMatrix &m)
 RMS of matrix entries.
 
diagMatrix diag (const matrix &m)
 obtain the real diagonal part of a hermitian matrix
 
diagMatrix eye (int N)
 identity
 
matrix zeroes (int nRows, int nCols)
 a dense-matrix of zeroes
 
matrix operator* (const matrix &m, const tiledBlockMatrix &tbm)
 multiply dense matrix by block matrix
 
template<typename scalar >
__hostanddev__ matrix3< scalar > operator* (scalar s, const matrix3< scalar > &m)
 
template<typename scalar >
__hostanddev__ matrix3< scalar > outer (const vector3< scalar > &a, const vector3< scalar > &b)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator* (const matrix3< scalar > &m, const vector3< scalar > &v)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator* (const matrix3< scalar > &m, const vector3< int > &v)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator* (const matrix3< int > &m, const vector3< scalar > &v)
 
__hostanddev__ vector3< int > operator* (const matrix3< int > &m, const vector3< int > &v)
 
__hostanddev__ vector3< complexoperator* (const matrix3< complex > &m, const vector3<> &v)
 
__hostanddev__ vector3< complexoperator* (const matrix3<> &m, const vector3< complex > &v)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator* (const vector3< scalar > &v, const matrix3< scalar > &m)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator* (const vector3< int > &v, const matrix3< scalar > &m)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator* (const vector3< scalar > &v, const matrix3< int > &m)
 
__hostanddev__ vector3< int > operator* (const vector3< int > &v, const matrix3< int > &m)
 
__hostanddev__ vector3< complexoperator* (const vector3<> &v, const matrix3< complex > &m)
 
__hostanddev__ vector3< complexoperator* (const vector3< complex > &v, const matrix3<> &m)
 
template<typename scalar >
__hostanddev__ matrix3< scalar > operator* (const matrix3< scalar > &m, const matrix3< scalar > &n)
 
template<typename scalar >
__hostanddev__ matrix3< scalar > operator* (const matrix3< scalar > &m, const matrix3< int > &n)
 
template<typename scalar >
__hostanddev__ matrix3< scalar > operator* (const matrix3< int > &m, const matrix3< scalar > &n)
 
__hostanddev__ matrix3< int > operator* (const matrix3< int > &m, const matrix3< int > &n)
 
template<typename scalar >
__hostanddev__ matrix3< scalar > & operator*= (matrix3< scalar > &m, const matrix3< scalar > &n)
 
template<typename scalar >
__hostanddev__ matrix3< scalar > Diag (vector3< scalar > v)
 
template<typename scalar >
__hostanddev__ scalar trace (const matrix3< scalar > &m)
 Trace of matrix.
 
template<typename scalar >
__hostanddev__ scalar det (const matrix3< scalar > &m)
 Determinant of matrix.
 
template<typename scalar >
__hostanddev__ matrix3< scalar > adjugate (const matrix3< scalar > &m)
 
__hostanddev__ matrix3 inv (const matrix3<> &m)
 
template<typename scalar >
__hostanddev__ scalar nrm2sq (const matrix3< scalar > &m)
 square of 2-norm of matrix
 
template<typename scalar >
__hostanddev__ double nrm2 (const matrix3< scalar > &m)
 2-norm of matrix
 
__hostanddev__ matrix3 rotation (double theta, int axis)
 Create a rotation matrix.
 
template<class T >
ceildiv (T num, T den)
 Ceiling of a positive integer division, templated over int types.
 
template<class T >
floorMultiple (T num, T den)
 Return largest multiple of den smaller than num, templated over int types.
 
__hostanddev__ double real (const complex &c)
 real part
 
__hostanddev__ double imag (const complex &c)
 imaginary part
 
__hostanddev__ double norm (const complex &c)
 absolute value squared
 
__hostanddev__ double abs (const complex &c)
 absolute value
 
__hostanddev__ double arg (const complex &c)
 argument (phase angle)
 
__hostanddev__ complex conj (const complex &c)
 complex conjugate
 
__hostanddev__ double conj (const double &c)
 identity operation provided to ease templating over complex and double
 
__hostanddev__ complex operator+ (double r, const complex &c)
 
__hostanddev__ complex operator- (double r, const complex &c)
 
__hostanddev__ complex operator* (double r, const complex &c)
 
__hostanddev__ complex cis (double x)
 Compute cis(x) = exp(iota x) = cos(x) + iota sin(x)
 
template<typename T >
std::vector< typename
T::DataType * > 
dataPref (TptrCollection &x)
 Extract a std::vector of data pointers from a ScalarFieldArray.
 
template<typename T >
std::vector< const typename
T::DataType * > 
constDataPref (const TptrCollection &x)
 Extract a std::vector of const data pointers from a const ScalarFieldArray.
 
template<typename T >
TptrCollection clone (const TptrCollection &x)
 Create a copy of the data (note operator= references same data since Tptr's are pointers!)
 
template<typename T >
TptrCollectionoperator*= (TptrCollection &x, double alpha)
 Scale.
 
template<class T >
TptrCollection operator* (const TptrCollection &in, double scaleFac)
 
template<class T >
TptrCollection operator* (double scaleFac, const TptrCollection &in)
 
template<class T >
TptrCollection operator* (TptrCollection &&in, double scaleFac)
 Add (destructible input)
 
template<class T >
TptrCollection operator* (double scaleFac, TptrCollection &&in)
 Add (destructible input)
 
template<typename T >
void axpy (double alpha, const TptrCollection &x, TptrCollection &y)
 y += alpha x
 
ScalarFieldArray operator* (const ScalarFieldArray &x, const ScalarFieldArray &y)
 elementise multiply
 
ScalarFieldArray operator* (ScalarFieldArray &&x, const ScalarFieldArray &y)
 
ScalarFieldArray operator* (const ScalarFieldArray &x, ScalarFieldArray &&y)
 
ScalarFieldArrayoperator*= (ScalarFieldArray &y, const ScalarField &x)
 
ScalarFieldArray operator* (const ScalarField &x, ScalarFieldArray &&y)
 
ScalarFieldArray operator* (ScalarFieldArray &&y, const ScalarField &x)
 
ScalarFieldArray operator* (const ScalarField &x, const ScalarFieldArray &y)
 
ScalarFieldArray operator* (const ScalarFieldArray &y, const ScalarField &x)
 
template<class T >
TptrCollectionoperator+= (TptrCollection &in, const TptrCollection &other)
 Increment.
 
template<class T >
TptrCollectionoperator-= (TptrCollection &in, const TptrCollection &other)
 Decrement.
 
template<class T >
TptrCollection operator+ (const TptrCollection &in1, const TptrCollection &in2)
 
template<class T >
TptrCollection operator+ (const TptrCollection &in1, TptrCollection &&in2)
 Add (destructible input)
 
template<class T >
TptrCollection operator+ (TptrCollection &&in1, const TptrCollection &in2)
 Add (destructible input)
 
template<class T >
TptrCollection operator+ (TptrCollection &&in1, TptrCollection &&in2)
 Add (destructible inputs)
 
template<class T >
TptrCollection operator- (const TptrCollection &in1, const TptrCollection &in2)
 
template<class T >
TptrCollection operator- (const TptrCollection &in1, TptrCollection &&in2)
 Add (destructible inputs)
 
template<class T >
TptrCollection operator- (TptrCollection &&in1, const TptrCollection &in2)
 Add (destructible inputs)
 
template<class T >
TptrCollection operator- (TptrCollection &&in1, TptrCollection &&in2)
 Add (destructible inputs)
 
template<typename T >
double dot (const TptrCollection &x, const TptrCollection &y)
 Inner product.
 
ScalarFieldTildeArray lGradient (const ScalarFieldTilde &, int l)
 spherical tensor gradient of order l (2l+1 outputs, multiplied by Ylm(Ghat) (iG)^l)
 
ScalarFieldTilde lDivergence (const ScalarFieldTildeArray &, int l)
 spherical tensor divergence of order l (2l+1 inputs, multiplied by Ylm(Ghat) (iG)^l, and summed)
 
template<typename T >
void initZero (TptrCollection &x)
 Initialize (non-null) data to zero.
 
template<typename T >
void nullToZero (TptrCollection &x, const GridInfo &gInfo, int N=0)
 
template<typename T >
void initRandomFlat (TptrCollection &x)
 Initialize to random numbers (uniform on 0 to 1)
 
template<typename T >
void randomize (TptrCollection &x)
 Initialize to normal random numbers:
 
template<typename T >
void loadFromFile (TptrCollection &x, const char *filename)
 Load an array of scalar fields from a file.
 
template<typename T >
void saveToFile (const TptrCollection &x, const char *filename)
 Save an array of scalar fields to a file.
 
ScalarFieldArray I (ScalarFieldTildeArray &&X)
 
ScalarFieldArray I (const ScalarFieldTildeArray &X)
 Reciprocal to real space.
 
ScalarFieldTildeArray J (const ScalarFieldArray &X)
 
ScalarFieldTildeArray Idag (const ScalarFieldArray &X)
 
ScalarFieldArray Jdag (ScalarFieldTildeArray &&X)
 
ScalarFieldArray Jdag (const ScalarFieldTildeArray &X)
 Hermitian conjugate of J.
 
template<typename T >
T & operator+= (T &y, const scaled< T > &x)
 
template<typename T >
T & operator-= (T &y, const scaled< T > &x)
 
template<typename T >
operator+ (const scaled< T > &x, const scaled< T > &y)
 
template<typename T >
operator- (const scaled< T > &x, const scaled< T > &y)
 
template<typename T >
scaled< T > operator- (const scaled< T > &x)
 
template<typename T >
scaled< T > operator* (double s, const scaled< T > &x)
 
template<typename T >
scaled< T > operator* (const scaled< T > &x, double s)
 
template<typename scalar >
__hostanddev__ tensor3< scalar > loadTensor (const tensor3< const scalar * > &tArr, int i)
 Load tensor from a constant tensor field.
 
template<typename scalar >
__hostanddev__ tensor3< scalar > loadTensor (const tensor3< scalar * > &tArr, int i)
 Load tensor from a tensor field.
 
template<typename scalar >
__hostanddev__ void storeTensor (const tensor3< scalar > &t, tensor3< scalar * > &tArr, int i)
 Store tensor to a tensor field.
 
template<typename scalar >
__hostanddev__ void accumTensor (const tensor3< scalar > &t, tensor3< scalar * > &tArr, int i)
 Accumulate tensor onto a tensor field.
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator+ (scalar s, const vector3< scalar > &a)
 
__hostanddev__ vector3< int > operator+ (int s, const vector3< int > &a)
 
__hostanddev__ vector3< int > operator+ (const vector3< int > &a, int s)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator+ (const vector3< scalar > &a, int s)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator+ (int s, const vector3< scalar > &a)
 
__hostanddev__ vector3 operator+ (vector3<> a, vector3< int > b)
 
__hostanddev__ vector3 operator+ (vector3< int > a, vector3<> b)
 
template<typename scalar >
__hostanddev__ vector3< scalar > & operator*= (vector3< scalar > &a, scalar s)
 
template<typename scalar >
__hostanddev__ vector3< scalar > & operator*= (vector3< scalar > &a, double s)
 
template<typename scalar >
__hostanddev__ vector3< scalar > & operator*= (vector3< scalar > &a, int s)
 
__hostanddev__ vector3operator*= (vector3<> &a, double s)
 
__hostanddev__ vector3operator*= (vector3<> &a, int s)
 
__hostanddev__ vector3< int > & operator*= (vector3< int > &a, int s)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator* (scalar s, const vector3< scalar > &a)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator* (double s, const vector3< scalar > &a)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator* (const vector3< scalar > &a, double s)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator* (int s, const vector3< scalar > &a)
 
template<typename scalar >
__hostanddev__ vector3< scalar > operator* (const vector3< scalar > &a, int s)
 
__hostanddev__ vector3 operator* (double s, const vector3<> &a)
 
__hostanddev__ vector3 operator* (const vector3<> &a, double s)
 
__hostanddev__ vector3 operator* (double s, const vector3< int > &a)
 
__hostanddev__ vector3 operator* (const vector3< int > &a, double s)
 
__hostanddev__ vector3< complexoperator* (complex s, const vector3<> &a)
 
__hostanddev__ vector3< complexoperator* (const vector3<> &a, complex s)
 
__hostanddev__ vector3< complexoperator* (complex s, const vector3< int > &a)
 
__hostanddev__ vector3< complexoperator* (const vector3< int > &a, complex s)
 
template<typename scalar >
__hostanddev__ scalar dot (const vector3< scalar > &a, const vector3< scalar > &b)
 
template<typename scalar >
__hostanddev__ scalar dot (const vector3< double > &a, const vector3< scalar > &b)
 
template<typename scalar >
__hostanddev__ scalar dot (const vector3< int > &a, const vector3< scalar > &b)
 
template<typename scalar >
__hostanddev__ scalar dot (const vector3< scalar > &a, const vector3< double > &b)
 
template<typename scalar >
__hostanddev__ scalar dot (const vector3< scalar > &a, const vector3< int > &b)
 
__hostanddev__ double dot (const vector3< double > &a, const vector3< double > &b)
 
__hostanddev__ double dot (const vector3< int > &a, const vector3< double > &b)
 
__hostanddev__ double dot (const vector3< double > &a, const vector3< int > &b)
 
__hostanddev__ int dot (const vector3< int > &a, const vector3< int > &b)
 
template<typename scalar >
__hostanddev__ vector3< scalar > cross (const vector3< scalar > &a, const vector3< scalar > &b)
 cross product
 
template<typename scalar >
__hostanddev__ scalar box (const vector3< scalar > &a, const vector3< scalar > &b, const vector3< scalar > &c)
 box product / triple product
 
__hostanddev__ vector3< double > normalize (const vector3<> &v)
 normalise vector
 
__hostanddev__ vector3< int > round (const vector3<> &v, double *err=0)
 Round vector3<> to vector3<int> (and optionally retrieve error)
 
__hostanddev__ double circDistanceSquared (const vector3<> &a, const vector3<> &b)
 
template<typename T >
gcd (T x, T y)
 GCD of integers, templated over integer types.
 
template<typename T >
vector3< T > gcdReduce (const vector3< T > &d)
 Reduce an integer vector by its gcd.
 
template<typename scalar >
__hostanddev__ vector3< scalar > loadVector (const vector3< const scalar * > &vArr, int i)
 Load vector from a constant vector field.
 
template<typename scalar >
__hostanddev__ vector3< scalar > loadVector (const vector3< scalar * > &vArr, int i)
 Load vector from a vector field.
 
template<typename scalar >
__hostanddev__ void storeVector (const vector3< scalar > &v, vector3< scalar * > &vArr, int i)
 Store vector to a vector field.
 
template<typename scalar >
__hostanddev__ void accumVector (const vector3< scalar > &v, vector3< scalar * > &vArr, int i)
 Accumulate vector onto a vector field.
 
template<class T , int N>
TptrMul clone (const TptrMul &X)
 Clone (NOTE: operator= is by reference for ScalarField multiplets)
 
template<class T , int N>
void initZero (TptrMul &X)
 Initialize data to 0 and scale factors to 1.
 
template<class T , int N>
void nullToZero (TptrMul &X, const GridInfo &gInfo)
 Allocate and initialize each component of X to 0 if null.
 
template<int N>
void initRandom (RptrMul &X, double cap=0.0)
 initialize element-wise with a unit-normal random number (with a cap if cap>0)
 
template<int N>
void initRandomFlat (RptrMul &X)
 initialize element-wise with a unit-flat [0:1) random number
 
template<int N>
void randomize (RptrMul &X)
 alternate interface required by Minimizable
 
template<class T , int N>
TptrMuloperator*= (TptrMul &in, const TptrMul &other)
 Elementwise multiply each component.
 
template<class T , int N>
TptrMul operator* (const TptrMul &in1, const TptrMul &in2)
 Elementwise multiply each component (preserve inputs)
 
template<class T , int N>
TptrMul operator* (TptrMul &&in1, const TptrMul &in2)
 Elementwise multiply each component (destructible input)
 
template<class T , int N>
TptrMul operator* (const TptrMul &in1, TptrMul &&in2)
 Elementwise multiply each component (destructible input)
 
template<class T , int N>
TptrMul operator* (TptrMul &&in1, TptrMul &&in2)
 Elementwise multiply each component (destructible inputs)
 
template<class T , int N>
TptrMuloperator*= (TptrMul &inM, const Tptr &inS)
 Elementwise multiply each component.
 
template<class T , int N>
TptrMul operator* (const TptrMul &inM, const Tptr &inS)
 Elementwise multiply each component (preserve inputs)
 
template<class T , int N>
TptrMul operator* (const Tptr &inS, const TptrMul &inM)
 Elementwise multiply each component (preserve inputs)
 
template<class T , int N>
TptrMul operator* (TptrMul &&inM, const Tptr &inS)
 Elementwise multiply each component (destructible input)
 
template<class T , int N>
TptrMul operator* (const Tptr &inS, TptrMul &&inM)
 Elementwise multiply each component (destructible input)
 
template<class T , int N>
TptrMuloperator*= (TptrMul &in, double scaleFac)
 Scale.
 
template<class T , int N>
TptrMul operator* (const TptrMul &in, double scaleFac)
 Scalar multiply (preserve input)
 
template<class T , int N>
TptrMul operator* (double scaleFac, const TptrMul &in)
 Scalar multiply (preserve input)
 
template<class T , int N>
TptrMul operator* (TptrMul &&in, double scaleFac)
 Scalar multiply (destructible input)
 
template<class T , int N>
TptrMul operator* (double scaleFac, TptrMul &&in)
 Scalar multiply (destructible input)
 
template<class T >
ScalarFieldMultiplet< T, 3 > operator* (vector3<> v, const Tptr &in)
 3-vector multiply
 
template<class T >
Tptr dot (vector3<> v, const ScalarFieldMultiplet< T, 3 > &in)
 3-vector multiply
 
template<class T , int N>
void axpy (double alpha, const TptrMul &X, TptrMul &Y)
 Linear combine Y += alpha * X.
 
template<class T , int N>
TptrMuloperator+= (TptrMul &in, const TptrMul &other)
 Increment.
 
template<class T , int N>
TptrMuloperator-= (TptrMul &in, const TptrMul &other)
 Decrement.
 
template<class T , int N>
TptrMul operator+ (const TptrMul &in1, const TptrMul &in2)
 Add (preserve inputs)
 
template<class T , int N>
TptrMul operator+ (const TptrMul &in1, TptrMul &&in2)
 Add (destructible input)
 
template<class T , int N>
TptrMul operator+ (TptrMul &&in1, const TptrMul &in2)
 Add (destructible input)
 
template<class T , int N>
TptrMul operator+ (TptrMul &&in1, TptrMul &&in2)
 Add (destructible inputs)
 
template<class T , int N>
TptrMul operator- (const TptrMul &in1, const TptrMul &in2)
 Subtract (preserve input)
 
template<class T , int N>
TptrMul operator- (const TptrMul &in1, TptrMul &&in2)
 Subtract (destructible input)
 
template<class T , int N>
TptrMul operator- (TptrMul &&in1, const TptrMul &in2)
 Subtract (destructible input)
 
template<class T , int N>
TptrMul operator- (TptrMul &&in1, TptrMul &&in2)
 Subtract (destructible inputs)
 
template<class T , int N>
TptrMul operator- (const TptrMul &in)
 Negate.
 
template<class T , int N>
TptrMul operator- (TptrMul &&in)
 Negate.
 
template<class T , int N>
void axpy (double alpha, const Tptr &X, TptrMul &Y)
 Linear combine Y += alpha * X.
 
template<class T , int N>
TptrMuloperator+= (TptrMul &in, const Tptr &other)
 Increment.
 
template<class T , int N>
TptrMuloperator-= (TptrMul &in, const Tptr &other)
 Decrement.
 
template<class T , int N>
TptrMul operator+ (const TptrMul &in1, const Tptr &in2)
 Add (preserve inputs)
 
template<class T , int N>
TptrMul operator+ (const Tptr &in1, const TptrMul &in2)
 Add (preserve inputs)
 
template<class T , int N>
TptrMul operator+ (const Tptr &in1, TptrMul &&in2)
 Add (destructible input)
 
template<class T , int N>
TptrMul operator+ (TptrMul &&in1, const Tptr &in2)
 Add (destructible input)
 
template<class T , int N>
TptrMul operator- (const TptrMul &in1, const Tptr &in2)
 Subtract (preserve input)
 
template<class T , int N>
TptrMul operator- (const Tptr &in1, const TptrMul &in2)
 Subtract (preserve input)
 
template<class T , int N>
TptrMul operator- (TptrMul &&in1, const Tptr &in2)
 Subtract (destructible input)
 
template<class T , int N>
TptrMul operator- (const Tptr &in1, TptrMul &&in2)
 Subtract (destructible input)
 
template<class T , int N>
double dot (const TptrMul &X, const TptrMul &Y)
 Inner product.
 
template<class T , int N>
double nrm2 (const TptrMul &X)
 2-norm
 
template<class T , int N>
double sum (const TptrMul &X)
 Sum of elements.
 
vector3 getGzero (const VectorFieldTilde &X)
 return G=0 components
 
void setGzero (const VectorFieldTilde &X, vector3<> v)
 set G=0 components
 
vector3 sumComponents (const VectorField &X)
 Sum of elements (component-wise)
 
ScalarField lengthSquared (const VectorField &X)
 Elementwise length squared.
 
ScalarField lengthSquaredWeighted (const vector3<> &w, const VectorField &X)
 Elementwise length squared (weighted)
 
ScalarField dotElemwise (const VectorField &X, const VectorField &Y)
 Elementwise dot.
 
template<int N>
RptrMuloperator+= (RptrMul &in, double scalar)
 Increment by scalar.
 
template<int N>
RptrMul operator+ (double scalar, const RptrMul &in)
 Add scalar (preserve input)
 
template<int N>
RptrMul operator+ (const RptrMul &in, double scalar)
 Add scalar (preserve input)
 
template<int N>
RptrMul operator+ (double scalar, RptrMul &&in)
 Add scalar (destructible input)
 
template<int N>
RptrMul operator+ (RptrMul &&in, double scalar)
 Add scalar (destructible input)
 
template<int N>
GptrMuloperator*= (GptrMul &in, const RealKernel &kernel)
 Multiply by kernel.
 
template<int N>
GptrMul operator* (const RealKernel &kernel, const GptrMul &in)
 Multiply by kernel (preserve input)
 
template<int N>
GptrMul operator* (const GptrMul &in, const RealKernel &kernel)
 Multiply by kernel (preserve input)
 
template<int N>
GptrMul operator* (const RealKernel &kernel, GptrMul &&in)
 Multiply by kernel (destructible input)
 
template<int N>
GptrMul operator* (GptrMul &&in, const RealKernel &kernel)
 Multiply by kernel (destructible input)
 
template<int N>
GptrMul O (GptrMul &&X)
 Inner product operator (diagonal in PW basis)
 
template<int N>
GptrMul O (const GptrMul &X)
 Inner product operator (diagonal in PW basis)
 
template<int N>
RptrMul I (GptrMul &&X)
 Forward transform: PW basis -> real space (destructible input)
 
template<int N>
GptrMul J (const RptrMul &X)
 Inverse transform: Real space -> PW basis.
 
template<int N>
GptrMul Idag (const RptrMul &X)
 Forward transform transpose: Real space -> PW basis.
 
template<int N>
RptrMul Jdag (GptrMul &&X)
 Inverse transform transpose: PW basis -> real space (destructible input)
 
template<int N>
RptrMul Jdag (const GptrMul &X)
 Inverse transform transpose: PW basis -> real space (preserve input)
 
template<int N>
RptrMul I (const GptrMul &X)
 Forward transform: PW basis -> real space (preserve input)
 
VectorFieldTilde gradient (const ScalarFieldTilde &)
 compute the gradient of a complex field, returns cartesian components
 
VectorField gradient (const ScalarField &)
 compute the gradient of a complex field, returns cartesian components
 
ScalarFieldTilde divergence (const VectorFieldTilde &)
 compute the divergence of a vector field specified in cartesian components
 
ScalarField divergence (const VectorField &)
 compute the divergence of a vector field specified in cartesian components
 
VectorFieldTilde operator* (const RadialFunctionG &, const VectorFieldTilde &)
 Convolve a vector field by a radial function (preserve input)
 
VectorFieldTilde operator* (const RadialFunctionG &, VectorFieldTilde &&)
 Convolve a vector field by a radial function (destructible input)
 
TensorFieldTilde tensorGradient (const ScalarFieldTilde &)
 symmetric traceless tensor second derivative of a scalar field
 
ScalarFieldTilde tensorDivergence (const TensorFieldTilde &)
 second derivative contraction of a symmetric traceless tensor field
 
template<int N>
void printStats (const RptrMul &, const char *name, FILE *fpLog=stdout)
 Print mean and standard deviation of each component array with specified name (debug utility)
 
void init (std::vector< ColumnBundle > &, int nbundles, int ncols=0, const Basis *basis=0, const ElecInfo *eInfo=0)
 Initialize an array of column bundles (with appropriate wavefunction sizes if ncols, basis, qnum and eInfo are all non-zero)
 
void randomize (std::vector< ColumnBundle > &, const ElecInfo &eInfo)
 randomize an array of columnbundles
 
ColumnBundle clone (const ColumnBundle &)
 
void randomize (ColumnBundle &x)
 Copies the input. More...
 
double dot (const ColumnBundle &x, const ColumnBundle &y)
 inner product
 
ColumnBundleoperator+= (ColumnBundle &Y, const scaled< ColumnBundle > &X)
 
ColumnBundleoperator-= (ColumnBundle &Y, const scaled< ColumnBundle > &X)
 
ColumnBundle operator+ (const scaled< ColumnBundle > &Y1, const scaled< ColumnBundle > &Y2)
 
ColumnBundle operator- (const scaled< ColumnBundle > &Y1, const scaled< ColumnBundle > &Y2)
 
ColumnBundleoperator*= (ColumnBundle &X, double s)
 
ColumnBundle operator* (double s, ColumnBundle &&Y)
 
ColumnBundle operator* (ColumnBundle &&Y, double s)
 
scaled< ColumnBundleoperator* (double s, const ColumnBundle &Y)
 
scaled< ColumnBundleoperator* (const ColumnBundle &Y, double s)
 
scaled< ColumnBundleoperator- (const ColumnBundle &Y)
 
ColumnBundleoperator*= (ColumnBundle &X, complex s)
 
ColumnBundle operator* (complex s, const ColumnBundle &Y)
 
ColumnBundle operator* (const ColumnBundle &Y, complex s)
 
ColumnBundleMatrixProduct operator* (const scaled< ColumnBundle > &sY, const matrixScaledTransOp &Mst)
 
ColumnBundleoperator+= (ColumnBundle &Y, const ColumnBundleMatrixProduct &XM)
 
ColumnBundleoperator-= (ColumnBundle &Y, const ColumnBundleMatrixProduct &XM)
 
ColumnBundle operator+ (const ColumnBundleMatrixProduct &XM1, const ColumnBundleMatrixProduct &XM2)
 
ColumnBundle operator- (const ColumnBundleMatrixProduct &XM1, const ColumnBundleMatrixProduct &XM2)
 
ColumnBundle operator* (const scaled< ColumnBundle > &, const diagMatrix &)
 
matrix operator^ (const scaled< ColumnBundle > &, const scaled< ColumnBundle > &)
 inner product
 
vector3< matrixspinOverlap (const scaled< ColumnBundle > &sY)
 spin-resolved inner product for a spinorial ColumnBundle
 
ColumnBundle Idag_DiagV_I (const ColumnBundle &C, const ScalarFieldArray &V)
 
ColumnBundle L (const ColumnBundle &Y)
 Apply Laplacian.
 
ColumnBundle Linv (const ColumnBundle &Y)
 Apply Laplacian inverse.
 
ColumnBundle O (const ColumnBundle &Y, std::vector< matrix > *VdagY=0)
 Apply overlap (and optionally retrieve pseudopotential projections for later reuse)
 
ColumnBundle D (const ColumnBundle &Y, int iDir)
 Compute the cartesian gradient of a column bundle in direction# iDir.
 
ColumnBundle DD (const ColumnBundle &Y, int iDir, int jDir)
 Compute second spatial derivative of a column bundle along directions# iDir, jDir.
 
void precond_inv_kinetic (ColumnBundle &Y, double KErollover)
 Apply inverse kinetic preconditioner (Roughly inv((k+G)^2/2)) in-place.
 
diagMatrix diagDot (const ColumnBundle &X, const ColumnBundle &Y)
 compute diag(X^Y) efficiently (avoid the off-diagonals)
 
void precond_inv_kinetic_band (ColumnBundle &Y, const diagMatrix &KEref)
 In-place inverse kinetic preconditioner with band-by-band KE reference (Used by BandDavidson)
 
ColumnBundle translate (ColumnBundle &&, vector3<> dr)
 translate a column-bundle by dr in lattice coordinates (destructible input)
 
ColumnBundle translate (const ColumnBundle &, vector3<> dr)
 translate a column-bundle by dr in lattice coordinates (preserve input)
 
void translateColumns (ColumnBundle &, const vector3<> *dr)
 translate each column of a column bundle by a different dr (in-place)
 
ColumnBundle switchBasis (const ColumnBundle &, const Basis &)
 return wavefunction projected to a different basis
 
complex traceinner (const diagMatrix &F, const ColumnBundle &X, const ColumnBundle &Y)
 Return trace(F*X^Y)
 
ScalarFieldArray diagouterI (const diagMatrix &F, const ColumnBundle &X, int nDensities, const GridInfo *gInfoOut=0)
 

Detailed Description

Macro Definition Documentation

#define MUL_MAT_MAT (   retType)
Value:
for(int i=0; i<3; i++) \
for(int j=0; j<3; j++) \
for(int k=0; k<3; k++) \
ret(i,j) += m(i,k) * n(k,j); \
return ret;
3x3 matrix
Definition: matrix3.h:31
#define MUL_MAT_VEC (   retType)
Value:
for(int i=0; i<3; i++) \
for(int j=0; j<3; j++) \
ret[i] += m(i,j) * v[j]; \
return ret;
Generic 3-vector.
Definition: vector3.h:34
#define MUL_VEC_MAT (   retType)
Value:
for(int i=0; i<3; i++) \
for(int j=0; j<3; j++) \
ret[j] += v[i] * m(i,j); \
return ret;
Generic 3-vector.
Definition: vector3.h:34

Function Documentation

template<typename scalar >
__hostanddev__ matrix3<scalar> adjugate ( const matrix3< scalar > &  m)
Parameters
mCalculate adjugate (matrix of signed cofactors)
void axpy ( complex  alpha,
const ManagedMemory< complex > &  x,
ManagedMemory< complex > &  y 
)

scale y *= x elementwise

standard blas y += alpha*x

__hostanddev__ double circDistanceSquared ( const vector3<> &  a,
const vector3<> &  b 
)

Return squared distance between a and b, interpreted as coordinates on a unit S1xS1xS1 embedded in R^6 i.e. sum of squares of distances between each coordinate put on a circle. This is useful for checking distances in a periodic cell safely.

matrix cisInv ( const matrix A,
matrix Bevecs = 0,
diagMatrix Beigs = 0 
)

Inverse function of cis: find matrix B such that A = exp(iota B). A must be unitary. Optionally retrieve the eigensystem of the result (note the real eigenvalues are of the input, not the result)

complex det ( const matrix A)

Compute the determinant of an arbitrary matrix A (via LU decomposition) If skipZeros is true, skips diagonal entries in the diagonal of the LU decomposition that are below a tolerance

template<typename scalar >
__hostanddev__ matrix3<scalar> Diag ( vector3< scalar >  v)
Parameters
vConstruct diagonal matrix
ScalarFieldArray diagouterI ( const diagMatrix F,
const ColumnBundle X,
int  nDensities,
const GridInfo gInfoOut = 0 
)

Returns diag((I*X)*F*(I*X)^) (Compute density from an orthonormal wavefunction ColumnBundle with some fillings F). nDensities is the number of scalar field in the output and controls how spin/spinors are handled: 1: return total density regardless of spin / spinor nature 2: return spin density in X.qnum->index()'th component of the output (valid for non-spinor X only) 4: return spin density-matrix (valid for spinor X only) If gInfoOut is specified, function ensures that the output is changed to that grid (in case tighter wfns grid is in use)

Parameters
XReciprocal to real space
ScalarFieldTildeArray Idag ( const ScalarFieldArray X)
inline
Parameters
XHermitian conjugate of I
ColumnBundle Idag_DiagV_I ( const ColumnBundle C,
const ScalarFieldArray V 
)

Return Idag V .* I C (evaluated columnwise) The handling of the spin structure of V parallels that of diagouterI, with V.size() taking the role of nDensities

matrix inv ( const matrix A)

Compute inverse of an arbitrary matrix A (via LU decomposition)

inverse of matrix

__hostanddev__ matrix3 inv ( const matrix3<> &  m)
Parameters
mMatrix inverse
matrix invsqrt ( const matrix A,
matrix Aevecs = 0,
diagMatrix Aeigs = 0,
bool *  isSingular = 0 
)

Compute matrix A^-0.5 and optionally the eigensystem of A (if non-null). If isSingular is provided, function will set it to true and return rather than stack-tracing in singular cases.

ScalarFieldTildeArray J ( const ScalarFieldArray X)
inline
Parameters
XReal to reciprocal space
ScalarFieldArray Jdag ( ScalarFieldTildeArray &&  X)
inline
Parameters
XHermitian conjugate of J
double nrm2 ( const matrixScaledTransOp A)

return transpose symmetric part: (A + AT)/2

return 2-norm (avoid explicit transpose or sub-matrix operations)

template<typename T >
void nullToZero ( TptrCollection x,
const GridInfo gInfo,
int  N = 0 
)

Allocate all null components and initialize them to zero If N is non-zero, resize collection to N scalar fields

template<class T >
TptrCollection operator+ ( const TptrCollection in1,
const TptrCollection in2 
)

< Add (destructible inputs)

template<class T >
TptrCollection operator- ( const TptrCollection in1,
const TptrCollection in2 
)

< Add (destructible inputs)

template<typename scalar >
__hostanddev__ matrix3<scalar> outer ( const vector3< scalar > &  a,
const vector3< scalar > &  b 
)
Parameters
bouter product
matrix pow ( const matrix A,
double  exponent,
matrix Aevecs = 0,
diagMatrix Aeigs = 0,
bool *  isSingular = 0 
)

Compute matrix A^exponent, and optionally the eigensystem of A (if non-null). If isSingular is provided, function will set it to true and return rather than stack-tracing in singular cases.

void randomize ( ColumnBundle x)

Copies the input.

Initialize to random numbers