JDFTx  1.3.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Simulation grid and data management

Files

file  Operators.h
 Operators on ScalarField's and ScalarFieldTilde's.
 
file  ScalarField.h
 Real and complex scalar fields in real and reciprocal space.
 
file  ScalarFieldIO.h
 I/O utilities for the data arrays.
 
file  VectorField.h
 Generic multiplet of data arrays (and specialized to triplets for vector fields in real/reciprocal space)
 

Classes

struct  FieldData
 Base class for ScalarFieldData and ScalarFieldTildeData. 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  ScalarFieldMultiplet< T, N >
 Generic multiplet object with overloaded arithmetic. More...
 

Macros

#define Tptr   std::shared_ptr<T>
 shorthand for writing the template operators (undef'd at end of header)
 
#define DECLARE_DATA_PREF_ACCESS
 
#define DECLARE_DATA_ACCESS
 
#define Tptr   std::shared_ptr<T>
 
#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 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

ScalarField Real (const complexScalarField &)
 real part of a complex scalar field (real-space)
 
ScalarFieldTilde Real (const complexScalarFieldTilde &)
 real part of a complex scalar field (reciprocal space)
 
ScalarField Imag (const complexScalarField &)
 imaginary part of a complex scalar field (real-space)
 
ScalarFieldTilde Imag (const complexScalarFieldTilde &)
 imaginary part of a complex scalar field (reciprocal space)
 
complexScalarField Complex (const ScalarField &)
 convert real to complex scalar field with zero imaginary part (real-space)
 
complexScalarField Complex (const ScalarField &re, const ScalarField &im)
 construct complex scalar field fromr eal and imaginary parts (real-space)
 
complexScalarFieldTilde Complex (const ScalarFieldTilde &)
 convert real to complex scalar field with zero imaginary part (reciprocal-space)
 
ScalarFieldTilde O (const ScalarFieldTilde &)
 Inner product operator (diagonal in PW basis)
 
ScalarFieldTilde O (ScalarFieldTilde &&)
 Inner product operator (diagonal in PW basis)
 
complexScalarFieldTilde O (const complexScalarFieldTilde &)
 Inner product operator (diagonal in PW basis)
 
complexScalarFieldTilde O (complexScalarFieldTilde &&)
 Inner product operator (diagonal in PW basis)
 
ScalarField I (const ScalarFieldTilde &, int nThreads=0)
 Forward transform: PW basis -> real space (preserve input)
 
ScalarField I (ScalarFieldTilde &&, int nThreads=0)
 Forward transform: PW basis -> real space (destructible input)
 
complexScalarField I (const complexScalarFieldTilde &, int nThreads=0)
 Forward transform: PW basis -> real space (preserve input)
 
complexScalarField I (complexScalarFieldTilde &&, int nThreads=0)
 Forward transform: PW basis -> real space (destructible input)
 
ScalarFieldTilde J (const ScalarField &, int nThreads=0)
 Inverse transform: Real space -> PW basis.
 
complexScalarFieldTilde J (const complexScalarField &, int nThreads=0)
 Inverse transform: Real space -> PW basis (preserve input)
 
complexScalarFieldTilde J (complexScalarField &&, int nThreads=0)
 Inverse transform: Real space -> PW basis (destructible input)
 
ScalarFieldTilde Idag (const ScalarField &, int nThreads=0)
 Forward transform transpose: Real space -> PW basis.
 
complexScalarFieldTilde Idag (const complexScalarField &, int nThreads=0)
 Forward transform transpose: Real space -> PW basis (preserve input)
 
complexScalarFieldTilde Idag (complexScalarField &&, int nThreads=0)
 Forward transform transpose: Real space -> PW basis (destructible input)
 
ScalarField Jdag (const ScalarFieldTilde &, int nThreads=0)
 Inverse transform transpose: PW basis -> real space (preserve input)
 
ScalarField Jdag (ScalarFieldTilde &&, int nThreads=0)
 Inverse transform transpose: PW basis -> real space (destructible input)
 
complexScalarField Jdag (const complexScalarFieldTilde &, int nThreads=0)
 Inverse transform transpose: PW basis -> real space (preserve input)
 
complexScalarField Jdag (complexScalarFieldTilde &&, int nThreads=0)
 Inverse transform transpose: PW basis -> real space (destructible input)
 
ScalarField JdagOJ (const ScalarField &)
 Evaluate Jdag(O(J())), which avoids 2 fourier transforms in PW basis (preserve input)
 
ScalarField JdagOJ (ScalarField &&)
 Evaluate Jdag(O(J())), which avoids 2 fourier transforms in PW basis (destructible input)
 
complexScalarField JdagOJ (const complexScalarField &)
 Evaluate Jdag(O(J())), which avoids 2 fourier transforms in PW basis (preserve input)
 
complexScalarField JdagOJ (complexScalarField &&)
 Evaluate Jdag(O(J())), which avoids 2 fourier transforms in PW basis (destructible input)
 
ScalarFieldTilde L (const ScalarFieldTilde &)
 Laplacian.
 
ScalarFieldTilde L (ScalarFieldTilde &&)
 Laplacian.
 
complexScalarFieldTilde L (const complexScalarFieldTilde &)
 Laplacian.
 
complexScalarFieldTilde L (complexScalarFieldTilde &&)
 Laplacian.
 
ScalarFieldTilde Linv (const ScalarFieldTilde &)
 Inverse Laplacian.
 
ScalarFieldTilde Linv (ScalarFieldTilde &&)
 Inverse Laplacian.
 
complexScalarFieldTilde Linv (const complexScalarFieldTilde &)
 Inverse Laplacian.
 
complexScalarFieldTilde Linv (complexScalarFieldTilde &&)
 Inverse Laplacian.
 
void zeroNyquist (RealKernel &Gdata)
 zeros out all the nyquist components of a real G-kernel
 
void zeroNyquist (ScalarFieldTilde &Gptr)
 zeros out all the nyquist components of a G-space data array
 
void zeroNyquist (ScalarField &Rptr)
 zeros out all the nyquist components of an R-space data array
 
ScalarField exp (const ScalarField &)
 Elementwise exponential (preserve input)
 
ScalarField exp (ScalarField &&)
 Elementwise exponential (destructible input)
 
ScalarField log (const ScalarField &)
 Elementwise logarithm (preserve input)
 
ScalarField log (ScalarField &&)
 Elementwise logarithm (destructible input)
 
ScalarField sqrt (const ScalarField &)
 Elementwise square root (preserve input)
 
ScalarField sqrt (ScalarField &&)
 Elementwise square root (destructible input)
 
ScalarField inv (const ScalarField &)
 Elementwise reciprocal (preserve input)
 
ScalarField inv (ScalarField &&)
 Elementwise reciprocal (destructible input)
 
ScalarField pow (const ScalarField &, double alpha)
 Elementwise power (preserve input)
 
ScalarField pow (ScalarField &&, double alpha)
 Elementwise power (destructible input)
 
template<class T >
Tptr clone (const Tptr &X)
 Clone (NOTE: operator= is by reference for the ScalarField classes)
 
template<class T >
Tptroperator*= (Tptr &in, double scaleFac)
 Scale.
 
template<class T >
Tptr operator* (const Tptr &in, double scaleFac)
 Scalar multiply (preserve input)
 
template<class T >
Tptr operator* (double scaleFac, const Tptr &in)
 Scalar multiply (preserve input)
 
template<class T >
Tptr operator* (Tptr &&in, double scaleFac)
 Scalar multiply (destructible input)
 
template<class T >
Tptr operator* (double scaleFac, Tptr &&in)
 Scalar multiply (destructible input)
 
template<class T >
Tptr conj (Tptr &&in)
 Generic elementwise conjugate for complex data:
 
template<class T >
Tptr conj (const Tptr &in)
 
template<class T >
Tptroperator*= (Tptr &in, const Tptr &other)
 Generic elementwise multiply for complex data:
 
ScalarFieldoperator*= (ScalarField &in, const ScalarField &other)
 Elementwise multiply for real data.
 
template<class T >
Tptr operator* (const Tptr &in1, const Tptr &in2)
 Elementwise multiply (preserve inputs)
 
template<class T >
Tptr operator* (const Tptr &in1, Tptr &&in2)
 Elementwise multiply (destructible input)
 
template<class T >
Tptr operator* (Tptr &&in1, const Tptr &in2)
 Elementwise multiply (destructible input)
 
template<class T >
Tptr operator* (Tptr &&in1, Tptr &&in2)
 Elementwise multiply (destructible inputs)
 
complexScalarFieldoperator*= (complexScalarField &, const ScalarField &)
 elementwise multiply
 
complexScalarField operator* (const complexScalarField &, const ScalarField &)
 elementwise multiply (preserve inputs)
 
complexScalarField operator* (const ScalarField &, const complexScalarField &)
 elementwise multiply (preserve inputs)
 
complexScalarField operator* (complexScalarField &&, const ScalarField &)
 elementwise multiply (destructible inputs)
 
complexScalarField operator* (const ScalarField &, complexScalarField &&)
 elementwise multiply (destructible inputs)
 
ScalarFieldTildeoperator*= (ScalarFieldTilde &, const RealKernel &)
 Elementwise multiply.
 
ScalarFieldTilde operator* (const RealKernel &, const ScalarFieldTilde &)
 Elementwise multiply (preserve inputs)
 
ScalarFieldTilde operator* (const ScalarFieldTilde &, const RealKernel &)
 Elementwise multiply (preserve inputs)
 
ScalarFieldTilde operator* (const RealKernel &, ScalarFieldTilde &&)
 Elementwise multiply (destructible input)
 
ScalarFieldTilde operator* (ScalarFieldTilde &&, const RealKernel &)
 Elementwise multiply (destructible input)
 
template<typename T >
void axpy (double alpha, const Tptr &X, Tptr &Y)
 Generic axpy for complex data types (Note: null pointers are treated as zero)
 
void axpy (double alpha, const ScalarField &X, ScalarField &Y)
 Real data Linear combine: Y += alpha * X (Note: null pointers are treated as zero)
 
template<class T >
Tptroperator+= (Tptr &in, const Tptr &other)
 Increment.
 
template<class T >
Tptroperator-= (Tptr &in, const Tptr &other)
 Decrement.
 
template<class T >
Tptr operator+ (const Tptr &in1, const Tptr &in2)
 Add (preserve inputs)
 
template<class T >
Tptr operator+ (const Tptr &in1, Tptr &&in2)
 Add (destructible input)
 
template<class T >
Tptr operator+ (Tptr &&in1, const Tptr &in2)
 Add (destructible input)
 
template<class T >
Tptr operator+ (Tptr &&in1, Tptr &&in2)
 Add (destructible inputs)
 
template<class T >
Tptr operator- (const Tptr &in1, const Tptr &in2)
 Subtract (preserve inputs)
 
template<class T >
Tptr operator- (const Tptr &in1, Tptr &&in2)
 Subtract (destructible input)
 
template<class T >
Tptr operator- (Tptr &&in1, const Tptr &in2)
 Subtract (destructible input)
 
template<class T >
Tptr operator- (Tptr &&in1, Tptr &&in2)
 Subtract (destructible inputs)
 
template<class T >
Tptr operator- (const Tptr &in)
 Negate.
 
template<class T >
Tptr operator- (Tptr &&in)
 Negate.
 
ScalarFieldoperator+= (ScalarField &, double)
 Increment by scalar.
 
ScalarField operator+ (double, const ScalarField &)
 Add scalar (preserve inputs)
 
ScalarField operator+ (const ScalarField &, double)
 Add scalar (preserve inputs)
 
ScalarField operator+ (double, ScalarField &&)
 Add scalar (destructible input)
 
ScalarField operator+ (ScalarField &&, double)
 Add scalar (destructible input)
 
ScalarFieldoperator-= (ScalarField &, double)
 Decrement by scalar.
 
ScalarField operator- (double, const ScalarField &)
 Subtract from scalar (preserve inputs)
 
ScalarField operator- (const ScalarField &, double)
 Subtract scalar (preserve inputs)
 
ScalarField operator- (double, ScalarField &&)
 Subtract from scalar (destructible input)
 
ScalarField operator- (ScalarField &&, double)
 Subtract scalar (destructible input)
 
template<typename T >
complex dot (const Tptr &X, const Tptr &Y)
 
template<typename T >
double nrm2 (const Tptr &X)
 
template<typename T >
complex sum (const Tptr &X)
 
double dot (const ScalarField &, const ScalarField &)
 Inner product.
 
double dot (const ScalarFieldTilde &, const ScalarFieldTilde &)
 Inner product.
 
double nrm2 (const ScalarField &)
 2-norm
 
double nrm2 (const ScalarFieldTilde &)
 2-norm
 
double sum (const ScalarField &)
 Sum of elements.
 
double sum (const ScalarFieldTilde &)
 Equivalent to dot() with a ScalarFieldTilde of all 1s (NOTE: sum(X) != sum(I(X)))
 
double integral (const ScalarField &)
 Integral in the unit cell (dV times sum())
 
double integral (const ScalarFieldTilde &)
 Integral in the unit cell (just fetches the G=0 component with correct prefactor)
 
complex integral (const complexScalarField &)
 Integral in the unit cell (dV times sum())
 
complex integral (const complexScalarFieldTilde &)
 Integral in the unit cell (just fetches the G=0 component with correct prefactor)
 
ScalarFieldTilde changeGrid (const ScalarFieldTilde &, const GridInfo &gInfoNew)
 
ScalarField changeGrid (const ScalarField &, const GridInfo &gInfoNew)
 
complexScalarFieldTilde changeGrid (const complexScalarFieldTilde &, const GridInfo &gInfoNew)
 
complexScalarField changeGrid (const complexScalarField &, const GridInfo &gInfoNew)
 
template<typename T >
void initZero (Tptr &X)
 
template<typename T >
void initZero (Tptr &X, const GridInfo &gInfo)
 
template<typename T >
void nullToZero (Tptr &X, const GridInfo &gInfo)
 If X is null, allocate and initialize to 0.
 
void initRandom (ScalarField &, double cap=0.0)
 initialize element-wise with a unit-normal random number (with a cap if cap>0)
 
void initRandomFlat (ScalarField &)
 initialize element-wise with a unit-flat [0:1) random number
 
void initGaussianKernel (RealKernel &, double x0)
 initialize to gaussian kernel exp(-(G x0/2)^2)
 
void initTranslation (ScalarFieldTilde &, const vector3<> &r)
 initialize to translation operator exp(-i G.r)
 
ScalarFieldTilde gaussConvolve (const ScalarFieldTilde &, double sigma)
 convolve with a gaussian
 
ScalarFieldTilde gaussConvolve (ScalarFieldTilde &&, double sigma)
 convolve with a gaussian (destructible input)
 
template<typename Func , typename... Args>
void applyFuncGsq (const GridInfo &gInfo, const Func &f, Args...args)
 Evaluate a function f(i, Gsq, args...) at each point in reciprocal space indexed by i.
 
template<typename Func , typename... Args>
void applyFunc_r (const GridInfo &gInfo, const Func &f, Args...args)
 Evaluate a function f(i, r, args...) at each point in real space index by i.
 
void printStats (const ScalarField &X, const char *name, FILE *fp=stdout)
 Print mean and standard deviation of array with specified name (debug utility)
 
template<typename Callable , typename Vec >
void checkSymmetry (Callable *func, const Vec &v1, const Vec &v2, const char *funcName)
 
template<typename T >
void saveRawBinary (const Tptr &X, FILE *fp)
 Save the data in raw binary format to stream.
 
template<typename T >
void saveRawBinary (const Tptr &X, const char *filename)
 Save the data in raw binary format to file.
 
template<typename T >
void loadRawBinary (Tptr &X, FILE *fp)
 Load the data in raw binary format from stream.
 
template<typename T >
void loadRawBinary (Tptr &X, const char *filename)
 Load the data in raw binary format from file.
 
void saveDX (const ScalarField &, const char *filenamePrefix)
 
std::vector< std::vector
< double > > 
sphericalize (const ScalarField *dataR, int nColumns, double drFac=1.0, vector3<> *center=0)
 
void saveSphericalized (const ScalarField *dataR, int nColumns, const char *filename, double drFac=1.0, vector3<> *center=0)
 
void saveSphericalized (const ScalarFieldTilde *dataG, int nColumns, const char *filename, double dGFac=1.0)
 
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 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
 
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)
 

