QMCPACK
LinearGrid< T, CT > Struct Template Reference

One-Dimensional linear-grid. More...

+ Inheritance diagram for LinearGrid< T, CT >:
+ Collaboration diagram for LinearGrid< T, CT >:

Public Member Functions

std::unique_ptr< OneDimGridBase< T, CT > > makeClone () const override
 
int locate (T r) const override
 evaluate the index of r More...
 
void set (T ri, T rf, int n) override
 Set the grid given the parameters. More...
 
Ugrid einspline_grid ()
 
- Public Member Functions inherited from OneDimGridBase< T, CT >
 OneDimGridBase ()
 
virtual ~OneDimGridBase ()=default
 
int getGridTag () const
 
int getIndex (T r) const
 
T & operator[] (int i)
 assign a value More...
 
T & operator() (int i)
 assign a value More...
 
operator[] (int i) const
 return a value More...
 
operator() (int i) const
 return a value More...
 
const T * data () const
 
T * data ()
 
dh () const
 return the differential spacing of the grid More...
 
r (int i) const
 returns $ r(i)$ More...
 
dr (int i) const
 returns $ r(i+1)-r(i)$ More...
 
int size () const
 returns the size of the grid More...
 
rmin () const
 return the first grid point More...
 
rmax () const
 return the last grid point More...
 
template<typename T1 >
int getIndexAndDistanceFromGridPoint (T r, T1 &dist) const
 

Additional Inherited Members

- Public Types inherited from OneDimGridBase< T, CT >
using value_type = T
 
using Array_t = CT
 
- Public Attributes inherited from OneDimGridBase< T, CT >
int GridTag
 
int num_points
 
value_type lower_bound
 
value_type upper_bound
 
value_type Delta
 differential spacing of the grid More...
 
double DeltaInv
 
Array_t X
 array to store the radial grid data More...
 

Detailed Description

template<class T, class CT = Vector<T>>
struct qmcplusplus::LinearGrid< T, CT >

One-Dimensional linear-grid.

The analytic form $ r_i = r_0 + i\left( \frac{r_f - r_0}{N-1} \right), $ where $ N $ is the number of points and the index $ i $ runs from 0 to $ N-1 $.

Definition at line 134 of file OneDimGridBase.h.

Member Function Documentation

◆ einspline_grid()

Ugrid einspline_grid ( )
inline

Definition at line 166 of file OneDimGridBase.h.

Referenced by QMCFiniteSize::getSkSpline().

167  {
168  Ugrid grid;
169  grid.start = lower_bound;
170  grid.end = upper_bound;
171  grid.num = num_points;
172  grid.delta = Delta;
173  grid.delta_inv = DeltaInv;
174 
175  return grid;
176  }
value_type Delta
differential spacing of the grid

◆ locate()

int locate ( r) const
inlineoverridevirtual

evaluate the index of r

Parameters
rcurrent position

The grid index satisfies $ X[Loc] \ge r < X[Loc+1]$.

Implements OneDimGridBase< T, CT >.

Definition at line 149 of file OneDimGridBase.h.

149 { return static_cast<int>((static_cast<double>(r) - X[0]) * DeltaInv); }
Array_t X
array to store the radial grid data
T r(int i) const
returns

◆ makeClone()

std::unique_ptr<OneDimGridBase<T, CT> > makeClone ( ) const
inlineoverridevirtual

Implements OneDimGridBase< T, CT >.

Definition at line 144 of file OneDimGridBase.h.

Referenced by qmcplusplus::createSpline4RbyVs_temp(), qmcplusplus::createSpline4RbyVsDeriv_temp(), and qmcplusplus::TEST_CASE().

145  {
146  return std::make_unique<LinearGrid<T, CT>>(*this);
147  }

◆ set()

void set ( ri,
rf,
int  n 
)
inlineoverridevirtual

Set the grid given the parameters.

Parameters
riinitial grid point
rffinal grid point
nnumber of grid points

Implements OneDimGridBase< T, CT >.

Definition at line 151 of file OneDimGridBase.h.

Referenced by CoulombPBCAB::add(), SkParserBase::compute_grid(), RadialJastrowBuilder::createJ1(), SkParserBase::get_grid(), CoulombPBCAA::initBreakup(), CoulombPBCAB::initBreakup(), BackflowBuilder::makeShortRange_twoBody(), QMCFiniteSize::spline_clamped(), and qmcplusplus::TEST_CASE().

152  {
154  lower_bound = ri;
155  upper_bound = rf;
156  num_points = n;
157  // Delta is the differential spacing
158  X.resize(n);
159  double Delta_DP = (static_cast<double>(rf) - static_cast<double>(ri)) / static_cast<double>(n - 1);
160  Delta = Delta_DP;
161  DeltaInv = 1.0 / Delta_DP;
162  for (int i = 0; i < n; i++)
163  X[i] = ri + Delta_DP * i;
164  }
Array_t X
array to store the radial grid data
value_type Delta
differential spacing of the grid

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