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 29 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.

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

◆ 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 173 of file OneBodyDensityMatrices.cpp.

References OperatorEstBase::data_locality_.

174  : OneBodyDensityMatrices(obdm)
175 {
176  data_locality_ = dl;
177 }
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 528 of file OneBodyDensityMatrices.cpp.

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

533 {
534  implAccumulate(walkers, psets, wfns, rng);
535 }
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 456 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().

457 {
459  dens = 0.0;
460  for (int i = 0; i < basis_size_; ++i)
461  {
462  Value b = basis_values_[i];
463  dens += std::abs(qmcplusplus::conj(b) * b);
464  }
465  dens /= basis_size_;
466 }
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
MakeReturn< UnaryNode< FnFabs, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t abs(const Vector< T1, C1 > &l)

◆ 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 483 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().

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

◆ 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 501 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().

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

◆ 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 468 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().

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

◆ calcFullDataSize()

size_t calcFullDataSize ( size_t  basis_size,
int  num_species 
)
private

Definition at line 197 of file OneBodyDensityMatrices.cpp.

Referenced by OneBodyDensityMatrices::OneBodyDensityMatrices().

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

◆ diffuse()

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

produce a position difference vector from timestep

Definition at line 440 of file OneBodyDensityMatrices.cpp.

References qmcplusplus::assignGaussRand(), and OHMMS_DIM.

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

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

◆ diffuseSpin()

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

spin diffusion

Definition at line 448 of file OneBodyDensityMatrices.cpp.

References qmcplusplus::assignGaussRand().

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

449 {
450  Real diff;
451  assignGaussRand(&diff, 1, rng);
452  diff *= sqt;
453  return diff;
454 }
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 549 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().

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

◆ 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 310 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().

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

◆ 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 369 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().

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

◆ 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 603 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().

604 {
606  int p = 0;
607  for (int s = 0; s < species_.size(); ++s)
608  {
609  int nb = 0;
610  Matrix<Value>& P_nb = phi_nb[s];
611  for (int n = 0; n < species_sizes_[s]; ++n, ++p)
612  {
613  if (is_spinor_)
615  else
617  for (int b = 0; b < basis_size_; ++b, ++nb)
618  P_nb(nb) = qmcplusplus::conj(basis_values_[b]);
619  }
620  }
621 }
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
float conj(const float &c)
Workaround to allow conj on scalar to return real instead of complex.
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
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 623 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().

626 {
628  int mb = 0;
629  for (int m = 0; m < samples_; ++m)
630  {
631  if (is_spinor_)
633  else
635  for (int b = 0; b < basis_size_; ++b, ++mb)
636  Phi_mb(mb) = basis_values_[b];
637  }
638 }
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 640 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().

643 {
645  for (int m = 0; m < samples_; ++m)
646  {
647  // get N ratios for the current sample point
648  if (!is_spinor_)
649  {
650  pset_target.makeVirtualMoves(rsamples_[m]);
651  psi_target.evaluateRatiosAlltoOne(pset_target, psi_ratios_);
652  }
653  else
654  {
655  //note: makeVirtualMoves updates distance tables
656  //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
657  int p = 0;
658  for (int s = 0; s < species_.size(); ++s)
659  for (int n = 0; n < species_sizes_[s]; ++n, ++p)
660  {
661  pset_target.makeMoveWithSpin(p, rsamples_[m] - pset_target.R[p], ssamples_[m] - pset_target.spins[p]);
662  psi_ratios_[p] = psi_target.calcRatio(pset_target, p);
663  pset_target.rejectMove(p);
664  psi_target.resetPhaseDiff();
665  }
666  }
667 
668  // collect ratios into per-species matrices
669  int p = 0;
670  for (int s = 0; s < species_.size(); ++s)
671  {
672  Matrix<Value>& P_nm = psi_nm[s];
673  for (int n = 0; n < species_sizes_[s]; ++n, ++p)
674  {
675  P_nm(n, m) = qmcplusplus::conj(psi_ratios_[p]);
676  }
677  }
678  }
679 }
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
float conj(const float &c)
Workaround to allow conj on scalar to return real instead of complex.
int size() const
return the number of species
Definition: SpeciesSet.h:53
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 207 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().

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

◆ generateUniformGrid()

void generateUniformGrid ( RandomBase< FullPrecReal > &  rng)
inlineprivate

Definition at line 278 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().

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

◆ generateUniformSamples()

void generateUniformSamples ( RandomBase< FullPrecReal > &  rng)
inlineprivate

Definition at line 299 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().

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

◆ 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 537 of file OneBodyDensityMatrices.cpp.

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

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

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

◆ normalizeBasis()

void normalizeBasis ( ParticleSet pset_target)
inlineprivate