Detailed Description

Macro Definition Documentation

#define DECLARE_DATA_ACCESS
Value:
DataType* data(bool shouldAbsorbScale=true) { return (DataType*)FieldData::data(shouldAbsorbScale); } \
const DataType* data(bool shouldAbsorbScale=true) const { return (const DataType*)FieldData::data(shouldAbsorbScale); } \
DataType* dataGpu(bool shouldAbsorbScale=true) { return (DataType*)FieldData::dataGpu(shouldAbsorbScale); } \
const DataType* dataGpu(bool shouldAbsorbScale=true) const { return (const DataType*)FieldData::dataGpu(shouldAbsorbScale); }
complex * dataGpu()
Get a gpu data pointer (must be called from GPU owner thread)
complex * data()
Return a pointer to the actual data Return a CPU pointer to the actual data, will move data from GPU ...
#define DECLARE_DATA_PREF_ACCESS
Value:
DataType* dataPref(bool shouldAbsorbScale=true) { return dataGpu(shouldAbsorbScale); } \
const DataType* dataPref(bool shouldAbsorbScale=true) const { return dataGpu(shouldAbsorbScale); }
std::vector< typename T::DataType * > dataPref(TptrCollection &x)
Extract a std::vector of data pointers from a ScalarFieldArray.
Definition: ScalarFieldArray.h:37

Function Documentation

template<typename Callable , typename Vec >
void checkSymmetry ( Callable *  func,
const Vec &  v1,
const Vec &  v2,
const char *  funcName 
)

Check the symmetry of a linear operator

Template Parameters
CallableAn operator with function call signature Vec Callable(const Vec&)
VecAny operand type representing an element of a vector space
template<typename T >
complex dot ( const Tptr X,
const Tptr Y 
)
Parameters
YGeneric inner product for complex types
template<typename T >
double nrm2 ( const Tptr X)
Parameters
XGeneric 2-norm for complex types
void saveDX ( const ScalarField ,
const char *  filenamePrefix 
)
Save data to a raw binary along with a DataExplorer header
Parameters
filenamePrefixBinary data is saved to filenamePrefix.bin with DataExplorer header filenamePrefix.dx
void saveSphericalized ( const ScalarField dataR,
int  nColumns,
const char *  filename,
double  drFac = 1.0,
vector3<> *  center = 0 
)
Saves an array of real space data pointers to a multicolumn 1D 'sphericalized' file (for gnuplot)
Parameters
dataRThe data to sphericalize and save
nColumnsNumber of ScalarField's in dataR[]
filenameOutput file in which column 1 will be the radius, column 2 to nColumns+1 would be the sphericalized versions of dataR[0 to nColumns-1]
drFacis the spacing in radius as a fraction of the diameter of the sample box (R ./ S) (drFac << 1 is likely to give noisy results, particularly close to r=0)
centerThe origin for spherical coordinates [default = center of box (if null pointer is passed)]
void saveSphericalized ( const ScalarFieldTilde dataG,
int  nColumns,
const char *  filename,
double  dGFac = 1.0 
)
Saves an array of reciprocal space data pointers to a multicolumn 1D 'sphericalized' file (for gnuplot)
Parameters
dataGThe data to sphericalize (about G=0) and save
nColumnsNumber of ScalarFieldTilde's in dataG[]
filenameOutput file in which column 1 will be the radius, column 2 to nColumns+1 would be the sphericalized versions of dataG[0 to nColumns-1]
dGFacis the spacing in radius as a fraction of the diameter of the Brillouin zone (dGFac << 1 is likely to give noisy results, particularly close to G=0)
std::vector< std::vector<double> > sphericalize ( const ScalarField dataR,
int  nColumns,
double  drFac = 1.0,
vector3<> *  center = 0 
)
Spherically average scalar fields about an arbitrary center (with Wigner-Seitz wrapping)
Parameters
dataRThe data to sphericalize and save
nColumnsNumber of ScalarField's in dataR[]
drFacis the spacing in radius as a fraction of the diameter of the sample box (R ./ S) (drFac << 1 is likely to give noisy results, particularly close to r=0)
centerThe origin for spherical coordinates [default = center of box (if null pointer is passed)]
Returns
The first array contains the radial grid, and the subsequent ones the spherically-averaged results, one for each dataR, and the last column contains the weight of the radial grid point
template<typename T >
complex sum ( const Tptr X)
Parameters
XGeneric sum for complex types