QMCPACK
SpinDensity Class Reference
+ Inheritance diagram for SpinDensity:
+ Collaboration diagram for SpinDensity:

Public Types

using Lattice_t = ParticleSet::ParticleLayout
 
using dens_t = std::vector< RealType >
 
using pts_t = std::vector< PosType >
 
- Public Types inherited from OperatorBase
enum  EnergyDomains { KINETIC = 0, POTENTIAL, NO_ENERGY_DOMAIN }
 enum to denote energy domain of operators More...
 
enum  QuantumDomains {
  NO_QUANTUM_DOMAIN = 0, CLASSICAL, QUANTUM, CLASSICAL_CLASSICAL,
  QUANTUM_CLASSICAL, QUANTUM_QUANTUM
}
 
enum  {
  PRIMARY = 0, OPTIMIZABLE = 1, RATIOUPDATE = 2, PHYSICAL = 3,
  COLLECTABLE = 4, NONLOCAL = 5
}
 enum for update_mode More...
 
using Return_t = FullPrecRealType
 type of return value of evaluate More...
 
using ValueMatrix = SPOSet::ValueMatrix
 For fast derivative evaluation. More...
 
using GradMatrix = SPOSet::GradMatrix
 
using BufferType = ParticleSet::Buffer_t
 typedef for the serialized buffer More...
 
using Walker_t = ParticleSet::Walker_t
 typedef for the walker More...
 
using ParticleScalar = ParticleSet::Scalar_t
 typedef for the ParticleScalar More...
 
using SPOMap = SPOSet::SPOMap
 typedef for SPOMap More...
 
- Public Types inherited from QMCTraits
enum  { DIM = OHMMS_DIM, DIM_VGL = OHMMS_DIM + 2 }
 
using QTBase = QMCTypes< OHMMS_PRECISION, DIM >
 
using QTFull = QMCTypes< OHMMS_PRECISION_FULL, DIM >
 
using RealType = QTBase::RealType
 
using ComplexType = QTBase::ComplexType
 
using ValueType = QTBase::ValueType
 
using PosType = QTBase::PosType
 
using GradType = QTBase::GradType
 
using TensorType = QTBase::TensorType
 
using IndexType = OHMMS_INDEXTYPE
 define other types More...
 
using FullPrecRealType = QTFull::RealType
 
using FullPrecValueType = QTFull::ValueType
 
using PropertySetType = RecordNamedProperty< FullPrecRealType >
 define PropertyList_t More...
 
using PtclGrpIndexes = std::vector< std::pair< int, int > >
 

Public Member Functions

 SpinDensity (ParticleSet &P)
 
 ~SpinDensity () override
 
std::string getClassName () const override
 return class name More...
 
std::unique_ptr< OperatorBasemakeClone (ParticleSet &P, TrialWaveFunction &psi) final
 
bool put (xmlNodePtr cur) override
 Read the input parameter. More...
 
Return_t evaluate (ParticleSet &P) override
 Evaluate the local energy contribution of this component. More...
 
void addObservables (PropertySetType &plist, BufferType &olist) override
 named values to the property list Default implementaton uses addValue(plist_) More...
 
void registerCollectables (std::vector< ObservableHelper > &h5desc, hdf_archive &file) const override
 
void resetTargetParticleSet (ParticleSet &P) override
 Reset the data with the target ParticleSet. More...
 
void setObservables (PropertySetType &plist) override
 Set the values evaluated by this object to plist Default implementation is to assign Value which is updated by evaluate function using my_index_. More...
 
void setParticlePropertyList (PropertySetType &plist, int offset) override
 
void checkout_scalar_arrays (TraceManager &tm)
 
void collect_scalar_samples ()
 
void delete_scalar_arrays ()
 
bool get (std::ostream &os) const override
 write about the class More...
 
void reset ()
 
void report (const std::string &pad)
 
void test (int moves, ParticleSet &P)
 
Return_t test_evaluate (ParticleSet &P, int &pmin, int &pmax)
 
- Public Member Functions inherited from OperatorBase
 OperatorBase ()
 Construct a new Operator Base object Default and unique empty constructor. More...
 
virtual ~OperatorBase ()=default
 
virtual bool dependsOnWaveFunction () const
 return true if this operator depends on a wavefunction More...
 
std::bitset< 8 > & getUpdateMode () noexcept
 get update_mode_ reference More...
 
Return_t getValue () const noexcept
 get a copy of value_ More...
 
std::string getName () const noexcept
 getter a copy of my_name_, rvalue small string optimization More...
 
void setName (const std::string name) noexcept
 Set my_name member, uses small string optimization (pass by value) More...
 
TraceRequestgetRequest () noexcept
 Get request_ member. More...
 
virtual void registerObservables (std::vector< ObservableHelper > &h5desc, hdf_archive &file) const
 add to observable descriptor for hdf5 The default implementation is to register a scalar for this->value_ More...
 
virtual void setHistories (Walker_t &ThisWalker)
 
virtual Return_t evaluateDeterministic (ParticleSet &P)
 Evaluate the local energy contribution of this component, deterministically based on current state. More...
 
virtual void mw_evaluate (const RefVectorWithLeader< OperatorBase > &o_list, const RefVectorWithLeader< TrialWaveFunction > &wf_list, const RefVectorWithLeader< ParticleSet > &p_list) const
 Evaluate the contribution of this component of multiple walkers. More...
 
virtual void mw_evaluatePerParticle (const RefVectorWithLeader< OperatorBase > &o_list, const RefVectorWithLeader< TrialWaveFunction > &wf_list, const RefVectorWithLeader< ParticleSet > &p_list, const std::vector< ListenerVector< RealType >> &listeners, const std::vector< ListenerVector< RealType >> &listeners_ions) const
 Evaluate the contribution of this component of multiple walkers per particle and report to registerd listeners from objects in Estimators. More...
 
virtual void mw_evaluateWithParameterDerivatives (const RefVectorWithLeader< OperatorBase > &o_list, const RefVectorWithLeader< ParticleSet > &p_list, const opt_variables_type &optvars, const RecordArray< ValueType > &dlogpsi, RecordArray< ValueType > &dhpsioverpsi) const
 TODO: add docs. More...
 
virtual Return_t rejectedMove (ParticleSet &P)
 TODO: add docs. More...
 
virtual Return_t evaluateWithToperator (ParticleSet &P)
 Evaluate the local energy contribution of this component with Toperators updated if requested. More...
 
virtual void mw_evaluateWithToperator (const RefVectorWithLeader< OperatorBase > &o_list, const RefVectorWithLeader< TrialWaveFunction > &wf_list, const RefVectorWithLeader< ParticleSet > &p_list) const
 Evaluate the contribution of this component of multiple walkers. More...
 
virtual void mw_evaluatePerParticleWithToperator (const RefVectorWithLeader< OperatorBase > &o_list, const RefVectorWithLeader< TrialWaveFunction > &wf_list, const RefVectorWithLeader< ParticleSet > &p_list, const std::vector< ListenerVector< RealType >> &listeners, const std::vector< ListenerVector< RealType >> &listeners_ions) const
 Evaluate the contribution of this component of multiple walkers per particle and report to registerd listeners from objects in Estimators. More...
 
virtual Return_t evaluateValueAndDerivatives (ParticleSet &P, const opt_variables_type &optvars, const Vector< ValueType > &dlogpsi, Vector< ValueType > &dhpsioverpsi)
 Evaluate value and derivatives wrt the optimizables. More...
 
virtual Return_t evaluateWithIonDerivs (ParticleSet &P, ParticleSet &ions, TrialWaveFunction &psi, ParticleSet::ParticlePos &hf_term, ParticleSet::ParticlePos &pulay_term)
 Evaluate contribution to local energy and derivatives w.r.t ionic coordinates from OperatorBase. More...
 
virtual Return_t evaluateWithIonDerivsDeterministic (ParticleSet &P, ParticleSet &ions, TrialWaveFunction &psi, ParticleSet::ParticlePos &hf_term, ParticleSet::ParticlePos &pulay_term)
 Evaluate contribution to local energy and derivatives w.r.t ionic coordinates from OperatorBase. More...
 
virtual void evaluateOneBodyOpMatrix (ParticleSet &P, const TWFFastDerivWrapper &psi, std::vector< ValueMatrix > &B)
 Evaluate "B" matrix for observable. More...
 
virtual void evaluateOneBodyOpMatrixForceDeriv (ParticleSet &P, ParticleSet &source, const TWFFastDerivWrapper &psi, const int iat, std::vector< std::vector< ValueMatrix >> &Bforce)
 Evaluate "dB/dR" matrices for observable. More...
 
virtual void updateSource (ParticleSet &s)
 Update data associated with a particleset. More...
 
virtual Return_t getEnsembleAverage ()
 Return an average value by collective operation. More...
 
virtual void createResource (ResourceCollection &collection) const
 Initialize a shared resource and hand it to a collection. More...
 
virtual void acquireResource (ResourceCollection &collection, const RefVectorWithLeader< OperatorBase > &o_list) const
 Acquire a shared resource from a collection. More...
 
virtual void releaseResource (ResourceCollection &collection, const RefVectorWithLeader< OperatorBase > &o_list) const
 Return a shared resource to a collection. More...
 
virtual void setRandomGenerator (RandomBase< FullPrecRealType > *rng)
 Set the Random Generator object TODO: add docs. More...
 
virtual void add2Hamiltonian (ParticleSet &qp, TrialWaveFunction &psi, QMCHamiltonian &targetH)
 TODO: add docs. More...
 
virtual void getRequiredTraces (TraceManager &tm)
 TODO: add docs. More...
 
virtual void informOfPerParticleListener ()
 
bool isClassical () const noexcept
 
bool isQuantum () const noexcept
 
bool isClassicalClassical () const noexcept
 
bool isQuantumClassical () const noexcept
 
bool isQuantumQuantum () const noexcept
 
bool getMode (const int i) const noexcept
 Return the mode i. More...
 
bool isNonLocal () const noexcept
 TODO: add docs. More...
 
bool hasListener () const noexcept
 
void contributeTraceQuantities ()
 Make trace quantities available. More...
 
void checkoutTraceQuantities (TraceManager &tm)
 Checkout trace arrays Derived classes must guard individual checkouts using request info. More...
 
void collectScalarTraces ()
 Collect scalar trace data. More...
 
void deleteTraceQuantities ()
 delete trace arrays More...
 

Public Attributes

ParticleSetPtmp
 
int nspecies
 
std::vector< int > species_size
 
std::vector< std::string > species_name
 
Lattice_t cell
 
PosType corner
 
TinyVector< int, DIMgrid
 
TinyVector< int, DIMgdims
 
int npoints
 

Additional Inherited Members

- Protected Member Functions inherited from OperatorBase
virtual void contributeScalarQuantities ()
 
virtual void checkoutScalarQuantities (TraceManager &tm)
 
virtual void collectScalarQuantities ()
 
virtual void deleteScalarQuantities ()
 
virtual void contributeParticleQuantities ()
 
virtual void checkoutParticleQuantities (TraceManager &tm)
 
virtual void deleteParticleQuantities ()
 
virtual void setComputeForces (bool compute)
 
void setEnergyDomain (EnergyDomains edomain)
 Set the Energy Domain. More...
 
void setQuantumDomain (QuantumDomains qdomain)
 set quantum domain More...
 
void oneBodyQuantumDomain (const ParticleSet &P)
 set quantum domain for one-body operator More...
 
void twoBodyQuantumDomain (const ParticleSet &P)
 set quantum domain for two-body operator More...
 
void twoBodyQuantumDomain (const ParticleSet &P1, const ParticleSet &P2)
 set quantum domain for two-body operator More...
 
void addValue (PropertySetType &plist)
 named values to the property list More...
 
- Protected Attributes inherited from OperatorBase
std::bitset< 8 > update_mode_
 set the current update mode More...
 
Return_t value_
 current value More...
 
std::string name_
 name of this object More...
 
TraceRequest request_
 whether traces are being collected More...
 
int my_index_
 starting index of this object More...
 
Return_t new_value_
 a new value for a proposed move More...
 
Walker_tt_walker_
 reference to the current walker More...
 
bool streaming_particles_
 
bool have_required_traces_
 

Detailed Description

Definition at line 22 of file SpinDensity.h.

Member Typedef Documentation

◆ dens_t

using dens_t = std::vector<RealType>

Definition at line 26 of file SpinDensity.h.

◆ Lattice_t

Definition at line 25 of file SpinDensity.h.

◆ pts_t

using pts_t = std::vector<PosType>

Definition at line 27 of file SpinDensity.h.

Constructor & Destructor Documentation

◆ SpinDensity()

Definition at line 21 of file SpinDensity.cpp.

References APP_ABORT, ParticleSet::getLattice(), ParticleSet::getSpeciesSet(), ParticleSet::groupsize(), SpinDensity::nspecies, SpinDensity::Ptmp, SpinDensity::reset(), qmcplusplus::Units::time::s, SpeciesSet::size(), SpinDensity::species_name, SpinDensity::species_size, SpeciesSet::speciesName, and qmcplusplus::SUPERCELL_OPEN.

22 {
23  // get particle information
24  SpeciesSet& species = P.getSpeciesSet();
25  nspecies = species.size();
26  for (int s = 0; s < nspecies; ++s)
27  species_size.push_back(P.groupsize(s));
28  for (int s = 0; s < nspecies; ++s)
29  species_name.push_back(species.speciesName[s]);
30  reset();
31 
32  //jtk: spin density only works for periodic bc's for now
33  // abort if using open boundary conditions
34  bool open_bcs = (P.getLattice().SuperCellEnum == SUPERCELL_OPEN);
35  if (open_bcs)
36  {
37  APP_ABORT("SpinDensity is not implemented for open boundary conditions at present\n please contact the developers "
38  "if you need this feature");
39  }
40 
41  Ptmp = &P;
42 }
std::vector< int > species_size
Definition: SpinDensity.h:33
#define APP_ABORT(msg)
Widely used but deprecated fatal error macros from legacy code.
Definition: AppAbort.h:27
std::vector< std::string > species_name
Definition: SpinDensity.h:34

◆ ~SpinDensity()

~SpinDensity ( )
inlineoverride

Definition at line 43 of file SpinDensity.h.

43 {}

Member Function Documentation

◆ addObservables()

void addObservables ( PropertySetType plist,
BufferType collectables 
)
overridevirtual

named values to the property list Default implementaton uses addValue(plist_)

Parameters
plistRecordNameProperty
collectablesObservables that are accumulated by evaluate

Reimplemented from OperatorBase.

Definition at line 182 of file SpinDensity.cpp.

References PooledData< T >::add(), PooledData< T >::current(), OperatorBase::my_index_, SpinDensity::npoints, and SpinDensity::nspecies.

183 {
184  my_index_ = collectables.current();
185  std::vector<RealType> tmp(nspecies * npoints);
186  collectables.add(tmp.begin(), tmp.end());
187 }
int my_index_
starting index of this object
Definition: OperatorBase.h:535

◆ checkout_scalar_arrays()

void checkout_scalar_arrays ( TraceManager tm)
inline

Definition at line 60 of file SpinDensity.h.

60 {}

◆ collect_scalar_samples()

void collect_scalar_samples ( )
inline

Definition at line 61 of file SpinDensity.h.

61 {}

◆ delete_scalar_arrays()

void delete_scalar_arrays ( )
inline

Definition at line 62 of file SpinDensity.h.

62 {}

◆ evaluate()

SpinDensity::Return_t evaluate ( ParticleSet P)
overridevirtual

Evaluate the local energy contribution of this component.

Parameters
Pinput configuration containing N particles
Returns
the value of the Hamiltonian component

Implements OperatorBase.

Definition at line 205 of file SpinDensity.cpp.

References SpinDensity::cell, ParticleSet::Collectables, SpinDensity::corner, QMCTraits::DIM, qmcplusplus::floor(), SpinDensity::gdims, SpinDensity::grid, OperatorBase::my_index_, SpinDensity::npoints, SpinDensity::nspecies, qmcplusplus::Units::time::ps, ParticleSet::R, qmcplusplus::Units::time::s, SpinDensity::species_size, OperatorBase::t_walker_, CrystalLattice< T, D >::toUnit(), and Walker< t_traits, p_traits >::Weight.

206 {
208  int p = 0;
209  int offset = my_index_;
210  for (int s = 0; s < nspecies; ++s, offset += npoints)
211  for (int ps = 0; ps < species_size[s]; ++ps, ++p)
212  {
213  PosType u = cell.toUnit(P.R[p] - corner);
214  //bool inside = true;
215  //for(int d=0;d<DIM;++d)
216  // inside &= u[d]>0.0 && u[d]<1.0;
217  //if(inside)
218  //{
219  int point = offset;
220  for (int d = 0; d < DIM; ++d)
221  point += gdims[d] * ((int)(grid[d] * (u[d] - std::floor(u[d])))); //periodic only
222  P.Collectables[point] += w;
223  //}
224  }
225  return 0.0;
226 }
std::vector< int > species_size
Definition: SpinDensity.h:33
int my_index_
starting index of this object
Definition: OperatorBase.h:535
TinyVector< int, DIM > gdims
Definition: SpinDensity.h:38
QMCTraits::PosType PosType
SingleParticlePos toUnit(const TinyVector< T1, D > &r) const
Convert a cartesian vector to a unit vector.
Walker_t * t_walker_
reference to the current walker
Definition: OperatorBase.h:541
QMCTraits::RealType RealType
FullPrecRealType Weight
Weight of the walker.
Definition: Walker.h:102
MakeReturn< UnaryNode< FnFloor, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t floor(const Vector< T1, C1 > &l)
TinyVector< int, DIM > grid
Definition: SpinDensity.h:37

◆ get()

bool get ( std::ostream &  os) const
inlineoverridevirtual

write about the class

Implements OperatorBase.

Definition at line 66 of file SpinDensity.h.

66 { return false; }

◆ getClassName()

std::string getClassName ( ) const
inlineoverridevirtual

return class name

Implements OperatorBase.

Definition at line 46 of file SpinDensity.h.

46 { return "SpinDensity"; }

◆ makeClone()

std::unique_ptr< OperatorBase > makeClone ( ParticleSet P,
TrialWaveFunction psi 
)
finalvirtual

Implements OperatorBase.

Definition at line 53 of file SpinDensity.cpp.

54 {
55  return std::make_unique<SpinDensity>(*this);
56 }

◆ put()

bool put ( xmlNodePtr  cur)
overridevirtual

Read the input parameter.

Parameters
curxml node for a OperatorBase object

Implements OperatorBase.

Definition at line 59 of file SpinDensity.cpp.

References OhmmsAttributeSet::add(), APP_ABORT, qmcplusplus::ceil(), SpinDensity::cell, CrystalLattice< T, D >::Center, SpinDensity::corner, QMCTraits::DIM, qmcplusplus::dot(), SpinDensity::gdims, ParticleSet::getLattice(), getXMLAttributeValue(), SpinDensity::grid, OperatorBase::name_, SpinDensity::npoints, SpinDensity::Ptmp, OhmmsAttributeSet::put(), putContent(), SpinDensity::report(), SpinDensity::reset(), CrystalLattice< T, D >::Rv, CrystalLattice< T, D >::set(), qmcplusplus::sqrt(), and SpinDensity::test().

60 {
61  using std::ceil;
62  using std::sqrt;
63  reset();
64  std::string write_report = "no";
65  OhmmsAttributeSet attrib;
66  attrib.add(name_, "name");
67  attrib.add(write_report, "report");
68  attrib.put(cur);
69 
70  bool have_dr = false;
71  bool have_grid = false;
72  bool have_center = false;
73  bool have_corner = false;
74  bool have_cell = false;
75 
76  PosType dr;
77  PosType center;
78  Tensor<RealType, DIM> axes;
79 
80  int test_moves = 0;
81 
82  xmlNodePtr element = cur->xmlChildrenNode;
83  while (element != NULL)
84  {
85  std::string ename((const char*)element->name);
86  if (ename == "parameter")
87  {
88  const std::string name(getXMLAttributeValue(element, "name"));
89  if (name == "dr")
90  {
91  have_dr = true;
92  putContent(dr, element);
93  }
94  else if (name == "grid")
95  {
96  have_grid = true;
97  putContent(grid, element);
98  }
99  else if (name == "corner")
100  {
101  have_corner = true;
102  putContent(corner, element);
103  }
104  else if (name == "center")
105  {
106  have_center = true;
107  putContent(center, element);
108  }
109  else if (name == "cell")
110  {
111  have_cell = true;
112  putContent(axes, element);
113  }
114  else if (name == "test_moves")
115  putContent(test_moves, element);
116  }
117  element = element->next;
118  }
119 
120  if (have_dr && have_grid)
121  {
122  APP_ABORT("SpinDensity::put dr and grid are provided, this is ambiguous");
123  }
124  else if (!have_dr && !have_grid)
125  APP_ABORT("SpinDensity::put must provide dr or grid");
126 
127  if (have_corner && have_center)
128  APP_ABORT("SpinDensity::put corner and center are provided, this is ambiguous");
129  if (have_cell)
130  {
131  cell.set(axes);
132  if (!have_corner && !have_center)
133  APP_ABORT("SpinDensity::put must provide corner or center");
134  }
135  else
136  cell = Ptmp->getLattice();
137 
138  if (have_center)
139  corner = center - cell.Center;
140 
141  if (have_dr)
142  for (int d = 0; d < DIM; ++d)
143  grid[d] = (int)ceil(sqrt(dot(cell.Rv[d], cell.Rv[d])) / dr[d]);
144 
145  npoints = 1;
146  for (int d = 0; d < DIM; ++d)
147  npoints *= grid[d];
148  gdims[0] = npoints / grid[0];
149  for (int d = 1; d < DIM; ++d)
150  gdims[d] = gdims[d - 1] / grid[d];
151 
152  if (write_report == "yes")
153  report(" ");
154  if (test_moves > 0)
155  test(test_moves, *Ptmp);
156 
157  return true;
158 }
void set(const Tensor< TT, D > &lat)
set the lattice vector from the command-line options
SingleParticlePos Center
Center of the cell sum(Rv[i])/2.
bool put(xmlNodePtr cur)
assign attributes to the set
Definition: AttributeSet.h:55
TinyVector< SingleParticlePos, D > Rv
Real-space unit vectors.
TinyVector< int, DIM > gdims
Definition: SpinDensity.h:38
QMCTraits::PosType PosType
std::string name_
name of this object
Definition: OperatorBase.h:527
MakeReturn< UnaryNode< FnCeil, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t ceil(const Vector< T1, C1 > &l)
class to handle a set of attributes of an xmlNode
Definition: AttributeSet.h:24
void report(const std::string &pad)
#define APP_ABORT(msg)
Widely used but deprecated fatal error macros from legacy code.
Definition: AppAbort.h:27
std::string getXMLAttributeValue(const xmlNodePtr cur, const std::string_view name)
get the value string for attribute name if name is unfound in cur you get an empty string back this i...
bool putContent(T &a, xmlNodePtr cur)
replaces a&#39;s value with the first "element" in the "string" returned by XMLNodeString{cur}.
Definition: libxmldefs.h:88
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
void test(int moves, ParticleSet &P)
const auto & getLattice() const
Definition: ParticleSet.h:251
Tensor< typename BinaryReturn< T1, T2, OpMultiply >::Type_t, D > dot(const AntiSymTensor< T1, D > &lhs, const AntiSymTensor< T2, D > &rhs)
void add(PDT &aparam, const std::string &aname, std::vector< PDT > candidate_values={}, TagStatus status=TagStatus::OPTIONAL)
add a new attribute
Definition: AttributeSet.h:42
TinyVector< int, DIM > grid
Definition: SpinDensity.h:37

◆ registerCollectables()

void registerCollectables ( std::vector< ObservableHelper > &  h5desc,
hdf_archive file 
) const
overridevirtual

Reimplemented from OperatorBase.

Definition at line 190 of file SpinDensity.cpp.

References OperatorBase::my_index_, OperatorBase::name_, SpinDensity::npoints, SpinDensity::nspecies, qmcplusplus::oh, qmcplusplus::Units::time::s, ObservableHelper::set_dimensions(), and SpinDensity::species_name.

191 {
192  std::vector<int> ng(1);
193  ng[0] = npoints;
194 
195  hdf_path hdf_name{name_};
196  for (int s = 0; s < nspecies; ++s)
197  {
198  h5desc.emplace_back(hdf_name / species_name[s]);
199  auto& oh = h5desc.back();
201  }
202 }
int my_index_
starting index of this object
Definition: OperatorBase.h:535
std::string name_
name of this object
Definition: OperatorBase.h:527
void set_dimensions(const std::vector< int > &dims, int first)
set the shape of this observable
std::vector< std::string > species_name
Definition: SpinDensity.h:34
ObservableHelper oh

◆ report()

void report ( const std::string &  pad)

Definition at line 161 of file SpinDensity.cpp.

References qmcplusplus::app_log(), SpinDensity::cell, CrystalLattice< T, D >::Center, SpinDensity::corner, QMCTraits::DIM, SpinDensity::gdims, SpinDensity::grid, SpinDensity::npoints, SpinDensity::nspecies, CrystalLattice< T, D >::Rv, qmcplusplus::Units::time::s, SpinDensity::species_name, and SpinDensity::species_size.

Referenced by SpinDensity::put().

162 {
163  app_log() << pad << "SpinDensity report" << std::endl;
164  app_log() << pad << " dim = " << DIM << std::endl;
165  app_log() << pad << " npoints = " << npoints << std::endl;
166  app_log() << pad << " grid = " << grid << std::endl;
167  app_log() << pad << " gdims = " << gdims << std::endl;
168  app_log() << pad << " corner = " << corner << std::endl;
169  app_log() << pad << " center = " << corner + cell.Center << std::endl;
170  app_log() << pad << " cell " << std::endl;
171  for (int d = 0; d < DIM; ++d)
172  app_log() << pad << " " << d << " " << cell.Rv[d] << std::endl;
173  app_log() << pad << " end cell " << std::endl;
174  app_log() << pad << " nspecies = " << nspecies << std::endl;
175  for (int s = 0; s < nspecies; ++s)
176  app_log() << pad << " species[" << s << "]"
177  << " = " << species_name[s] << " " << species_size[s] << std::endl;
178  app_log() << pad << "end SpinDensity report" << std::endl;
179 }
SingleParticlePos Center
Center of the cell sum(Rv[i])/2.
std::vector< int > species_size
Definition: SpinDensity.h:33
std::ostream & app_log()
Definition: OutputManager.h:65
TinyVector< SingleParticlePos, D > Rv
Real-space unit vectors.
TinyVector< int, DIM > gdims
Definition: SpinDensity.h:38
std::vector< std::string > species_name
Definition: SpinDensity.h:34
TinyVector< int, DIM > grid
Definition: SpinDensity.h:37

◆ reset()

void reset ( )

Definition at line 45 of file SpinDensity.cpp.

References OperatorBase::COLLECTABLE, SpinDensity::corner, OperatorBase::name_, and OperatorBase::update_mode_.

Referenced by SpinDensity::put(), and SpinDensity::SpinDensity().

46 {
47  name_ = "SpinDensity";
48  update_mode_.set(COLLECTABLE, 1);
49  corner = 0.0;
50 }
std::string name_
name of this object
Definition: OperatorBase.h:527
std::bitset< 8 > update_mode_
set the current update mode
Definition: OperatorBase.h:521

◆ resetTargetParticleSet()

void resetTargetParticleSet ( ParticleSet P)
inlineoverridevirtual

Reset the data with the target ParticleSet.

Parameters
Pnew target ParticleSet

Implements OperatorBase.

Definition at line 56 of file SpinDensity.h.

56 {}

◆ setObservables()

void setObservables ( PropertySetType plist)
inlineoverridevirtual

Set the values evaluated by this object to plist Default implementation is to assign Value which is updated by evaluate function using my_index_.

Parameters
plistRecordNameProperty

Reimplemented from OperatorBase.

Definition at line 57 of file SpinDensity.h.

57 {}

◆ setParticlePropertyList()

