QMCPACK
OneBodyDensityMatrices Class Reference

Per crowd Estimator for OneBodyDensityMatrices aka 1RDM DensityMatrices1B. More...

+ Inheritance diagram for OneBodyDensityMatrices:
+ Collaboration diagram for OneBodyDensityMatrices:

Classes

struct  OneBodyDensityMatrixTimers
 

Public Types

enum  Sampling { VOLUME_BASED, METROPOLIS, NO_SAMPLING }
 
using Value = QMCTraits::ValueType
 
using FullPrecValue = QMCTraits::FullPrecValueType
 
using Real = RealAlias< Value >
 
using FullPrecReal = RealAlias< FullPrecValue >
 
using Grad = TinyVector< Value, OHMMS_DIM >
 
using Lattice = PtclOnLatticeTraits::ParticleLayout
 
using Position = QMCTraits::PosType
 
using Evaluator = OneBodyDensityMatricesInput::Evaluator
 
using Integrator = OneBodyDensityMatricesInput::Integrator
 
using SPOMap = SPOSet::SPOMap
 
- Public Types inherited from OperatorEstBase
using QMCT = QMCTraits
 
using FullPrecRealType = QMCT::FullPrecRealType
 
using MCPWalker = Walker< QMCTraits, PtclOnLatticeTraits >
 
using Data = std::vector< QMCT::RealType >
 

Private Attributes

OneBodyDensityMatricesInput input_
 
Lattice lattice_
 
SpeciesSet species_
 
CompositeSPOSet basis_functions_
 
Vector< Valuebasis_values_
 
Vector< Valuebasis_norms_
 
Vector< Gradbasis_gradients_
 
Vector< Valuebasis_laplacians_
 
Vector< Valuebasis_spin_gradients_
 
std::vector< Positionrsamples_
 
std::vector< Realssamples_
 
Vector< Realsamples_weights_
 
int basis_size_
 
std::vector< int > species_sizes_
 
std::vector< std::string > species_names_
 
size_t samples_
 

More...
 
Sampling sampling_
 Sampling method, this derived values in input_. More...
 
Position center_
 If not defined in OBDMI taken from lattice_. More...
 
Position rcorner_
 with respect to center_ using lattice_; More...
 
Real volume_
 
bool periodic_
 
std::vector< Valuepsi_ratios_
 , ' . More...
 
std::vector< Matrix< Value > > Phi_NB_
 row major per sample workspaces More...
 
std::vector< Matrix< Value > > Psi_NM_
 ratio per particle per sample size: particles * samples vector is over species each matrix row: particle col: sample More...
 
std::vector< Matrix< Value > > Phi_Psi_NB_
 
std::vector< Matrix< Value > > N_BB_
 
Matrix< ValuePhi_MB_
 basis_values_ at each r of rsamples_ row: sample col: basis_value size: samples * basis_size More...
 
int nmoves_ = 0
 

More...
 
int naccepted_ = 0
 number of accepted samples More...
 
Real acceptance_ratio_ = 0.0
 running acceptance ratio over all samples More...
 
int ind_dims_ [OHMMS_DIM]
 
bool warmed_up_ = false
 
Real metric_ = 1.0
 }@ More...
 
Position rpcur_
 current position – As long Positions are TinyVectors they are intialized to zero vectors More...
 
Position dpcur_
 current drift More...
 
Real rhocur_ = 0.0
 current density More...
 
Real spcur_ = 0.0
 spin related variables More...
 
Real dspcur_ = 0.0
 
const bool is_spinor_
 
OneBodyDensityMatrixTimers timers_
 
template<typename T >
class testing::OneBodyDensityMatricesTests
 
 OneBodyDensityMatrices (OneBodyDensityMatricesInput &&obdmi, const Lattice &lattice, const SpeciesSet &species, const SPOMap &spomap, const ParticleSet &pset_target)
 Standard Constructor Call this to make a new OBDM this is what you should be calling. More...
 
 OneBodyDensityMatrices (const OneBodyDensityMatrices &obdm, DataLocality dl)
 Constructor used when spawing crowd clones needs to be public so std::make_unique can call it. More...
 
std::unique_ptr< OperatorEstBasespawnCrowdClone () const override
 
void accumulate (const RefVector< MCPWalker > &walkers, const RefVector< ParticleSet > &psets, const RefVector< TrialWaveFunction > &wfns, const RefVector< QMCHamiltonian > &hams, RandomBase< FullPrecReal > &rng) override
 
void startBlock (int steps) override
 
void registerOperatorEstimator (hdf_archive &file) override
 create and tie OperatorEstimator's observable_helper hdf5 wrapper to stat.h5 file More...
 
 OneBodyDensityMatrices (const OneBodyDensityMatrices &obdm)=default
 Default copy constructor. More...
 
void implAccumulate (const RefVector< MCPWalker > &walkers, const RefVector< ParticleSet > &psets, const RefVector< TrialWaveFunction > &wfns, RandomBase< FullPrecReal > &rng)
 Unfortunate design RandomGenerator type aliasing and virtual inheritance requires this for testing. More...
 
size_t calcFullDataSize (size_t basis_size, int num_species)
 
void normalizeBasis (ParticleSet &pset_target)
 
void report (const std::string &pad="")
 
void evaluateMatrix (ParticleSet &pset_target, TrialWaveFunction &psi_target, const MCPWalker &walker, RandomBase< FullPrecReal > &rng)
 
void generateSamples (const Real weight, ParticleSet &pset_target, RandomBase< FullPrecReal > &rng, int steps=0)
 Dispatch method to difference methods of generating samples. More...
 
void generateUniformGrid (RandomBase< FullPrecReal > &rng)
 
void generateUniformSamples (RandomBase< FullPrecReal > &rng)
 
void generateDensitySamples (bool save, int steps, RandomBase< FullPrecReal > &rng, ParticleSet &pset_target)
 generate samples for density integration More...
 
void generateDensitySamplesWithSpin (bool save, int steps, RandomBase< FullPrecReal > &rng, ParticleSet &pset_target)
 same as above, but with spin variables included More...
 
void generateSampleRatios (ParticleSet &pset_target, TrialWaveFunction &psi_target, std::vector< Matrix< Value >> &Psi_nm)
 
Position diffuse (const Real sqt, RandomBase< FullPrecReal > &rng)
 produce a position difference vector from timestep More...
 
Real diffuseSpin (const Real sqt, RandomBase< FullPrecReal > &rng)
 spin diffusion More...
 
void calcDensity (const Position &r, Real &dens, ParticleSet &pset_target)
 calculate density based on r More...
 
void calcDensityWithSpin (const Position &r, const Real &s, Real &dens, ParticleSet &pset_target)
 same as above, but with spin move More...
 
void calcDensityDrift (const Position &r, Real &dens, Position &drift, ParticleSet &pset_target)
 calculate density and drift bashed on r More...
 
void calcDensityDriftWithSpin (const Position &r, const Real &s, Real &dens, Position &drift, Real &sdrift, ParticleSet &pset_target)
 same as above, but with spin move More...
 
void generateSampleBasis (Matrix< Value > &Phi_mb, ParticleSet &pset_target, TrialWaveFunction &psi_target)
 set Phi_mp to basis vaules per sample sideeffects: More...
 
void generateParticleBasis (ParticleSet &pset_target, std::vector< Matrix< Value >> &phi_nb)
 set phi_nb to basis values per target particleset particle sideeffects: More...
 
void updateBasis (const Position &r, ParticleSet &pset_target)
 basis set updates More...
 
void updateBasisWithSpin (const Position &r, const Real &s, ParticleSet &pset_target)
 basis set updates with spin More...
 
void updateBasisD012 (const Position &r, ParticleSet &pset_target)
 evaluates vgl on basis_functions_ for r sideeffects: More...
 
void updateBasisD012WithSpin (const Position &r, const Real &s, ParticleSet &pset_target)
 same as above, but includes spin gradients More...
 
void warmupSampling (ParticleSet &pset_target, RandomBase< FullPrecReal > &rng)
 does some warmup sampling i.e. More...
 

Additional Inherited Members

- Public Member Functions inherited from OperatorEstBase
 OperatorEstBase (DataLocality dl)
 constructor More...
 
 OperatorEstBase (const OperatorEstBase &oth)
 Shallow copy constructor! This alows us to keep the default copy constructors for derived classes which is quite useful to the spawnCrowdClone design. More...
 
virtual ~OperatorEstBase ()=default
 virtual destructor More...
 
virtual void accumulate (const RefVector< MCPWalker > &walkers, const RefVector< ParticleSet > &psets, const RefVector< TrialWaveFunction > &wfns, const RefVector< QMCHamiltonian > &hams, RandomBase< FullPrecRealType > &rng)=0
 Accumulate whatever it is you are accumulating with respect to walkers. More...
 
virtual void collect (const RefVector< OperatorEstBase > &oebs)
 Reduce estimator result data from crowds to rank. More...
 
virtual void normalize (QMCT::RealType invToWgt)
 
std::vector< QMCT::RealType > & get_data ()
 
void write (hdf_archive &file)
 Write to previously registered observable_helper hdf5 wrapper. More...
 
void zero ()
 zero data appropriately for the DataLocality More...
 
QMCT::FullPrecRealType get_walkers_weight () const
 Return the total walker weight for this block. More...
 
const std::string & get_my_name () const
 
virtual void registerListeners (QMCHamiltonian &ham_leader)
 Register 0-many listeners with a leading QMCHamiltonian instance i.e. More...
 
bool isListenerRequired ()
 
DataLocality get_data_locality () const
 
- Protected Attributes inherited from OperatorEstBase
DataLocality data_locality_
 locality for accumulation of estimator data. More...
 
std::string my_name_
 name of this object – only used for debugging and h5 output More...
 
QMCT::FullPrecRealType walkers_weight_
 
std::vector< ObservableHelperh5desc_
 
Data data_
 
bool requires_listener_ = false
 

Detailed Description

Per crowd Estimator for OneBodyDensityMatrices aka 1RDM DensityMatrices1B.

Todo:

most matrices are written to by incrementing a single vector index into their memory. This isn't compatible with aligned rows and ignores much less error prone accessing that Matrix supplys. Fix it.

functions favor output arguments or state updates over return values, simplify.

Definition at line 45 of file OneBodyDensityMatrices.h.

Member Typedef Documentation

◆ Evaluator

◆ FullPrecReal

Definition at line 51 of file OneBodyDensityMatrices.h.

◆ FullPrecValue

◆ Grad

Definition at line 52 of file OneBodyDensityMatrices.h.

◆ Integrator

◆ Lattice

◆ Position

Definition at line 54 of file OneBodyDensityMatrices.h.

◆ Real

using Real = RealAlias<Value>

Definition at line 50 of file OneBodyDensityMatrices.h.

◆ SPOMap

Definition at line 58 of file OneBodyDensityMatrices.h.

◆ Value

Definition at line 48 of file OneBodyDensityMatrices.h.

Member Enumeration Documentation

◆ Sampling

enum Sampling
strong
Enumerator
VOLUME_BASED 
METROPOLIS 
NO_SAMPLING 

Definition at line 60 of file OneBodyDensityMatrices.h.

61  {
62  VOLUME_BASED,
63  METROPOLIS,
64  NO_SAMPLING
65  };

Constructor & Destructor Documentation

◆ OneBodyDensityMatrices() [1/3]

OneBodyDensityMatrices ( OneBodyDensityMatricesInput &&  obdmi,
const Lattice lattice,
const SpeciesSet species,
const SPOMap spomap,
const ParticleSet pset_target 
)

Standard Constructor Call this to make a new OBDM this is what you should be calling.

Definition at line 30 of file OneBodyDensityMatrices.cpp.

References CompositeSPOSet::add(), OneBodyDensityMatrices::basis_functions_, OneBodyDensityMatrices::basis_gradients_, OneBodyDensityMatrices::basis_laplacians_, OneBodyDensityMatrices::basis_norms_, OneBodyDensityMatrices::basis_size_, OneBodyDensityMatrices::basis_spin_gradients_, OneBodyDensityMatrices::basis_values_, OneBodyDensityMatrices::calcFullDataSize(), CrystalLattice< T, D >::Center, OneBodyDensityMatrices::center_, qmcplusplus::crowd, OperatorEstBase::data_, OneBodyDensityMatricesInput::DENSITY, qmcplusplus::exp(), OneBodyDensityMatricesInput::get_basis_sets(), OneBodyDensityMatricesInput::get_center(), OneBodyDensityMatricesInput::get_center_defined(), OneBodyDensityMatricesInput::get_corner(), OneBodyDensityMatricesInput::get_corner_defined(), OneBodyDensityMatricesInput::get_evaluator(), OneBodyDensityMatricesInput::get_integrator(), OneBodyDensityMatricesInput::get_normalized(), OneBodyDensityMatricesInput::get_points(), OneBodyDensityMatricesInput::get_samples(), OneBodyDensityMatricesInput::get_scale(), OneBodyDensityMatricesInput::get_volume_normalized(), SpeciesSet::getAttribute(), SpeciesSet::hasAttribute(), OneBodyDensityMatrices::ind_dims_, OneBodyDensityMatrices::input_, OneBodyDensityMatrices::is_spinor_, OneBodyDensityMatrices::lattice_, qmcplusplus::log(), OneBodyDensityMatricesInput::MATRIX, OneBodyDensityMatrices::metric_, OneBodyDensityMatrices::METROPOLIS, OperatorEstBase::my_name_, OneBodyDensityMatrices::N_BB_, OneBodyDensityMatrices::normalizeBasis(), OHMMS_DIM, OneBodyDensityMatrices::periodic_, OneBodyDensityMatrices::Phi_MB_, OneBodyDensityMatrices::Phi_NB_, OneBodyDensityMatrices::Phi_Psi_NB_, qmcplusplus::pow(), qmcplusplus::pset_target, OneBodyDensityMatrices::Psi_NM_, OneBodyDensityMatrices::psi_ratios_, OneBodyDensityMatrices::rcorner_, CrystalLattice< T, D >::reset(), Matrix< T, Alloc >::resize(), Vector< T, Alloc >::resize(), OneBodyDensityMatrices::rsamples_, qmcplusplus::Units::time::s, OneBodyDensityMatrices::samples_, OneBodyDensityMatrices::samples_weights_, OneBodyDensityMatrices::sampling_, SpeciesSet::size(), SPOSet::size(), OneBodyDensityMatrices::species_, OneBodyDensityMatrices::species_names_, OneBodyDensityMatrices::species_sizes_, SpeciesSet::speciesName, qmcplusplus::spomap, qmcplusplus::sqrt(), OneBodyDensityMatrices::ssamples_, qmcplusplus::SUPERCELL_OPEN, CrystalLattice< T, D >::SuperCellEnum, OneBodyDensityMatricesInput::UNIFORM, OneBodyDensityMatricesInput::UNIFORM_GRID, CrystalLattice< T, D >::Volume, OneBodyDensityMatrices::volume_, and OneBodyDensityMatrices::VOLUME_BASED.

36  input_(obdmi),
38  species_(species),
39  basis_functions_("OneBodyDensityMatrices::basis"),
40  is_spinor_(pset_target.isSpinor()),
41  timers_("OneBodyDensityMatrix")
42 {
43  my_name_ = "OneBodyDensityMatrices";
44  lattice_.reset();
45 
47  {
50  }
51  else
52  {
55  else
58  }
59 
62 
63  // Here we discover sampling is derived (this may belong in input class)
64  switch (input_.get_integrator())
65  {
71  for (int d = 1; d < OHMMS_DIM; ++d)
72  ind_dims_[d] = ind_dims_[d - 1] / input_.get_points();
73  break;
78  break;
82  metric_ = 1.0 / samples_;
83  break;
84  }
85  rsamples_.resize(samples_);
86  if (is_spinor_)
87  ssamples_.resize(samples_);
88 
90  throw UniformCommunicateError("OneBodyDensityMatrices::OneBodyDensityMatrices only density sampling implemented "
91  "for calculations using spinors");
92 
93 
94  // get the sposets that form the basis
95  auto& sposets = input_.get_basis_sets();
96 
97  for (int i = 0; i < sposets.size(); ++i)
98  {
99  auto spo_it = spomap.find(sposets[i]);
100  if (spo_it == spomap.end())
101  throw UniformCommunicateError("OneBodyDensityMatrices::OneBodyDensityMatrices sposet " + sposets[i] +
102  " does not exist.");
103  basis_functions_.add(spo_it->second->makeClone());
104  }
106 
107  if (basis_size_ < 1)
108  throw UniformCommunicateError("OneBodyDensityMatrices::OneBodyDensityMatrices basis_size must be greater than one");
109 
110  int nspecies = species.size();
111  if (!species_.hasAttribute("membersize"))
112  throw UniformCommunicateError("OneBodyDensityMatrices::OneBodyDensityMatrices error: Species set does not have the "
113  "required attribute 'membersize'");
114  int isize = species.getAttribute("membersize");
115  // We have the count per species at least a fundamental as the total particles.
116  // the sum of species membersize and total particles should be an invariant.
117  int nparticles = 0;
118  for (int s = 0; s < nspecies; ++s)
119  nparticles += species(isize, s);
120  for (int s = 0; s < nspecies; ++s)
121  species_sizes_.push_back(species(isize, s));
122  for (int s = 0; s < nspecies; ++s)
123  species_names_.push_back(species.speciesName[s]);
124 
127 
128  Real bn_standard = 1.0;
130  bn_standard = 1.0 / std::sqrt(volume_);
131  for (int i = 0; i < basis_size_; ++i)
132  basis_norms_[i] = bn_standard;
133 
135  psi_ratios_.resize(nparticles);
136 
138  {
140  Phi_NB_.reserve(nspecies);
141  Psi_NM_.reserve(nspecies);
142  Phi_Psi_NB_.reserve(nspecies);
143  N_BB_.reserve(nspecies);
144  for (int s = 0; s < nspecies; ++s)
145  {
146  int specs_size = species_sizes_[s];
147  Phi_NB_.emplace_back(specs_size, basis_size_);
148  Psi_NM_.emplace_back(specs_size, samples_);
149  Phi_Psi_NB_.emplace_back(specs_size, basis_size_);
150  N_BB_.emplace_back(basis_size_, basis_size_);
151  }
152  }
153 
155  {
158  if (is_spinor_)
160  }
161 
162  // so if the input is not normalized, normalize it.
163  // with respect to what?
164  if (!input_.get_normalized())
165  {
166  //Since the following is not a const method we copy particle set
167  ParticleSet pset_temp(pset_target);
168  normalizeBasis(pset_temp);
169  }
170 
172 }
void resize(size_type n, Type_t val=Type_t())
Resize the container.
Definition: OhmmsVector.h:166
int SuperCellEnum
supercell enumeration
SingleParticlePos Center
Center of the cell sum(Rv[i])/2.
std::vector< std::string > species_names_
void reset()
Evaluate the reciprocal vectors, volume and metric tensor.
bool hasAttribute(const std::string &aname) const
Check for attribute presence This replaces code that gets numAttributes then tries to addAttribute fo...
Definition: SpeciesSet.cpp:70
Position center_
If not defined in OBDMI taken from lattice_.
int size() const
return the size of the orbital set Ye: this needs to be replaced by getOrbitalSetSize(); ...
Definition: SPOSet.h:75
Position rcorner_
with respect to center_ using lattice_;
void resize(size_type n, size_type m)
Resize the container.
Definition: OhmmsMatrix.h:99
std::vector< Matrix< Value > > Phi_NB_
row major per sample workspaces
size_t calcFullDataSize(size_t basis_size, int num_species)
Matrix< Value > Phi_MB_
basis_values_ at each r of rsamples_ row: sample col: basis_value size: samples * basis_size ...
void normalizeBasis(ParticleSet &pset_target)
int size() const
return the number of species
Definition: SpeciesSet.h:53
#define OHMMS_DIM
Definition: config.h:64
MakeReturn< BinaryNode< FnPow, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t, typename CreateLeaf< Vector< T2, C2 > >::Leaf_t > >::Expression_t pow(const Vector< T1, C1 > &l, const Vector< T2, C2 > &r)
CrystalLattice< OHMMS_PRECISION, OHMMS_DIM > lattice
Scalar_t Volume
Physical properties of a supercell.
OperatorEstBase(DataLocality dl)
constructor
OneBodyDensityMatricesInput obdmi(node)
MakeReturn< UnaryNode< FnExp, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t exp(const Vector< T1, C1 > &l)
void add(std::unique_ptr< SPOSet > component)
add a sposet component to this composite sposet
Sampling sampling_
Sampling method, this derived values in input_.
MakeReturn< UnaryNode< FnLog, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t log(const Vector< T1, C1 > &l)
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
const std::vector< std::string > & get_basis_sets() const
std::vector< Matrix< Value > > Phi_Psi_NB_
std::vector< Matrix< Value > > N_BB_
std::vector< Value > psi_ratios_
, &#39; .
std::string my_name_
name of this object – only used for debugging and h5 output
std::vector< Matrix< Value > > Psi_NM_
ratio per particle per sample size: particles * samples vector is over species each matrix row: parti...

◆ OneBodyDensityMatrices() [2/3]

Constructor used when spawing crowd clones needs to be public so std::make_unique can call it.

Do not use directly unless you've really thought it through.

Definition at line 174 of file OneBodyDensityMatrices.cpp.

References OperatorEstBase::data_locality_.

175  : OneBodyDensityMatrices(obdm)
176 {
177  data_locality_ = dl;
178 }
DataLocality data_locality_
locality for accumulation of estimator data.
OneBodyDensityMatrices(OneBodyDensityMatricesInput &&obdmi, const Lattice &lattice, const SpeciesSet &species, const SPOMap &spomap, const ParticleSet &pset_target)
Standard Constructor Call this to make a new OBDM this is what you should be calling.

◆ OneBodyDensityMatrices() [3/3]

OneBodyDensityMatrices ( const OneBodyDensityMatrices obdm)
privatedefault

Default copy constructor.

Instances of this estimator is assume to be thread scope, i.e. never called by more than one thread at a time. note the OperatorEstBase copy constructor does not copy or even allocate data_

Member Function Documentation

◆ accumulate()

void accumulate ( const RefVector< MCPWalker > &  walkers,
const RefVector< ParticleSet > &  psets,
const RefVector< TrialWaveFunction > &  wfns,
const RefVector< QMCHamiltonian > &  hams,
RandomBase< FullPrecReal > &  rng 
)
override

Definition at line 529 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::implAccumulate(), and qmcplusplus::hdf::walkers.

534 {
535  implAccumulate(walkers, psets, wfns, rng);
536 }
const char walkers[]
Definition: HDFVersion.h:36
void implAccumulate(const RefVector< MCPWalker > &walkers, const RefVector< ParticleSet > &psets, const RefVector< TrialWaveFunction > &wfns, RandomBase< FullPrecReal > &rng)
Unfortunate design RandomGenerator type aliasing and virtual inheritance requires this for testing...

◆ calcDensity()

void calcDensity ( const Position r,
Real dens,
ParticleSet pset_target 
)
inlineprivate

calculate density based on r

Parameters
[in]rposition
[out]densdensity

called by generateDensitySamples to get trial dens. also called by test_derivatives. sideeffects:

  • updateBasis is called

Definition at line 457 of file OneBodyDensityMatrices.cpp.

References qmcplusplus::abs(), OneBodyDensityMatrices::basis_size_, OneBodyDensityMatrices::basis_values_, qmcplusplus::conj(), qmcplusplus::pset_target, and OneBodyDensityMatrices::updateBasis().

Referenced by OneBodyDensityMatrices::generateDensitySamples().

458 {
460  dens = 0.0;
461  for (int i = 0; i < basis_size_; ++i)
462  {
463  Value b = basis_values_[i];
464  dens += std::abs(qmcplusplus::conj(b) * b);
465  }
466  dens /= basis_size_;
467 }
MakeReturn< UnaryNode< FnFabs, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t abs(const Vector< T1, C1 > &l)
float conj(const float &c)
Workaround to allow conj on scalar to return real instead of complex.
void updateBasis(const Position &r, ParticleSet &pset_target)
basis set updates

◆ calcDensityDrift()

void calcDensityDrift ( const Position r,
Real dens,
Position drift,
ParticleSet pset_target 
)
private

calculate density and drift bashed on r

Parameters
[in]rposition
[out]densdensity
[out]driftdrift

called by warmupSamples to get an initial drift and density. called by generateDensitySamples to get trial drift and trial dens. also called by test_derivatives. sideeffects:

  • updateBasisD012 is called

Definition at line 484 of file OneBodyDensityMatrices.cpp.

References qmcplusplus::abs(), OneBodyDensityMatrices::basis_gradients_, OneBodyDensityMatrices::basis_size_, OneBodyDensityMatrices::basis_values_, qmcplusplus::conj(), OneBodyDensityMatricesInput::get_timestep(), OneBodyDensityMatrices::input_, OHMMS_DIM, qmcplusplus::pset_target, and OneBodyDensityMatrices::updateBasisD012().

Referenced by OneBodyDensityMatrices::generateDensitySamples(), and OneBodyDensityMatrices::warmupSampling().

485 {
487  dens = 0.0;
488  drift = 0.0;
489  for (int i = 0; i < basis_size_; ++i)
490  {
491  const Grad& bg = basis_gradients_[i];
492  Value b = basis_values_[i];
493  Value bc = qmcplusplus::conj(b);
494  dens += std::abs(bc * b);
495  for (int d = 0; d < OHMMS_DIM; ++d)
496  drift[d] += std::real(bc * bg[d]);
497  }
498  drift *= input_.get_timestep() / dens;
499  dens /= basis_size_;
500 }
QMCTraits::RealType real
void updateBasisD012(const Position &r, ParticleSet &pset_target)
evaluates vgl on basis_functions_ for r sideeffects:
MakeReturn< UnaryNode< FnFabs, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t abs(const Vector< T1, C1 > &l)
#define OHMMS_DIM
Definition: config.h:64
TinyVector< Value, OHMMS_DIM > Grad
float conj(const float &c)
Workaround to allow conj on scalar to return real instead of complex.

◆ calcDensityDriftWithSpin()

void calcDensityDriftWithSpin ( const Position r,
const Real s,
Real dens,
Position drift,
Real sdrift,
ParticleSet pset_target 
)
private

same as above, but with spin move

Definition at line 502 of file OneBodyDensityMatrices.cpp.

References qmcplusplus::abs(), OneBodyDensityMatrices::basis_gradients_, OneBodyDensityMatrices::basis_size_, OneBodyDensityMatrices::basis_spin_gradients_, OneBodyDensityMatrices::basis_values_, qmcplusplus::conj(), OneBodyDensityMatricesInput::get_timestep(), OneBodyDensityMatrices::input_, OHMMS_DIM, qmcplusplus::pset_target, qmcplusplus::Units::time::s, and OneBodyDensityMatrices::updateBasisD012WithSpin().

Referenced by OneBodyDensityMatrices::generateDensitySamplesWithSpin(), and OneBodyDensityMatrices::warmupSampling().

508 {
510  dens = 0.0;
511  drift = 0.0;
512  sdrift = 0.0;
513  for (int i = 0; i < basis_size_; ++i)
514  {
515  const Grad& bg = basis_gradients_[i];
516  const Value& bsg = basis_spin_gradients_[i];
517  Value b = basis_values_[i];
518  Value bc = qmcplusplus::conj(b);
519  dens += std::abs(bc * b);
520  for (int d = 0; d < OHMMS_DIM; ++d)
521  drift[d] += std::real(bc * bg[d]);
522  sdrift += std::real(bc * bsg);
523  }
524  drift *= input_.get_timestep() / dens;
525  sdrift *= input_.get_timestep() / dens;
526  dens /= basis_size_;
527 }
QMCTraits::RealType real
MakeReturn< UnaryNode< FnFabs, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t abs(const Vector< T1, C1 > &l)
#define OHMMS_DIM
Definition: config.h:64
TinyVector< Value, OHMMS_DIM > Grad
float conj(const float &c)
Workaround to allow conj on scalar to return real instead of complex.
void updateBasisD012WithSpin(const Position &r, const Real &s, ParticleSet &pset_target)
same as above, but includes spin gradients

◆ calcDensityWithSpin()

void calcDensityWithSpin ( const Position r,
const Real s,
Real dens,
ParticleSet pset_target 
)
inlineprivate

same as above, but with spin move

Definition at line 469 of file OneBodyDensityMatrices.cpp.

References qmcplusplus::abs(), OneBodyDensityMatrices::basis_size_, OneBodyDensityMatrices::basis_values_, qmcplusplus::conj(), qmcplusplus::pset_target, qmcplusplus::Units::time::s, and OneBodyDensityMatrices::updateBasisWithSpin().

Referenced by OneBodyDensityMatrices::generateDensitySamplesWithSpin().

473 {
475  dens = 0.0;
476  for (int i = 0; i < basis_size_; ++i)
477  {
478  Value b = basis_values_[i];
479  dens += std::abs(qmcplusplus::conj(b) * b);
480  }
481  dens /= basis_size_;
482 }
MakeReturn< UnaryNode< FnFabs, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t abs(const Vector< T1, C1 > &l)
void updateBasisWithSpin(const Position &r, const Real &s, ParticleSet &pset_target)
basis set updates with spin
float conj(const float &c)
Workaround to allow conj on scalar to return real instead of complex.

◆ calcFullDataSize()

size_t calcFullDataSize ( size_t  basis_size,
int  num_species 
)
private

Definition at line 198 of file OneBodyDensityMatrices.cpp.

Referenced by OneBodyDensityMatrices::OneBodyDensityMatrices().

199 {
200  if constexpr (IsComplex_t<Value>::value)
201  return 2 * basis_size * basis_size * nspecies;
202  else
203  return basis_size * basis_size * nspecies;
204 }

◆ diffuse()

OneBodyDensityMatrices::Position diffuse ( const Real  sqt,
RandomBase< FullPrecReal > &  rng 
)
private

produce a position difference vector from timestep

Definition at line 441 of file OneBodyDensityMatrices.cpp.

References qmcplusplus::assignGaussRand(), and OHMMS_DIM.

Referenced by OneBodyDensityMatrices::generateDensitySamples(), OneBodyDensityMatrices::generateDensitySamplesWithSpin(), and OneBodyDensityMatrices::warmupSampling().

442 {
443  Position diff;
444  assignGaussRand(&diff[0], OHMMS_DIM, rng);
445  diff *= sqt;
446  return diff;
447 }
#define OHMMS_DIM
Definition: config.h:64
void assignGaussRand(T *restrict a, unsigned n, RG &rng)

◆ diffuseSpin()

OneBodyDensityMatrices::Real diffuseSpin ( const Real  sqt,
RandomBase< FullPrecReal > &  rng 
)
private

spin diffusion

Definition at line 449 of file OneBodyDensityMatrices.cpp.

References qmcplusplus::assignGaussRand().

Referenced by OneBodyDensityMatrices::generateDensitySamplesWithSpin(), and OneBodyDensityMatrices::warmupSampling().

450 {
451  Real diff;
452  assignGaussRand(&diff, 1, rng);
453  diff *= sqt;
454  return diff;
455 }
void assignGaussRand(T *restrict a, unsigned n, RG &rng)

◆ evaluateMatrix()

void evaluateMatrix ( ParticleSet pset_target,
TrialWaveFunction psi_target,
const MCPWalker walker,
RandomBase< FullPrecReal > &  rng 
)
private

Definition at line 550 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::OneBodyDensityMatrixTimers::accumulate_timer, OneBodyDensityMatrices::basis_size_, OperatorEstBase::data_, qmcplusplus::MatrixOperators::diag_product(), OneBodyDensityMatrices::generateParticleBasis(), OneBodyDensityMatrices::generateSampleBasis(), OneBodyDensityMatrices::generateSampleRatios(), OneBodyDensityMatrices::generateSamples(), qmcplusplus::imag(), OneBodyDensityMatrices::OneBodyDensityMatrixTimers::matrix_products_timer, OneBodyDensityMatrices::metric_, qmcplusplus::n, OneBodyDensityMatrices::N_BB_, OneBodyDensityMatrices::Phi_MB_, OneBodyDensityMatrices::Phi_NB_, OneBodyDensityMatrices::Phi_Psi_NB_, qmcplusplus::MatrixOperators::product(), qmcplusplus::MatrixOperators::product_AtB(), qmcplusplus::pset_target, OneBodyDensityMatrices::Psi_NM_, qmcplusplus::real(), qmcplusplus::Units::time::s, OneBodyDensityMatrices::samples_weights_, SpeciesSet::size(), OneBodyDensityMatrices::species_, OneBodyDensityMatrices::timers_, qmcplusplus::walker, and OneBodyDensityMatrices::warmupSampling().

Referenced by OneBodyDensityMatrices::implAccumulate(), and OneBodyDensityMatricesTests< T >::testEvaluateMatrix().

554 {
555  //perform warmup sampling the first time
557  // get weight and single particle energy trace data
558  Real weight = walker.Weight * metric_;
559 
560  // compute sample positions (monte carlo or deterministic)
561  generateSamples(weight, pset_target, rng);
562  // compute basis and wavefunction ratio values in matrix form
563  generateSampleBasis(Phi_MB_, pset_target, psi_target); // basis : samples x basis_size
564  generateSampleRatios(pset_target, psi_target, Psi_NM_); // conj(Psi ratio) : particles x samples
565  generateParticleBasis(pset_target, Phi_NB_); // conj(basis) : particles x basis_size
566 
567  // \todo separate testable and optimizable block, should be function
568  // perform integration via matrix products
569  {
571  for (int s = 0; s < species_.size(); ++s)
572  {
573  Matrix<Value>& Psi_nm = Psi_NM_[s];
574  Matrix<Value>& Phi_Psi_nb = Phi_Psi_NB_[s];
575  Matrix<Value>& Phi_nb = Phi_NB_[s];
576  diag_product(Psi_nm, samples_weights_, Psi_nm);
577  product(Psi_nm, Phi_MB_, Phi_Psi_nb); // ratio*basis : particles x basis_size
578  product_AtB(Phi_nb, Phi_Psi_nb, N_BB_[s]); // conj(basis)^T*ratio*basis : basis_size^2
579  }
580  }
581  // accumulate data for this walker
582  {
583  ScopedTimer local_timer(timers_.accumulate_timer);
584  const int basis_size_sq = basis_size_ * basis_size_;
585  int ij = 0;
586  for (int s = 0; s < species_.size(); ++s)
587  {
588  //int ij=nindex; // for testing
589  const Matrix<Value>& NDM = N_BB_[s];
590  for (int n = 0; n < basis_size_sq; ++n)
591  {
592  Value val = NDM(n);
593  data_[ij] += real(val);
594  ij++;
595 #if defined(QMC_COMPLEX)
596  data_[ij] += imag(val);
597  ij++;
598 #endif
599  }
600  }
601  }
602 }
void generateSampleRatios(ParticleSet &pset_target, TrialWaveFunction &psi_target, std::vector< Matrix< Value >> &Psi_nm)
float real(const float &c)
real part of a scalar. Cannot be replaced by std::real due to AFQMC specific needs.
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
std::vector< Matrix< Value > > Phi_NB_
row major per sample workspaces
void generateSampleBasis(Matrix< Value > &Phi_mb, ParticleSet &pset_target, TrialWaveFunction &psi_target)
set Phi_mp to basis vaules per sample sideeffects:
Matrix< Value > Phi_MB_
basis_values_ at each r of rsamples_ row: sample col: basis_value size: samples * basis_size ...
int size() const
return the number of species
Definition: SpeciesSet.h:53
void generateSamples(const Real weight, ParticleSet &pset_target, RandomBase< FullPrecReal > &rng, int steps=0)
Dispatch method to difference methods of generating samples.
float imag(const float &c)
imaginary part of a scalar. Cannot be replaced by std::imag due to AFQMC specific needs...
void warmupSampling(ParticleSet &pset_target, RandomBase< FullPrecReal > &rng)
does some warmup sampling i.e.
void product(const Matrix< T > &A, const Matrix< T > &B, Matrix< T > &C)
static function to perform C=AB for real matrices
void generateParticleBasis(ParticleSet &pset_target, std::vector< Matrix< Value >> &phi_nb)
set phi_nb to basis values per target particleset particle sideeffects:
std::vector< Matrix< Value > > Phi_Psi_NB_
std::vector< Matrix< Value > > N_BB_
void product_AtB(const Matrix< T > &A, const Matrix< T > &B, Matrix< T > &C)
void diag_product(const Matrix< T1 > &A, const Vector< T2 > &B, Matrix< T3 > &C)
C = A*diag(B)
std::vector< Matrix< Value > > Psi_NM_
ratio per particle per sample size: particles * samples vector is over species each matrix row: parti...

◆ generateDensitySamples()

void generateDensitySamples ( bool  save,
int  steps,
RandomBase< FullPrecReal > &  rng,
ParticleSet pset_target 
)
inlineprivate

generate samples for density integration

Parameters
[in]saveif false throw out the samples
[in]stepsactually the number of samples which are basically steps.
[in]rngrandom generator
[in]pset_targetwill be returned to its initial state but is mutated.

sideeffects:

Definition at line 311 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::acceptance_ratio_, qmcplusplus::app_log(), OneBodyDensityMatrices::calcDensity(), OneBodyDensityMatrices::calcDensityDrift(), OneBodyDensityMatrices::diffuse(), qmcplusplus::dot(), OneBodyDensityMatrices::dpcur_, qmcplusplus::exp(), OneBodyDensityMatricesInput::get_timestep(), OneBodyDensityMatricesInput::get_use_drift(), OneBodyDensityMatricesInput::get_write_acceptance_ratio(), OneBodyDensityMatrices::input_, OneBodyDensityMatrices::naccepted_, OneBodyDensityMatrices::nmoves_, omp_get_thread_num(), qmcplusplus::pset_target, OneBodyDensityMatrices::rhocur_, OneBodyDensityMatrices::rpcur_, OneBodyDensityMatrices::rsamples_, qmcplusplus::Units::time::s, OneBodyDensityMatrices::samples_weights_, and qmcplusplus::sqrt().

Referenced by OneBodyDensityMatrices::generateSamples().

315 {
316  const auto timestep = input_.get_timestep();
317  Real sqt = std::sqrt(timestep);
318  Real ot = 1.0 / timestep;
319  Position r = rpcur_; //current position
320  Position d = dpcur_; //current drift
321  Real rho = rhocur_; //current density
322  for (int s = 0; s < steps; ++s)
323  {
324  nmoves_++;
325  Position rp; // trial pos
326  Position dp; // trial drift
327  Position ds; // drift sum
328  Real rhop; // trial density
329  Real ratio; // dens ratio
330  Real Pacc; // acc prob
331  Position diff = diffuse(sqt, rng); // get diffusion
332  if (input_.get_use_drift())
333  {
334  rp = r + diff + d; //update trial position
335  calcDensityDrift(rp, rhop, dp, pset_target); //get trial drift and density
336  ratio = rhop / rho; //density ratio
337  ds = dp + d; //drift sum
338  Pacc = ratio * std::exp(-ot * (dot(diff, ds) + .5 * dot(ds, ds))); //acceptance probability
339  }
340  else
341  {
342  rp = r + diff; //update trial position
343  calcDensity(rp, rhop, pset_target); //get trial density
344  ratio = rhop / rho; //density ratio
345  Pacc = ratio; //acceptance probability
346  }
347  if (rng() < Pacc)
348  { //accept move
349  r = rp;
350  d = dp;
351  rho = rhop;
352  naccepted_++;
353  }
354  if (save)
355  {
356  rsamples_[s] = r;
357  samples_weights_[s] = 1.0 / rho;
358  }
359  }
361 
363  app_log() << "dm1b acceptance_ratio = " << acceptance_ratio_ << std::endl;
364 
365  rpcur_ = r;
366  dpcur_ = d;
367  rhocur_ = rho;
368 }
Real acceptance_ratio_
running acceptance ratio over all samples
std::ostream & app_log()
Definition: OutputManager.h:65
void calcDensity(const Position &r, Real &dens, ParticleSet &pset_target)
calculate density based on r
void calcDensityDrift(const Position &r, Real &dens, Position &drift, ParticleSet &pset_target)
calculate density and drift bashed on r
omp_int_t omp_get_thread_num()
Definition: OpenMP.h:25
Position diffuse(const Real sqt, RandomBase< FullPrecReal > &rng)
produce a position difference vector from timestep
MakeReturn< UnaryNode< FnExp, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t exp(const Vector< T1, C1 > &l)
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
Tensor< typename BinaryReturn< T1, T2, OpMultiply >::Type_t, D > dot(const AntiSymTensor< T1, D > &lhs, const AntiSymTensor< T2, D > &rhs)
Position rpcur_
current position – As long Positions are TinyVectors they are intialized to zero vectors ...
int naccepted_
number of accepted samples

◆ generateDensitySamplesWithSpin()

void generateDensitySamplesWithSpin ( bool  save,
int  steps,
RandomBase< FullPrecReal > &  rng,
ParticleSet pset_target 
)
inlineprivate

same as above, but with spin variables included

Definition at line 370 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::acceptance_ratio_, qmcplusplus::app_log(), OneBodyDensityMatrices::calcDensityDriftWithSpin(), OneBodyDensityMatrices::calcDensityWithSpin(), OneBodyDensityMatrices::diffuse(), OneBodyDensityMatrices::diffuseSpin(), qmcplusplus::dot(), OneBodyDensityMatrices::dpcur_, OneBodyDensityMatrices::dspcur_, qmcplusplus::exp(), OneBodyDensityMatricesInput::get_timestep(), OneBodyDensityMatricesInput::get_use_drift(), OneBodyDensityMatricesInput::get_write_acceptance_ratio(), OneBodyDensityMatrices::input_, OneBodyDensityMatrices::naccepted_, OneBodyDensityMatrices::nmoves_, omp_get_thread_num(), qmcplusplus::pset_target, OneBodyDensityMatrices::rhocur_, OneBodyDensityMatrices::rpcur_, OneBodyDensityMatrices::rsamples_, qmcplusplus::Units::time::s, OneBodyDensityMatrices::samples_weights_, OneBodyDensityMatrices::spcur_, qmcplusplus::sqrt(), and OneBodyDensityMatrices::ssamples_.

Referenced by OneBodyDensityMatrices::generateSamples().

374 {
375  const auto timestep = input_.get_timestep();
376  Real sqt = std::sqrt(timestep);
377  Real ot = 1.0 / timestep;
378  Position r = rpcur_; //current position
379  Position d = dpcur_; //current drift
380  Real rho = rhocur_; //current density
381  for (int s = 0; s < steps; ++s)
382  {
383  nmoves_++;
384  Position rp; // trial pos
385  Position dp; // trial drift
386  Position ds; // drift sum
387  Real rhop; // trial density
388  Real ratio; // dens ratio
389  Real Pacc; // acc prob
390  Position diff = diffuse(sqt, rng); // get diffusion
391 
392  //now do spin variables
393  Real spinp; // trial spin
394  Real dspinp = 0; // trial spindrifty
395  Real sdiff = diffuseSpin(sqt, rng); //spin diffusion
396  if (input_.get_use_drift())
397  {
398  rp = r + diff + d; //update trial position
399  spinp = spcur_ + sdiff + dspcur_; //update trial spin
400  calcDensityDriftWithSpin(rp, spinp, rhop, dp, dspinp, pset_target); //get trial drift and density
401  ratio = rhop / rho; //density ratio
402  ds = dp + d; //drift sum
403  auto spin_ds = dspinp + dspcur_; //spin drift sum
404  Pacc = ratio * std::exp(-ot * (dot(diff, ds) + .5 * dot(ds, ds))) *
405  std::exp(-ot * (sdiff * spin_ds) + 0.5 * spin_ds * spin_ds); //acceptance probability
406  }
407  else
408  {
409  rp = r + diff; //update trial position
410  spinp = spcur_ + sdiff; //update trial spin
411  calcDensityWithSpin(rp, spinp, rhop, pset_target); //get trial density
412  ratio = rhop / rho; //density ratio
413  Pacc = ratio; //acceptance probability
414  }
415  if (rng() < Pacc)
416  { //accept move
417  r = rp;
418  d = dp;
419  spcur_ = spinp;
420  dspcur_ = dspinp;
421  rho = rhop;
422  naccepted_++;
423  }
424  if (save)
425  {
426  rsamples_[s] = r;
427  samples_weights_[s] = 1.0 / rho;
428  ssamples_[s] = spcur_;
429  }
430  }
432 
434  app_log() << "dm1b acceptance_ratio = " << acceptance_ratio_ << std::endl;
435 
436  rpcur_ = r;
437  dpcur_ = d;
438  rhocur_ = rho;
439 }
Real acceptance_ratio_
running acceptance ratio over all samples
std::ostream & app_log()
Definition: OutputManager.h:65
omp_int_t omp_get_thread_num()
Definition: OpenMP.h:25
Position diffuse(const Real sqt, RandomBase< FullPrecReal > &rng)
produce a position difference vector from timestep
MakeReturn< UnaryNode< FnExp, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t exp(const Vector< T1, C1 > &l)
void calcDensityDriftWithSpin(const Position &r, const Real &s, Real &dens, Position &drift, Real &sdrift, ParticleSet &pset_target)
same as above, but with spin move
void calcDensityWithSpin(const Position &r, const Real &s, Real &dens, ParticleSet &pset_target)
same as above, but with spin move
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
Real diffuseSpin(const Real sqt, RandomBase< FullPrecReal > &rng)
spin diffusion
Tensor< typename BinaryReturn< T1, T2, OpMultiply >::Type_t, D > dot(const AntiSymTensor< T1, D > &lhs, const AntiSymTensor< T2, D > &rhs)
Position rpcur_
current position – As long Positions are TinyVectors they are intialized to zero vectors ...
int naccepted_
number of accepted samples

◆ generateParticleBasis()

void generateParticleBasis ( ParticleSet pset_target,
std::vector< Matrix< Value >> &  phi_nb 
)
private

set phi_nb to basis values per target particleset particle sideeffects:

  • updates basis_values_ to last rsample

Definition at line 604 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::basis_size_, OneBodyDensityMatrices::basis_values_, qmcplusplus::conj(), OneBodyDensityMatrices::OneBodyDensityMatrixTimers::gen_particle_basis_timer, OneBodyDensityMatrices::is_spinor_, qmcplusplus::n, qmcplusplus::pset_target, qmcplusplus::Units::time::s, SpeciesSet::size(), OneBodyDensityMatrices::species_, OneBodyDensityMatrices::species_sizes_, OneBodyDensityMatrices::timers_, OneBodyDensityMatrices::updateBasis(), and OneBodyDensityMatrices::updateBasisWithSpin().

Referenced by OneBodyDensityMatrices::evaluateMatrix().

605 {
607  int p = 0;
608  for (int s = 0; s < species_.size(); ++s)
609  {
610  int nb = 0;
611  Matrix<Value>& P_nb = phi_nb[s];
612  for (int n = 0; n < species_sizes_[s]; ++n, ++p)
613  {
614  if (is_spinor_)
616  else
618  for (int b = 0; b < basis_size_; ++b, ++nb)
619  P_nb(nb) = qmcplusplus::conj(basis_values_[b]);
620  }
621  }
622 }
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
void updateBasisWithSpin(const Position &r, const Real &s, ParticleSet &pset_target)
basis set updates with spin
int size() const
return the number of species
Definition: SpeciesSet.h:53
float conj(const float &c)
Workaround to allow conj on scalar to return real instead of complex.
void updateBasis(const Position &r, ParticleSet &pset_target)
basis set updates

◆ generateSampleBasis()

void generateSampleBasis ( Matrix< Value > &  Phi_mb,
ParticleSet pset_target,
TrialWaveFunction psi_target 
)
private

set Phi_mp to basis vaules per sample sideeffects:

  • updates basis_values_ to last rsample

Definition at line 624 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::basis_size_, OneBodyDensityMatrices::basis_values_, OneBodyDensityMatrices::OneBodyDensityMatrixTimers::gen_sample_basis_timer, OneBodyDensityMatrices::is_spinor_, qmcplusplus::Units::distance::m, qmcplusplus::pset_target, OneBodyDensityMatrices::rsamples_, OneBodyDensityMatrices::samples_, OneBodyDensityMatrices::ssamples_, OneBodyDensityMatrices::timers_, OneBodyDensityMatrices::updateBasis(), and OneBodyDensityMatrices::updateBasisWithSpin().

Referenced by OneBodyDensityMatrices::evaluateMatrix().

