QMCPACK
ShortRangeCuspFunctor< T > Struct Template Reference
+ Inheritance diagram for ShortRangeCuspFunctor< T >:
+ Collaboration diagram for ShortRangeCuspFunctor< T >:

Public Member Functions

 ShortRangeCuspFunctor (const std::string &my_name)
 default constructor More...
 
void setCusp (real_type cusp) override
 sets the cusp condition and disables optimization of the cusp-determining parameter More...
 
OptimizableFunctorBasemakeClone () const override
 clone the functor More...
 
void reset () override
 Implement the reset function, which was pure virtual in OptimizableFunctorBase, even though we don't need it. More...
 
real_type evaluate (real_type r) const
 compute U(r) at a particular value of r More...
 
real_type evaluate (real_type r, real_type &dudr, real_type &d2udr2) const
 compute U(r), dU/dr, and d^2U/dr^2 at a particular value of r More...
 
real_type evaluate (real_type r, real_type &dudr, real_type &d2udr2, real_type &d3udr3) const
 compute U(r), dU/dr, d^2U/dr^2, and d^3U/dr^3 More...
 
real_type evaluateV (const int iat, const int iStart, const int iEnd, const T *restrict _distArray, T *restrict distArrayCompressed) const
 compute U(r) at multiple distances and sum the results More...
 
void evaluateVGL (const int iat, const int iStart, const int iEnd, const T *distArray, T *restrict valArray, T *restrict gradArray, T *restrict laplArray, T *restrict distArrayCompressed, int *restrict distIndices) const
 compute U(r), dU/dr, and d^2U/dr^2 at multiple values of r More...
 
real_type f (real_type r) override
 compute U(r) at a particular distance or return zero if beyond the cutoff More...
 
real_type df (real_type r) override
 compute dU/dr at a particular distance or return zero if beyond the cutoff More...
 
bool evaluateDerivatives (real_type r, std::vector< TinyVector< real_type, 3 >> &derivs) override
 compute derivatives of U(r), dU/dr, and d^2U/dr^2 with respect to the variational parameters More...
 
bool evaluateDerivatives (real_type r, std::vector< real_type > &derivs) override
 compute derivatives of U(r) with respect to the variational parameters More...
 
template<class U >
void set_variable_from_xml (ReportEngine &PRE, xmlNodePtr cur, U &variable_to_set, std::string &id_to_set, bool &opt_to_set)
 set up a variational parameter using the supplied xml node More...
 
bool put (xmlNodePtr cur) override
 read in information about the functor from an xml node More...
 
void checkInVariablesExclusive (opt_variables_type &active) override
 check in variational parameters to the global list of parameters used by the optimizer. More...
 
void checkOutVariables (const opt_variables_type &active) override
 check out variational optimizable variables More...
 
void resetParametersExclusive (const opt_variables_type &active) override
 reset the parameters during optimizations. More...
 
- Public Member Functions inherited from OptimizableFunctorBase
 OptimizableFunctorBase (const std::string &name="")
 default constructor More...
 
virtual ~OptimizableFunctorBase ()=default
 virtual destrutor More...
 
void getIndex (const opt_variables_type &active)
 
virtual real_type f (real_type r)=0
 evaluate the value at r More...
 
virtual real_type df (real_type r)=0
 evaluate the first derivative More...
 
virtual void setDensity (real_type n)
 empty virtual function to help builder classes More...
 
virtual void setCusp (real_type cusp)
 empty virtual function to help builder classes More...
 
virtual void setPeriodic (bool periodic)
 empty virtual function to help builder classes More...
 
virtual bool evaluateDerivatives (real_type r, std::vector< qmcplusplus::TinyVector< real_type, 3 >> &derivs)
 
virtual bool evaluateDerivatives (real_type r, std::vector< real_type > &derivs)
 
virtual void setGridManager (bool willmanage)
 
virtual void checkInVariablesExclusive (opt_variables_type &active)=0
 check in variational parameters to the global list of parameters used by the optimizer. More...
 
virtual void resetParametersExclusive (const opt_variables_type &active)=0
 reset the parameters during optimizations More...
 
- Public Member Functions inherited from OptimizableObject
 OptimizableObject (const std::string &name)
 
const std::string & getName () const
 
bool isOptimized () const
 
virtual void reportStatus (std::ostream &os)
 print the state, e.g., optimizables More...
 
void setOptimization (bool state)
 
virtual void writeVariationalParameters (hdf_archive &hout)
 Write the variational parameters for this object to the VP HDF file. More...
 
virtual void readVariationalParameters (hdf_archive &hin)
 Read the variational parameters for this object from the VP HDF file. More...
 

Static Public Member Functions

static void mw_evaluateV (const int num_groups, const ShortRangeCuspFunctor *const functors[], const int n_src, const int *grp_ids, const int num_pairs, const int *ref_at, const T *mw_dist, const int dist_stride, T *mw_vals, Vector< char, OffloadPinnedAllocator< char >> &transfer_buffer)
 evaluate sum of the pair potentials FIXME More...
 

Public Attributes

bool Opt_A
 true, if A is optimizable More...
 
bool Opt_R0
 true, if R0 is optimizable More...
 
bool Opt_B
 true, if B variables are optimizable More...
 
real_type A
 variable that controls the cusp More...
 
real_type R0
 the soft cutoff distance that controls how short ranged the functor is More...
 
std::vector< real_typeB
 variables that add detail through an expansion in sigmoidal functions More...
 
std::string ID_A
 id of A More...
 
std::string ID_R0
 id of R0 More...
 
std::string ID_B
 id of B More...
 
- Public Attributes inherited from OptimizableFunctorBase
real_type cutoff_radius = 0.0
 maximum cutoff More...
 
opt_variables_type myVars
 set of variables to be optimized More...
 

Additional Inherited Members

- Public Types inherited from OptimizableFunctorBase
using real_type = optimize::VariableSet::real_type
 typedef for real values More...
 
using opt_variables_type = optimize::VariableSet
 typedef for variableset: this is going to be replaced More...
 

Detailed Description

template<class T>
struct qmcplusplus::ShortRangeCuspFunctor< T >

Definition at line 57 of file ShortRangeCuspFunctor.h.

Constructor & Destructor Documentation

◆ ShortRangeCuspFunctor()

ShortRangeCuspFunctor ( const std::string &  my_name)
inline

default constructor

Definition at line 87 of file ShortRangeCuspFunctor.h.

References OptimizableFunctorBase::cutoff_radius, and ShortRangeCuspFunctor< T >::reset().

Referenced by ShortRangeCuspFunctor< T >::makeClone().

88  : OptimizableFunctorBase(my_name),
89  Opt_A(false),
90  Opt_R0(true),
91  Opt_B(true),
92  A(1.0),
93  R0(0.06),
94  ID_A("string_not_set"),
95  ID_R0("string_not_set"),
96  ID_B("string_not_set")
97  {
98  cutoff_radius = 1.0e4; //some big range
99  reset();
100  }
bool Opt_R0
true, if R0 is optimizable
void reset() override
Implement the reset function, which was pure virtual in OptimizableFunctorBase, even though we don&#39;t ...
OptimizableFunctorBase(const std::string &name="")
default constructor
real_type A
variable that controls the cusp
bool Opt_A
true, if A is optimizable
bool Opt_B
true, if B variables are optimizable
real_type R0
the soft cutoff distance that controls how short ranged the functor is

Member Function Documentation

◆ checkInVariablesExclusive()

void checkInVariablesExclusive ( opt_variables_type active)
inlineoverridevirtual

check in variational parameters to the global list of parameters used by the optimizer.

Parameters
activea super set of optimizable variables

The existing checkInVariables implementation in WFC/SPO/.. are inclusive and it calls checkInVariables of its members class A: public SPOSet {} class B: public WFC { A objA; checkInVariables() { objA.checkInVariables(); } };

With OptimizableObject, class A: public OptimizableObject {} class B: public OptimizableObject { A objA; checkInVariablesExclusive() { // should not call objA.checkInVariablesExclusive() if objA has been extracted; } }; A vector of OptimizableObject, will be created by calling extractOptimizableObjects(). All the checkInVariablesExclusive() will be called through this vector and thus checkInVariablesExclusive implementation should only handle non-OptimizableObject members.

Implements OptimizableObject.

Definition at line 636 of file ShortRangeCuspFunctor.h.

References VariableSet::insertFrom(), and OptimizableFunctorBase::myVars.

Referenced by qmcplusplus::TEST_CASE().

637  {
638  active.insertFrom(myVars);
639  }
opt_variables_type myVars
set of variables to be optimized

◆ checkOutVariables()

void checkOutVariables ( const opt_variables_type active)
inlineoverridevirtual

check out variational optimizable variables

Parameters
activea super set of optimizable variables

Implements OptimizableFunctorBase.

Definition at line 641 of file ShortRangeCuspFunctor.h.

References VariableSet::getIndex(), and OptimizableFunctorBase::myVars.

642  {
643  myVars.getIndex(active);
644  }
int getIndex(const std::string &vname) const
return the Index vaule for the named parameter
opt_variables_type myVars
set of variables to be optimized

◆ df()

real_type df ( real_type  r)
inlineoverride

compute dU/dr at a particular distance or return zero if beyond the cutoff

Definition at line 275 of file ShortRangeCuspFunctor.h.

References OptimizableFunctorBase::cutoff_radius, and ShortRangeCuspFunctor< T >::evaluate().

276  {
277  if (r >= cutoff_radius)
278  return 0.0;
279  real_type du, d2u;
280  evaluate(r, du, d2u);
281  return du;
282  }
real_type evaluate(real_type r) const
compute U(r) at a particular value of r
OHMMS_PRECISION real_type

◆ evaluate() [1/3]

real_type evaluate ( real_type  r) const
inline

compute U(r) at a particular value of r

Definition at line 121 of file ShortRangeCuspFunctor.h.

References ShortRangeCuspFunctor< T >::A, ShortRangeCuspFunctor< T >::B, OptimizableFunctorBase::cutoff_radius, qmcplusplus::exp(), ShortRangeCuspFunctor< T >::R0, and qmcplusplus::Units::time::s.

Referenced by ShortRangeCuspFunctor< T >::df(), ShortRangeCuspFunctor< T >::evaluateV(), ShortRangeCuspFunctor< T >::evaluateVGL(), ShortRangeCuspFunctor< T >::f(), and qmcplusplus::TEST_CASE().

122  {
123  // the functor will be almost exactly zero at long range, so just return that if r is beyond the hard cutoff
124  if (r >= cutoff_radius)
125  return 0.0;
126 
127  // get the ratio of the distance and the soft cutoff distance
128  const real_type s = r / R0;
129 
130  // sum up the sigmoidal function expansion
131  real_type sig_sum = 0.0;
132  real_type sn = s * s; // s^n
133  for (int i = 0; i < B.size(); i++)
134  {
135  sig_sum += B[i] * sn / (1.0 + sn);
136  sn *= s; // update s^n
137  }
138 
139  // return U(r)
140  return -1.0 * std::exp(-s) * (A * R0 + sig_sum);
141  }
real_type A
variable that controls the cusp
MakeReturn< UnaryNode< FnExp, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t exp(const Vector< T1, C1 > &l)
OHMMS_PRECISION real_type
std::vector< real_type > B
variables that add detail through an expansion in sigmoidal functions
real_type R0
the soft cutoff distance that controls how short ranged the functor is

◆ evaluate() [2/3]

real_type evaluate ( real_type  r,
real_type dudr,
real_type d2udr2 
) const
inline

compute U(r), dU/dr, and d^2U/dr^2 at a particular value of r

Definition at line 144 of file ShortRangeCuspFunctor.h.

References ShortRangeCuspFunctor< T >::A, ShortRangeCuspFunctor< T >::B, OptimizableFunctorBase::cutoff_radius, qmcplusplus::exp(), qmcplusplus::n, ShortRangeCuspFunctor< T >::R0, and qmcplusplus::Units::time::s.

145  {
146  // the functor will be almost exactly zero at long range, so just return that if r is beyond the hard cutoff
147  if (r >= cutoff_radius)
148  {
149  dudr = 0.0;
150  d2udr2 = 0.0;
151  return 0.0;
152  }
153 
154  // get the ratio of the distance and the soft cutoff distance
155  const real_type s = r / R0;
156 
157  // get the exponential factor
158  const real_type ex = std::exp(-s);
159 
160  // sum the terms related to the sigmoidal functions that are needed for U, dU, and d^2U
161  real_type sig_sum_0 = 0.0; // contributes to U(r)
162  real_type sig_sum_1 = 0.0; // contributes to dU/dr
163  real_type sig_sum_2 = 0.0; // contributes to d2U/dr2
164  real_type n = 2.0;
165  real_type sn = s * s; // s^{n }
166  real_type snm1 = s; // s^{n-1}
167  real_type snm2 = 1.0; // s^{n-2}
168  for (int i = 0; i < B.size(); i++)
169  {
170  const real_type d = 1.0 + sn;
171  const real_type d2 = d * d;
172  const real_type soverd = s / d;
173  const real_type noverd2 = n / d2;
174  sig_sum_0 += B[i] * sn / d;
175  sig_sum_1 += B[i] * (sn * sn + sn - n * snm1) / d2;
176  sig_sum_2 += B[i] * snm2 * (d * noverd2 * noverd2 * (sn - 1.0) + noverd2 * (1.0 + 2.0 * s) - soverd * s);
177  snm2 = snm1; // update s^{n-2}
178  snm1 = sn; // update s^{n-1}
179  sn *= s; // update s^{n }
180  n += 1.0; // update n
181  }
182 
183  // get some useful ratios ( ex / R0 and ex / R0^2 )
184  const real_type exoverR0 = ex / R0;
185  const real_type exoverR02 = exoverR0 / R0;
186 
187  // evaluate dU/dr
188  dudr = A * ex + exoverR0 * sig_sum_1;
189 
190  // evaluate d^2U/dr^2
191  d2udr2 = -A * exoverR0 + exoverR02 * sig_sum_2;
192 
193  // return U(r)
194  return -1.0 * ex * (A * R0 + sig_sum_0);
195  }
real_type A
variable that controls the cusp
MakeReturn< UnaryNode< FnExp, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t exp(const Vector< T1, C1 > &l)
OHMMS_PRECISION real_type
std::vector< real_type > B
variables that add detail through an expansion in sigmoidal functions
real_type R0
the soft cutoff distance that controls how short ranged the functor is

◆ evaluate() [3/3]

real_type evaluate ( real_type  r,
real_type dudr,
real_type d2udr2,
real_type d3udr3 
) const
inline

compute U(r), dU/dr, d^2U/dr^2, and d^3U/dr^3

Definition at line 198 of file ShortRangeCuspFunctor.h.

References ReportEngine::error().

199  {
200  ReportEngine PRE("ShortRangeCuspFunctor", "evaluate(r, dudr, d2udr2, d3udr3)");
201  PRE.error("evaluate(r, dudr, d2udr2, d3udr3) not implemented for ShortRangeCuspFunctor", true);
202  return 0.0;
203  }

◆ evaluateDerivatives() [1/2]

bool evaluateDerivatives ( real_type  r,
std::vector< TinyVector< real_type, 3 >> &  derivs 
)
inlineoverride

compute derivatives of U(r), dU/dr, and d^2U/dr^2 with respect to the variational parameters

Definition at line 285 of file ShortRangeCuspFunctor.h.

References ShortRangeCuspFunctor< T >::A, ShortRangeCuspFunctor< T >::B, qmcplusplus::exp(), qmcplusplus::n, ShortRangeCuspFunctor< T >::Opt_A, ShortRangeCuspFunctor< T >::Opt_B, ShortRangeCuspFunctor< T >::Opt_R0, ShortRangeCuspFunctor< T >::R0, and qmcplusplus::Units::time::s.

Referenced by qmcplusplus::TEST_CASE().

286  {
287  // get the ratio of the distance and the soft cutoff distance
288  const real_type s = r / R0;
289 
290  // get the exponential factor
291  const real_type ex = std::exp(-s);
292 
293  // get some useful ratios
294  const real_type exoverR0 = ex / R0;
295  const real_type exoverR02 = exoverR0 / R0;
296  const real_type exoverR03 = exoverR02 / R0;
297 
298  // initialize the index that will track where to put different parameters' derivatives
299  int i = 0;
300 
301  // evaluate derivatives with respect to the cusp condition variable A
302  if (Opt_A)
303  {
304  derivs[i][0] = -ex * R0; //dU/da
305  derivs[i][1] = ex; //d(dU/da)/dr
306  derivs[i][2] = -exoverR0; //d^2 (dU/da)/dr^2
307  ++i;
308  }
309 
310  // evaluate derivatives with respect to the soft cutoff radius
311  if (Opt_R0)
312  {
313  // sum up terms related to the sigmoidal exansion
314  real_type n = 2.0;
315  real_type sn = s * s; // s^{n }
316  real_type snm1 = s; // s^{n-1}
317  real_type snm2 = 1.0; // s^{n-2}
318  real_type sig_sum_0 = 0.0; // contributes to dUdR0
319  real_type sig_sum_1 = 0.0; // contributes to d(dU/dR0)/dr
320  real_type sig_sum_2 = 0.0; // contributes to d^2(dU/dR0)/dr^2
321  for (int j = 0; j < B.size(); j++)
322  {
323  const real_type d = 1.0 + sn;
324  const real_type d2 = d * d;
325  const real_type soverd = s / d;
326  const real_type noverd2 = n / d2;
327  sig_sum_0 += B[j] * sn * (noverd2 - soverd);
328  sig_sum_1 += B[j] * snm1 * (d * noverd2 * noverd2 * (1.0 - sn) - 2.0 * noverd2 * s + soverd * (s - 1.0));
329  sig_sum_2 += B[j] * snm2 *
330  (s * (3.0 * s - 1.0) * noverd2 - s * (s - 2.0) * soverd +
331  (n + n * sn * (sn - 4.0) + (3.0 * s + 1.0) * (sn * sn - 1.0)) * noverd2 * noverd2);
332  snm2 = snm1; // update s^{n-2}
333  snm1 = sn; // update s^{n-1}
334  sn *= s; // update s^{n }
335  n += 1.0; // update n
336  }
337 
338  // compute terms related to the cusp-inducing term
339  const real_type cusp_part_0 = -A * ex * (s + 1.0); // contributes to dUdR0
340  const real_type cusp_part_1 = A * exoverR0 * s; // contributes to d(dU/dR0)/dr
341  const real_type cusp_part_2 = -A * exoverR02 * (s - 1.0); // contributes to d^2(dU/dR0)/dr^2
342 
343  // compute the desired derivatives
344  derivs[i][0] = cusp_part_0 + exoverR0 * sig_sum_0; // dUdR0
345  derivs[i][1] = cusp_part_1 + exoverR02 * sig_sum_1; // d(dU/dR0)/dr
346  derivs[i][2] = cusp_part_2 + exoverR03 * sig_sum_2; // d^2(dU/dR0)/dr^2
347 
348  // increment the index tracking where to put derivatives
349  ++i;
350  }
351 
352  // evaluate derivatives with respect to the sigmoidal expansion coefficients
353  if (Opt_B)
354  {
355  // loop over the terms in the expansion over sigmoidal functions
356  real_type n = 2.0;
357  real_type sn = s * s; // s^{n }
358  real_type snm1 = s; // s^{n-1}
359  real_type snm2 = 1.0; // s^{n-2}
360  for (int j = 0; j < B.size(); j++)
361  {
362  // compute some useful values
363  const real_type d = 1.0 + sn;
364  const real_type d2 = d * d;
365  const real_type soverd = s / d;
366  const real_type noverd2 = n / d2;
367 
368  // dU/dB_j
369  derivs[i][0] = -ex * sn / d;
370 
371  // d(dU/dB_j)/dr
372  derivs[i][1] = exoverR0 * (sn * sn + sn - n * snm1) / d2;
373 
374  // d^2(dU/dB_j)/dr^2
375  derivs[i][2] = exoverR02 * snm2 * (d * noverd2 * noverd2 * (sn - 1.0) + noverd2 * (1.0 + 2.0 * s) - soverd * s);
376 
377  snm2 = snm1; // update s^{n-2}
378  snm1 = sn; // update s^{n-1}
379  sn *= s; // update s^{n }
380  n += 1.0; // update n
381 
382  // increment the index tracking where to put derivatives
383  ++i;
384  }
385  }
386 
387  return true;
388  }
bool Opt_R0
true, if R0 is optimizable
real_type A
variable that controls the cusp
MakeReturn< UnaryNode< FnExp, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t exp(const Vector< T1, C1 > &l)
bool Opt_A
true, if A is optimizable
bool Opt_B
true, if B variables are optimizable
OHMMS_PRECISION real_type
std::vector< real_type > B
variables that add detail through an expansion in sigmoidal functions
real_type R0
the soft cutoff distance that controls how short ranged the functor is

◆ evaluateDerivatives() [2/2]

bool evaluateDerivatives ( real_type  r,
std::vector< real_type > &  derivs 
)
inlineoverride

compute derivatives of U(r) with respect to the variational parameters

Definition at line 391 of file ShortRangeCuspFunctor.h.

References ShortRangeCuspFunctor< T >::A, ShortRangeCuspFunctor< T >::B, qmcplusplus::exp(), qmcplusplus::n, ShortRangeCuspFunctor< T >::Opt_A, ShortRangeCuspFunctor< T >::Opt_B, ShortRangeCuspFunctor< T >::Opt_R0, ShortRangeCuspFunctor< T >::R0, and qmcplusplus::Units::time::s.

392  {
393  // get the ratio of the distance and the soft cutoff distance
394  const real_type s = r / R0;
395 
396  // get the exponential factor
397  const real_type ex = std::exp(-s);
398 
399  // get some useful ratios
400  const real_type exoverR0 = ex / R0;
401 
402  // initialize the index that will track where to put different parameters' derivatives
403  int i = 0;
404 
405  // evaluate derivatives with respect to the cusp condition variable A
406  if (Opt_A)
407  {
408  derivs[i] = -ex * R0; //dU/da
409  ++i;
410  }
411 
412  // evaluate derivatives with respect to the soft cutoff radius
413  if (Opt_R0)
414  {
415  // sum up terms related to the sigmoidal exansion
416  real_type n = 2.0;
417  real_type sn = s * s; // s^n
418  real_type sig_sum = 0.0;
419  for (int j = 0; j < B.size(); j++)
420  {
421  const real_type d = 1.0 + sn;
422  const real_type d2 = d * d;
423  const real_type soverd = s / d;
424  const real_type noverd2 = n / d2;
425  sig_sum += B[j] * sn * (noverd2 - soverd);
426  sn *= s; // update s^n
427  n += 1.0; // update n
428  }
429 
430  // compute term related to the cusp-inducing term
431  const real_type cusp_part = -A * ex * (s + 1.0);
432 
433  // compute dU/dR0
434  derivs[i] = cusp_part + exoverR0 * sig_sum;
435 
436  // increment the index tracking where to put derivatives
437  ++i;
438  }
439 
440  // evaluate derivatives with respect to the sigmoidal expansion coefficients
441  if (Opt_B)
442  {
443  real_type sn = s * s; // s^n
444  for (int j = 0; j < B.size(); j++)
445  {
446  derivs[i] = -ex * sn / (1.0 + sn); // dU/dB_j
447  sn *= s; // update s^n
448  ++i; // increment the index tracking where to put derivatives
449  }
450  }
451 
452  return true;
453  }
bool Opt_R0
true, if R0 is optimizable
real_type A
variable that controls the cusp
MakeReturn< UnaryNode< FnExp, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t exp(const Vector< T1, C1 > &l)
bool Opt_A
true, if A is optimizable
bool Opt_B
true, if B variables are optimizable
OHMMS_PRECISION real_type
std::vector< real_type > B
variables that add detail through an expansion in sigmoidal functions
real_type R0
the soft cutoff distance that controls how short ranged the functor is

◆ evaluateV()

real_type evaluateV ( const int  iat,
const int  iStart,
const int  iEnd,
const T *restrict  _distArray,
T *restrict  distArrayCompressed 
) const
inline

compute U(r) at multiple distances and sum the results

Definition at line 206 of file ShortRangeCuspFunctor.h.

References ShortRangeCuspFunctor< T >::evaluate().

211  {
212  real_type sum(0);
213  for (int idx = iStart; idx < iEnd; idx++)
214  if (idx != iat)
215  sum += evaluate(_distArray[idx]);
216  return sum;
217  }
real_type evaluate(real_type r) const
compute U(r) at a particular value of r
OHMMS_PRECISION real_type

◆ evaluateVGL()

void evaluateVGL ( const int  iat,
const int  iStart,
const int  iEnd,
const T *  distArray,
T *restrict  valArray,
T *restrict  gradArray,
T *restrict  laplArray,
T *restrict  distArrayCompressed,
int *restrict  distIndices 
) const
inline

compute U(r), dU/dr, and d^2U/dr^2 at multiple values of r

Definition at line 247 of file ShortRangeCuspFunctor.h.

References ShortRangeCuspFunctor< T >::evaluate().

256  {
257  for (int idx = iStart; idx < iEnd; idx++)
258  {
259  valArray[idx] = evaluate(distArray[idx], gradArray[idx], laplArray[idx]);
260  gradArray[idx] /= distArray[idx];
261  }
262  if (iat >= iStart && iat < iEnd)
263  valArray[iat] = gradArray[iat] = laplArray[iat] = T(0);
264  }
real_type evaluate(real_type r) const
compute U(r) at a particular value of r

◆ f()

real_type f ( real_type  r)
inlineoverride

compute U(r) at a particular distance or return zero if beyond the cutoff

Definition at line 267 of file ShortRangeCuspFunctor.h.

References OptimizableFunctorBase::cutoff_radius, and ShortRangeCuspFunctor< T >::evaluate().

268  {
269  if (r >= cutoff_radius)
270  return 0.0;
271  return evaluate(r);
272  }
real_type evaluate(real_type r) const
compute U(r) at a particular value of r

◆ makeClone()

OptimizableFunctorBase* makeClone ( ) const
inlineoverridevirtual

clone the functor

Implements OptimizableFunctorBase.

Definition at line 112 of file ShortRangeCuspFunctor.h.

References ShortRangeCuspFunctor< T >::ShortRangeCuspFunctor().

112 { return new ShortRangeCuspFunctor(*this); }
ShortRangeCuspFunctor(const std::string &my_name)
default constructor

◆ mw_evaluateV()

static void mw_evaluateV ( const int  num_groups,
const ShortRangeCuspFunctor< T > *const  functors[],
const int  n_src,
const int *  grp_ids,
const int  num_pairs,
const int *  ref_at,
const T *  mw_dist,
const int  dist_stride,
T *  mw_vals,
Vector< char, OffloadPinnedAllocator< char >> &  transfer_buffer 
)
inlinestatic

evaluate sum of the pair potentials FIXME

Returns
$\sum u(r_j)$ for r_j < cutoff_radius

Definition at line 222 of file ShortRangeCuspFunctor.h.

232  {
233  for (int ip = 0; ip < num_pairs; ip++)
234  {
235  mw_vals[ip] = 0;
236  for (int j = 0; j < n_src; j++)
237  {
238  const int ig = grp_ids[j];
239  auto& functor(*functors[ig]);
240  if (j != ref_at[ip])
241  mw_vals[ip] += functor.evaluate(mw_dist[ip * dist_stride + j]);
242  }
243  }
244  }

◆ put()

bool put ( xmlNodePtr  cur)
inlineoverridevirtual

read in information about the functor from an xml node

Implements OptimizableFunctorBase.

Definition at line 500 of file ShortRangeCuspFunctor.h.

References ShortRangeCuspFunctor< T >::A, OhmmsAttributeSet::add(), qmcplusplus::app_summary(), ShortRangeCuspFunctor< T >::B, castXMLCharToChar(), OptimizableFunctorBase::cutoff_radius, ReportEngine::error(), ShortRangeCuspFunctor< T >::ID_A, ShortRangeCuspFunctor< T >::ID_B, ShortRangeCuspFunctor< T >::ID_R0, VariableSet::insert(), qmcplusplus::lowerCase(), OptimizableFunctorBase::myVars, ShortRangeCuspFunctor< T >::Opt_A, ShortRangeCuspFunctor< T >::Opt_B, ShortRangeCuspFunctor< T >::Opt_R0, optimize::OTHER_P, OhmmsAttributeSet::put(), putContent(), ShortRangeCuspFunctor< T >::R0, ShortRangeCuspFunctor< T >::reset(), ShortRangeCuspFunctor< T >::set_variable_from_xml(), and ShortRangeCuspFunctor< T >::setCusp().

Referenced by qmcplusplus::TEST_CASE().

501  {
502  // set up an object for info / warning / error reporting
503  ReportEngine PRE("ShortRangeCuspFunctor", "put(xmlNodePtr)");
504 
505  // create some variables to read information in to
506  real_type radius = -1.0; // hard cutoff radius
507  real_type cusp_in = -1.0; // cusp condition
508 
509  // read from the xml node
510  OhmmsAttributeSet rAttrib;
511  rAttrib.add(radius, "rcut");
512  rAttrib.add(radius, "cutoff");
513  rAttrib.add(cusp_in, "cusp");
514  rAttrib.put(cur);
515 
516  // set the hard cutoff radius if we have it
517  if (radius > 0.0)
518  cutoff_radius = radius;
519 
520  // set the cusp if we have it
521  if (cusp_in > 0.0)
522  setCusp(cusp_in);
523 
524  // loop over this node's children to read in variational parameter information
525  xmlNodePtr childPtr = cur->xmlChildrenNode;
526  while (childPtr != NULL)
527  {
528  // skip blank nodes
529  if (xmlIsBlankNode(childPtr))
530  {
531  childPtr = childPtr->next;
532  continue;
533  }
534 
535  // get the name of the child node
536  std::string cname(lowerCase(castXMLCharToChar(childPtr->name)));
537  // read in a variable
538  if (cname == "var")
539  {
540  // read in the name of the variable
541  std::string v_name("string_not_set");
542  OhmmsAttributeSet att;
543  att.add(v_name, "name");
544  att.put(childPtr);
545 
546  // read in the variable's info
547  if (v_name == "A")
548  set_variable_from_xml(PRE, childPtr, A, ID_A, Opt_A);
549  else if (v_name == "R0")
550  set_variable_from_xml(PRE, childPtr, R0, ID_R0, Opt_R0);
551  else if (v_name == "string_not_set")
552  PRE.error("variable name not set", true);
553  else
554  PRE.error("unrecognized variable name: " + v_name, true);
555  }
556 
557  // read in the B coefficients
558  else if (cname == "coefficients")
559  {
560  // read in the id, whether to optimize, and the node type
561  std::string id("string_not_set");
562  std::string type("string_not_set");
563  std::string optimize("string_not_set");
564  OhmmsAttributeSet att;
565  att.add(id, "id");
566  att.add(type, "type");
567  att.add(optimize, "optimize");
568  att.put(childPtr);
569 
570  // sanity check that this node has been specified as an array type
571  if (type != "Array")
572  PRE.error("ShortRangeCuspFunctor expected coefficients parameter array type to be \"Array\"", true);
573 
574  // read in the vector of coefficients
575  putContent(B, childPtr);
576 
577  // set the id if we have it
578  if (id != "string_not_set")
579  ID_B = id;
580 
581  // if the coefficients are to be optimized, add them to the variable set
583  if (optimize == "yes")
584  {
585  if (id == "string_not_set")
586  PRE.error("\"id\" must be set if we are going to optimize B coefficients", true);
587  Opt_B = true;
588  for (int i = 0; i < B.size(); i++)
589  {
590  std::stringstream sstr;
591  sstr << ID_B << "_" << i;
592  myVars.insert(sstr.str(), B.at(i), Opt_B, optimize::OTHER_P);
593  }
594  }
595  else if (optimize == "no")
596  {
597  Opt_B = false;
598  }
599  else if (optimize != "string_not_set")
600  {
601  PRE.error("Unrecognized value for \"optimize\". Should be either yes or no", true);
602  }
603  }
604 
605  // error if cname is not recognized
606  else
607  {
608  PRE.error("\"" + cname +
609  "\" is not a recognized value for \"cname\". Allowed values are \"var\" and \"coefficients\"",
610  true);
611  }
612 
613  // go to the next node
614  childPtr = childPtr->next;
615  }
616 
617  // summarize what was read in
618  app_summary() << " ---------------------------------------------------------" << std::endl;
619  app_summary() << " cusp variable A: " << A << std::endl;
620  app_summary() << " soft cutoff R0: " << R0 << std::endl;
621  app_summary() << " sigmoidal expansion coefficients B:";
622  for (int i = 0; i < B.size(); i++)
623  app_summary() << " " << B.at(i);
624  app_summary() << std::endl;
625  app_summary() << " hard cutoff radius: " << cutoff_radius << std::endl;
626  app_summary() << " Opt_A: " << (Opt_A ? "yes" : "no") << std::endl;
627  app_summary() << " Opt_R0: " << (Opt_R0 ? "yes" : "no") << std::endl;
628  app_summary() << " Opt_B: " << (Opt_B ? "yes" : "no") << std::endl;
629 
630  // reset the functor now that we've read its info in (currently does nothing)
631  reset();
632 
633  return true;
634  }
void set_variable_from_xml(ReportEngine &PRE, xmlNodePtr cur, U &variable_to_set, std::string &id_to_set, bool &opt_to_set)
set up a variational parameter using the supplied xml node
bool Opt_R0
true, if R0 is optimizable
std::ostream & app_summary()
Definition: OutputManager.h:63
bool put(xmlNodePtr cur)
assign attributes to the set
Definition: AttributeSet.h:55
void reset() override
Implement the reset function, which was pure virtual in OptimizableFunctorBase, even though we don&#39;t ...
real_type A
variable that controls the cusp
void insert(const std::string &vname, real_type v, bool enable=true, int type=OTHER_P)
Definition: VariableSet.h:133
class to handle a set of attributes of an xmlNode
Definition: AttributeSet.h:24
std::string lowerCase(const std::string_view s)
++17
char * castXMLCharToChar(xmlChar *c)
assign a value from a node. Use specialization for classes.
Definition: libxmldefs.h:62
bool Opt_A
true, if A is optimizable
bool Opt_B
true, if B variables are optimizable
bool putContent(T &a, xmlNodePtr cur)
replaces a&#39;s value with the first "element" in the "string" returned by XMLNodeString{cur}.
Definition: libxmldefs.h:88
void setCusp(real_type cusp) override
sets the cusp condition and disables optimization of the cusp-determining parameter ...
OHMMS_PRECISION real_type
std::vector< real_type > B
variables that add detail through an expansion in sigmoidal functions
void add(PDT &aparam, const std::string &aname, std::vector< PDT > candidate_values={}, TagStatus status=TagStatus::OPTIONAL)
add a new attribute
Definition: AttributeSet.h:42
opt_variables_type myVars
set of variables to be optimized
real_type R0
the soft cutoff distance that controls how short ranged the functor is

◆ reset()

void reset ( )
inlineoverridevirtual

Implement the reset function, which was pure virtual in OptimizableFunctorBase, even though we don't need it.

Implements OptimizableFunctorBase.

Definition at line 115 of file ShortRangeCuspFunctor.h.

Referenced by ShortRangeCuspFunctor< T >::put(), ShortRangeCuspFunctor< T >::resetParametersExclusive(), ShortRangeCuspFunctor< T >::setCusp(), and ShortRangeCuspFunctor< T >::ShortRangeCuspFunctor().

116  {
117  //cutoff_radius = 1.0e4; //some big range
118  }

◆ resetParametersExclusive()

void resetParametersExclusive ( const opt_variables_type active)
inlineoverridevirtual

reset the parameters during optimizations.

Exclusive, see checkInVariablesExclusive

Implements OptimizableObject.

Definition at line 646 of file ShortRangeCuspFunctor.h.

References ShortRangeCuspFunctor< T >::A, ShortRangeCuspFunctor< T >::B, OptimizableFunctorBase::myVars, ShortRangeCuspFunctor< T >::Opt_A, ShortRangeCuspFunctor< T >::Opt_B, ShortRangeCuspFunctor< T >::Opt_R0, ShortRangeCuspFunctor< T >::R0, ShortRangeCuspFunctor< T >::reset(), VariableSet::size(), and VariableSet::where().

Referenced by qmcplusplus::TEST_CASE().

647  {
648  if (myVars.size())
649  {
650  int ia = myVars.where(0);
651  if (ia > -1)
652  {
653  int i = 0;
654  if (Opt_A)
655  A = std::real(myVars[i++] = active[ia++]);
656  if (Opt_R0)
657  R0 = std::real(myVars[i++] = active[ia++]);
658  for (int j = 0; Opt_B && j < B.size(); j++)
659  B.at(j) = std::real(myVars[i++] = active[ia++]);
660  }
661  reset();
662  }
663  }
QMCTraits::RealType real
bool Opt_R0
true, if R0 is optimizable
void reset() override
Implement the reset function, which was pure virtual in OptimizableFunctorBase, even though we don&#39;t ...
real_type A
variable that controls the cusp
int where(int i) const
return the locator of the i-th Index
Definition: VariableSet.h:90
bool Opt_A
true, if A is optimizable
size_type size() const
return the size
Definition: VariableSet.h:88
bool Opt_B
true, if B variables are optimizable
std::vector< real_type > B
variables that add detail through an expansion in sigmoidal functions
opt_variables_type myVars
set of variables to be optimized
real_type R0
the soft cutoff distance that controls how short ranged the functor is

◆ set_variable_from_xml()

void set_variable_from_xml ( ReportEngine PRE,
xmlNodePtr  cur,
U &  variable_to_set,
std::string &  id_to_set,
bool &  opt_to_set 
)
inline

set up a variational parameter using the supplied xml node

Definition at line 457 of file ShortRangeCuspFunctor.h.

References OhmmsAttributeSet::add(), ReportEngine::error(), VariableSet::insert(), qmcplusplus::lowerCase(), OptimizableFunctorBase::myVars, optimize::OTHER_P, OhmmsAttributeSet::put(), and putContent().

Referenced by ShortRangeCuspFunctor< T >::put().

462  {
463  // get id, name, and optimize info
464  std::string id("string_not_set");
465  std::string name("string_not_set");
466  std::string optimize("string_not_set");
467  OhmmsAttributeSet rAttrib;
468  rAttrib.add(id, "id");
469  rAttrib.add(name, "name");
470  rAttrib.add(optimize, "optimize");
471  rAttrib.put(cur);
472 
473  // read in the variable
474  putContent(variable_to_set, cur);
475 
476  // set the id if we have it
477  if (id != "string_not_set")
478  id_to_set = id;
479 
480  // if we are to optimize the variable, add it to the optimizable variable set
482  if (optimize == "yes")
483  {
484  if (id == "string_not_set")
485  PRE.error("\"id\" must be set if we are going to optimize variable " + name, true);
486  opt_to_set = true;
487  myVars.insert(id, variable_to_set, opt_to_set, optimize::OTHER_P);
488  }
489  else if (optimize == "no")
490  {
491  opt_to_set = false;
492  }
493  else if (optimize != "string_not_set")
494  {
495  PRE.error("Unrecognized value for \"optimize\". Should be either yes or no", true);
496  }
497  }
bool put(xmlNodePtr cur)
assign attributes to the set
Definition: AttributeSet.h:55
void insert(const std::string &vname, real_type v, bool enable=true, int type=OTHER_P)
Definition: VariableSet.h:133
class to handle a set of attributes of an xmlNode
Definition: AttributeSet.h:24
std::string lowerCase(const std::string_view s)
++17
bool putContent(T &a, xmlNodePtr cur)
replaces a&#39;s value with the first "element" in the "string" returned by XMLNodeString{cur}.
Definition: libxmldefs.h:88
void add(PDT &aparam, const std::string &aname, std::vector< PDT > candidate_values={}, TagStatus status=TagStatus::OPTIONAL)
add a new attribute
Definition: AttributeSet.h:42
opt_variables_type myVars
set of variables to be optimized

◆ setCusp()

void setCusp ( real_type  cusp)
inlineoverride

sets the cusp condition and disables optimization of the cusp-determining parameter

Definition at line 103 of file ShortRangeCuspFunctor.h.

References ShortRangeCuspFunctor< T >::A, ShortRangeCuspFunctor< T >::Opt_A, and ShortRangeCuspFunctor< T >::reset().

Referenced by ShortRangeCuspFunctor< T >::put().

104  {
105  //throw std::runtime_error("ShortRangeCuspFunctor::setCusp was called");
106  A = cusp;
107  Opt_A = false;
108  reset();
109  }
void reset() override
Implement the reset function, which was pure virtual in OptimizableFunctorBase, even though we don&#39;t ...
real_type A
variable that controls the cusp
bool Opt_A
true, if A is optimizable

Member Data Documentation

◆ A

◆ B

◆ ID_A

std::string ID_A

id of A

Definition at line 76 of file ShortRangeCuspFunctor.h.

Referenced by ShortRangeCuspFunctor< T >::put(), and qmcplusplus::TEST_CASE().

◆ ID_B

std::string ID_B

id of B

Definition at line 80 of file ShortRangeCuspFunctor.h.

Referenced by ShortRangeCuspFunctor< T >::put(), and qmcplusplus::TEST_CASE().

◆ ID_R0

std::string ID_R0

id of R0

Definition at line 78 of file ShortRangeCuspFunctor.h.

Referenced by ShortRangeCuspFunctor< T >::put(), and qmcplusplus::TEST_CASE().

◆ Opt_A

◆ Opt_B

◆ Opt_R0

◆ R0


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