void setParticlePropertyList ( PropertySetType plist,
int  offset 
)
inlineoverridevirtual

Reimplemented from OperatorBase.

Definition at line 58 of file SpinDensity.h.

58 {}

◆ test()

void test ( int  moves,
ParticleSet P 
)

Definition at line 229 of file SpinDensity.cpp.

References APP_ABORT, qmcplusplus::app_log(), QMCTraits::DIM, ParticleSet::getLattice(), ParticleSet::getTotalNum(), qmcplusplus::Units::distance::m, omptarget::min(), ParticleSet::R, and SpinDensity::test_evaluate().

Referenced by SpinDensity::put().

230 {
231  app_log() << " SpinDensity test" << std::endl;
232  RandomGenerator rng;
233  int particles = P.getTotalNum();
234  int pmin = std::numeric_limits<int>::max();
235  int pmax = std::numeric_limits<int>::min();
236  for (int m = 0; m < moves; ++m)
237  {
238  for (int p = 0; p < particles; ++p)
239  {
240  PosType u;
241  for (int d = 0; d < DIM; ++d)
242  u[d] = rng();
243  P.R[p] = P.getLattice().toCart(u);
244  }
245  test_evaluate(P, pmin, pmax);
246  }
247  app_log() << " end SpinDensity test" << std::endl;
248  APP_ABORT("SpinDensity::test test complete");
249 }
std::ostream & app_log()
Definition: OutputManager.h:65
StdRandom< OHMMS_PRECISION_FULL > RandomGenerator
Return_t test_evaluate(ParticleSet &P, int &pmin, int &pmax)
T min(T a, T b)
QMCTraits::PosType PosType
#define APP_ABORT(msg)
Widely used but deprecated fatal error macros from legacy code.
Definition: AppAbort.h:27

◆ test_evaluate()

SpinDensity::Return_t test_evaluate ( ParticleSet P,
int &  pmin,
int &  pmax 
)

Definition at line 252 of file SpinDensity.cpp.

References qmcplusplus::app_log(), SpinDensity::cell, SpinDensity::corner, QMCTraits::DIM, SpinDensity::gdims, SpinDensity::grid, omptarget::min(), SpinDensity::npoints, SpinDensity::nspecies, qmcplusplus::Units::time::ps, ParticleSet::R, qmcplusplus::Units::time::s, SpinDensity::species_size, and CrystalLattice< T, D >::toUnit().

Referenced by SpinDensity::test().

253 {
254  int p = 0;
255  int offset = 0;
256  for (int s = 0; s < nspecies; ++s, offset += npoints)
257  for (int ps = 0; ps < species_size[s]; ++ps, ++p)
258  {
259  PosType u = cell.toUnit(P.R[p] - corner);
260  bool inside = true;
261  for (int d = 0; d < DIM; ++d)
262  inside &= u[d] > 0.0 && u[d] < 1.0;
263  if (inside)
264  {
265  int point = offset;
266  for (int d = 0; d < DIM; ++d)
267  point += gdims[d] * ((int)(u[d] * grid[d]));
268  pmin = std::min(pmin, point - offset);
269  pmax = std::max(pmax, point - offset);
270  }
271  }
272  app_log() << " pmin = " << pmin << " pmax = " << pmax << " npoints = " << npoints << std::endl;
273  return 0.0;
274 }
std::vector< int > species_size
Definition: SpinDensity.h:33
std::ostream & app_log()
Definition: OutputManager.h:65
T min(T a, T b)
TinyVector< int, DIM > gdims
Definition: SpinDensity.h:38
QMCTraits::PosType PosType
SingleParticlePos toUnit(const TinyVector< T1, D > &r) const
Convert a cartesian vector to a unit vector.
TinyVector< int, DIM > grid
Definition: SpinDensity.h:37

Member Data Documentation

◆ cell

◆ corner

◆ gdims

◆ grid

◆ npoints

◆ nspecies

◆ Ptmp

ParticleSet* Ptmp

Definition at line 29 of file SpinDensity.h.

Referenced by SpinDensity::put(), and SpinDensity::SpinDensity().

◆ species_name

std::vector<std::string> species_name

◆ species_size

std::vector<int> species_size

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