627 {
629  int mb = 0;
630  for (int m = 0; m < samples_; ++m)
631  {
632  if (is_spinor_)
634  else
636  for (int b = 0; b < basis_size_; ++b, ++mb)
637  Phi_mb(mb) = basis_values_[b];
638  }
639 }
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
void updateBasisWithSpin(const Position &r, const Real &s, ParticleSet &pset_target)
basis set updates with spin
void updateBasis(const Position &r, ParticleSet &pset_target)
basis set updates

◆ generateSampleRatios()

void generateSampleRatios ( ParticleSet pset_target,
TrialWaveFunction psi_target,
std::vector< Matrix< Value >> &  Psi_nm 
)
private

Definition at line 641 of file OneBodyDensityMatrices.cpp.

References TrialWaveFunction::calcRatio(), qmcplusplus::conj(), TrialWaveFunction::evaluateRatiosAlltoOne(), OneBodyDensityMatrices::OneBodyDensityMatrixTimers::gen_sample_ratios_timer, OneBodyDensityMatrices::is_spinor_, qmcplusplus::Units::distance::m, qmcplusplus::n, qmcplusplus::pset_target, OneBodyDensityMatrices::psi_ratios_, TrialWaveFunction::resetPhaseDiff(), OneBodyDensityMatrices::rsamples_, qmcplusplus::Units::time::s, OneBodyDensityMatrices::samples_, SpeciesSet::size(), OneBodyDensityMatrices::species_, OneBodyDensityMatrices::species_sizes_, OneBodyDensityMatrices::ssamples_, and OneBodyDensityMatrices::timers_.

Referenced by OneBodyDensityMatrices::evaluateMatrix().

644 {
646  for (int m = 0; m < samples_; ++m)
647  {
648  // get N ratios for the current sample point
649  if (!is_spinor_)
650  {
651  pset_target.makeVirtualMoves(rsamples_[m]);
652  psi_target.evaluateRatiosAlltoOne(pset_target, psi_ratios_);
653  }
654  else
655  {
656  //note: makeVirtualMoves updates distance tables
657  //There is not a corresponding "distance table" for the spins, so we need to brute force this by moving each electron and using makeMoveWithSpin, calcRatio, and rejectMove, and resetPhaseDiff, similar to how NonLocalECP works for evaluating quadrature points without virtual particles
658  int p = 0;
659  for (int s = 0; s < species_.size(); ++s)
660  for (int n = 0; n < species_sizes_[s]; ++n, ++p)
661  {
662  pset_target.makeMoveWithSpin(p, rsamples_[m] - pset_target.R[p], ssamples_[m] - pset_target.spins[p]);
663  psi_ratios_[p] = psi_target.calcRatio(pset_target, p);
664  pset_target.rejectMove(p);
665  psi_target.resetPhaseDiff();
666  }
667  }
668 
669  // collect ratios into per-species matrices
670  int p = 0;
671  for (int s = 0; s < species_.size(); ++s)
672  {
673  Matrix<Value>& P_nm = psi_nm[s];
674  for (int n = 0; n < species_sizes_[s]; ++n, ++p)
675  {
676  P_nm(n, m) = qmcplusplus::conj(psi_ratios_[p]);
677  }
678  }
679  }
680 }
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
int size() const
return the number of species
Definition: SpeciesSet.h:53
float conj(const float &c)
Workaround to allow conj on scalar to return real instead of complex.
std::vector< Value > psi_ratios_
, &#39; .

◆ generateSamples()

void generateSamples ( const Real  weight,
ParticleSet pset_target,
RandomBase< FullPrecReal > &  rng,
int  steps = 0 
)
private

Dispatch method to difference methods of generating samples.

dispatch determined by Integrator.

Parameters
[in]weightof this walker's samples
[in]pset_targetwill be returned to its initial state but is mutated.
[in]rngrandom generator. templated for testing without dependency on app level rng.
[in]stepsIf integrator_ = Integrator::DENSITY steps is a key parameter otherwise ignored. when set to 0 it is reset to samples_ internally

sideeffects:

  • samples_weights_ are set.
  • rsamples_ are set. for Density
    • update basis_values_, basis_gradients_, basis_laplacians_

Definition at line 208 of file OneBodyDensityMatrices.cpp.

References qmcplusplus::app_log(), Vector< T, Alloc >::begin(), OneBodyDensityMatricesInput::DENSITY, Vector< T, Alloc >::end(), OneBodyDensityMatrices::OneBodyDensityMatrixTimers::gen_samples_timer, OneBodyDensityMatrices::generateDensitySamples(), OneBodyDensityMatrices::generateDensitySamplesWithSpin(), OneBodyDensityMatrices::generateUniformGrid(), OneBodyDensityMatrices::generateUniformSamples(), OneBodyDensityMatricesInput::get_integrator(), OneBodyDensityMatricesInput::get_rstats(), OneBodyDensityMatrices::input_, OneBodyDensityMatrices::is_spinor_, OneBodyDensityMatrices::METROPOLIS, omptarget::min(), OHMMS_DIM, omp_get_thread_num(), qmcplusplus::pset_target, OneBodyDensityMatrices::rsamples_, qmcplusplus::Units::time::s, OneBodyDensityMatrices::samples_, OneBodyDensityMatrices::samples_weights_, OneBodyDensityMatrices::sampling_, qmcplusplus::sqrt(), OneBodyDensityMatrices::timers_, OneBodyDensityMatricesInput::UNIFORM, and OneBodyDensityMatricesInput::UNIFORM_GRID.

Referenced by OneBodyDensityMatrices::evaluateMatrix(), OneBodyDensityMatricesTests< T >::testGenerateSamples(), and OneBodyDensityMatrices::warmupSampling().

212 {
214 
215  // Steps will always be 0 unless these are samples for warmup which is only for metropolis
216  // This is not a clear way to write this
217  // \todo rewrite to make algorithm more clears
218  bool save = false;
219  if (steps == 0)
220  {
221  save = true;
222  steps = samples_;
223  }
224 
225  switch (input_.get_integrator())
226  {
228  generateUniformGrid(rng);
229  break;
230  case Integrator::UNIFORM:
232  break;
233  case Integrator::DENSITY:
234  if (is_spinor_)
235  generateDensitySamplesWithSpin(save, steps, rng, pset_target);
236  else
237  generateDensitySamples(save, steps, rng, pset_target);
238  break;
239  }
240 
241  if (save)
242  {
244  samples_weights_ *= weight;
245  else
246  {
247  std::fill(samples_weights_.begin(), samples_weights_.end(), weight);
248  }
249  }
250 
251  // optional check
252  if (input_.get_rstats() && omp_get_thread_num() == 0)
253  {
254  Position rmin = std::numeric_limits<Real>::max();
255  Position rmax = -std::numeric_limits<Real>::max();
256  Position rmean = 0.0;
257  Position rstd = 0.0;
258  for (int s = 0; s < rsamples_.size(); ++s)
259  for (int d = 0; d < OHMMS_DIM; ++d)
260  {
261  Real rd = rsamples_[s][d];
262  rmin[d] = std::min(rmin[d], rd);
263  rmax[d] = std::max(rmax[d], rd);
264  rmean[d] += rd;
265  rstd[d] += rd * rd;
266  }
267  rmean /= rsamples_.size();
268  rstd /= rsamples_.size();
269  for (int d = 0; d < OHMMS_DIM; ++d)
270  rstd[d] = std::sqrt(rstd[d] - rmean[d] * rmean[d]);
271  app_log() << "\nrsamples properties:" << std::endl;
272  app_log() << " rmin = " << rmin << std::endl;
273  app_log() << " rmax = " << rmax << std::endl;
274  app_log() << " rmean = " << rmean << std::endl;
275  app_log() << " rstd = " << rstd << std::endl;
276  }
277 }
std::ostream & app_log()
Definition: OutputManager.h:65
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
#define OHMMS_DIM
Definition: config.h:64
void generateDensitySamplesWithSpin(bool save, int steps, RandomBase< FullPrecReal > &rng, ParticleSet &pset_target)
same as above, but with spin variables included
T min(T a, T b)
omp_int_t omp_get_thread_num()
Definition: OpenMP.h:25
void generateDensitySamples(bool save, int steps, RandomBase< FullPrecReal > &rng, ParticleSet &pset_target)
generate samples for density integration
void generateUniformGrid(RandomBase< FullPrecReal > &rng)
Sampling sampling_
Sampling method, this derived values in input_.
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
void generateUniformSamples(RandomBase< FullPrecReal > &rng)

◆ generateUniformGrid()

void generateUniformGrid ( RandomBase< FullPrecReal > &  rng)
inlineprivate

Definition at line 279 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatricesInput::get_points(), OneBodyDensityMatricesInput::get_scale(), OneBodyDensityMatrices::ind_dims_, OneBodyDensityMatrices::input_, OneBodyDensityMatrices::lattice_, OHMMS_DIM, OneBodyDensityMatrices::rcorner_, OneBodyDensityMatrices::rsamples_, qmcplusplus::Units::time::s, OneBodyDensityMatrices::samples_, and CrystalLattice< T, D >::toCart().

Referenced by OneBodyDensityMatrices::generateSamples().

280 {
281  Position rp;
282  Position ushift = 0.0;
284  for (int d = 0; d < OHMMS_DIM; ++d)
285  ushift[d] += rng() * du;
286  for (int s = 0; s < samples_; ++s)
287  {
288  int nrem = s;
289  for (int d = 0; d < OHMMS_DIM - 1; ++d)
290  {
291  int ind = nrem / ind_dims_[d];
292  rp[d] = ind * du + ushift[d];
293  nrem -= ind * ind_dims_[d];
294  }
295  rp[OHMMS_DIM - 1] = nrem * du + ushift[OHMMS_DIM - 1];
297  }
298 }
Position rcorner_
with respect to center_ using lattice_;
#define OHMMS_DIM
Definition: config.h:64
SingleParticlePos toCart(const TinyVector< T1, D > &c) const
Convert a unit vector to a cartesian vector.

◆ generateUniformSamples()

