QMCPACK
TraceManager.h
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: Jeongnim Kim, jeongnim.kim@gmail.com, University of Illinois at Urbana-Champaign
8 // Cynthia Gu, zg1@ornl.gov, Oak Ridge National Laboratory
9 // Norbert Podhorszki, pnorbert@ornl.gov, Oak Ridge National Laboratory
10 // Paul R. C. Kent, kentpr@ornl.gov, Oak Ridge National Laboratory
11 // Jaron T. Krogel, krogeljt@ornl.gov, Oak Ridge National Laboratory
12 // Mark Dewing, markdewing@gmail.com, University of Illinois at Urbana-Champaign
13 // Mark A. Berrill, berrillma@ornl.gov, Oak Ridge National Laboratory
14 //
15 // File created by: Jeongnim Kim, jeongnim.kim@gmail.com, University of Illinois at Urbana-Champaign
16 //////////////////////////////////////////////////////////////////////////////////////
17 
18 
19 #ifndef QMCPLUSPLUS_TRACEMANAGER_H
20 #define QMCPLUSPLUS_TRACEMANAGER_H
21 
22 
23 #if !defined(DISABLE_TRACEMANAGER)
24 
25 
26 #include <Configuration.h>
28 #include "OhmmsData/AttributeSet.h"
29 #include "OhmmsPETE/OhmmsArray.h"
30 #include "Particle/ParticleSet.h"
32 #include "ModernStringUtils.hpp"
33 #include "Message/Communicate.h"
34 #include "hdf/hdf_archive.h"
35 #include "Concurrency/OpenMP.h"
36 
37 #include <algorithm>
38 #include <array>
39 #include <map>
40 #include <memory>
41 #include <set>
42 
43 namespace qmcplusplus
44 {
45 //#define TRACE_CHECK
46 
47 const unsigned int DMAX = 4;
48 using TraceInt = long;
50 using TraceComp = std::complex<TraceReal>;
51 
52 
54 {
55  std::string name;
68 
69  inline TraceQuantity()
70  {
71  default_quantity = false;
72  combined_quantity = false;
73  scalar_available = false;
74  array_available = false;
76  array_stream_requested = false;
77  scalar_write_requested = false;
78  array_write_requested = false;
79  stream_scalar = false;
80  stream_array = false;
81  write_scalar = false;
82  write_array = false;
83  }
84 
85  inline void incorporate(const TraceQuantity& other)
86  {
87  if (name != other.name)
88  {
89  APP_ABORT("TraceQuantity::incorporate\n cannot merge quantities with differing names\n names: " + name + " " +
90  other.name);
91  }
98  }
99 };
100 
101 
102 //means of control of traces
103 // which quantities are available, streaming, writing, etc
104 // handles global (user + all observable) trace request
105 //medium of exchange of trace information between trace manager and observables
106 // handles individual trace requests from observables
108 {
109  //switch to allow or disallow streams/writes of requested quantities
112 
113  //switch to allow or disallow scalar/array streaming/writing
115  bool arrays_on;
116 
117  //all scalar/array quantities should be provided for streaming
120 
121  //all scalar/array quantities should be provided for writing
124 
125  //records whether default scalars/arrays should stream/write
130 
131 
132  //quantities made available or requested for streaming or writing
133  std::map<std::string, TraceQuantity> quantities;
134 
135  //dependency lists for combined quantities
136  std::map<std::string, std::set<std::string>> combined_dependencies;
137 
138  //used to screen checked out quantities for writing
139  std::string scalar_domain;
140 
141  inline TraceRequest() { reset(); }
142 
143  inline void reset()
144  {
145  allow_streams = false;
146  allow_writes = false;
147  scalars_on = false;
148  arrays_on = false;
149  stream_all_scalars = false;
150  stream_all_arrays = false;
151  write_all_scalars = false;
152  write_all_arrays = false;
154  streaming_default_arrays = false;
155  writing_default_scalars = false;
156  writing_default_arrays = false;
157 
158  quantities.clear();
159  combined_dependencies.clear();
160  }
161 
162  inline void set_scalar_domain(const std::string& domain) { scalar_domain = domain; }
163 
164  inline bool screen_sample(const std::string& domain, const std::string& name, bool& write)
165  {
166  bool scalar = domain == scalar_domain;
167  bool present = quantities.find(name) != quantities.end();
168  bool stream = false;
169  write = false;
170  if (present)
171  {
172  TraceQuantity& q = quantities[name];
173  if (scalar)
174  {
175  stream = q.stream_scalar;
176  write = q.write_scalar;
177  }
178  else
179  {
180  stream = q.stream_array;
181  write = q.write_array;
182  }
183  }
184  return stream;
185  }
186 
187  //Contributor API (OperatorBase and others)
188  //declare that scalars are available for a quantity
189  inline void contribute_scalar(const std::string& name, bool default_quantity = false)
190  {
191  guarantee_presence(name);
192  quantities[name].scalar_available = true;
193  if (default_quantity)
194  quantities[name].default_quantity = true;
195  }
196 
197  //declare that arrays are available for a quantity
198  inline void contribute_array(const std::string& name, bool default_quantity = false)
199  {
200  guarantee_presence(name);
201  quantities[name].array_available = true;
202  if (default_quantity)
203  quantities[name].default_quantity = true;
204  }
205 
206  //declare that a combined quantity could be made available
207  inline void contribute_combined(const std::string& name,
208  std::vector<std::string>& deps,
209  bool scalar = false,
210  bool array = false,
211  bool default_quantity = false)
212  {
213  guarantee_presence(name, true);
214  TraceQuantity& q = quantities[name];
215  q.combined_quantity = true;
217  q.array_available = array;
218  if (default_quantity)
219  q.default_quantity = true;
220  //combined_dependencies[name] = deps;
221  }
222 
223  //declare that a scalar quantity is desired for streaming
224  inline void request_scalar(const std::string& name, bool write = false)
225  {
226  guarantee_presence(name);
227  quantities[name].scalar_stream_requested = true;
228  if (write)
229  quantities[name].scalar_write_requested = true;
230  }
231 
232  //declare that a array quantity is desired for streaming
233  inline void request_array(const std::string& name, bool write = false)
234  {
235  guarantee_presence(name);
236  quantities[name].array_stream_requested = true;
237  if (write)
238  quantities[name].array_write_requested = true;
239  }
240 
241  //query whether a scalar quantity should/will be streaming
242  inline bool streaming_scalar(const std::string& name)
243  {
244  check_presence(name);
245  return quantities[name].stream_scalar;
246  }
247 
248  //query whether a array quantity should/will be streaming
249  inline bool streaming_array(const std::string& name)
250  {
251  check_presence(name);
252  return quantities[name].stream_array;
253  }
254 
255  //query whether a quantity will be streaming in any fashion
256  inline bool streaming(const std::string& name)
257  {
258  check_presence(name);
259  TraceQuantity& q = quantities[name];
260  return q.stream_scalar || q.stream_array;
261  }
262 
263 
264  //TraceManager API
265  //declare that scalar quantities are desired for streaming
266  inline void request_scalar(const std::set<std::string>& names, bool write = false)
267  {
268  std::set<std::string>::iterator name;
269  for (name = names.begin(); name != names.end(); ++name)
270  request_scalar(*name, write);
271  }
272 
273  //declare that array quantities are desired for streaming
274  inline void request_array(const std::set<std::string>& names, bool write = false)
275  {
276  std::set<std::string>::iterator name;
277  for (name = names.begin(); name != names.end(); ++name)
278  request_array(*name, write);
279  }
280 
281  //merge in all quantities from a contributor request
282  inline void incorporate(TraceRequest& other)
283  {
284  std::map<std::string, TraceQuantity>::iterator it;
285  for (it = other.quantities.begin(); it != other.quantities.end(); ++it)
286  {
287  const TraceQuantity& q = it->second;
288  if (quantity_present(q.name))
289  quantities[q.name].incorporate(q);
290  else
291  quantities[q.name] = q;
292  }
293  std::map<std::string, std::set<std::string>>::iterator d;
294  for (d = other.combined_dependencies.begin(); d != other.combined_dependencies.end(); ++d)
295  {
296  const std::string& name = d->first;
297  std::set<std::string>& deps = d->second;
298  if (combined_dependencies.find(name) != combined_dependencies.end())
299  combined_dependencies[name].insert(deps.begin(), deps.end());
300  else
301  combined_dependencies[name] = deps;
302  }
303  }
304 
305  //balance requests with availability and turn streaming/writing on/off
307  {
308  std::map<std::string, TraceQuantity>::iterator it;
309  for (it = quantities.begin(); it != quantities.end(); ++it)
310  {
311  TraceQuantity& q = it->second;
312 
313  q.stream_scalar =
315 
317 
318  q.stream_array =
320 
322  }
323  // default quantities stream and write if any others do
325  writing_default_scalars = false;
326  streaming_default_arrays = false;
327  writing_default_arrays = false;
328  for (it = quantities.begin(); it != quantities.end(); ++it)
329  {
330  TraceQuantity& q = it->second;
335  }
336  for (it = quantities.begin(); it != quantities.end(); ++it)
337  {
338  TraceQuantity& q = it->second;
339  if (q.default_quantity)
340  {
345  }
346  }
347  // if any combined quantities are streaming, their dependencies must also
348  for (it = quantities.begin(); it != quantities.end(); ++it)
349  {
350  TraceQuantity& q = it->second;
351  if (q.combined_quantity)
352  {
353  std::set<std::string>& deps = combined_dependencies[q.name];
354  std::set<std::string>::iterator name;
355  if (q.stream_scalar || q.stream_array)
356  {
357  for (name = deps.begin(); name != deps.end(); ++name)
358  {
359  check_presence(*name);
360  TraceQuantity& qd = quantities[*name];
362  qd.stream_array |= (qd.array_available && q.stream_array);
363  }
364  }
365  }
366  }
367  combined_dependencies.clear();
368  }
369 
370  //relay updated streaming information to contributor
371  inline void relay_stream_info(TraceRequest& other)
372  {
374  other.allow_writes = allow_writes;
375  other.scalars_on = scalars_on;
376  other.arrays_on = arrays_on;
385  std::map<std::string, TraceQuantity>::iterator it;
386  for (it = other.quantities.begin(); it != other.quantities.end(); ++it)
387  {
388  TraceQuantity& q = it->second;
389  check_presence(q.name);
390  q = quantities[q.name];
391  }
392  other.combined_dependencies.clear();
393  }
394 
395  inline void report()
396  {
397  //app_log()<<"\n TraceRequest"<< std::endl;
398  app_log() << " allow_streams = " << allow_streams << std::endl;
399  app_log() << " allow_writes = " << allow_writes << std::endl;
400  app_log() << " scalars_on = " << scalars_on << std::endl;
401  app_log() << " arrays_on = " << arrays_on << std::endl;
402  app_log() << " stream_all_scalars = " << stream_all_scalars << std::endl;
403  app_log() << " stream_all_arrays = " << stream_all_arrays << std::endl;
404  app_log() << " write_all_scalars = " << write_all_scalars << std::endl;
405  app_log() << " write_all_arrays = " << write_all_arrays << std::endl;
406  app_log() << " streaming_default_scalars = " << streaming_default_scalars << std::endl;
407  app_log() << " streaming_default_arrays = " << streaming_default_arrays << std::endl;
408  app_log() << " writing_default_scalars = " << writing_default_scalars << std::endl;
409  app_log() << " writing_default_arrays = " << writing_default_arrays << std::endl;
410 
411  write_selected("scalars available", "scalar_available");
412  write_selected("arrays available", "array_available");
413 
414  write_selected("scalar streams requested", "scalar_stream_requested");
415  write_selected("array streams requested", "array_stream_requested");
416 
417  write_selected("scalar writes requested", "scalar_write_requested");
418  write_selected("array writes requested", "array_write_requested");
419 
420  write_selected("scalar streams occurring", "stream_scalar");
421  write_selected("array streams occurring", "stream_array");
422 
423  write_selected("scalar writes occurring", "write_scalar");
424  write_selected("array writes occurring", "write_array");
425  }
426 
427  inline void write_selected(const std::string& header, const std::string& selector)
428  {
429  app_log() << " " << header << ":";
430  int n = 0;
431  std::map<std::string, TraceQuantity>::iterator it;
432  for (it = quantities.begin(); it != quantities.end(); ++it)
433  {
434  TraceQuantity& q = it->second;
435  bool selected = false;
436  if (selector == "scalar_available")
437  selected = q.scalar_available;
438  else if (selector == "array_available")
439  selected = q.array_available;
440  else if (selector == "scalar_stream_requested")
441  selected = q.scalar_stream_requested;
442  else if (selector == "array_stream_requested")
443  selected = q.array_stream_requested;
444  else if (selector == "scalar_write_requested")
445  selected = q.scalar_write_requested;
446  else if (selector == "array_write_requested")
447  selected = q.array_write_requested;
448  else if (selector == "stream_scalar")
449  selected = q.stream_scalar;
450  else if (selector == "stream_array")
451  selected = q.stream_array;
452  else if (selector == "write_scalar")
453  selected = q.write_scalar;
454  else if (selector == "write_array")
455  selected = q.write_array;
456  else
457  APP_ABORT("TraceRequest::write_selected unrecognized selector: " + selector);
458  if (selected)
459  {
460  if (n % 5 == 0)
461  app_log() << std::endl << " ";
462  n++;
463  app_log() << " " << q.name;
464  }
465  }
466  app_log() << std::endl;
467  }
468 
469 
470  //private (internal) API
471  //query whether a quantity is present
472  inline bool quantity_present(const std::string& name) { return quantities.find(name) != quantities.end(); }
473 
474  //create a quantity if it is not present
475  inline void guarantee_presence(const std::string& name, bool combined = false)
476  {
477  if (!quantity_present(name))
478  {
479  TraceQuantity q;
480  q.name = name;
481  quantities[name] = q;
482  }
483  if (combined)
484  if (combined_dependencies.find(name) == combined_dependencies.end())
485  {
486  std::set<std::string> stmp;
487  combined_dependencies[name] = stmp;
488  }
489  }
490 
491  //abort if a quantity is not present
492  inline void check_presence(const std::string& name)
493  {
494  if (!quantity_present(name))
495  {
496  APP_ABORT("TraceRequest::check_presence quantity " + name + " is not present");
497  }
498  }
499 
500  //query whether any quantities are streaming
502 
503  //query whether any quantities are writing
505 
506  //query whether any scalar quantities are streaming
508 
509  //query whether any array quantities are streaming
510  inline bool streaming_arrays() { return streaming_default_arrays; }
511 };
512 
513 
514 template<typename T>
516 {
517  std::string domain;
518  std::string name;
519  int index;
522  int size;
527  bool write;
529  std::map<std::string, TraceInt> meta_int;
530  std::map<std::string, TraceReal> meta_real;
531  std::map<std::string, std::string> meta_string;
532  bool verbose;
533 
534  inline TraceSample(const std::string& sdomain, const std::string& sname, int sindex, int sdim, Vector<T>& ssample)
535  : sample(ssample), verbose(false)
536  {
537  initialize(sdomain, sname, sindex, sdim);
538  }
539 
540 
541  inline TraceSample(const std::string& sdomain,
542  const std::string& sname,
543  int sindex,
544  int sdim,
545  TinyVector<int, DMAX> sshape,
546  Vector<T>& ssample)
547  : sample(ssample), verbose(false)
548  {
549  initialize(sdomain, sname, sindex, sdim);
550  shape = sshape;
551  size = sample.size();
552  check_shape();
553  }
554 
555  inline virtual ~TraceSample() = default;
556 
557  inline void initialize(const std::string& sdomain, const std::string& sname, int sindex, int sdim)
558  {
559  domain = sdomain, name = sname;
560  dimension = sdim;
561  index = sindex;
562  array_trace = false;
563  write = false;
564  buffer_start = -1;
565  buffer_end = -1;
566  }
567 
568 
569  inline void set_unit_size(int usize) { unit_size = usize; }
570 
571 
572  inline void set_data_size() { data_size = size * unit_size; }
573 
574 
575  inline void check_shape()
576  {
577  bool correct_shape = true;
578  bool correct_dimension = dimension <= DMAX;
579  if (correct_dimension)
580  {
581  int tsize = 1;
582  for (int d = 0; d < dimension; ++d)
583  {
584  tsize *= shape[d];
585  correct_dimension = correct_dimension && shape[d] > 0;
586  }
587  correct_shape = tsize == size;
588  }
589  if (!correct_dimension)
590  APP_ABORT("TraceSample::check_shape dimension of sample array is incorrect");
591  if (!correct_shape)
592  APP_ABORT("TraceSample::check_shape shape and size of sample array do not match");
593  }
594 
595 
596  inline bool same_shape(TraceSample<T>* other)
597  {
598  bool same = dimension == other->dimension && size == other->size;
599  if (same)
600  for (int d = 0; d < dimension; ++d)
601  same = same && shape[d] == other->shape[d];
602  return same;
603  }
604 
605  virtual bool is_combined() { return false; }
606 
607  inline void set_buffer_range(int& bstart)
608  {
609  set_data_size();
610  if (write)
611  {
612  buffer_start = bstart;
613  buffer_end = bstart + data_size;
614  bstart = buffer_end;
615  }
616  }
617 
618 
619  inline T sum()
620  {
621  T s(0);
622  for (int i = 0; i < sample.size(); ++i)
623  s += sample[i];
624  return s;
625  }
626 
627  inline void write_summary(int ind = -1, std::string pad = " ")
628  {
629  std::string pad2 = pad + " ";
630  if (ind == -1)
631  app_log() << pad << " TraceSample " << name << std::endl;
632  else
633  app_log() << pad << ind << " TraceSample " << name << std::endl;
634  app_log() << pad2 << "domain = " << domain << std::endl;
635  app_log() << pad2 << "name = " << name << std::endl;
636  app_log() << pad2 << "index = " << index << std::endl;
637  app_log() << pad2 << "array_trace = " << array_trace << std::endl;
638  app_log() << pad2 << "dimension = " << dimension << std::endl;
639  app_log() << pad2 << "size = " << size << std::endl;
640  app_log() << pad2 << "unit_size = " << unit_size << std::endl;
641  app_log() << pad2 << "data_size = " << data_size << std::endl;
642  app_log() << pad2 << "shape = " << shape << std::endl;
643  app_log() << pad2 << "write = " << write << std::endl;
644  app_log() << pad2 << "buffer range = [" << buffer_start << "," << buffer_end << ")" << std::endl;
645  }
646 };
647 
648 
649 template<typename T>
651 {
652  bool combined;
653  std::vector<TraceReal> weights;
654  std::vector<TraceSample<T>*> components;
655 
656 
657  inline CombinedTraceSample(const std::string& sdomain,
658  const std::string& sname,
659  int sindex,
660  int sdim,
661  Vector<T>& ssample)
662  : TraceSample<T>(sdomain, sname, sindex, sdim, ssample)
663  {
664  reset();
665  }
666 
667 
668  inline CombinedTraceSample(const std::string& sdomain,
669  const std::string& sname,
670  int sindex,
671  int sdim,
672  TinyVector<int, DMAX> sshape,
673  Vector<T>& ssample)
674  : TraceSample<T>(sdomain, sname, sindex, sdim, sshape, ssample)
675  {
676  reset();
677  }
678 
679  bool is_combined() override { return true; }
680 
681  inline void reset() { combined = false; }
682 
683 
684  inline void add_component(TraceSample<T>* component, TraceReal weight)
685  {
686  if (components.size() == 0)
687  {
688  this->dimension = component->dimension;
689  this->size = component->size;
690  this->shape = component->shape;
691  this->data_size = component->data_size;
692  this->array_trace = component->array_trace;
693  this->sample.resize(component->size);
694  }
695  else if (!this->same_shape(component))
696  {
697  APP_ABORT("CombinedTraceSample::add_component attempted to add a different shape component\n my domain: " +
698  this->domain + "\n my name: " + this->name + "\n component domain: " + component->domain +
699  "\n component name: " + component->name);
700  }
701  else if (this->domain != component->domain)
702  APP_ABORT("CombinedTraceSample::add_component attempted to add a different domain component\n my domain: " +
703  this->domain + "\n my name: " + this->name + "\n component domain: " + component->domain +
704  "\n component name: " + component->name);
705  weights.push_back(weight);
706  components.push_back(component);
707  }
708 
709 
710  inline void combine()
711  {
712  std::fill(this->sample.begin(), this->sample.end(), T(0));
713  for (int i = 0; i < components.size(); ++i)
714  {
715  T weight = weights[i];
716  auto& component = components[i]->sample;
717  for (int j = 0; j < this->sample.size(); ++j)
718  this->sample[j] += weight * component[j];
719  }
720  combined = true;
721  }
722 
723 
724  inline void write_summary_combined(int ind, std::string pad = " ")
725  {
726  std::string pad2 = pad + " ";
727  std::string pad3 = pad2 + " ";
728  app_log() << pad << ind << " CombinedTraceSample " << this->name << std::endl;
729  app_log() << pad2 << "domain = " << this->domain << std::endl;
730  app_log() << pad2 << "ncomponents = " << components.size() << std::endl;
731  app_log() << pad2 << "components" << std::endl;
732  for (int i = 0; i < components.size(); ++i)
733  {
734  TraceSample<T>& c = *components[i];
735  app_log() << pad3 << c.name << " " << c.index << " " << weights[i] << std::endl;
736  }
737  app_log() << pad2 << "end components" << std::endl;
738  app_log() << pad2 << "vector address = " << (size_t) & this->sample << std::endl;
739  }
740 };
741 
742 
743 template<typename T>
745 {
746  return left->data_size < right->data_size;
747 }
748 
749 
750 template<typename T>
752 {
753  std::vector<TraceSample<T>*> samples;
754  std::map<std::string, std::map<std::string, int>> sample_indices;
755  std::vector<TraceSample<T>*> ordered_samples;
756  std::vector<CombinedTraceSample<T>*> combined_samples;
757  std::vector<Vector<T>*> combined_sample_vectors;
758  bool verbose;
759 
760  inline TraceSamples() : verbose(false) {}
761 
762  inline ~TraceSamples() { finalize(); }
763 
764 
765  inline void set_verbose(bool v) { verbose = v; }
766 
767 
768  inline int size() { return samples.size(); }
769 
770 
771  inline void assign_sample_index(const std::string& domain, const std::string& name, int index, std::string label = "")
772  {
773  if (sample_indices.count(domain) > 0 && sample_indices[domain].count(name) > 0)
774  {
775  APP_ABORT("TraceSamples::checkout " + label + " variable " + name + " already exists in domain " + domain);
776  }
777  else
778  {
779  sample_indices[domain][name] = index;
780  }
781  }
782 
783  template<int D>
784  inline Array<T, D>* checkout_array(const std::string& domain, const std::string& name, TinyVector<int, DMAX> shape)
785  {
786  int index = samples.size();
787  assign_sample_index(domain, name, index, "array");
788  std::array<size_t, D> subshape;
789  for (int idim = 0; idim < D; idim++)
790  subshape[idim] = shape[idim];
791  Array<T, D>* a = new Array<T, D>(subshape);
792  TraceSample<T>* s = new TraceSample<T>(domain, name, index, D, shape, a->storage());
793  samples.push_back(s);
794  if (verbose)
795  app_log() << "TraceSamples::checkout_array " << domain << " " << name << " " << index << std::endl;
796  return a;
797  }
798 
799 
800  template<int D>
801  inline Array<T, D>* checkout_array(const ParticleSet& P, const std::string& name, TinyVector<int, DMAX> shape)
802  {
803  const std::string& domain = P.parentName();
804  int index = samples.size();
805  assign_sample_index(domain, name, index, "array");
806  std::array<size_t, D> subshape;
807  for (int idim = 0; idim < D; idim++)
808  subshape[idim] = shape[idim];
809  Array<T, D>* a = new Array<T, D>(subshape);
810  TraceSample<T>* s = new TraceSample<T>(domain, name, index, D, shape, a->storage());
811  samples.push_back(s);
812  s->array_trace = true;
813  if (verbose)
814  app_log() << "TraceSamples::checkout_array " << domain << " " << name << " " << index << std::endl;
815  return a;
816  }
817 
818 
819  inline TraceSample<T>* get_trace(const std::string& domain, const std::string& name)
820  {
821  TraceSample<T>* ts = NULL;
822  for (int i = 0; i < samples.size(); ++i)
823  {
824  TraceSample<T>& tsc = *samples[i];
825  if (tsc.domain == domain && tsc.name == name)
826  {
827  ts = &tsc;
828  break;
829  }
830  }
831  if (ts == NULL)
832  APP_ABORT("TraceSamples::get_trace failed to get trace for quantity " + name + " in domain " + domain);
833  return ts;
834  }
835 
836 
837  inline CombinedTraceSample<T>* get_combined_trace(const std::string& domain, const std::string& name)
838  {
839  CombinedTraceSample<T>* ts = NULL;
840  for (int i = 0; i < combined_samples.size(); ++i)
841  {
843  if (tsc.domain == domain && tsc.name == name)
844  {
845  ts = &tsc;
846  break;
847  }
848  }
849  return ts;
850  }
851 
852 
853  inline bool make_combined_trace(const std::string& name,
854  std::vector<std::string>& names,
855  std::vector<TraceReal>& weights)
856  {
857  bool created = false;
858  if (names.size() != weights.size())
859  APP_ABORT("TraceSamples::make_combined_trace names and weights must be the same size");
860  std::map<std::string, std::map<std::string, int>>::iterator it;
861  for (it = sample_indices.begin(); it != sample_indices.end(); it++)
862  {
863  std::string domain = it->first;
864  std::map<std::string, int>& indices = it->second;
865  bool any_present = false;
866  for (int i = 0; i < names.size(); ++i)
867  any_present = any_present || indices.count(names[i]) > 0;
868  if (any_present)
869  {
870  int index = samples.size();
871  auto* sample = new Vector<T>;
872  CombinedTraceSample<T>* combined = new CombinedTraceSample<T>(domain, name, index, 0, *sample);
873  for (int i = 0; i < names.size(); ++i)
874  {
875  if (indices.count(names[i]) > 0)
876  {
877  TraceSample<T>* component = samples[indices[names[i]]];
878  combined->add_component(component, weights[i]);
879  }
880  }
881  assign_sample_index(domain, name, index);
882  samples.push_back(combined);
883  combined_samples.push_back(combined);
884  combined_sample_vectors.push_back(sample);
885  created = true;
886  }
887  }
888  return created;
889  }
890 
891 
892  inline void set_unit_size(int usize)
893  {
894  for (int i = 0; i < samples.size(); i++)
895  samples[i]->set_unit_size(usize);
896  }
897 
898 
899  inline void screen_writes(TraceRequest& request)
900  {
901  for (int i = 0; i < samples.size(); i++)
902  {
903  TraceSample<T>& s = *samples[i];
904  bool stream = request.screen_sample(s.domain, s.name, s.write);
905  if (verbose)
906  app_log() << "TraceRequest screening " << s.name << " in domain " << s.domain << ". stream: " << stream
907  << " write: " << s.write << std::endl;
908  if (!stream && !s.is_combined())
909  app_log() << "warning: quantity " + s.name + " in domain " + s.domain +
910  " was not requested but is streaming anyway"
911  << std::endl;
912  }
913  }
914 
915 
916  inline void order_by_size()
917  {
918  for (int i = 0; i < samples.size(); i++)
919  samples[i]->set_data_size();
920  ordered_samples.resize(samples.size());
921  copy(samples.begin(), samples.end(), ordered_samples.begin());
922  sort(ordered_samples.begin(), ordered_samples.end(), TraceSample_comp<T>);
923  }
924 
925 
926  inline void set_buffer_ranges(int& starting_index)
927  {
928  for (int i = 0; i < ordered_samples.size(); i++)
929  {
930  TraceSample<T>& sample = *ordered_samples[i];
931  sample.set_buffer_range(starting_index);
932  }
933  }
934 
935 
936  inline int total_size()
937  {
938  int s = 0;
939  for (int i = 0; i < samples.size(); i++)
940  s += samples[i]->sample.size() * samples[i]->unit_size;
941  return s;
942  }
943 
944 
945  inline int min_buffer_index()
946  {
947  int min_index = 2000000000;
948  for (int i = 0; i < samples.size(); i++)
949  min_index = std::min(min_index, samples[i]->buffer_start);
950  return min_index;
951  }
952 
953 
954  inline int max_buffer_index()
955  {
956  int max_index = -1;
957  for (int i = 0; i < samples.size(); i++)
958  max_index = std::max(max_index, samples[i]->buffer_end);
959  return max_index;
960  }
961 
962 
963  inline void combine_samples()
964  {
965  for (int i = 0; i < combined_samples.size(); ++i)
966  combined_samples[i]->combine();
967  }
968 
969 
971  {
972  for (int i = 0; i < combined_samples.size(); ++i)
973  combined_samples[i]->reset();
974  }
975 
976 
977  inline void finalize()
978  {
979  //note combined_samples pointers are a subset of those in samples
980  // and so only samples needs to be deleted
981  delete_iter(samples.begin(), samples.end());
983  samples.resize(0);
984  ordered_samples.resize(0);
985  combined_samples.resize(0);
986  combined_sample_vectors.resize(0);
987  sample_indices.clear();
988  }
989 
990 
992  {
993  std::map<std::string, std::map<std::string, int>>::iterator it;
994  std::map<std::string, int>::iterator it2;
995  for (it = sample_indices.begin(); it != sample_indices.end(); it++)
996  {
997  const std::string& domain = it->first;
998  std::map<std::string, int>& indices = it->second;
999  f.push(domain);
1000  for (it2 = indices.begin(); it2 != indices.end(); ++it2)
1001  {
1002  const std::string& quantity = it2->first;
1003  TraceSample<T>& sample = *samples[it2->second];
1004  if (sample.write)
1005  {
1006  f.push(quantity);
1007  f.write(sample.dimension, "dimension");
1008  f.write(sample.shape, "shape");
1009  f.write(sample.size, "size");
1010  f.write(sample.unit_size, "unit_size");
1011  f.write(sample.buffer_start, "row_start");
1012  f.write(sample.buffer_end, "row_end");
1013  f.pop();
1014  }
1015  }
1016  f.pop();
1017  }
1018  }
1019 
1020 
1021  inline void write_summary(std::string type, std::string pad = " ")
1022  {
1023  std::string pad2 = pad + " ";
1024  std::string pad3 = pad2 + " ";
1025  std::string pad4 = pad3 + " ";
1026  app_log() << pad << "TraceSamples<" << type << ">" << std::endl;
1027  app_log() << pad2 << "nsamples = " << samples.size() << std::endl;
1028  app_log() << pad2 << "ncombined_samples = " << combined_samples.size() << std::endl;
1029  app_log() << pad2 << "sample_indices" << std::endl;
1030  std::map<std::string, std::map<std::string, int>>::iterator it;
1031  std::map<std::string, int>::iterator it2;
1032  for (it = sample_indices.begin(); it != sample_indices.end(); it++)
1033  {
1034  const std::string& domain = it->first;
1035  std::map<std::string, int>& indices = it->second;
1036  app_log() << pad3 << domain << std::endl;
1037  for (it2 = indices.begin(); it2 != indices.end(); ++it2)
1038  app_log() << pad4 << it2->first << " = " << it2->second << std::endl;
1039  }
1040  app_log() << pad2 << "end sample_indices" << std::endl;
1041  app_log() << pad2 << "combined_sample_vectors = ";
1042  for (int i = 0; i < combined_sample_vectors.size(); ++i)
1043  app_log() << (size_t)combined_sample_vectors[i] << " ";
1044  app_log() << std::endl;
1045  app_log() << pad2 << "combined_samples" << std::endl;
1046  for (int i = 0; i < combined_samples.size(); ++i)
1047  combined_samples[i]->write_summary_combined(i, pad3);
1048  app_log() << pad2 << "end combined_samples" << std::endl;
1049  app_log() << pad2 << "samples" << std::endl;
1050  for (int i = 0; i < ordered_samples.size(); ++i)
1051  ordered_samples[i]->write_summary(i, pad3);
1052  //for(int i=0; i<samples.size(); ++i)
1053  // samples[i]->write_summary(i,pad3);
1054  app_log() << pad2 << "end samples" << std::endl;
1055  app_log() << pad << "end TraceSamples<" << type << ">" << std::endl;
1056  }
1057 
1058 
1059  inline void user_report(const std::string& type, const std::string& pad = " ")
1060  {
1061  std::string pad2 = pad + " ";
1062  app_log() << pad << type << " traces provided by estimators" << std::endl;
1063  std::map<std::string, std::map<std::string, int>>::iterator it;
1064  std::map<std::string, int>::iterator it2;
1065  for (it = sample_indices.begin(); it != sample_indices.end(); it++)
1066  {
1067  const std::string& domain = it->first;
1068  std::map<std::string, int>& indices = it->second;
1069  app_log() << pad2 << "domain " << domain << ": ";
1070  int n = 0;
1071  for (it2 = indices.begin(); it2 != indices.end(); ++it2)
1072  {
1073  if (n % 5 == 0)
1074  app_log() << std::endl << pad2 << " ";
1075  n++;
1076  const std::string& quantity = it2->first;
1077  app_log() << quantity << " ";
1078  }
1079  app_log() << std::endl;
1080  }
1081  }
1082 };
1083 
1084 
1085 template<typename T>
1087 {
1091  std::string type;
1093  bool verbose;
1094 
1095  //hdf variables
1096  std::string top;
1097  hsize_t dims[2];
1099 
1100 
1102  {
1103  type = "?";
1104  has_complex = false;
1105  reset();
1106  }
1107 
1108 
1109  inline void set_verbose(bool v) { verbose = v; }
1110 
1111 
1112  inline void set_type(std::string stype)
1113  {
1114  type = stype;
1115  top = type + "_data";
1116  }
1117 
1118 
1119  inline void reset() { buffer.resize(0, buffer.size(1)); }
1120 
1121 
1122  inline void set_samples(TraceSamples<T>& s) { samples = &s; }
1123 
1124 
1125  inline void set_samples(TraceSamples<std::complex<T>>& s)
1126  {
1127  complex_samples = &s;
1128  has_complex = true;
1129  }
1130 
1131 
1132  inline void make_combined_trace(const std::string& name,
1133  std::vector<std::string>& names,
1134  std::vector<TraceReal>& weights)
1135  {
1136  bool created_real = samples->make_combined_trace(name, names, weights);
1137  if (has_complex)
1138  {
1139  bool created_complex = complex_samples->make_combined_trace(name, names, weights);
1140  if (created_real && created_complex)
1141  APP_ABORT("TraceBuffer<" + type +
1142  ">::make_combined_trace\n cannot create real and complex combined traces for the same quantity\n "
1143  "attempted for quantity " +
1144  name);
1145  }
1146  }
1147 
1148 
1149  inline void order_and_resize()
1150  {
1151  //put the sample data in size order
1152  samples->set_unit_size(1);
1153  samples->order_by_size();
1154  if (has_complex)
1155  {
1158  }
1159  //assign buffer ranges to each sample
1160  int sample_size = 0;
1161  samples->set_buffer_ranges(sample_size);
1162  if (has_complex)
1163  complex_samples->set_buffer_ranges(sample_size);
1164 #if defined(TRACE_CHECK)
1165  test_buffer_write(sample_size);
1166 #endif
1167  //resize the buffer
1168  int nsamples_init = 1;
1169  buffer.resize(nsamples_init, sample_size);
1170  }
1171 
1172 
1173  inline bool same_as(TraceBuffer<T>& ref) { return buffer.size(1) == ref.buffer.size(1); }
1174 
1175 
1176  inline void collect_sample()
1177  {
1178  if (verbose)
1179  app_log() << " TraceBuffer<" << type << ">::collect_sample()" << std::endl;
1180  //make more room, if necessary
1181  int nrows = buffer.size(0);
1182  int row_size = buffer.size(1);
1183  if (row_size > 0)
1184  {
1185  //make space for the row, if necessary
1186  int current_row = nrows;
1187  nrows++;
1188  // resizing buffer(type Array) doesn't preserve data. Thus keep old data and copy over
1189  auto buffer_old(buffer);
1190  buffer.resize(nrows, row_size);
1191  std::copy_n(buffer_old.data(), buffer_old.size(), buffer.data());
1192  if (verbose)
1193  app_log() << " increasing # of rows to " << nrows << std::endl;
1194  //combine samples
1195  samples->combine_samples();
1196  if (has_complex)
1198  //collect data from all samples into the buffer row
1199  {
1200  std::vector<TraceSample<T>*>& ordered_samples = samples->ordered_samples;
1201  for (int s = 0; s < ordered_samples.size(); s++)
1202  {
1203  TraceSample<T>& tsample = *ordered_samples[s];
1204  if (tsample.write)
1205  {
1206  auto& sample = tsample.sample;
1207  for (int i = 0; i < sample.size(); ++i)
1208  buffer(current_row, tsample.buffer_start + i) = sample[i];
1209  }
1210  }
1211  }
1212  if (has_complex)
1213  {
1214  std::vector<TraceSample<std::complex<T>>*>& ordered_samples = complex_samples->ordered_samples;
1215  for (int s = 0; s < ordered_samples.size(); s++)
1216  {
1217  TraceSample<std::complex<T>>& tsample = *ordered_samples[s];
1218  if (tsample.write)
1219  {
1220  auto& sample = tsample.sample;
1221  for (int i = 0, ib = 0; i < sample.size(); ++i, ib += 2)
1222  {
1223  buffer(current_row, tsample.buffer_start + ib) = sample[i].real();
1224  buffer(current_row, tsample.buffer_start + ib + 1) = sample[i].imag();
1225  }
1226  }
1227  }
1228  }
1229  //reset combined samples so they can be recombined on the next step
1230  samples->reset_combined_samples();
1231  if (has_complex)
1233 #if defined(TRACE_CHECK)
1234  test_buffer_collect(current_row);
1235 #endif
1236  }
1237  }
1238 
1239 
1240  inline void write() { APP_ABORT("TraceBuffer::write has not yet been implemented"); }
1241 
1242 
1243  inline void write_summary(std::string pad = " ")
1244  {
1245  std::string pad2 = pad + " ";
1246  app_log() << pad << "TraceBuffer<" << type << ">" << std::endl;
1247  app_log() << pad2 << "nrows = " << buffer.size(0) << std::endl;
1248  app_log() << pad2 << "row_size = " << buffer.size(1) << std::endl;
1249  app_log() << pad2 << "has_complex = " << has_complex << std::endl;
1250  samples->write_summary(type, pad2);
1251  if (has_complex)
1252  complex_samples->write_summary("complex " + type, pad2);
1253  app_log() << pad << "end TraceBuffer<" << type << ">" << std::endl;
1254  }
1255 
1256 
1257  inline void user_report(const std::string& pad = " ")
1258  {
1259  samples->user_report(type, pad);
1260  if (has_complex)
1261  complex_samples->user_report("complex " + type, pad);
1262  }
1263 
1265  {
1266  f.push(top);
1267  f.push("layout");
1268  samples->register_hdf_data(f);
1269  if (has_complex)
1271  f.pop();
1272  f.pop();
1273  if (!f.open_groups())
1274  APP_ABORT("TraceBuffer<" + type +
1275  ">::register_hdf_data() some hdf groups are still open at the end of registration");
1276  hdf_file_pointer = 0;
1277  }
1278 
1279 
1281 
1282 
1283  inline void write_hdf(hdf_archive& f, hsize_t& file_pointer)
1284  {
1285  if (verbose)
1286  app_log() << "TraceBuffer<" << type << ">::write_hdf() " << file_pointer << " " << buffer.size(0) << " "
1287  << buffer.size(1) << std::endl;
1288  dims[0] = buffer.size(0);
1289  dims[1] = buffer.size(1);
1290  if (dims[0] > 0)
1291  {
1292  f.push(top);
1293  h5d_append(f.top(), "traces", file_pointer, buffer.dim(), dims, buffer.data());
1294  f.pop();
1295  }
1296  f.flush();
1297  }
1298 
1299 
1300  inline void test_buffer_write(int sample_size)
1301  {
1302  //check that the size is correct
1303  int ssize = samples->total_size();
1304  if (has_complex)
1305  ssize += complex_samples->total_size();
1306  if (sample_size != ssize)
1307  {
1308  app_log() << "sample_size = " << sample_size << "\ntotal_size = " << ssize << std::endl;
1309  APP_ABORT("TraceBuffer::test_buffer_write sample_size and total_size do not match");
1310  }
1311  //check that buffer indices fall in the expected range
1312  int nsamples = samples->size();
1313  int min_index = samples->min_buffer_index();
1314  int max_index = samples->max_buffer_index();
1315  if (has_complex)
1316  {
1317  nsamples += complex_samples->size();
1318  min_index = std::min(min_index, complex_samples->min_buffer_index());
1319  max_index = std::max(max_index, complex_samples->max_buffer_index());
1320  }
1321  if (nsamples > 0)
1322  {
1323  if (min_index != 0)
1324  APP_ABORT("TraceBuffer::test_buffer_write min_index!=0\n min_index=" << min_index);
1325  if (max_index != sample_size)
1326  APP_ABORT("TraceBuffer::test_buffer_write max_index!=sample_size");
1327  //check that no overlap exists in writes to buffer
1328  Array<int, 2> test_buffer;
1329  test_buffer.resize(1, sample_size);
1330  std::fill(test_buffer.begin(), test_buffer.end(), 0);
1331  int row = 0;
1332  int row_size = test_buffer.size(1);
1333  int offset = row * row_size;
1334  int* loc1 = &test_buffer(offset);
1335  int* loc2 = &test_buffer(row, 0);
1336  if (loc1 != loc2)
1337  APP_ABORT("TraceBuffer::test_buffer_write serialized buffer offset improperly computed");
1338  {
1339  int boffset;
1340  std::vector<TraceSample<T>*>& ordered_samples = samples->ordered_samples;
1341  for (int s = 0; s < ordered_samples.size(); s++)
1342  {
1343  TraceSample<T>& tsample = *ordered_samples[s];
1344  std::vector<T>& sample = tsample.sample;
1345  boffset = offset + tsample.buffer_start;
1346  for (int i = 0; i < sample.size(); ++i)
1347  test_buffer(boffset + i) = 1;
1348  }
1349  }
1350  if (has_complex)
1351  {
1352  int boffset;
1353  std::vector<TraceSample<std::complex<T>>*>& ordered_samples = complex_samples->ordered_samples;
1354  for (int s = 0; s < ordered_samples.size(); s++)
1355  {
1356  TraceSample<std::complex<T>>& tsample = *ordered_samples[s];
1357  std::vector<std::complex<T>>& sample = tsample.sample;
1358  boffset = offset + tsample.buffer_start;
1359  for (int i = 0, ib = 0; i < sample.size(); ++i, ib += tsample.unit_size)
1360  {
1361  test_buffer(boffset + ib) = 1;
1362  test_buffer(boffset + ib + 1) = 1;
1363  }
1364  }
1365  }
1366  //app_log()<<"test_buffer:"<< std::endl;
1367  //for(int i=0;i<row_size;++i)
1368  // app_log()<<" "<<i<<" "<<test_buffer(row,i)<< std::endl;
1369  for (int i = 0; i < row_size; ++i)
1370  if (!test_buffer(row, i))
1371  APP_ABORT("TraceBuffer::test_buffer_write write to row is not contiguous");
1372  }
1373  }
1374 
1375 
1376  inline void test_buffer_collect(int current_row)
1377  {
1378  if (verbose)
1379  app_log() << "TraceBuffer::test_buffer_collect" << std::endl;
1380  std::string scalars = "scalars";
1381  std::map<std::string, std::map<std::string, int>>::iterator dom;
1382  std::map<std::string, int>::iterator var;
1383  std::map<std::string, std::map<std::string, int>>& domains = samples->sample_indices;
1384  std::vector<TraceSample<T>*>& tsamples = samples->samples;
1385  std::map<std::string, int>& scalar_vars = domains[scalars];
1386  for (var = scalar_vars.begin(); var != scalar_vars.end(); var++)
1387  {
1388  const std::string& name = var->first;
1389  TraceSample<T>& sample = *tsamples[var->second];
1390  T value = buffer(current_row, sample.buffer_start);
1391  T svalue = 0;
1392  bool any_present = false;
1393  for (dom = domains.begin(); dom != domains.end(); dom++)
1394  {
1395  const std::string& domain = dom->first;
1396  if (domain != scalars)
1397  {
1398  std::map<std::string, int>& vars = dom->second;
1399  if (vars.count(name) > 0)
1400  {
1401  any_present = true;
1402  TraceSample<T>& ssample = *tsamples[vars[name]];
1403  int start = ssample.buffer_start;
1404  int end = ssample.buffer_end;
1405  for (int i = start; i < end; i++)
1406  svalue += buffer(current_row, i);
1407  }
1408  }
1409  }
1410  if (any_present)
1411  {
1412  if (verbose)
1413  app_log() << " " << name << " " << value << " " << svalue << std::endl;
1414  }
1415  }
1416  }
1417 };
1418 
1419 
1421 {
1422 private:
1423  //collections of samples for a single walker step
1424  // the associated arrays will be updated following evaluate
1428 
1429  //buffers for storing samples
1430  // single row of buffer is a single sample from one walker
1431  // number of rows adjusts to accommodate walker samples
1434 
1435 public:
1436  static double trace_tol;
1437 
1439 
1441  std::string default_domain;
1446  bool verbose;
1447  std::string format;
1449  std::string file_root;
1451  std::unique_ptr<hdf_archive> hdf_file;
1452 
1454  {
1456  master_copy = true;
1457  communicator = comm;
1458  throttle = 1;
1459  format = "hdf";
1460  default_domain = "scalars";
1462  int_buffer.set_type("int");
1463  real_buffer.set_type("real");
1467  }
1468 
1469 
1471  {
1472  if (verbose)
1473  app_log() << "TraceManager::makeClone " << master_copy << std::endl;
1474  if (!master_copy)
1475  APP_ABORT("TraceManager::makeClone only the master copy should call this function");
1476  TraceManager* tm = new TraceManager();
1477  tm->master_copy = false;
1478  tm->transfer_state_from(*this);
1479  tm->distribute();
1480  return tm;
1481  }
1482 
1483 
1484  inline void transfer_state_from(const TraceManager& tm)
1485  {
1487  request = tm.request;
1490  throttle = tm.throttle;
1491  verbose = tm.verbose;
1492  format = tm.format;
1493  hdf_format = tm.hdf_format;
1495  }
1496 
1497 
1498  inline void distribute()
1499  {
1505  }
1506 
1507 
1508  inline void reset_permissions()
1509  {
1510  method_allows_traces = false;
1511  streaming_traces = false;
1512  writing_traces = false;
1513  verbose = false;
1514  hdf_format = false;
1515  request.reset();
1516  }
1517 
1518 
1519  inline void put(xmlNodePtr cur, bool allow_traces, std::string series_root)
1520  {
1522  method_allows_traces = allow_traces;
1523  file_root = series_root;
1524  bool traces_requested = cur != NULL;
1525  streaming_traces = traces_requested && method_allows_traces;
1526  if (streaming_traces)
1527  {
1528  if (omp_get_thread_num() == 0)
1529  {
1530  app_log() << "\n TraceManager::put() " << master_copy << std::endl;
1531  app_log() << " traces requested : " << traces_requested << std::endl;
1532  app_log() << " method allows traces : " << method_allows_traces << std::endl;
1533  app_log() << " streaming traces : " << streaming_traces << std::endl;
1534  app_log() << std::endl;
1535  }
1536  //read trace attributes
1537  std::string writing = "yes";
1538  std::string scalar = "yes";
1539  std::string array = "yes";
1540  std::string scalar_defaults = "yes";
1541  std::string array_defaults = "yes";
1542  std::string verbose_write = "no";
1543  OhmmsAttributeSet attrib;
1544  attrib.add(writing, "write");
1545  attrib.add(scalar, "scalar");
1546  attrib.add(array, "array");
1547  attrib.add(scalar_defaults, "scalar_defaults");
1548  attrib.add(array_defaults, "array_defaults");
1549  attrib.add(format, "format");
1550  attrib.add(throttle, "throttle");
1551  attrib.add(verbose_write, "verbose");
1552  attrib.add(array, "particle"); //legacy
1553  attrib.add(array_defaults, "particle_defaults"); //legacy
1554  attrib.put(cur);
1555  writing_traces = writing == "yes";
1556  bool scalars_on = scalar == "yes";
1557  bool arrays_on = array == "yes";
1558  bool use_scalar_defaults = scalar_defaults == "yes";
1559  bool use_array_defaults = array_defaults == "yes";
1560  verbose = verbose_write == "yes";
1561  format = lowerCase(format);
1562  if (format == "hdf")
1563  {
1564  hdf_format = true;
1565  }
1566  else
1567  {
1568  APP_ABORT("TraceManager::put " + format + " is not a valid file format for traces\n valid options is: hdf");
1569  }
1570 
1571  //read scalar and array elements
1572  // each requests that certain traces be computed
1573  std::set<std::string> scalar_requests;
1574  std::set<std::string> array_requests;
1575  xmlNodePtr element = cur->children;
1576  while (element != NULL)
1577  {
1578  std::string name((const char*)element->name);
1579  if (name == "scalar_traces")
1580  {
1581  std::string defaults = "no";
1582  OhmmsAttributeSet eattrib;
1583  eattrib.add(defaults, "defaults");
1584  eattrib.put(element);
1585  use_scalar_defaults = use_scalar_defaults && defaults == "yes";
1586  if (!use_scalar_defaults)
1587  {
1588  std::vector<std::string> scalar_list;
1589  putContent(scalar_list, element);
1590  scalar_requests.insert(scalar_list.begin(), scalar_list.end());
1591  }
1592  }
1593  else if (name == "array_traces" || name == "particle_traces")
1594  {
1595  std::string defaults = "no";
1596  OhmmsAttributeSet eattrib;
1597  eattrib.add(defaults, "defaults");
1598  eattrib.put(element);
1599  use_array_defaults = use_array_defaults && defaults == "yes";
1600  if (!use_array_defaults)
1601  {
1602  std::vector<std::string> array_list;
1603  putContent(array_list, element);
1604  array_requests.insert(array_list.begin(), array_list.end());
1605  }
1606  }
1607  else if (name != "text")
1608  {
1609  APP_ABORT("TraceManager::put " + name +
1610  " is not a valid sub-element of <trace/>\n valid options are: scalar_traces, array_traces");
1611  }
1612  element = element->next;
1613  }
1614 
1616 
1617  //input user quantity requests into the traces request
1620  request.scalars_on = scalars_on;
1621  request.arrays_on = arrays_on;
1622  request.stream_all_scalars = use_scalar_defaults;
1623  request.stream_all_arrays = use_array_defaults;
1626  request.request_scalar(scalar_requests, writing_traces);
1627  request.request_array(array_requests, writing_traces);
1628 
1629  //distribute verbosity level to buffer and sample objects
1630  distribute();
1631  }
1632 
1633  //streaming_traces = false;
1634  //writing_traces = false;
1635  }
1636 
1637 
1638  inline void update_status()
1639  {
1642  }
1643 
1644 
1645  inline void screen_writes()
1646  {
1650  }
1651 
1652  inline void initialize_traces()
1653  {
1654  if (streaming_traces)
1655  {
1656  if (verbose)
1657  app_log() << "TraceManager::initialize_traces " << master_copy << std::endl;
1658  //organize trace samples and initialize buffers
1659  if (writing_traces)
1660  {
1663  }
1664  }
1665  }
1666 
1667 
1668  inline void finalize_traces()
1669  {
1670  if (verbose)
1671  app_log() << "TraceManager::finalize_traces " << master_copy << std::endl;
1675  }
1676 
1677 
1678  //checkout functions to be used by any OperatorBase or Estimator
1679  // the array checked out should be updated during evaluate
1680  // object calling checkout is responsible for deleting the new array
1681 
1682  // checkout integer arrays
1683  template<int D>
1684  inline Array<TraceInt, D>* checkout_int(const std::string& name, int n1 = 1, int n2 = 0, int n3 = 0, int n4 = 0)
1685  {
1686  return checkout_int<D>(default_domain, name, n1, n2, n3, n4);
1687  }
1688  template<int D>
1689  inline Array<TraceInt, D>* checkout_int(const std::string& domain,
1690  const std::string& name,
1691  int n1 = 1,
1692  int n2 = 0,
1693  int n3 = 0,
1694  int n4 = 0)
1695  {
1696  TinyVector<int, DMAX> shape(n1, n2, n3, n4);
1697  return int_samples.checkout_array<D>(domain, name, shape);
1698  }
1699  template<int D>
1700  inline Array<TraceInt, D>* checkout_int(const std::string& name,
1701  const ParticleSet& P,
1702  int n2 = 0,
1703  int n3 = 0,
1704  int n4 = 0)
1705  {
1706  TinyVector<int, DMAX> shape(P.getTotalNum(), n2, n3, n4);
1707  return int_samples.checkout_array<D>(P, name, shape);
1708  }
1709 
1710 
1711  // checkout real arrays
1712  template<int D>
1713  inline Array<TraceReal, D>* checkout_real(const std::string& name, int n1 = 1, int n2 = 0, int n3 = 0, int n4 = 0)
1714  {
1715  return checkout_real<D>(default_domain, name, n1, n2, n3, n4);
1716  }
1717  template<int D>
1718  inline Array<TraceReal, D>* checkout_real(const std::string& domain,
1719  const std::string& name,
1720  int n1 = 1,
1721  int n2 = 0,
1722  int n3 = 0,
1723  int n4 = 0)
1724  {
1725  TinyVector<int, DMAX> shape(n1, n2, n3, n4);
1726  return real_samples.checkout_array<D>(domain, name, shape);
1727  }
1728  template<int D>
1729  inline Array<TraceReal, D>* checkout_real(const std::string& name,
1730  const ParticleSet& P,
1731  int n2 = 0,
1732  int n3 = 0,
1733  int n4 = 0)
1734  {
1735  TinyVector<int, DMAX> shape(P.getTotalNum(), n2, n3, n4);
1736  return real_samples.checkout_array<D>(P, name, shape);
1737  }
1738 
1739 
1740  // checkout complex arrays
1741  template<int D>
1742  inline Array<TraceComp, D>* checkout_complex(const std::string& name, int n1 = 1, int n2 = 0, int n3 = 0, int n4 = 0)
1743  {
1744  return checkout_complex<D>(default_domain, name, n1, n2, n3, n4);
1745  }
1746  template<int D>
1747  inline Array<TraceComp, D>* checkout_complex(const std::string& domain,
1748  const std::string& name,
1749  int n1 = 1,
1750  int n2 = 0,
1751  int n3 = 0,
1752  int n4 = 0)
1753  {
1754  TinyVector<int, DMAX> shape(n1, n2, n3, n4);
1755  return comp_samples.checkout_array<D>(domain, name, shape);
1756  }
1757  template<int D>
1758  inline Array<TraceComp, D>* checkout_complex(const std::string& name,
1759  const ParticleSet& P,
1760  int n2 = 0,
1761  int n3 = 0,
1762  int n4 = 0)
1763  {
1764  TinyVector<int, DMAX> shape(P.getTotalNum(), n2, n3, n4);
1765  return comp_samples.checkout_array<D>(P, name, shape);
1766  }
1767 
1768 
1769  //get trace functions
1770  // used by estimators that require trace information
1771  inline TraceSample<TraceInt>* get_int_trace(const std::string& name) { return get_int_trace(default_domain, name); }
1772  inline TraceSample<TraceInt>* get_int_trace(const std::string& domain, const std::string& name)
1773  {
1774  return int_samples.get_trace(domain, name);
1775  }
1776  inline TraceSample<TraceInt>* get_int_trace(const ParticleSet& P, const std::string& name)
1777  {
1778  return int_samples.get_trace(P.parentName(), name);
1779  }
1780 
1781  inline TraceSample<TraceReal>* get_real_trace(const std::string& name)
1782  {
1783  return get_real_trace(default_domain, name);
1784  }
1785  inline TraceSample<TraceReal>* get_real_trace(const std::string& domain, const std::string& name)
1786  {
1787  return real_samples.get_trace(domain, name);
1788  }
1789  inline TraceSample<TraceReal>* get_real_trace(const ParticleSet& P, const std::string& name)
1790  {
1791  return real_samples.get_trace(P.parentName(), name);
1792  }
1793 
1794  inline TraceSample<TraceComp>* get_complex_trace(const std::string& name)
1795  {
1796  return get_complex_trace(default_domain, name);
1797  }
1798  inline TraceSample<TraceComp>* get_complex_trace(const std::string& domain, const std::string& name)
1799  {
1800  return comp_samples.get_trace(domain, name);
1801  }
1802  inline TraceSample<TraceComp>* get_complex_trace(const ParticleSet& P, const std::string& name)
1803  {
1804  return comp_samples.get_trace(P.parentName(), name);
1805  }
1806 
1807 
1809  {
1810  return get_int_combined_trace(default_domain, name);
1811  }
1812  inline CombinedTraceSample<TraceInt>* get_int_combined_trace(const std::string& domain, const std::string& name)
1813  {
1814  return int_samples.get_combined_trace(domain, name);
1815  }
1816  inline CombinedTraceSample<TraceInt>* get_int_combined_trace(const ParticleSet& P, const std::string& name)
1817  {
1818  return int_samples.get_combined_trace(P.parentName(), name);
1819  }
1820 
1822  {
1824  }
1825  inline CombinedTraceSample<TraceReal>* get_real_combined_trace(const std::string& domain, const std::string& name)
1826  {
1827  return real_samples.get_combined_trace(domain, name);
1828  }
1829  inline CombinedTraceSample<TraceReal>* get_real_combined_trace(const ParticleSet& P, const std::string& name)
1830  {
1831  return real_samples.get_combined_trace(P.parentName(), name);
1832  }
1833 
1835  {
1837  }
1838  inline CombinedTraceSample<TraceComp>* get_complex_combined_trace(const std::string& domain, const std::string& name)
1839  {
1840  return comp_samples.get_combined_trace(domain, name);
1841  }
1842  inline CombinedTraceSample<TraceComp>* get_complex_combined_trace(const ParticleSet& P, const std::string& name)
1843  {
1844  return comp_samples.get_combined_trace(P.parentName(), name);
1845  }
1846 
1847 
1848  //create combined trace out of existing traces
1849  inline void make_combined_trace(const std::string& name, std::vector<std::string>& names)
1850  {
1851  std::vector<TraceReal> weights;
1852  weights.resize(names.size());
1853  std::fill(weights.begin(), weights.end(), 1.0);
1854  make_combined_trace(name, names, weights);
1855  }
1856 
1857  inline void make_combined_trace(const std::string& name,
1858  std::vector<std::string>& names,
1859  std::vector<TraceReal>& weights)
1860  {
1861  if (streaming_traces)
1862  {
1863  if (verbose)
1864  app_log() << "TraceManager::make_combined_trace " << master_copy << " " << name << std::endl;
1865  real_buffer.make_combined_trace(name, names, weights);
1866  }
1867  }
1868 
1869 
1870  inline void check_clones(std::vector<TraceManager*>& clones)
1871  {
1872  if (writing_traces && clones.size() > 0)
1873  {
1874  if (verbose)
1875  app_log() << "TraceManager::check_clones" << std::endl;
1876  bool all_same = true;
1877  bool int_same;
1878  bool real_same;
1879  TraceManager& ref = *clones[0];
1880  for (int i = 0; i < clones.size(); ++i)
1881  {
1882  TraceManager& tm = *clones[i];
1883  int_same = tm.int_buffer.same_as(ref.int_buffer);
1884  real_same = tm.real_buffer.same_as(ref.real_buffer);
1885  all_same = all_same && int_same && real_same;
1886  }
1887  if (!all_same)
1888  {
1889  for (int i = 0; i < clones.size(); ++i)
1890  clones[i]->write_summary();
1891  APP_ABORT("TraceManager::check_clones trace buffer widths of clones do not match\n contiguous write is "
1892  "impossible\n this was first caused by clones contributing array traces from identical, but "
1893  "differently named, particlesets such as e, e2, e3 ... (fixed)\n please check the TraceManager "
1894  "summaries printed above");
1895  }
1896  }
1897  }
1898 
1899 
1900  inline void reset_buffers()
1901  {
1902  if (writing_traces)
1903  {
1904  if (verbose)
1905  app_log() << "TraceManager::reset_buffers " << master_copy << std::endl;
1906  int_buffer.reset();
1907  real_buffer.reset();
1908  }
1909  }
1910 
1911 
1912  //store the full sample from a single walker step in buffers
1913  inline void buffer_sample(int current_step)
1914  {
1915  if (writing_traces && current_step % throttle == 0)
1916  {
1917  if (verbose)
1918  app_log() << " TraceManager::buffer_sample() " << master_copy << std::endl;
1921  }
1922  }
1923 
1924 
1925  //write buffered trace data to file
1926  inline void write_buffers(std::vector<TraceManager*>& clones, int block)
1927  {
1928  if (master_copy)
1929  {
1930  if (writing_traces)
1931  {
1932  //double tstart = MPI_Wtime();
1933  if (verbose)
1934  app_log() << "TraceManager::write_buffers " << master_copy << std::endl;
1935  if (hdf_format)
1936  {
1937  write_buffers_hdf(clones);
1938  }
1939  }
1940  }
1941  else
1942  APP_ABORT("TraceManager::write_buffers should not be called from non-master copy");
1943  }
1944 
1945 
1946  inline void open_file(std::vector<TraceManager*>& clones)
1947  {
1948  if (master_copy)
1949  {
1950  if (writing_traces)
1951  {
1952  if (verbose)
1953  app_log() << "TraceManager::open_file " << master_copy << std::endl;
1954  if (verbose)
1955  clones[0]->write_summary();
1956  if (hdf_format)
1957  {
1958  open_hdf_file(clones);
1959  }
1960  }
1961  }
1962  else
1963  APP_ABORT("TraceManager::open_file should not be called from non-master copy");
1964  }
1965 
1966 
1967  inline void close_file()
1968  {
1969  if (master_copy)
1970  {
1971  if (writing_traces)
1972  {
1973  if (verbose)
1974  app_log() << "TraceManager::close_file " << master_copy << std::endl;
1975  if (hdf_format)
1976  {
1977  close_hdf_file();
1978  }
1979  }
1980  }
1981  else
1982  APP_ABORT("TraceManager::close_file should not be called from non-master copy");
1983  }
1984 
1985 
1986  inline void startRun(int blocks, std::vector<TraceManager*>& clones)
1987  {
1988  if (verbose)
1989  app_log() << "TraceManager::startRun " << master_copy << std::endl;
1990  if (master_copy)
1991  {
1993  check_clones(clones);
1994  open_file(clones);
1995  }
1996  else
1997  APP_ABORT("TraceManager::startRun should not be called from non-master copy");
1998  }
1999 
2000 
2001  inline void stopRun()
2002  {
2003  if (verbose)
2004  app_log() << "TraceManager::stopRun " << master_copy << std::endl;
2005  if (master_copy)
2006  {
2007  close_file();
2008  finalize_traces();
2009  }
2010  else
2011  APP_ABORT("TraceManager::stopRun should not be called from non-master copy");
2012  }
2013 
2014 
2015  inline void startBlock(int nsteps)
2016  {
2017  if (verbose)
2018  app_log() << "TraceManager::startBlock " << master_copy << std::endl;
2019  reset_buffers();
2020  }
2021 
2022 
2023  inline void stopBlock()
2024  {
2025  if (verbose)
2026  app_log() << "TraceManager::stopBlock " << master_copy << std::endl;
2027  }
2028 
2029 
2030  inline void write_summary(std::string pad = " ")
2031  {
2032  std::string pad2 = pad + " ";
2033  app_log() << std::endl;
2034  app_log() << pad << "TraceManager (detailed summary)" << std::endl;
2035  app_log() << pad2 << "master_copy = " << master_copy << std::endl;
2036  app_log() << pad2 << "method_allows_traces = " << method_allows_traces << std::endl;
2037  app_log() << pad2 << "streaming_traces = " << streaming_traces << std::endl;
2038  app_log() << pad2 << "writing_traces = " << writing_traces << std::endl;
2039  app_log() << pad2 << "format = " << format << std::endl;
2040  app_log() << pad2 << "hdf format = " << hdf_format << std::endl;
2041  app_log() << pad2 << "default_domain = " << default_domain << std::endl;
2042  int_buffer.write_summary(pad2);
2043  real_buffer.write_summary(pad2);
2044  app_log() << pad << "end TraceManager" << std::endl;
2045  }
2046 
2047 
2048  inline void user_report(std::string pad = " ")
2049  {
2050  std::string pad2 = pad + " ";
2051  std::string pad3 = pad2 + " ";
2052  app_log() << std::endl;
2053  app_log() << pad << "Traces report" << std::endl;
2054  request.report();
2055  app_log() << pad2 << "Type and domain breakdown of streaming quantities:" << std::endl;
2056  std::set<std::string>::iterator req;
2057  int_buffer.user_report(pad3);
2058  real_buffer.user_report(pad3);
2059  app_log() << std::endl;
2060  //if(verbose)
2061  // write_summary(pad);
2062  }
2063 
2064  //hdf file operations
2065  inline void open_hdf_file(std::vector<TraceManager*>& clones)
2066  {
2067  if (clones.size() == 0)
2068  APP_ABORT("TraceManager::open_hdf_file no trace clones exist, cannot open file");
2069  int nprocs = communicator->size();
2070  int rank = communicator->rank();
2071  std::array<char, 32> ptoken;
2072  std::string file_name = file_root;
2073  if (nprocs > 1)
2074  {
2075  int length{0};
2076  if (nprocs > 10000)
2077  length = std::snprintf(ptoken.data(), ptoken.size(), ".p%05d", rank);
2078  else if (nprocs > 1000)
2079  length = std::snprintf(ptoken.data(), ptoken.size(), ".p%04d", rank);
2080  else
2081  length = std::snprintf(ptoken.data(), ptoken.size(), ".p%03d", rank);
2082  if (length < 0)
2083  throw std::runtime_error("Error generating filename");
2084  file_name.append(ptoken.data(), length);
2085  }
2086  file_name += ".traces.h5";
2087  if (verbose)
2088  app_log() << "TraceManager::open_hdf_file opening traces hdf file " << file_name << std::endl;
2089  hdf_file = std::make_unique<hdf_archive>();
2090  bool successful = hdf_file->create(file_name);
2091  if (!successful)
2092  APP_ABORT("TraceManager::open_hdf_file failed to open hdf file " + file_name);
2093  // only clones have active buffers and associated data
2094  TraceManager& tm = *clones[0];
2095  //tm.write_summary();
2098  }
2099 
2100 
2101  inline void write_buffers_hdf(std::vector<TraceManager*>& clones)
2102  {
2103  if (verbose)
2104  app_log() << "TraceManager::write_buffers_hdf " << master_copy << std::endl;
2105  for (int ip = 0; ip < clones.size(); ++ip)
2106  {
2107  TraceManager& tm = *clones[ip];
2110  }
2111  }
2112 
2113  inline void close_hdf_file() { hdf_file.reset(); }
2114 };
2115 
2116 
2117 } // namespace qmcplusplus
2118 
2119 
2120 #else
2121 
2122 
2123 // make a vacuous class for TraceManager for lighter compilation
2124 // disabling TraceManager should not affect other runtime behavior
2125 
2126 #include "Particle/ParticleSet.h"
2127 
2128 namespace qmcplusplus
2129 {
2130 using TraceInt = long;
2131 using TraceReal = double;
2132 using TraceComp = std::complex<TraceReal>;
2133 
2134 struct TraceRequest
2135 {
2137 
2138  inline TraceRequest() { streaming_default_scalars = false; }
2139 
2140  //inline bool screen_sample(const std::string& domain,const std::string& name,bool& write) { return false; }
2141  inline bool streaming_scalar(const std::string& name) { return false; }
2142  inline bool streaming_array(const std::string& name) { return false; }
2143  inline bool streaming(const std::string& name) { return false; }
2144  //inline bool quantity_present(const std::string& name) { return false; }
2145  inline bool streaming() { return false; }
2146  //inline bool writing() { return false; }
2147  //inline bool streaming_scalars() { return false; }
2148  //inline bool streaming_arrays() { return false; }
2149 
2150 
2151  //inline void set_scalar_domain(const std::string& domain) { }
2152  inline void request_scalar(const std::string& name, bool write = false) {}
2153  inline void request_array(const std::string& name, bool write = false) {}
2154  //inline void request_scalar(const std::set<std::string>& names,bool write=false) { }
2155  //inline void request_array(const std::set<std::string>& names,bool write=false) { }
2156  inline void incorporate(TraceRequest& other) {}
2157  inline void determine_stream_write() {}
2158  inline void relay_stream_info(TraceRequest& other) {}
2159  //inline void report() { }
2160  //inline void write_selected(const std::string& header,const std::string& selector) { }
2161  //inline void guarantee_presence(const std::string& name,bool combined=false) { }
2162  //inline void check_presence(const std::string& name) { }
2163  inline void contribute_scalar(const std::string& name, bool default_quantity = false) {}
2164  inline void contribute_array(const std::string& name, bool default_quantity = false) {}
2165  inline void contribute_combined(const std::string& name,
2166  std::vector<std::string>& deps,
2167  bool scalar = false,
2168  bool array = false,
2169  bool default_quantity = false)
2170  {}
2171 };
2172 
2173 
2174 template<typename T>
2175 struct TraceSample
2176 {
2177  std::vector<T>& sample;
2178 };
2179 
2180 
2181 template<typename T>
2182 struct CombinedTraceSample : public TraceSample<T>
2183 {
2184  inline void combine() {}
2185 };
2186 
2187 
2188 struct TraceManager
2189 {
2190  TraceRequest request;
2191  bool streaming_traces;
2192 
2193 
2194  TraceManager(Communicate* comm = 0) { streaming_traces = false; }
2195 
2196  inline TraceManager* makeClone() { return new TraceManager(); }
2197 
2198  inline void transfer_state_from(const TraceManager& tm) {}
2199  //inline void distribute() { }
2200  //inline void reset_permissions() { }
2201  inline void put(xmlNodePtr cur, bool allow_traces, std::string series_root) {}
2202  inline void update_status() {}
2203  inline void screen_writes() {}
2204  inline void initialize_traces() {}
2205  inline void finalize_traces() {}
2206 
2207  template<int D>
2208  inline Array<TraceInt, D>* checkout_int(const std::string& name, int n1 = 1, int n2 = 0, int n3 = 0, int n4 = 0)
2209  {
2210  return 0;
2211  }
2212  //template<int D> inline Array<TraceInt,D>* checkout_int( const std::string& domain,const std::string& name, int n1=1,int n2=0,int n3=0,int n4=0) { return 0; }
2213  template<int D>
2214  inline Array<TraceInt, D>* checkout_int(const std::string& name,
2215  const ParticleSet& P,
2216  int n2 = 0,
2217  int n3 = 0,
2218  int n4 = 0)
2219  {
2220  return 0;
2221  }
2222  template<int D>
2223  inline Array<TraceReal, D>* checkout_real(const std::string& name, int n1 = 1, int n2 = 0, int n3 = 0, int n4 = 0)
2224  {
2225  return 0;
2226  }
2227  //template<int D> inline Array<TraceReal,D>* checkout_real( const std::string& domain,const std::string& name, int n1=1,int n2=0,int n3=0,int n4=0) { return 0; }
2228  template<int D>
2229  inline Array<TraceReal, D>* checkout_real(const std::string& name,
2230  const ParticleSet& P,
2231  int n2 = 0,
2232  int n3 = 0,
2233  int n4 = 0)
2234  {
2235  return 0;
2236  }
2237  //template<int D> inline Array<TraceComp,D>* checkout_complex(const std::string& name, int n1=1,int n2=0,int n3=0,int n4=0) { return 0; }
2238  //template<int D> inline Array<TraceComp,D>* checkout_complex(const std::string& domain,const std::string& name, int n1=1,int n2=0,int n3=0,int n4=0) { return 0; }
2239  template<int D>
2240  inline Array<TraceComp, D>* checkout_complex(const std::string& name,
2241  const ParticleSet& P,
2242  int n2 = 0,
2243  int n3 = 0,
2244  int n4 = 0)
2245  {
2246  return 0;
2247  }
2248 
2249  //inline TraceSample<TraceInt>* get_int_trace(const std::string& name) { return 0; }
2250  //inline TraceSample<TraceInt>* get_int_trace(const std::string& domain, const std::string& name) { return 0; }
2251  //inline TraceSample<TraceInt>* get_int_trace(const ParticleSet& P, const std::string& name) { return 0; }
2252  inline TraceSample<TraceReal>* get_real_trace(const std::string& name) { return 0; }
2253  //inline TraceSample<TraceReal>* get_real_trace(const std::string& domain, const std::string& name) { return 0; }
2254  inline TraceSample<TraceReal>* get_real_trace(const ParticleSet& P, const std::string& name) { return 0; }
2255  //inline TraceSample<TraceComp>* get_complex_trace(const std::string& name) { return 0; }
2256  //inline TraceSample<TraceComp>* get_complex_trace(const std::string& domain, const std::string& name) { return 0; }
2257  inline TraceSample<TraceComp>* get_complex_trace(const ParticleSet& P, const std::string& name) { return 0; }
2258 
2259  //inline CombinedTraceSample<TraceInt>* get_int_combined_trace(const std::string& name) { return 0; }
2260  //inline CombinedTraceSample<TraceInt>* get_int_combined_trace(const std::string& domain, const std::string& name) { return 0; }
2261  //inline CombinedTraceSample<TraceInt>* get_int_combined_trace(const ParticleSet& P, const std::string& name) { return 0; }
2262  //inline CombinedTraceSample<TraceReal>* get_real_combined_trace(const std::string& name) { return 0; }
2263  //inline CombinedTraceSample<TraceReal>* get_real_combined_trace(const std::string& domain, const std::string& name) { return 0; }
2264  inline CombinedTraceSample<TraceReal>* get_real_combined_trace(const ParticleSet& P, const std::string& name)
2265  {
2266  return 0;
2267  }
2268  //inline CombinedTraceSample<TraceComp>* get_complex_combined_trace(const std::string& name) { return 0; }
2269  //inline CombinedTraceSample<TraceComp>* get_complex_combined_trace(const std::string& domain, const std::string& name) { return 0; }
2270  //inline CombinedTraceSample<TraceComp>* get_complex_combined_trace(const ParticleSet& P, const std::string& name) { return 0; }
2271 
2272  inline void make_combined_trace(const std::string& name, std::vector<std::string>& names) {}
2273  //inline void make_combined_trace(const std::string& name,std::vector<std::string>& names,std::vector<TraceReal>& weights) { }
2274  //inline void check_clones(std::vector<TraceManager*>& clones) { }
2275  //inline void reset_buffers() { }
2276  inline void buffer_sample(int current_step) {}
2277  inline void write_buffers(std::vector<TraceManager*>& clones, int block) {}
2278  //inline void open_file(std::vector<TraceManager*>& clones) { }
2279  //inline void close_file() { }
2280  inline void startRun(int blocks, std::vector<TraceManager*>& clones) {}
2281  inline void stopRun() {}
2282  inline void startBlock(int nsteps) {}
2283  inline void stopBlock() {}
2284  //inline void write_summary( std::string pad=" ") { }
2285  inline void user_report(std::string pad = " ") {}
2286 };
2287 
2288 
2289 } // namespace qmcplusplus
2290 
2291 #endif
2292 
2293 
2294 #endif
TinyVector< int, DMAX > shape
Definition: TraceManager.h:525
void resize(size_type n, Type_t val=Type_t())
Resize the container.
Definition: OhmmsVector.h:166
TraceSample< TraceComp > * get_complex_trace(const std::string &domain, const std::string &name)
Array< TraceInt, D > * checkout_int(const std::string &name, const ParticleSet &P, int n2=0, int n3=0, int n4=0)
CombinedTraceSample< T > * get_combined_trace(const std::string &domain, const std::string &name)
Definition: TraceManager.h:837
void set_samples(TraceSamples< T > &s)
TraceBuffer< TraceInt > int_buffer
void put(xmlNodePtr cur, bool allow_traces, std::string series_root)
Array< TraceComp, D > * checkout_complex(const std::string &name, int n1=1, int n2=0, int n3=0, int n4=0)
bool quantity_present(const std::string &name)
Definition: TraceManager.h:472
void write(T &data, const std::string &aname)
write the data to the group aname and check status runtime error is issued on I/O error ...
Definition: hdf_archive.h:259
void initialize(const std::string &sdomain, const std::string &sname, int sindex, int sdim)
Definition: TraceManager.h:557
bool h5d_append(hid_t grp, const std::string &aname, hsize_t &current, hsize_t ndims, const hsize_t *const dims, const T *const first, hsize_t chunk_size=1, hid_t xfer_plist=H5P_DEFAULT)
TraceSamples< TraceReal > real_samples
void delete_iter(IT first, IT last)
delete the pointers in [first,last)
std::map< std::string, TraceQuantity > quantities
Definition: TraceManager.h:133
helper functions for EinsplineSetBuilder
Definition: Configuration.h:43
void open_hdf_file(std::vector< TraceManager *> &clones)
hid_t top() const
return the top of the group stack
Definition: hdf_archive.h:192
int rank() const
return the rank
Definition: Communicate.h:116
void incorporate(const TraceQuantity &other)
Definition: TraceManager.h:85
void write_hdf(hdf_archive &f)
bool streaming_array(const std::string &name)
Definition: TraceManager.h:249
bool open_groups()
check if any groups are open group stack will have entries if so
Definition: hdf_archive.h:198
void write_buffers(std::vector< TraceManager *> &clones, int block)
void register_hdf_data(hdf_archive &f)
Definition: TraceManager.h:991
TraceSample< TraceReal > * get_real_trace(const ParticleSet &P, const std::string &name)
size_t getTotalNum() const
Definition: ParticleSet.h:493
void startRun(int blocks, std::vector< TraceManager *> &clones)
#define OHMMS_PRECISION
Definition: config.h:70
std::ostream & app_log()
Definition: OutputManager.h:65
bool put(xmlNodePtr cur)
assign attributes to the set
Definition: AttributeSet.h:55
Type_t * data()
Definition: OhmmsArray.h:87
CombinedTraceSample< TraceComp > * get_complex_combined_trace(const std::string &name)
CombinedTraceSample< TraceComp > * get_complex_combined_trace(const std::string &domain, const std::string &name)
std::vector< Vector< T > * > combined_sample_vectors
Definition: TraceManager.h:757
void user_report(const std::string &pad=" ")
Array< TraceComp, D > * checkout_complex(const std::string &name, const ParticleSet &P, int n2=0, int n3=0, int n4=0)
void buffer_sample(int current_step)
Array< TraceComp, D > * checkout_complex(const std::string &domain, const std::string &name, int n1=1, int n2=0, int n3=0, int n4=0)
void set_type(std::string stype)
class to handle hdf file
Definition: hdf_archive.h:51
bool streaming_scalar(const std::string &name)
Definition: TraceManager.h:242
void write_selected(const std::string &header, const std::string &selector)
Definition: TraceManager.h:427
void write_hdf(hdf_archive &f, hsize_t &file_pointer)
CombinedTraceSample< TraceInt > * get_int_combined_trace(const ParticleSet &P, const std::string &name)
std::vector< TraceSample< T > * > ordered_samples
Definition: TraceManager.h:755
virtual bool is_combined()
Definition: TraceManager.h:605
Array< TraceInt, D > * checkout_int(const std::string &name, int n1=1, int n2=0, int n3=0, int n4=0)
void resize(const std::array< SIZET, D > &dims)
Resize the container.
Definition: OhmmsArray.h:65
int size() const
return the number of tasks
Definition: Communicate.h:118
void check_presence(const std::string &name)
Definition: TraceManager.h:492
bool same_as(TraceBuffer< T > &ref)
void startBlock(int nsteps)
CombinedTraceSample< TraceInt > * get_int_combined_trace(const std::string &domain, const std::string &name)
void write_summary(int ind=-1, std::string pad=" ")
Definition: TraceManager.h:627
T min(T a, T b)
Array< TraceReal, D > * checkout_real(const std::string &domain, const std::string &name, int n1=1, int n2=0, int n3=0, int n4=0)
CombinedTraceSample< TraceReal > * get_real_combined_trace(const std::string &domain, const std::string &name)
Declaration of OhmmsElementBase and define xml-related macros.
bool streaming(const std::string &name)
Definition: TraceManager.h:256
void copy(const Array< T1, 3 > &src, Array< T2, 3 > &dest)
Definition: Blitz.h:639
std::vector< TraceSample< T > * > components
Definition: TraceManager.h:654
Wrapping information on parallelism.
Definition: Communicate.h:68
void write_summary(std::string type, std::string pad=" ")
omp_int_t omp_get_thread_num()
Definition: OpenMP.h:25
TraceSample< TraceReal > * get_real_trace(const std::string &domain, const std::string &name)
Specialized paritlce class for atomistic simulations.
Definition: ParticleSet.h:55
Container_t & storage()
Definition: OhmmsArray.h:60
TraceManager(Communicate *comm=0)
void set_unit_size(int usize)
Definition: TraceManager.h:569
size_type size() const
return the current size
Definition: OhmmsVector.h:162
bool TraceSample_comp(TraceSample< T > *left, TraceSample< T > *right)
Definition: TraceManager.h:744
std::map< std::string, TraceInt > meta_int
Definition: TraceManager.h:529
void make_combined_trace(const std::string &name, std::vector< std::string > &names)
void write_summary(std::string pad=" ")
void make_combined_trace(const std::string &name, std::vector< std::string > &names, std::vector< TraceReal > &weights)
TraceSample< TraceInt > * get_int_trace(const ParticleSet &P, const std::string &name)
CombinedTraceSample(const std::string &sdomain, const std::string &sname, int sindex, int sdim, TinyVector< int, DMAX > sshape, Vector< T > &ssample)
Definition: TraceManager.h:668
TraceSamples< TraceInt > int_samples
void transfer_state_from(const TraceManager &tm)
void open_file(std::vector< TraceManager *> &clones)
class to handle a set of attributes of an xmlNode
Definition: AttributeSet.h:24
Array< T, D > * checkout_array(const ParticleSet &P, const std::string &name, TinyVector< int, DMAX > shape)
Definition: TraceManager.h:801
void test_buffer_collect(int current_row)
void contribute_scalar(const std::string &name, bool default_quantity=false)
Definition: TraceManager.h:189
std::vector< TraceReal > weights
Definition: TraceManager.h:653
std::string lowerCase(const std::string_view s)
++17
const std::string & parentName() const
return parent&#39;s name
Definition: ParticleSet.h:236
bool same_shape(TraceSample< T > *other)
Definition: TraceManager.h:596
size_t size() const
Definition: OhmmsArray.h:57
TraceSample< TraceInt > * get_int_trace(const std::string &name)
#define APP_ABORT(msg)
Widely used but deprecated fatal error macros from legacy code.
Definition: AppAbort.h:27
void user_report(std::string pad=" ")
void make_combined_trace(const std::string &name, std::vector< std::string > &names, std::vector< TraceReal > &weights)
TraceSample< TraceComp > * get_complex_trace(const ParticleSet &P, const std::string &name)
TraceSample(const std::string &sdomain, const std::string &sname, int sindex, int sdim, Vector< T > &ssample)
Definition: TraceManager.h:534
void incorporate(TraceRequest &other)
Definition: TraceManager.h:282
void write_buffers_hdf(std::vector< TraceManager *> &clones)
void user_report(const std::string &type, const std::string &pad=" ")
TraceManager * makeClone()
void set_unit_size(int usize)
Definition: TraceManager.h:892
std::vector< CombinedTraceSample< T > * > combined_samples
Definition: TraceManager.h:756
TraceSamples< TraceComp > comp_samples
TraceSamples< std::complex< T > > * complex_samples
TraceSample< TraceInt > * get_int_trace(const std::string &domain, const std::string &name)
void set_buffer_range(int &bstart)
Definition: TraceManager.h:607
CombinedTraceSample< TraceReal > * get_real_combined_trace(const ParticleSet &P, const std::string &name)
void request_array(const std::string &name, bool write=false)
Definition: TraceManager.h:233
Array< TraceReal, D > * checkout_real(const std::string &name, int n1=1, int n2=0, int n3=0, int n4=0)
void push(const std::string &gname, bool createit=true)
push a group to the group stack
double scalar
Definition: Blitz.h:190
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
TraceSample(const std::string &sdomain, const std::string &sname, int sindex, int sdim, TinyVector< int, DMAX > sshape, Vector< T > &ssample)
Definition: TraceManager.h:541
void contribute_combined(const std::string &name, std::vector< std::string > &deps, bool scalar=false, bool array=false, bool default_quantity=false)
Definition: TraceManager.h:207
void request_scalar(const std::set< std::string > &names, bool write=false)
Definition: TraceManager.h:266
std::vector< TraceSample< T > * > samples
Definition: TraceManager.h:753
Array< TraceInt, D > * checkout_int(const std::string &domain, const std::string &name, int n1=1, int n2=0, int n3=0, int n4=0)
void relay_stream_info(TraceRequest &other)
Definition: TraceManager.h:371
Array< TraceReal, D > * checkout_real(const std::string &name, const ParticleSet &P, int n2=0, int n3=0, int n4=0)
sycl::event copy_n(sycl::queue &aq, const T1 *restrict VA, size_t array_size, T2 *restrict VC, const std::vector< sycl::event > &events)
Definition: syclBLAS.cpp:548
void assign_sample_index(const std::string &domain, const std::string &name, int index, std::string label="")
Definition: TraceManager.h:771
TraceSample< TraceComp > * get_complex_trace(const std::string &name)
void set_buffer_ranges(int &starting_index)
Definition: TraceManager.h:926
std::map< std::string, TraceReal > meta_real
Definition: TraceManager.h:530
CombinedTraceSample< TraceInt > * get_int_combined_trace(const std::string &name)
void request_array(const std::set< std::string > &names, bool write=false)
Definition: TraceManager.h:274
void contribute_array(const std::string &name, bool default_quantity=false)
Definition: TraceManager.h:198
void register_hdf_data(hdf_archive &f)
std::complex< TraceReal > TraceComp
Definition: TraceManager.h:50
Container_t::iterator begin()
Definition: OhmmsArray.h:80
CombinedTraceSample(const std::string &sdomain, const std::string &sname, int sindex, int sdim, Vector< T > &ssample)
Definition: TraceManager.h:657
void set_samples(TraceSamples< std::complex< T >> &s)
void write_summary(std::string pad=" ")
unsigned dim() const
Definition: OhmmsArray.h:55
void test_buffer_write(int sample_size)
void add_component(TraceSample< T > *component, TraceReal weight)
Definition: TraceManager.h:684
void write_summary_combined(int ind, std::string pad=" ")
Definition: TraceManager.h:724
std::map< std::string, std::map< std::string, int > > sample_indices
Definition: TraceManager.h:754
std::map< std::string, std::set< std::string > > combined_dependencies
Definition: TraceManager.h:136
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
TraceSamples< T > * samples
bool screen_sample(const std::string &domain, const std::string &name, bool &write)
Definition: TraceManager.h:164
const unsigned int DMAX
Definition: TraceManager.h:47
void set_scalar_domain(const std::string &domain)
Definition: TraceManager.h:162
std::unique_ptr< hdf_archive > hdf_file
TraceSample< TraceReal > * get_real_trace(const std::string &name)
TraceBuffer< TraceReal > real_buffer
void screen_writes(TraceRequest &request)
Definition: TraceManager.h:899
void request_scalar(const std::string &name, bool write=false)
Definition: TraceManager.h:224
bool make_combined_trace(const std::string &name, std::vector< std::string > &names, std::vector< TraceReal > &weights)
Definition: TraceManager.h:853
TraceSample< T > * get_trace(const std::string &domain, const std::string &name)
Definition: TraceManager.h:819
CombinedTraceSample< TraceComp > * get_complex_combined_trace(const ParticleSet &P, const std::string &name)
void check_clones(std::vector< TraceManager *> &clones)
Array< T, D > * checkout_array(const std::string &domain, const std::string &name, TinyVector< int, DMAX > shape)
Definition: TraceManager.h:784
void flush()
flush a file
Definition: hdf_archive.h:146
void guarantee_presence(const std::string &name, bool combined=false)
Definition: TraceManager.h:475
std::map< std::string, std::string > meta_string
Definition: TraceManager.h:531
virtual ~TraceSample()=default
Container_t::iterator end()
Definition: OhmmsArray.h:81
OHMMS_PRECISION TraceReal
Definition: TraceManager.h:49