QMCPACK
EnergyDensityEstimator.cpp
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////////////
2 // This file is distributed under the University of Illinois/NCSA Open Source License.
3 // See LICENSE file in top directory for details.
4 //
5 // Copyright (c) 2016 Jeongnim Kim and QMCPACK developers.
6 //
7 // File developed by: Ye Luo, yeluo@anl.gov, Argonne National Laboratory
8 // Jaron T. Krogel, krogeljt@ornl.gov, Oak Ridge National Laboratory
9 // Mark A. Berrill, berrillma@ornl.gov, Oak Ridge National Laboratory
10 //
11 // File created by: Ye Luo, yeluo@anl.gov, Argonne National Laboratory
12 //////////////////////////////////////////////////////////////////////////////////////
13 
14 
15 #include "EnergyDensityEstimator.h"
16 #include "OhmmsData/AttributeSet.h"
18 #include "Particle/DistanceTable.h"
20 #include "Utilities/string_utils.h"
21 #include <string>
22 #include <vector>
23 
24 namespace qmcplusplus
25 {
26 EnergyDensityEstimator::EnergyDensityEstimator(const PSPool& PSP, const std::string& defaultKE)
27  : psetpool(PSP), Pdynamic(0), Pstatic(0), w_trace(0), Td_trace(0), Vd_trace(0), Vs_trace(0)
28 {
29  update_mode_.set(COLLECTABLE, 1);
30  defKE = defaultKE;
31  nsamples = 0;
32  ion_points = false;
33  nions = -1;
34  request_.request_scalar("weight");
35  request_.request_array("Kinetic");
36  request_.request_array("LocalPotential");
37 }
38 
39 
41 
42 
43 bool EnergyDensityEstimator::put(xmlNodePtr cur, ParticleSet& Pdyn)
44 {
45  Pdynamic = &Pdyn;
46  return put(cur);
47 }
48 
49 
50 /** check xml elements
51  */
52 bool EnergyDensityEstimator::put(xmlNodePtr cur)
53 {
54  input_xml = cur;
55  //initialize simple xml attributes
56  name_ = "EnergyDensity";
57  std::string dyn, stat = "";
58  ion_points = false;
59  OhmmsAttributeSet attrib;
60  attrib.add(name_, "name");
61  attrib.add(dyn, "dynamic");
62  attrib.add(stat, "static");
63  attrib.add(ion_points, "ion_points");
64  attrib.put(cur);
65  //collect particle sets
66  if (!Pdynamic)
68  if (Pdynamic->hasSK())
71  std::vector<ParticleSet*> Pref;
72  if (stat == "")
73  {
74  Pstatic = 0;
75  dtable_index = -1;
76  }
77  else
78  {
79  Pstatic = get_particleset(stat);
80  if (Pstatic->hasSK())
83  Pref.resize(1);
84  Pref[0] = Pstatic;
85  if (!ion_points)
87  else
89  }
90  //size arrays
93  if (ion_points)
94  {
96  Rion.resize(nions, DIM);
97  for (int i = 0; i < nions; i++)
98  for (int d = 0; d < DIM; d++)
99  Rion(i, d) = Pstatic->R[i][d];
100  }
102  fill(particles_outside.begin(), particles_outside.end(), true);
103  //read xml element contents
104  xmlNodePtr element;
105  bool stop = false;
106  //initialize reference points
107  app_log() << "Initializing reference points" << std::endl;
108  bool has_ref = false;
109  element = cur->children;
110  while (element != NULL)
111  {
112  std::string name((const char*)element->name);
113  if (name == "reference_points")
114  {
115  if (has_ref)
116  {
117  APP_ABORT("EnergyDensityEstimator::put: EDE can only have one instance of reference_points.");
118  }
119  else
120  {
121  bool ref_succeeded = ref.put(element, *Pdynamic, Pref);
122  stop = stop || !ref_succeeded;
123  has_ref = true;
124  }
125  }
126  element = element->next;
127  }
128  if (!has_ref)
129  {
130  bool ref_succeeded = ref.put(*Pdynamic, Pref);
131  stop = stop || !ref_succeeded;
132  }
133  //initialize grids or other cell partitions
134  bool periodic = Pdynamic->getLattice().SuperCellEnum != SUPERCELL_OPEN;
135  bool grid_succeeded;
136  element = cur->children;
137  int nvalues = (int)nEDValues;
138  while (element != NULL)
139  {
140  std::string name = (const char*)element->name;
141  if (name == "spacegrid")
142  {
143  SpaceGrid* sg = new SpaceGrid(nvalues);
144  spacegrids.push_back(sg);
145  if (Pstatic)
146  {
147  set_ptcl();
148  grid_succeeded = sg->put(element, ref.points, Rptcl, Zptcl, Pdynamic->getTotalNum(), periodic, false);
149  unset_ptcl();
150  }
151  else
152  grid_succeeded = sg->put(element, ref.points, periodic, false);
153  stop = stop || !grid_succeeded;
154  }
155  element = element->next;
156  }
157  if (stop == true)
158  {
159  APP_ABORT("EnergyDensityEstimator::put");
160  }
161  return true;
162 }
163 
164 
166 {
167  ParticleSet& P = *Pstatic;
168  SpeciesSet& species(P.getSpeciesSet());
169  int ChargeAttribIndx = species.addAttribute("charge");
170  int nspecies = species.TotalNum;
171  int nps = P.getTotalNum();
172  std::vector<RealType> Zspec;
173  Zspec.resize(nspecies);
174  Zptcl.resize(nps);
175  for (int spec = 0; spec < nspecies; spec++)
176  Zspec[spec] = species(ChargeAttribIndx, spec);
177  for (int i = 0; i < nps; i++)
178  Zptcl[i] = Zspec[P.GroupID[i]];
179  Rptcl.resize(P.R.size());
180  for (int i = 0; i < P.R.size(); i++)
181  Rptcl[i] = P.R[i];
182  if (P.getLattice().SuperCellEnum != SUPERCELL_OPEN)
184 }
185 
187 {
188  Zptcl.clear();
189  Rptcl.clear();
190 }
191 
192 
194 {
195  auto pit(psetpool.find(psname));
196  if (pit == psetpool.end())
197  {
198  app_log() << " ParticleSet " << psname << " does not exist" << std::endl;
199  APP_ABORT("EnergyDensityEstimator::put");
200  }
201  return pit->second.get();
202 }
203 
204 
206 {
207  bool write = omp_get_thread_num() == 0;
208  w_trace = tm.get_real_trace("weight");
209  Td_trace = tm.get_real_trace(*Pdynamic, "Kinetic");
210  Vd_trace = tm.get_real_combined_trace(*Pdynamic, "LocalPotential");
211  if (Pstatic)
212  Vs_trace = tm.get_real_combined_trace(*Pstatic, "LocalPotential");
213  have_required_traces_ = true;
214 }
215 
216 
218 {
219  os << "EnergyDensityEstimator::write_description" << std::endl;
220  os << std::endl;
221  os << " EnergyDensityEstimator details" << std::endl;
222  os << std::endl;
223  std::string indent = " ";
224  os << indent + "nparticles = " << nparticles << std::endl;
225  os << indent + "nspacegrids = " << spacegrids.size() << std::endl;
226  os << std::endl;
227  ref.write_description(os, indent);
228  os << std::endl;
229  for (int i = 0; i < spacegrids.size(); i++)
230  {
231  spacegrids[i]->write_description(os, indent);
232  }
233  os << std::endl;
234  os << " end EnergyDensityEstimator details" << std::endl;
235  os << std::endl;
236  os << "end EnergyDensityEstimator::write_description" << std::endl;
237  return;
238 }
239 
240 
241 bool EnergyDensityEstimator::get(std::ostream& os) const
242 {
243  os << "EDM replace this " << std::endl;
244  APP_ABORT("EnergyDensityEstimator::get");
245  return true;
246 }
247 
248 
250 {
251  //remains empty
252 }
253 
254 
255 //#define ENERGYDENSITY_CHECK
256 
258 {
260  {
261  Pdynamic = &P;
262  //Collect positions from ParticleSets
263  int p = 0;
264  {
265  const ParticlePos& Rs = Pdynamic->R;
266  for (int i = 0; i < Rs.size(); i++)
267  {
268  R[p] = Rs[i];
269  p++;
270  }
271  }
272  if (Pstatic && !ion_points)
273  {
274  const ParticlePos& Rs = Pstatic->R;
275  for (int i = 0; i < Rs.size(); i++)
276  {
277  R[p] = Rs[i];
278  p++;
279  }
280  }
281  if (P.getLattice().SuperCellEnum != SUPERCELL_OPEN)
282  P.applyMinimumImage(R);
283  //Convert information accumulated in ParticleSets into EnergyDensity quantities
284  RealType w = w_trace->sample[0];
285  p = 0;
286  {
287  Vd_trace->combine();
288  const ParticleSet& Ps = *Pdynamic;
289  const auto& Ts = Td_trace->sample;
290  const auto& Vs = Vd_trace->sample;
291  for (int i = 0; i < Ps.getTotalNum(); i++)
292  {
293  EDValues(p, W) = w;
294  EDValues(p, T) = w * Ts[i];
295  EDValues(p, V) = w * Vs[i];
296  p++;
297  }
298  }
299  if (Pstatic)
300  {
301  Vs_trace->combine();
302  const ParticleSet& Ps = *Pstatic;
303  const auto& Vs = Vs_trace->sample;
304  if (!ion_points)
305  for (int i = 0; i < Ps.getTotalNum(); i++)
306  {
307  EDValues(p, W) = w;
308  EDValues(p, T) = 0.0;
309  EDValues(p, V) = w * Vs[i];
310  p++;
311  }
312  else
313  for (int i = 0; i < Ps.getTotalNum(); i++)
314  {
315  EDIonValues(i, W) = w;
316  EDIonValues(i, T) = 0.0;
317  EDIonValues(i, V) = w * Vs[i];
318  }
319  }
320  //Accumulate energy density in spacegrids
321  const auto& dtab(P.getDistTableAB(dtable_index));
322  fill(particles_outside.begin(), particles_outside.end(), true);
323  for (int i = 0; i < spacegrids.size(); i++)
324  {
325  SpaceGrid& sg = *spacegrids[i];
327  }
328  //Accumulate energy density of particles outside any spacegrid
329  int bi, v;
330  const int bimax = outside_buffer_offset + (int)nEDValues;
331  for (int p = 0; p < particles_outside.size(); p++)
332  {
333  if (particles_outside[p])
334  {
335  for (bi = outside_buffer_offset, v = 0; bi < bimax; bi++, v++)
336  {
337  P.Collectables[bi] += EDValues(p, v);
338  }
339  }
340  }
341  if (ion_points)
342  {
343  // Accumulate energy density for ions at a point field
344  bi = ion_buffer_offset;
345  for (int i = 0; i < nions; i++)
346  for (v = 0; v < (int)nEDValues; v++, bi++)
347  {
348  P.Collectables[bi] += EDIonValues(i, v);
349  }
350  }
351  nsamples++;
352 #if defined(ENERGYDENSITY_CHECK)
353  int thread = omp_get_thread_num();
355  RealType Eref = P.PropertyList[WP::LOCALENERGY];
356  RealType Vref = P.PropertyList[WP::LOCALPOTENTIAL];
357  RealType Tref = Eref - Vref;
358 #pragma omp critical(edcheck)
359  {
360  RealType Dsum = 0.0;
361  RealType Tsum = 0.0;
362  RealType Vsum = 0.0;
363  RealType Esum = 0.0;
364  for (int p = 0; p < nparticles; p++)
365  {
366  Dsum += EDValues(p, W);
367  Tsum += EDValues(p, T);
368  Vsum += EDValues(p, V);
369  }
370  if (ion_points)
371  for (int i = 0; i < nions; i++)
372  {
373  Dsum += EDIonValues(i, W);
374  Tsum += EDIonValues(i, T);
375  Vsum += EDIonValues(i, V);
376  }
377  Esum = Tsum + Vsum;
378  static int cnt = 0;
379  //app_log()<<"eval ED Dsum"<<cnt<<" "<<Dsum<< std::endl;
380  app_log() << thread << " eval ED " << cnt << " " << Tsum << " " << Vsum << " " << Esum << std::endl;
381  int nvals = (int)nEDValues;
382  RealType edvals[nvals];
383  RealType edtmp[nvals];
384  for (int v = 0; v < nvals; v++)
385  edvals[v] = 0.0;
386  for (int i = 0; i < spacegrids.size(); i++)
387  {
388  SpaceGrid& sg = *spacegrids[i];
389  sg.sum(P.Collectables, edtmp);
390  for (int v = 0; v < nvals; v++)
391  edvals[v] += edtmp[v];
392  }
393  for (int v = 0; v < nvals; v++)
394  edvals[v] += P.Collectables[outside_buffer_offset + v];
395  if (ion_points)
396  {
397  for (int i = 0; i < nions; i++)
398  for (int v = 0; v < nvals; v++)
399  edvals[v] += P.Collectables[ion_buffer_offset + i * nvals + v];
400  }
401  //app_log()<<"eval ES Dsum"<<cnt<<" "<<edvals[W]<< std::endl;
402  app_log() << thread << " eval ES " << cnt << " " << edvals[T] << " " << edvals[V] << " " << edvals[T] + edvals[V]
403  << std::endl;
404  app_log() << thread << " ref E " << cnt << " " << Tref << " " << Vref << " " << Eref << std::endl;
405  cnt++;
406  }
407 #endif
408  //APP_ABORT("EnergyDensityEstimator::evaluate");
409  }
410 
411  return 0.0;
412 }
413 
414 
415 void EnergyDensityEstimator::write_Collectables(std::string& label, int& cnt, ParticleSet& P)
416 {
417  //for(int v=0;v<nEDValues;v++){
418  int ii = spacegrids[0]->buffer_offset;
419  int io = outside_buffer_offset;
420  double Ti = P.Collectables[ii + 1] / P.Collectables[ii] * 12.0;
421  double To = P.Collectables[io + 1] / P.Collectables[io] * 12.0;
422  app_log() << "EDcoll " << label << cnt << " " << Ti << " " << To << std::endl;
423  //}
424 }
425 
426 
428 {
429  app_log() << "EDValues" << std::endl;
430  for (int p = 0; p < nparticles; p++)
431  fprintf(stdout, " %d %e %e %e\n", p, EDValues(p, 0), EDValues(p, 1), EDValues(p, 2));
432 }
433 
435 {
436  app_log() << "Nonzero domains" << std::endl;
437  int nd = 1;
438  for (int i = 0; i < spacegrids.size(); i++)
439  nd += spacegrids[i]->nDomains();
440  for (int i = 0; i < nd; i++)
441  {
442  bool nonzero = false;
443  int n = outside_buffer_offset + i * nEDValues;
444  for (int v = 0; v < nEDValues; v++)
445  {
446  nonzero = nonzero || std::abs(P.Collectables[n + v]) > 1e-8;
447  }
448  if (nonzero)
449  {
450  // fprintf(stdout," %d %e %e %e %e %e %e\n",i,P.Collectables[n],
451  fprintf(stdout, " %d %e %e %e \n", i, P.Collectables[n], P.Collectables[n + 1], P.Collectables[n + 2]);
452  }
453  }
454 }
455 
456 
458 {
459  my_index_ = collectables.size();
460  //allocate space for energy density outside of any spacegrid
461  outside_buffer_offset = collectables.size();
462  int nvalues = (int)nEDValues;
463  std::vector<RealType> tmp(nvalues);
464  collectables.add(tmp.begin(), tmp.end());
465  //allocate space for spacegrids
466  for (int i = 0; i < spacegrids.size(); i++)
467  {
468  spacegrids[i]->allocate_buffer_space(collectables);
469  }
470  if (ion_points)
471  {
472  ion_buffer_offset = collectables.size();
473  nvalues = nions * ((int)nEDValues);
474  std::vector<RealType> tmp2(nvalues);
475  collectables.add(tmp2.begin(), tmp2.end());
476  }
477 }
478 
479 
480 void EnergyDensityEstimator::registerCollectables(std::vector<ObservableHelper>& h5desc, hdf_archive& file) const
481 {
482  hdf_path hdf_name{name_};
483  h5desc.emplace_back(hdf_name / "variables");
484  auto& oh = h5desc.back();
485  oh.addProperty(const_cast<int&>(nparticles), "nparticles", file);
486  int nspacegrids = spacegrids.size();
487  oh.addProperty(const_cast<int&>(nspacegrids), "nspacegrids", file);
488  oh.addProperty(const_cast<int&>(nsamples), "nsamples", file);
489  if (ion_points)
490  {
491  oh.addProperty(const_cast<int&>(nions), "nions", file);
492  oh.addProperty(const_cast<Matrix<RealType>&>(Rion), "ion_positions", file);
493  }
494 
495  ref.save(h5desc, file);
496  h5desc.emplace_back(hdf_name / "outside");
497  auto& ohOutside = h5desc.back();
498  std::vector<int> ng(1);
499  ng[0] = (int)nEDValues;
500  ohOutside.set_dimensions(ng, outside_buffer_offset);
501  for (int i = 0; i < spacegrids.size(); i++)
502  {
503  SpaceGrid& sg = *spacegrids[i];
504  sg.registerCollectables(h5desc, file, i);
505  }
506  if (ion_points)
507  {
508  std::vector<int> ng2(2);
509  ng2[0] = nions;
510  ng2[1] = (int)nEDValues;
511 
512  h5desc.emplace_back(hdf_name / "ions");
513  auto& ohIons = h5desc.back();
514  ohIons.set_dimensions(ng2, ion_buffer_offset);
515  }
516 }
517 
519 {
520  //remains empty
521 }
522 
524 {
525  //remains empty
526 }
527 
528 
529 std::unique_ptr<OperatorBase> EnergyDensityEstimator::makeClone(ParticleSet& qp, TrialWaveFunction& psi)
530 {
531  bool write = omp_get_thread_num() == 0;
532  if (write)
533  app_log() << "EnergyDensityEstimator::makeClone" << std::endl;
534 
535  std::unique_ptr<EnergyDensityEstimator> edclone = std::make_unique<EnergyDensityEstimator>(psetpool, defKE);
536  edclone->put(input_xml, qp);
537  //int thread = omp_get_thread_num();
538  //app_log()<<thread<<"make edclone"<< std::endl;
539  //edclone->Pdynamic = Pdynamic->get_clone(thread);
540  //edclone->Pstatic = Pstatic->get_clone(thread);
541  return edclone;
542 }
543 
544 
545 } // namespace qmcplusplus
void resize(size_type n, Type_t val=Type_t())
Resize the container.
Definition: OhmmsVector.h:166
size_type size() const
return the size of the data
Definition: PooledData.h:48
Return_t evaluate(ParticleSet &P) override
Evaluate the local energy contribution of this component.
void setObservables(PropertySetType &plist) override
Set the values evaluated by this object to plist Default implementation is to assign Value which is u...
void delete_iter(IT first, IT last)
delete the pointers in [first,last)
helper functions for EinsplineSetBuilder
Definition: Configuration.h:43
void applyMinimumImage(ParticlePos &pinout) const
bool put(xmlNodePtr cur) override
check xml elements
MakeReturn< UnaryNode< FnFabs, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t abs(const Vector< T1, C1 > &l)
size_t getTotalNum() const
Definition: ParticleSet.h:493
void setParticlePropertyList(PropertySetType &plist, int offset) override
int my_index_
starting index of this object
Definition: OperatorBase.h:535
std::ostream & app_log()
Definition: OutputManager.h:65
if(c->rank()==0)
bool put(xmlNodePtr cur, ParticleSet &P, std::vector< ParticleSet *> &Pref)
bool put(xmlNodePtr cur)
assign attributes to the set
Definition: AttributeSet.h:55
void registerCollectables(std::vector< ObservableHelper > &h5desc, hdf_archive &file) const override
class to handle hdf file
Definition: hdf_archive.h:51
Vectorized record engine for scalar properties.
void addProperty(T &p, const std::string &pname, hdf_archive &file)
add named property to describe the collectable this helper class handles
void resize(size_type n, size_type m)
Resize the container.
Definition: OhmmsMatrix.h:99
int addAttribute(const std::string &aname)
for a new attribute, allocate the data, !More often used to get the index of a species ...
Definition: SpeciesSet.cpp:45
ParticleSet * get_particleset(std::string &psname)
TraceRequest request_
whether traces are being collected
Definition: OperatorBase.h:531
ParticleIndex GroupID
Species ID.
Definition: ParticleSet.h:77
void sum(const BufferType &buf, RealType *vals)
Definition: SpaceGrid.cpp:1071
void evaluate(const ParticlePos &R, const Matrix< RealType > &values, BufferType &buf, std::vector< bool > &particles_outside, const DistanceTableAB &dtab)
Definition: SpaceGrid.cpp:821
PropertySetType PropertyList
name-value map of Walker Properties
Definition: ParticleSet.h:112
void write_nonzero_domains(const ParticleSet &P)
const DistanceTableAB & getDistTableAB(int table_ID) const
get a distance table by table_ID and dyanmic_cast to DistanceTableAB
void turnOnPerParticleSK()
Turn on per particle storage in Structure Factor.
int addTable(const ParticleSet &psrc, DTModes modes=DTModes::ALL_OFF)
add a distance table
omp_int_t omp_get_thread_num()
Definition: OpenMP.h:25
std::string name_
name of this object
Definition: OperatorBase.h:527
Specialized paritlce class for atomistic simulations.
Definition: ParticleSet.h:55
void getRequiredTraces(TraceManager &tm) override
TODO: add docs.
size_type size() const
return the current size
Definition: OhmmsVector.h:162
void add(T &x)
Definition: PooledData.h:88
void addObservables(PropertySetType &plist)
class to handle a set of attributes of an xmlNode
Definition: AttributeSet.h:24
std::map< std::string, Point > points
void write_Collectables(std::string &label, int &cnt, ParticleSet &P)
void resetTargetParticleSet(ParticleSet &P) override
Reset the data with the target ParticleSet.
std::unique_ptr< OperatorBase > makeClone(ParticleSet &qp, TrialWaveFunction &psi) final
EnergyDensityEstimator(const PSPool &PSP, const std::string &defaultKE)
ParticlePos R
Position.
Definition: ParticleSet.h:79
CombinedTraceSample< TraceReal > * Vs_trace
#define APP_ABORT(msg)
Widely used but deprecated fatal error macros from legacy code.
Definition: AppAbort.h:27
void registerCollectables(std::vector< ObservableHelper > &h5desc, hdf_archive &file, int grid_index) const
Definition: SpaceGrid.cpp:704
void request_array(const std::string &name, bool write=false)
Definition: TraceManager.h:233
Buffer_t Collectables
observables in addition to those registered in Properties/PropertyList
Definition: ParticleSet.h:126
SpeciesSet & getSpeciesSet()
retrun the SpeciesSet of this particle set
Definition: ParticleSet.h:231
FullPrecRealType Return_t
type of return value of evaluate
Definition: OperatorBase.h:64
Class to represent a many-body trial wave function.
Indexes
an enum denoting index of physical properties
void write_description(std::ostream &os, std::string &indent)
std::map< std::string, const std::unique_ptr< ParticleSet > > PSPool
ParticleAttrib< SingleParticlePos > ParticlePos
Definition: Configuration.h:92
void clear()
clear
Definition: OhmmsVector.h:191
const auto & getLattice() const
Definition: ParticleSet.h:251
Define a LRHandler with two template parameters.
Custom container for set of attributes for a set of species.
Definition: SpeciesSet.h:33
CombinedTraceSample< TraceReal > * Vd_trace
CombinedTraceSample< TraceReal > * get_real_combined_trace(const std::string &name)
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
bool get(std::ostream &os) const override
write about the class
Declaration of a MCWalkerConfiguration.
TraceSample< TraceReal > * get_real_trace(const std::string &name)
void request_scalar(const std::string &name, bool write=false)
Definition: TraceManager.h:224
std::bitset< 8 > update_mode_
set the current update mode
Definition: OperatorBase.h:521
ObservableHelper oh
bool put(xmlNodePtr cur, std::map< std::string, Point > &points, ParticlePos &R, std::vector< RealType > &Z, int ndp, bool is_periodic, bool abort_on_fail=true)
Definition: SpaceGrid.h:34
void save(std::vector< ObservableHelper > &h5desc, hdf_archive &file) const