void generateUniformSamples ( RandomBase< FullPrecReal > &  rng)
inlineprivate

Definition at line 300 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatricesInput::get_scale(), OneBodyDensityMatrices::input_, OneBodyDensityMatrices::lattice_, OHMMS_DIM, OneBodyDensityMatrices::rcorner_, OneBodyDensityMatrices::rsamples_, qmcplusplus::Units::time::s, OneBodyDensityMatrices::samples_, and CrystalLattice< T, D >::toCart().

Referenced by OneBodyDensityMatrices::generateSamples().

301 {
302  Position rp;
303  for (int s = 0; s < samples_; ++s)
304  {
305  for (int d = 0; d < OHMMS_DIM; ++d)
306  rp[d] = input_.get_scale() * rng();
308  }
309 }
Position rcorner_
with respect to center_ using lattice_;
#define OHMMS_DIM
Definition: config.h:64
SingleParticlePos toCart(const TinyVector< T1, D > &c) const
Convert a unit vector to a cartesian vector.

◆ implAccumulate()

void implAccumulate ( const RefVector< MCPWalker > &  walkers,
const RefVector< ParticleSet > &  psets,
const RefVector< TrialWaveFunction > &  wfns,
RandomBase< FullPrecReal > &  rng 
)
private

Unfortunate design RandomGenerator type aliasing and virtual inheritance requires this for testing.

Definition at line 538 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::evaluateMatrix(), qmcplusplus::hdf::walkers, and OperatorEstBase::walkers_weight_.

Referenced by OneBodyDensityMatrices::accumulate(), and OneBodyDensityMatricesTests< T >::testAccumulate().

542 {
543  for (int iw = 0; iw < walkers.size(); ++iw)
544  {
545  walkers_weight_ += walkers[iw].get().Weight;
546  evaluateMatrix(psets[iw], wfns[iw], walkers[iw], rng);
547  }
548 }
const char walkers[]
Definition: HDFVersion.h:36
void evaluateMatrix(ParticleSet &pset_target, TrialWaveFunction &psi_target, const MCPWalker &walker, RandomBase< FullPrecReal > &rng)
QMCT::FullPrecRealType walkers_weight_

◆ normalizeBasis()

void normalizeBasis ( ParticleSet pset_target)
inlineprivate

Definition at line 751 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::basis_norms_, OneBodyDensityMatrices::basis_size_, OneBodyDensityMatrices::basis_values_, Vector< T, Alloc >::begin(), qmcplusplus::conj(), Vector< T, Alloc >::end(), OneBodyDensityMatricesInput::get_points(), OneBodyDensityMatricesInput::get_scale(), OneBodyDensityMatrices::input_, OneBodyDensityMatrices::lattice_, OHMMS_DIM, qmcplusplus::pow(), qmcplusplus::pset_target, OneBodyDensityMatrices::rcorner_, qmcplusplus::real(), Vector< T, Alloc >::resize(), qmcplusplus::sqrt(), CrystalLattice< T, D >::toCart(), OneBodyDensityMatrices::updateBasis(), and OneBodyDensityMatrices::volume_.

Referenced by OneBodyDensityMatrices::OneBodyDensityMatrices().

752 {
753  int ngrid = std::max(200, input_.get_points());
754  int ngtot = pow(ngrid, OHMMS_DIM);
755  Real du = input_.get_scale() / ngrid;
756  Real dV = volume_ / ngtot;
757  Position rp;
758  Vector<Value> bnorms;
759  int gdims[OHMMS_DIM];
760  gdims[0] = pow(ngrid, OHMMS_DIM - 1);
761  for (int d = 1; d < OHMMS_DIM; ++d)
762  gdims[d] = gdims[d - 1] / ngrid;
763  bnorms.resize(basis_size_);
764  for (int i = 0; i < basis_size_; ++i)
765  bnorms[i] = 0.0;
766  std::fill(basis_norms_.begin(), basis_norms_.end(), 1.0);
767  for (int p = 0; p < ngtot; ++p)
768  {
769  int nrem = p;
770  for (int d = 0; d < OHMMS_DIM - 1; ++d)
771  {
772  int ind = nrem / gdims[d];
773  rp[d] = ind * du + du / 2;
774  nrem -= ind * gdims[d];
775  }
776  rp[OHMMS_DIM - 1] = nrem * du + du / 2;
777  rp = lattice_.toCart(rp) + rcorner_;
779  for (int i = 0; i < basis_size_; ++i)
780  bnorms[i] += qmcplusplus::conj(basis_values_[i]) * basis_values_[i] * dV;
781  }
782  for (int i = 0; i < basis_size_; ++i)
783  basis_norms_[i] = 1.0 / std::sqrt(real(bnorms[i]));
784 }
float real(const float &c)
real part of a scalar. Cannot be replaced by std::real due to AFQMC specific needs.
Position rcorner_
with respect to center_ using lattice_;
#define OHMMS_DIM
Definition: config.h:64
MakeReturn< BinaryNode< FnPow, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t, typename CreateLeaf< Vector< T2, C2 > >::Leaf_t > >::Expression_t pow(const Vector< T1, C1 > &l, const Vector< T2, C2 > &r)
float conj(const float &c)
Workaround to allow conj on scalar to return real instead of complex.
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
void updateBasis(const Position &r, ParticleSet &pset_target)
basis set updates
SingleParticlePos toCart(const TinyVector< T1, D > &c) const
Convert a unit vector to a cartesian vector.

◆ registerOperatorEstimator()

void registerOperatorEstimator ( hdf_archive file)
overridevirtual

create and tie OperatorEstimator's observable_helper hdf5 wrapper to stat.h5 file

Parameters
gidhdf5 group to which the observables belong

The default implementation does nothing. The derived classes which compute big data, e.g. density, should overwrite this function.

Reimplemented from OperatorEstBase.

Definition at line 786 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::basis_size_, OperatorEstBase::h5desc_, OperatorEstBase::my_name_, qmcplusplus::oh, qmcplusplus::Units::time::s, ObservableHelper::set_dimensions(), SpeciesSet::size(), OneBodyDensityMatrices::species_, and SpeciesSet::speciesName.

Referenced by OneBodyDensityMatricesTests< T >::testRegisterAndWrite().

787 {
788  std::vector<int> my_indexes(2, basis_size_);
789  if constexpr (IsComplex_t<Value>::value)
790  {
791  my_indexes.push_back(2);
792  }
793  int nentries = std::accumulate(my_indexes.begin(), my_indexes.end(), 1);
794 
795  int spin_data_size = 0;
796  if constexpr (IsComplex_t<Value>::value)
797  spin_data_size = 2 * basis_size_ * basis_size_;
798  else
799  spin_data_size = basis_size_ * basis_size_;
800 
801  hdf_path hdf_name{my_name_};
802  hdf_name /= "number_matrix";
803  for (int s = 0; s < species_.size(); ++s)
804  {
805  h5desc_.emplace_back(hdf_name / species_.speciesName[s]);
806  auto& oh = h5desc_.back();
807  oh.set_dimensions(my_indexes, s * spin_data_size);
808  }
809 }
std::vector< ObservableHelper > h5desc_
int size() const
return the number of species
Definition: SpeciesSet.h:53
std::vector< std::string > speciesName
Species name list.
Definition: SpeciesSet.h:44
void set_dimensions(const std::vector< int > &dims, int first)
set the shape of this observable
ObservableHelper oh
std::string my_name_
name of this object – only used for debugging and h5 output

◆ report()

void report ( const std::string &  pad = "")
private

◆ spawnCrowdClone()

std::unique_ptr< OperatorEstBase > spawnCrowdClone ( ) const
overridevirtual

Implements OperatorEstBase.

Definition at line 180 of file OneBodyDensityMatrices.cpp.

References OperatorEstBase::data_, OperatorEstBase::data_locality_, qmcplusplus::queue, and qmcplusplus::rank.

181 {
182  std::size_t data_size = data_.size();
183  auto spawn_data_locality = data_locality_;
184 
186  {
187  // This is just a stub until a memory saving optimization is deemed necessary
188  spawn_data_locality = DataLocality::queue;
189  data_size = 0;
190  throw std::runtime_error("There is no memory savings implementation for OneBodyDensityMatrices");
191  }
192 
193  auto spawn = std::make_unique<OneBodyDensityMatrices>(*this, spawn_data_locality);
194  spawn->get_data().resize(data_size, 0.0);
195  return spawn;
196 }
DataLocality data_locality_
locality for accumulation of estimator data.

◆ startBlock()

void startBlock ( int  steps)
overridevirtual

Implements OperatorEstBase.

Definition at line 206 of file OneBodyDensityMatrices.cpp.

206 {}

◆ updateBasis()

void updateBasis ( const Position r,
ParticleSet pset_target 
)
inlineprivate

basis set updates

Definition at line 682 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::basis_functions_, OneBodyDensityMatrices::basis_norms_, OneBodyDensityMatrices::basis_size_, OneBodyDensityMatrices::basis_values_, CompositeSPOSet::evaluateValue(), and qmcplusplus::pset_target.

Referenced by OneBodyDensityMatrices::calcDensity(), OneBodyDensityMatrices::generateParticleBasis(), OneBodyDensityMatrices::generateSampleBasis(), and OneBodyDensityMatrices::normalizeBasis().

683 {
684  // This is ridiculous in the case of splines, still necessary for hybrid/LCAO
685  pset_target.makeMove(0, r - pset_target.R[0]);
687  pset_target.rejectMove(0);
688  for (int i = 0; i < basis_size_; ++i)
689  basis_values_[i] *= basis_norms_[i];
690 }
void evaluateValue(const ParticleSet &P, int iat, ValueVector &psi) override
evaluate the values of this single-particle orbital set

◆ updateBasisD012()

void updateBasisD012 ( const Position r,
ParticleSet pset_target 
)
inlineprivate

evaluates vgl on basis_functions_ for r sideeffects:

  • sets basis_values_, basis_gradients_, basis_laplacians_ all are normalized by basis norms_

Definition at line 702 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::basis_functions_, OneBodyDensityMatrices::basis_gradients_, OneBodyDensityMatrices::basis_laplacians_, OneBodyDensityMatrices::basis_norms_, OneBodyDensityMatrices::basis_size_, OneBodyDensityMatrices::basis_values_, CompositeSPOSet::evaluateVGL(), and qmcplusplus::pset_target.

Referenced by OneBodyDensityMatrices::calcDensityDrift().

703 {
704  pset_target.makeMove(0, r - pset_target.R[0]);
706  pset_target.rejectMove(0);
707  for (int i = 0; i < basis_size_; ++i)
708  basis_values_[i] *= basis_norms_[i];
709  for (int i = 0; i < basis_size_; ++i)
711  for (int i = 0; i < basis_size_; ++i)
713 }
void evaluateVGL(const ParticleSet &P, int iat, ValueVector &psi, GradVector &dpsi, ValueVector &d2psi) override
evaluate the values, gradients and laplacians of this single-particle orbital set ...

◆ updateBasisD012WithSpin()

void updateBasisD012WithSpin ( const Position r,
const Real s,
ParticleSet pset_target 
)
inlineprivate

same as above, but includes spin gradients

Definition at line 715 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::basis_functions_, OneBodyDensityMatrices::basis_gradients_, OneBodyDensityMatrices::basis_laplacians_, OneBodyDensityMatrices::basis_norms_, OneBodyDensityMatrices::basis_size_, OneBodyDensityMatrices::basis_spin_gradients_, OneBodyDensityMatrices::basis_values_, CompositeSPOSet::evaluateVGL_spin(), qmcplusplus::pset_target, and qmcplusplus::Units::time::s.

Referenced by OneBodyDensityMatrices::calcDensityDriftWithSpin().

716 {
717  pset_target.makeMoveWithSpin(0, r - pset_target.R[0], s - pset_target.spins[0]);
720  pset_target.rejectMove(0);
721  for (int i = 0; i < basis_size_; ++i)
722  {
723  basis_values_[i] *= basis_norms_[i];
727  }
728 }
void evaluateVGL_spin(const ParticleSet &P, int iat, ValueVector &psi, GradVector &dpsi, ValueVector &d2psi, ValueVector &dspin_psi) override
evaluate the values, gradients and laplacians and spin gradient of this single-particle orbital set ...

◆ updateBasisWithSpin()

void updateBasisWithSpin ( const Position r,
const Real s,
ParticleSet pset_target 
)
inlineprivate

basis set updates with spin

Definition at line 692 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::basis_functions_, OneBodyDensityMatrices::basis_norms_, OneBodyDensityMatrices::basis_size_, OneBodyDensityMatrices::basis_values_, CompositeSPOSet::evaluateValue(), qmcplusplus::pset_target, and qmcplusplus::Units::time::s.

Referenced by OneBodyDensityMatrices::calcDensityWithSpin(), OneBodyDensityMatrices::generateParticleBasis(), and OneBodyDensityMatrices::generateSampleBasis().

693 {
694  // This is ridiculous in the case of splines, still necessary for hybrid/LCAO
695  pset_target.makeMoveWithSpin(0, r - pset_target.R[0], s - pset_target.spins[0]);
697  pset_target.rejectMove(0);
698  for (int i = 0; i < basis_size_; ++i)
699  basis_values_[i] *= basis_norms_[i];
700 }
void evaluateValue(const ParticleSet &P, int iat, ValueVector &psi) override
evaluate the values of this single-particle orbital set

◆ warmupSampling()

void warmupSampling ( ParticleSet pset_target,
RandomBase< FullPrecReal > &  rng 
)
private

does some warmup sampling i.e.

samples but throws away the results Only when integrator_ = Integrator::DENSITY sets rpcur_ initial rpcur + one diffusion step sets initial rhocur_ and dpcur_ Then calls generateSamples with number of input warmup samples.

Definition at line 730 of file OneBodyDensityMatrices.cpp.

References OneBodyDensityMatrices::calcDensityDrift(), OneBodyDensityMatrices::calcDensityDriftWithSpin(), OneBodyDensityMatrices::center_, OneBodyDensityMatrices::diffuse(), OneBodyDensityMatrices::diffuseSpin(), OneBodyDensityMatrices::dpcur_, OneBodyDensityMatrices::dspcur_, OneBodyDensityMatrices::generateSamples(), OneBodyDensityMatricesInput::get_timestep(), OneBodyDensityMatricesInput::get_warmup_samples(), OneBodyDensityMatrices::input_, OneBodyDensityMatrices::is_spinor_, OneBodyDensityMatrices::METROPOLIS, qmcplusplus::pset_target, OneBodyDensityMatrices::rhocur_, OneBodyDensityMatrices::rpcur_, OneBodyDensityMatrices::sampling_, OneBodyDensityMatrices::spcur_, qmcplusplus::sqrt(), and OneBodyDensityMatrices::warmed_up_.

Referenced by OneBodyDensityMatrices::evaluateMatrix().

731 {
733  {
734  if (!warmed_up_)
735  {
737  rpcur_ += center_;
738  if (!is_spinor_)
740  else
741  {
744  }
745  }
747  warmed_up_ = true;
748  }
749 }
Position center_
If not defined in OBDMI taken from lattice_.
void calcDensityDrift(const Position &r, Real &dens, Position &drift, ParticleSet &pset_target)
calculate density and drift bashed on r
void generateSamples(const Real weight, ParticleSet &pset_target, RandomBase< FullPrecReal > &rng, int steps=0)
Dispatch method to difference methods of generating samples.
Position diffuse(const Real sqt, RandomBase< FullPrecReal > &rng)
produce a position difference vector from timestep
void calcDensityDriftWithSpin(const Position &r, const Real &s, Real &dens, Position &drift, Real &sdrift, ParticleSet &pset_target)
same as above, but with spin move
Sampling sampling_
Sampling method, this derived values in input_.
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
Real diffuseSpin(const Real sqt, RandomBase< FullPrecReal > &rng)
spin diffusion
Position rpcur_
current position – As long Positions are TinyVectors they are intialized to zero vectors ...

Friends And Related Function Documentation

◆ testing::OneBodyDensityMatricesTests

Definition at line 345 of file OneBodyDensityMatrices.h.

Member Data Documentation

◆ acceptance_ratio_

Real acceptance_ratio_ = 0.0
private

running acceptance ratio over all samples

Definition at line 135 of file OneBodyDensityMatrices.h.

Referenced by OneBodyDensityMatrices::generateDensitySamples(), and OneBodyDensityMatrices::generateDensitySamplesWithSpin().

◆ basis_functions_

◆ basis_gradients_

◆ basis_laplacians_

◆ basis_norms_

◆ basis_size_

◆ basis_spin_gradients_

◆ basis_values_

◆ center_

Position center_
private

If not defined in OBDMI taken from lattice_.

Definition at line 80 of file OneBodyDensityMatrices.h.

Referenced by OneBodyDensityMatrices::OneBodyDensityMatrices(), and OneBodyDensityMatrices::warmupSampling().

◆ dpcur_

◆ dspcur_

◆ ind_dims_

◆ input_

◆ is_spinor_

◆ lattice_

◆ metric_

◆ N_BB_

std::vector<Matrix<Value> > N_BB_
private

◆ naccepted_

int naccepted_ = 0
private

◆ nmoves_

◆ periodic_

bool periodic_
private

◆ Phi_MB_

Matrix<Value> Phi_MB_
private

basis_values_ at each r of rsamples_ row: sample col: basis_value size: samples * basis_size

Definition at line 124 of file OneBodyDensityMatrices.h.

Referenced by OneBodyDensityMatrices::evaluateMatrix(), and OneBodyDensityMatrices::OneBodyDensityMatrices().

◆ Phi_NB_

std::vector<Matrix<Value> > Phi_NB_
private

row major per sample workspaces

conj(basis_values) for each particle size: samples * basis_size vector is over species each matrix row: particle column: basis_value

Definition at line 113 of file OneBodyDensityMatrices.h.

Referenced by OneBodyDensityMatrices::evaluateMatrix(), and OneBodyDensityMatrices::OneBodyDensityMatrices().

◆ Phi_Psi_NB_

std::vector<Matrix<Value> > Phi_Psi_NB_
private

◆ Psi_NM_

std::vector<Matrix<Value> > Psi_NM_
private

ratio per particle per sample size: particles * samples vector is over species each matrix row: particle col: sample

Definition at line 119 of file OneBodyDensityMatrices.h.

Referenced by OneBodyDensityMatrices::evaluateMatrix(), and OneBodyDensityMatrices::OneBodyDensityMatrices().

◆ psi_ratios_

std::vector<Value> psi_ratios_
private

, ' .

per particle ratios size: particles

Definition at line 105 of file OneBodyDensityMatrices.h.

Referenced by OneBodyDensityMatrices::generateSampleRatios(), and OneBodyDensityMatrices::OneBodyDensityMatrices().

◆ rcorner_

◆ rhocur_

◆ rpcur_

Position rpcur_
private

current position – As long Positions are TinyVectors they are intialized to zero vectors

Definition at line 144 of file OneBodyDensityMatrices.h.

Referenced by OneBodyDensityMatrices::generateDensitySamples(), OneBodyDensityMatrices::generateDensitySamplesWithSpin(), and OneBodyDensityMatrices::warmupSampling().

◆ rsamples_

◆ samples_

◆ samples_weights_

◆ sampling_

◆ spcur_

◆ species_

◆ species_names_

std::vector<std::string> species_names_
private

◆ species_sizes_

◆ ssamples_

◆ timers_

◆ volume_

◆ warmed_up_

bool warmed_up_ = false
private

Definition at line 137 of file OneBodyDensityMatrices.h.

Referenced by OneBodyDensityMatrices::warmupSampling().


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