Definition at line 750 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().

751 {
752  int ngrid = std::max(200, input_.get_points());
753  int ngtot = pow(ngrid, OHMMS_DIM);
754  Real du = input_.get_scale() / ngrid;
755  Real dV = volume_ / ngtot;
756  Position rp;
757  Vector<Value> bnorms;
758  int gdims[OHMMS_DIM];
759  gdims[0] = pow(ngrid, OHMMS_DIM - 1);
760  for (int d = 1; d < OHMMS_DIM; ++d)
761  gdims[d] = gdims[d - 1] / ngrid;
762  bnorms.resize(basis_size_);
763  for (int i = 0; i < basis_size_; ++i)
764  bnorms[i] = 0.0;
765  std::fill(basis_norms_.begin(), basis_norms_.end(), 1.0);
766  for (int p = 0; p < ngtot; ++p)
767  {
768  int nrem = p;
769  for (int d = 0; d < OHMMS_DIM - 1; ++d)
770  {
771  int ind = nrem / gdims[d];
772  rp[d] = ind * du + du / 2;
773  nrem -= ind * gdims[d];
774  }
775  rp[OHMMS_DIM - 1] = nrem * du + du / 2;
776  rp = lattice_.toCart(rp) + rcorner_;
778  for (int i = 0; i < basis_size_; ++i)
779  bnorms[i] += qmcplusplus::conj(basis_values_[i]) * basis_values_[i] * dV;
780  }
781  for (int i = 0; i < basis_size_; ++i)
782  basis_norms_[i] = 1.0 / std::sqrt(real(bnorms[i]));
783 }
SingleParticlePos toCart(const TinyVector< T1, D > &c) const
Convert a unit vector to a cartesian vector.
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
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
#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 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_;

◆ 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 785 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().

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

◆ report()

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

◆ spawnCrowdClone()

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

Implements OperatorEstBase.

Definition at line 179 of file OneBodyDensityMatrices.cpp.

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

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

◆ startBlock()

void startBlock ( int  steps)
overridevirtual

Implements OperatorEstBase.

Definition at line 205 of file OneBodyDensityMatrices.cpp.

205 {}

◆ updateBasis()

void updateBasis ( const Position r,
ParticleSet pset_target 
)
inlineprivate

basis set updates

Definition at line 681 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().

682 {
683  // This is ridiculous in the case of splines, still necessary for hybrid/LCAO
684  pset_target.makeMove(0, r - pset_target.R[0]);
686  pset_target.rejectMove(0);
687  for (int i = 0; i < basis_size_; ++i)
688  basis_values_[i] *= basis_norms_[i];
689 }
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 701 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().

702 {
703  pset_target.makeMove(0, r - pset_target.R[0]);
705  pset_target.rejectMove(0);
706  for (int i = 0; i < basis_size_; ++i)
707  basis_values_[i] *= basis_norms_[i];
708  for (int i = 0; i < basis_size_; ++i)
710  for (int i = 0; i < basis_size_; ++i)
712 }
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 714 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().

715 {
716  pset_target.makeMoveWithSpin(0, r - pset_target.R[0], s - pset_target.spins[0]);
719  pset_target.rejectMove(0);
720  for (int i = 0; i < basis_size_; ++i)
721  {
722  basis_values_[i] *= basis_norms_[i];
726  }
727 }
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 691 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().

692 {
693  // This is ridiculous in the case of splines, still necessary for hybrid/LCAO
694  pset_target.makeMoveWithSpin(0, r - pset_target.R[0], s - pset_target.spins[0]);
696  pset_target.rejectMove(0);
697  for (int i = 0; i < basis_size_; ++i)
698  basis_values_[i] *= basis_norms_[i];
699 }
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 729 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().

730 {
732  {
733  if (!warmed_up_)
734  {
736  rpcur_ += center_;
737  if (!is_spinor_)
739  else
740  {
743  }
744  }
746  warmed_up_ = true;
747  }
748 }
Position center_
If not defined in OBDMI taken from lattice_.
Position diffuse(const Real sqt, RandomBase< FullPrecReal > &rng)
produce a position difference vector from timestep
void calcDensityDrift(const Position &r, Real &dens, Position &drift, ParticleSet &pset_target)
calculate density and drift bashed on r
Sampling sampling_
Sampling method, this derived values in input_.
Real diffuseSpin(const Real sqt, RandomBase< FullPrecReal > &rng)
spin diffusion
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(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
Position rpcur_
current position – As long Positions are TinyVectors they are intialized to zero vectors ...
void generateSamples(const Real weight, ParticleSet &pset_target, RandomBase< FullPrecReal > &rng, int steps=0)
Dispatch method to difference methods of generating samples.

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: