MADNESS 0.10.1
Public Member Functions | Public Attributes | Static Private Member Functions | List of all members
madness::GaussianConvolution1D< Q > Class Template Reference

1D convolution with (derivative) Gaussian; coeff and expnt given in simulation coordinates [0,1] More...

#include <convolution1d.h>

Inheritance diagram for madness::GaussianConvolution1D< Q >:
Inheritance graph
[legend]
Collaboration diagram for madness::GaussianConvolution1D< Q >:
Collaboration graph
[legend]

Public Member Functions

 GaussianConvolution1D (int k, Q coeff, double expnt, int m, bool periodic, double bloch_k=0.0, KernelRange rng={})
 
virtual ~GaussianConvolution1D ()
 
bool issmall (Level n, Translation lx) const final
 
virtual Level natural_level () const final
 Returns the level for projection.
 
Tensor< Qrnlp (Level n, const Translation lx) const final
 Compute the projection of the operator onto the double order polynomials.
 
- Public Member Functions inherited from madness::Convolution1D< Q >
 Convolution1D (int k, int npt, int maxR, double bloch_k=0.0, KernelRange rng={})
 
virtual ~Convolution1D ()
 
bool get_issmall (Level n, Translation lx) const
 
const Tensor< Q > & get_rnlp (Level n, Translation lx) const
 
bool lattice_summed () const
 
const ConvolutionData1D< Q > * mod_nonstandard (const Key< 2 > &op_key) const
 Returns a pointer to the cached modified make_nonstandard form of the operator.
 
const ConvolutionData1D< Q > * nonstandard (Level n, Translation lx) const
 Returns a pointer to the cached make_nonstandard form of the operator.
 
Q phase (double R) const
 
bool range_restricted () const
 
const Tensor< Q > & rnlij (Level n, Translation lx, bool do_transpose=false) const
 Computes the transition matrix elements for the convolution for n,l.
 
bool rnlp_is_zero (Level n, Translation l) const
 

Public Attributes

const Q coeff
 Coefficient.
 
const double expnt
 Exponent.
 
const int m
 Order of derivative (0, 1, or 2 only)
 
const Level natlev
 Level to evaluate.
 
- Public Attributes inherited from madness::Convolution1D< Q >
double bloch_k
 k in exp(i k R) Bloch phase factor folded into lattice sum
 
Tensor< double > c
 
Tensor< double > hg
 
Tensor< double > hgT
 
Tensor< double > hgT2k
 
int k
 Wavelet order.
 
int maxR
 Number of lattice translations for sum.
 
SimpleCache< ConvolutionData1D< Q >, 2 > mod_ns_cache
 
int npt
 Number of quadrature points (is this used?)
 
SimpleCache< ConvolutionData1D< Q >, 1 > ns_cache
 
Tensor< double > quad_w
 
Tensor< double > quad_x
 
KernelRange range
 if range is nonnull, kernel range limited to to range (in simulation cell units), useful for finite-range convolutions with periodic functions
 
SimpleCache< Tensor< Q >, 1 > rnlij_cache
 
SimpleCache< Tensor< Q >, 1 > rnlp_cache
 

Static Private Member Functions

static int maxR (bool periodic, double expnt, const KernelRange &rng={})
 

Additional Inherited Members

- Public Types inherited from madness::Convolution1D< Q >
typedef Q opT
 The apply function uses this to infer resultT=opT*inputT.
 

Detailed Description

template<typename Q>
class madness::GaussianConvolution1D< Q >

1D convolution with (derivative) Gaussian; coeff and expnt given in simulation coordinates [0,1]

Note that the derivative is computed in simulation coordinates so you must be careful to scale the coefficients correctly.

Constructor & Destructor Documentation

◆ GaussianConvolution1D()

template<typename Q >
madness::GaussianConvolution1D< Q >::GaussianConvolution1D ( int  k,
Q  coeff,
double  expnt,
int  m,
bool  periodic,
double  bloch_k = 0.0,
KernelRange  rng = {} 
)
inlineexplicit

◆ ~GaussianConvolution1D()

template<typename Q >
virtual madness::GaussianConvolution1D< Q >::~GaussianConvolution1D ( )
inlinevirtual

Member Function Documentation

◆ issmall()

template<typename Q >
bool madness::GaussianConvolution1D< Q >::issmall ( Level  n,
Translation  lx 
) const
inlinefinalvirtual
Returns
true if the block of [r^n_l]_ij is expected to be small

Implements madness::Convolution1D< Q >.

References beta, madness::GaussianConvolution1D< Q >::expnt, and pow().

◆ maxR()

template<typename Q >
static int madness::GaussianConvolution1D< Q >::maxR ( bool  periodic,
double  expnt,
const KernelRange rng = {} 
)
inlinestaticprivate

◆ natural_level()

template<typename Q >
virtual Level madness::GaussianConvolution1D< Q >::natural_level ( ) const
inlinefinalvirtual

Returns the level for projection.

Reimplemented from madness::Convolution1D< Q >.

References madness::GaussianConvolution1D< Q >::natlev.

◆ rnlp()

template<typename Q >
Tensor< Q > madness::GaussianConvolution1D< Q >::rnlp ( Level  n,
const Translation  lx 
) const
inlinefinalvirtual

Compute the projection of the operator onto the double order polynomials.

The returned reference is to a cached tensor ... if you want to modify it, take a copy first.

Return in v[p] p=0..2*k-1

r(n,l,p) = 2^(-n) * int(K(2^(-n)*(z+l)) * phi(p,z), z=0..1)
char * p(char *buf, const char *name, int k, int initial_level, double thresh, int order)
Definition derivatives.cc:72
vecfuncT K(vecfuncT &ket, vecfuncT &bra, vecfuncT &vf)

The kernel is coeff*exp(-expnt*z^2)*z^m (with m>0). This is equivalent to

r(n,l,p) = 2^(-n*(m+1))*coeff * int( ((d/dz)^m exp(-beta*z^2)) * phi(p,z-l), z=l..l+1)
const int m
Order of derivative (0, 1, or 2 only)
Definition convolution1d.h:753
const Q coeff
Coefficient.
Definition convolution1d.h:750
const double beta
Definition gygi_soltion.cc:62
static const double d
Definition nonlinschro.cc:121

where

beta = alpha * 2^(-2*n)
static const double alpha
Definition testcosine.cc:10

Implements madness::Convolution1D< Q >.

References std::abs(), beta, madness::GaussianConvolution1D< Q >::coeff, e(), madness::GaussianConvolution1D< Q >::expnt, madness::KernelRange::finite(), madness::KernelRange::finite_hard(), madness::KernelRange::finite_soft(), h(), madness::KernelRange::iextent_x2(), madness::Convolution1D< Q >::k, L, madness::legendre_scaling_functions(), madness::GaussianConvolution1D< Q >::m, MADNESS_ASSERT, madness::KernelRange::N(), madness::Convolution1D< Q >::npt, p(), pow(), Q(), madness::Convolution1D< Q >::quad_w, madness::Convolution1D< Q >::quad_x, madness::Convolution1D< Q >::range, madness::Convolution1D< Q >::rnlp_is_zero(), madness::KernelRange::sigma(), v, madness::KernelRange::value(), and x0.

Referenced by main(), madness::test_rnlp_rangelimited(), and madness::test_rnlp_rangelimited_erf().

Member Data Documentation

◆ coeff

template<typename Q >
const Q madness::GaussianConvolution1D< Q >::coeff

Coefficient.

Referenced by madness::GaussianConvolution1D< Q >::rnlp().

◆ expnt

template<typename Q >
const double madness::GaussianConvolution1D< Q >::expnt

◆ m

template<typename Q >
const int madness::GaussianConvolution1D< Q >::m

Order of derivative (0, 1, or 2 only)

Referenced by madness::GaussianConvolution1D< Q >::rnlp().

◆ natlev

template<typename Q >
const Level madness::GaussianConvolution1D< Q >::natlev

Level to evaluate.

Referenced by madness::GaussianConvolution1D< Q >::natural_level().


The documentation for this class was generated from the following file: