QMCPACK
test_OneBodyDensityMatrices.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) 2024 QMCPACK developers.
6 //
7 // File developed by: Peter Doak, doakpw@ornl.gov, Oak Ridge National Lab
8 //
9 // File created by: Peter Doak, doakpw@ornl.gov, Oak Ridge National Lab
10 //////////////////////////////////////////////////////////////////////////////////////
11 
12 /** \file
13  * Tests for the OneBodyDensityMatrices.
14  * We can't reason about the state of the global Random in tests. A User can run only some tests,
15  * new tests will get added, other tests modified so global Random is called more times or fewer.
16  * Use of FakeRandom in unit tests in other tests of this executable can result in ambiguity
17  * about which global Random this test will have have access to as well.
18  * But several ParticleSet functions use the global Random so we have to avoid the normal
19  * sequence of particleset state transforms and set particle positions explicitly.
20  */
21 
22 #include "catch.hpp"
23 
24 #include <iostream>
25 #include "OneBodyDensityMatrices.h"
28 #include "EstimatorTesting.h"
30 #include "EstimatorInput.h"
31 #include "ParticleSet.h"
32 #include "TrialWaveFunction.h"
33 #include "OhmmsData/Libxml2Doc.h"
37 #include "Utilities/StdRandom.h"
39 #include "Utilities/ProjectData.h"
41 #include "CPU/math.hpp"
42 
43 namespace qmcplusplus
44 {
45 
46 /** set to true to generate new random R for particles.
47  * EXPECT all tests to fail in this case!
48  */
49 constexpr bool generate_test_data = false;
50 // set to true to dump obdm for same R but perhaps different RNG effecting the integration.
51 constexpr bool dump_obdm = false;
52 
53 namespace testing
54 {
56 
57 template<typename T>
59 {
60 public:
67 
69  {
70  OneBodyDensityMatrices obdm2(obdm);
71  CHECK(obdm.sampling_ == obdm2.sampling_);
72  CHECK(obdm.data_ != obdm2.data_);
73  }
74 
75  OneBodyDensityMatricesTests() = default;
79  StdRandom<T>& rng)
80  {
82  switch (input)
83  {
84  case (Input::valid::VANILLA):
85  obdm.generateSamples(1.0, pset_target, rng);
86  CHECK(obdm.nmoves_ == 64);
87  break;
88  case (Input::valid::SCALE):
89  obdm.generateSamples(1.0, pset_target, rng);
90  CHECK(obdm.nmoves_ == 0);
91  break;
92  case (Input::valid::GRID):
93  obdm.generateSamples(1.0, pset_target, rng);
94  CHECK(obdm.nmoves_ == 0);
95  CHECK(obdm.samples_ == pow(22, OHMMS_DIM));
96  break;
97  }
98  }
99 
103  StdRandom<T>& rng)
104  {
106  //confirm that the spins were sampled, will be zero if not
107  CHECK(obdm.is_spinor_ == true);
108  CHECK(std::abs(obdm.spcur_) > 1e-8);
109  }
110 
111  /** Checking approximate equality for complex valued data as
112  * two reals is not consistent with testing practices elsewhere in the code.
113  * Values that are slightly off may now fall in approximation limit properly.
114  *
115  * The smell from OperatorEstBase continuing the tradition of
116  * packing everything into a Real buffer is not lost on me.
117  */
118  void checkData(Real* ref_in, Real* test_in, size_t size)
119  {
121  {
122  size /= 2;
123  auto* ref_data = reinterpret_cast<std::complex<Real>*>(ref_in);
124  auto* test_data = reinterpret_cast<std::complex<Real>*>(test_in);
125  for (size_t id = 0; id < size; id += 2)
126 #if defined(MIXED_PRECISION)
127  CHECK(ref_data[id] == ComplexApprox(test_data[id]).epsilon(6e-4));
128 #else
129  CHECK(ref_data[id] == ComplexApprox(test_data[id]));
130 #endif
131  }
132  else
133  {
134  for (size_t id = 0; id < size; ++id)
135 #if defined(MIXED_PRECISION)
136  CHECK(ref_in[id] == Approx(test_in[id]).epsilon(1e-4));
137 #else
138  CHECK(ref_in[id] == Approx(test_in[id]));
139 #endif
140  }
141  }
142 
143  /** no change test for accumulate
144  */
147  RefVector<ParticleSet>& psets,
149  StdRandom<T>& rng)
150  {
151  obdm.implAccumulate(walkers, psets, twfcs, rng);
152  Data data(getAccumulateData());
153  auto& returned_data = obdm.data_;
154  if constexpr (generate_test_data)
155  FAIL_CHECK("Test always fails when generating new test reference data.");
156  else
157  checkData(data.data(), returned_data.data(), data.size());
158  }
159 
160  /** no change test for evaluateMatrix.
161  */
163  ParticleSet& pset,
164  TrialWaveFunction& trial_wavefunction,
165  MCPWalker& walker,
166  StdRandom<T>& rng)
167  {
168  obdm.evaluateMatrix(pset, trial_wavefunction, walker, rng);
170  auto& returned_data = obdm.data_;
171  if constexpr (generate_test_data)
172  FAIL_CHECK("Test always fails when generating new test reference data.");
173  else
174  checkData(returned_data.data(), data.data(), data.size());
175  }
176 
178  {
179  //this test is just going to set some arbitrary data, not actually calculate anything.
180  //then we will write this data to the hdf5
181  //then we will open and read the hdf5 and make sure the up and down data are properly written
182 
183  obdm.data_ = getUniqueSpinData();
184  hdf_archive hd;
185  std::string test_file{"1rdm_test.hdf"};
186  bool okay = hd.create(test_file);
187  REQUIRE(okay);
188 
189  obdm.registerOperatorEstimator(hd);
190  obdm.write(hd);
191  hd.close();
192 
193  hdf_archive hd_read;
194  okay = hd_read.open(test_file);
195  REQUIRE(okay);
196  Data up_data;
197  Data dn_data;
198 
199  hdf_path up_path = {"OneBodyDensityMatrices/number_matrix/u/value"};
200  hdf_path dn_path = {"OneBodyDensityMatrices/number_matrix/d/value"};
201 
202  size_t nb = obdm.basis_size_;
203  size_t down_offset = 0;
205  {
206  std::array<size_t, 4> shape = {1, nb, nb, 2};
207  hd_read.readSlabReshaped(up_data, shape, up_path.string());
208  hd_read.readSlabReshaped(dn_data, shape, dn_path.string());
209  down_offset = nb * nb * 2; //for real and imag
210  }
211  else if constexpr (std::is_floating_point<OneBodyDensityMatrices::Value>::value)
212  {
213  std::array<size_t, 3> shape = {1, nb, nb};
214  hd_read.readSlabReshaped(up_data, shape, up_path.string());
215  hd_read.readSlabReshaped(dn_data, shape, dn_path.string());
216  down_offset = nb * nb;
217  }
218  //The data in obdm is the reference data, we need to check that
219  //what was read from the hdf5 is consistent
220  checkData(obdm.data_.data(), up_data.data(), up_data.size());
221  checkData(obdm.data_.data() + down_offset, dn_data.data(), dn_data.size());
222  }
223 
225  {
226  std::cout << "Here is what is in your OneBodyDensityMatrices:\n" << NativePrint(obdm.data_) << '\n';
227  }
228 
229 private:
233 };
234 
235 } // namespace testing
236 
237 
238 TEST_CASE("OneBodyDensityMatrices::OneBodyDensityMatrices", "[estimators]")
239 {
243  bool okay = doc.parseFromString(Input::xml[Input::valid::VANILLA]);
244  if (!okay)
245  throw std::runtime_error("cannot parse OneBodyDensitMatricesInput section");
246  xmlNodePtr node = doc.getRoot();
249  auto species_set = testing::makeSpeciesSet(SpeciesCases::GOOD);
250 
252  Communicate* comm;
254 
256  auto wavefunction_pool =
258  auto& pset_target = *(particle_pool.getParticleSet("e"));
259  auto& spo_map = wavefunction_pool.getWaveFunction("wavefunction")->getSPOMap();
260 
261  // Good constructor
262  OneBodyDensityMatrices obdm(std::move(obdmi), lattice, species_set, spo_map, pset_target);
263  // make sure there is something in obdm's data
265  oeba[0] = 1.0;
267  obdmt.testCopyConstructor(obdm);
268 
269  species_set = testing::makeSpeciesSet(SpeciesCases::NO_MEMBERSIZE);
270  CHECK_THROWS_AS(OneBodyDensityMatrices(std::move(obdmi), lattice, species_set, spo_map, pset_target),
272 }
273 
274 TEST_CASE("OneBodyDensityMatrices::generateSamples", "[estimators]")
275 {
277 
279 
282 
284  auto wavefunction_pool =
286  auto& pset_target = *(particle_pool.getParticleSet("e"));
287  auto& species_set = pset_target.getSpeciesSet();
288  auto& spo_map = wavefunction_pool.getWaveFunction("wavefunction")->getSPOMap();
289 
290  auto samplingCaseRunner = [&pset_target, &species_set, &spo_map](Input::valid test_case) {
292 
293  bool okay = doc.parseFromString(Input::xml[test_case]);
294  if (!okay)
295  throw std::runtime_error("cannot parse OneBodyDensitMatricesInput section");
296  xmlNodePtr node = doc.getRoot();
298 
299  OneBodyDensityMatrices obDenMat(std::move(obdmi), pset_target.getLattice(), species_set, spo_map, pset_target);
300 
302  //Get control over which rng is used.
303  //we don't want FakeRandom.
305  obdmt.testGenerateSamples(test_case, obDenMat, pset_target, rng);
306  };
307 
308  samplingCaseRunner(Input::valid::VANILLA);
309  samplingCaseRunner(Input::valid::SCALE);
310  samplingCaseRunner(Input::valid::GRID);
311 }
312 
313 #ifdef QMC_COMPLEX
314 TEST_CASE("OneBodyDensityMatrices::generateSamplesForSpinor", "[estimators]")
315 {
316  using Input = testing::ValidOneBodyDensityMatricesInput;
318 
319  ProjectData test_project("test", ProjectData::DriverVersion::BATCH);
321 
323  auto wavefunction_pool =
325  auto& pset_target = *(particle_pool.getParticleSet("e"));
326  auto& species_set = pset_target.getSpeciesSet();
327  auto& spo_map = wavefunction_pool.getWaveFunction("wavefunction")->getSPOMap();
328 
329  auto samplingCaseRunner = [&pset_target, &species_set, &spo_map](Input::valid test_case) {
331 
332  bool okay = doc.parseFromString(Input::xml[test_case]);
333  if (!okay)
334  throw std::runtime_error("cannot parse OneBodyDensitMatricesInput section");
335  xmlNodePtr node = doc.getRoot();
336  OneBodyDensityMatricesInput obdmi(node);
337 
338  OneBodyDensityMatrices obDenMat(std::move(obdmi), pset_target.getLattice(), species_set, spo_map, pset_target);
339 
340  testing::OneBodyDensityMatricesTests<QMCTraits::FullPrecRealType> obdmt;
341  //Get control over which rng is used.
342  //we don't want FakeRandom.
343  StdRandom<OneBodyDensityMatrices::FullPrecRealType> rng;
344  obdmt.testGenerateSamplesForSpinor(test_case, obDenMat, pset_target, rng);
345  };
346 
347  //Spin sampling only added for density sampling
348  samplingCaseRunner(Input::valid::VANILLA);
349 }
350 #endif
351 
352 TEST_CASE("OneBodyDensityMatrices::spawnCrowdClone()", "[estimators]")
353 {
354  using Input = testing::ValidOneBodyDensityMatricesInput;
355 
357 
358  ProjectData test_project("test", ProjectData::DriverVersion::BATCH);
360 
362  auto wavefunction_pool =
364  auto& pset_target = *(particle_pool.getParticleSet("e"));
365  auto& species_set = pset_target.getSpeciesSet();
366  auto& spomap = wavefunction_pool.getWaveFunction("wavefunction")->getSPOMap();
367 
369  bool okay = doc.parseFromString(Input::xml[Input::valid::VANILLA]);
370  if (!okay)
371  throw std::runtime_error("cannot parse OneBodyDensitMatricesInput section");
372  xmlNodePtr node = doc.getRoot();
374 
377  REQUIRE(clone != nullptr);
378  REQUIRE(clone.get() != &original);
379  REQUIRE(dynamic_cast<decltype(&original)>(clone.get()) != nullptr);
380 }
381 
382 TEST_CASE("OneBodyDensityMatrices::accumulate", "[estimators]")
383 {
386 
389 
391  bool okay = doc.parseFromString(Input::xml[Input::valid::VANILLA]);
392  if (!okay)
393  throw std::runtime_error("cannot parse OneBodyDensitMatricesInput section");
394  xmlNodePtr node = doc.getRoot();
397  auto wavefunction_pool =
399  auto& spomap = wavefunction_pool.getWaveFunction("wavefunction")->getSPOMap();
400  auto& pset_target = *(particle_pool.getParticleSet("e"));
401  auto& pset_source = *(particle_pool.getParticleSet("ion"));
402  auto& species_set = pset_target.getSpeciesSet();
403  OneBodyDensityMatrices obdm(std::move(obdmi), pset_target.getLattice(), species_set, spomap, pset_target);
404 
405  std::vector<MCPWalker> walkers;
406  int nwalkers = 3;
407  for (int iw = 0; iw < nwalkers; ++iw)
408  walkers.emplace_back(8);
409 
410  std::vector<ParticleSet::ParticlePos> deterministic_rs = {{
411  {-0.6759092808, 0.835668385, 1.985307097},
412  {0.09710352868, -0.76751858, -1.89306891},
413  {-0.5605484247, -0.9578875303, 1.476860642},
414  {2.585144997, 1.862680197, 3.282609463},
415  {-0.1961335093, 1.111888766, -0.578481257},
416  {1.794641614, 1.6000278, -0.9474347234},
417  {2.157717228, 0.9254754186, 2.263158321},
418  {1.883366346, 2.136350632, 3.188981533},
419  },
420  {
421  {-0.2079261839, -0.2796236873, 0.5512072444},
422  {-0.2823159397, 0.7537326217, 0.01526880637},
423  {3.533515453, 2.433290243, 0.9281452894},
424  {2.051767349, 2.312927485, 0.7089259624},
425  {-1.043096781, 0.8190526962, -0.1958218962},
426  {0.9210210443, 0.7726522088, 0.3962054551},
427  {2.043324947, 0.3482068777, 3.39059639},
428  {0.9103830457, 2.167978764, 2.341906071},
429  },
430  {
431  {-0.466550231, 0.09173964709, -0.3779250085},
432  {-0.4211375415, -2.017466068, -1.691870451},
433  {2.090800285, 1.88529861, 2.152359247},
434  {2.973145723, 1.718174577, 3.822324753},
435  {-0.8552014828, -0.3484517336, -0.2870049179},
436  {0.2349359095, -0.5025780797, 0.2305756211},
437  {-0.03547382355, 2.279159069, 3.057915211},
438  {2.535993099, 1.637133598, 3.689830303},
439  }};
440  std::vector<ParticleSet> psets =
441  testing::generateRandomParticleSets<generate_test_data>(pset_target, pset_source, deterministic_rs, nwalkers);
442 
443  auto& trial_wavefunction = *(wavefunction_pool.getPrimary());
444  std::vector<UPtr<TrialWaveFunction>> twfcs(nwalkers);
445  for (int iw = 0; iw < nwalkers; ++iw)
446  twfcs[iw] = trial_wavefunction.makeClone(psets[iw]);
447 
449  rng.init(101);
450 
451  auto updateWalker = [](auto& walker, auto& pset_target, auto& trial_wavefunction) {
452  pset_target.update(true);
453  pset_target.donePbyP();
454  trial_wavefunction.evaluateLog(pset_target);
455  pset_target.saveWalker(walker);
456  };
457 
458  for (int iw = 0; iw < nwalkers; ++iw)
459  updateWalker(walkers[iw], psets[iw], *(twfcs[iw]));
460 
461  auto ref_walkers(makeRefVector<MCPWalker>(walkers));
462  auto ref_psets(makeRefVector<ParticleSet>(psets));
463  auto ref_twfcs(convertUPtrToRefVector(twfcs));
464 
466  obdmt.testAccumulate(obdm, ref_walkers, ref_psets, ref_twfcs, rng);
467 
468  if constexpr (dump_obdm)
469  obdmt.dumpData(obdm);
470 }
471 
472 TEST_CASE("OneBodyDensityMatrices::evaluateMatrix", "[estimators]")
473 {
476 
479 
480  for (auto valid_integrator :
481  std::vector<Input::valid>{Input::valid::VANILLA, Input::valid::SCALE, Input::valid::GRID})
482  {
484  bool okay = doc.parseFromString(Input::xml[valid_integrator]);
485  if (!okay)
486  throw std::runtime_error("cannot parse OneBodyDensitMatricesInput section");
487  xmlNodePtr node = doc.getRoot();
489 
490  std::string integrator_str =
492  std::cout << "Test evaluateMatrix for: " << integrator_str << '\n';
493 
495  auto wavefunction_pool =
497  auto& spomap = wavefunction_pool.getWaveFunction("wavefunction")->getSPOMap();
498  auto& pset_target = *(particle_pool.getParticleSet("e"));
499  auto& species_set = pset_target.getSpeciesSet();
500  OneBodyDensityMatrices obdm(std::move(obdmi), pset_target.getLattice(), species_set, spomap, pset_target);
501  auto& trial_wavefunction = *(wavefunction_pool.getPrimary());
502 
503  // Because we can't control or consistent know the global random state we must initialize particle positions to known values.
505  {4.120557308, 2.547962427, 2.11555481}, {2.545657158, 2.021627665, 3.17555666},
506  {1.251996636, 1.867651463, 0.7268046737}, {4.749059677, 5.845647812, 3.871560574},
507  {5.18129015, 4.168475151, 2.748870373}, {6.24560833, 4.087143421, 4.187825203},
508  {3.173382998, 3.651777267, 2.970916748}, {1.576967478, 2.874752045, 3.687536716},
509  };
510 
512  rng.init(101);
514  // Now we have to bring the pset, trial_wavefunction and walker to valid state.
515  //pset.loadWalker(walker, false);
516  pset_target.update(true);
517  pset_target.donePbyP();
518  trial_wavefunction.evaluateLog(pset_target);
519  pset_target.saveWalker(walker);
521  obdmt.testEvaluateMatrix(obdm, pset_target, trial_wavefunction, walker, rng);
522  if constexpr (dump_obdm)
523  obdmt.dumpData(obdm);
524  }
525 }
526 
527 TEST_CASE("OneBodyDensityMatrices::registerAndWrite", "[estimators]")
528 {
531 
534 
536  bool okay = doc.parseFromString(Input::xml[Input::valid::VANILLA]);
537  if (!okay)
538  throw std::runtime_error("cannot parse OneBodyDensitMatricesInput section");
539  xmlNodePtr node = doc.getRoot();
541 
542  std::string integrator_str =
544  std::cout << "Test registerAndWrite for: " << integrator_str << '\n';
545 
547  auto wavefunction_pool =
549  auto& spomap = wavefunction_pool.getWaveFunction("wavefunction")->getSPOMap();
550  auto& pset_target = *(particle_pool.getParticleSet("e"));
551  auto& species_set = pset_target.getSpeciesSet();
552  OneBodyDensityMatrices obdm(std::move(obdmi), pset_target.getLattice(), species_set, spomap, pset_target);
553 
555  obdmt.testRegisterAndWrite(obdm);
556 }
557 
558 namespace testing
559 {
560 // The test result data is defined down here for readability of the test code.
561 template<typename T>
563  OBDMI::Integrator integrator)
564 {
565  Data data;
566  switch (integrator)
567  {
570  {
571  data = {
572  0.8559739634, -2.775557562e-16, -0.0009312358165, -0.0002137464399, -0.00188481782, -0.00010560363,
573  -0.1690957318, 4.68545619e-11, 0.2135706414, 4.229001871e-11, -0.0003234768609, -0.0008471008737,
574  0.0004354250735, -0.000629313951, 0.2632686421, -4.718447855e-16, -0.0009312358165, 0.0002137464399,
575  0.6904253737, -5.551115123e-17, -2.29895663e-05, -6.759071163e-05, 8.571154293e-05, -2.029129991e-05,
576  -0.0002755957962, 6.15929991e-05, 0.01997752724, -0.03959823579, -0.03441868427, -0.03703684606,
577  -0.000228978625, 5.163016816e-05, -0.00188481782, 0.00010560363, -2.29895663e-05, 6.759071163e-05,
578  0.6910137193, 1.110223025e-15, 0.0001789239062, -9.720928471e-06, -0.0005431524679, 3.125021981e-05,
579  -0.0370646689, -0.03445096274, 0.03964615719, -0.02000771091, -0.0004552815242, 2.596515681e-05,
580  -0.1690957318, -4.685549171e-11, 8.571154293e-05, 2.029129991e-05, 0.0001789239062, 9.720928472e-06,
581  0.6510614284, 4.996003611e-16, -0.3264829535, -1.120145643e-11, -1.628936668e-05, -5.104526467e-05,
582  2.191572411e-05, -3.7934521e-05, -0.1113259098, -1.13588236e-11, 0.9949076208, 5.254666641e-08,
583  0.09716631627, -0.008698974489, -0.07670746649, 0.01101924465, -1.054191923, -1.009885029e-07,
584  0.6429960606, 5.478465426e-08, 0.006654585037, -0.003422468497, -0.008957381204, -0.002542568933,
585  0.3883611988, 2.492687778e-08, -0.2567403319, -2.274337275, -0.070057294, -2.130173974,
586  -0.185057423, -2.774684629, 0.494079399, 4.39071663, -0.2679771972, -2.380684367,
587  -0.2537442403, 0.1010428549, -0.2026493778, -0.04440749192, -0.1216357873, -1.078602297,
588  0.345585874, -1.689635575, 0.7414178343, -1.48274947, 0.742418168, -1.971367257,
589  -0.6650567573, 3.261922039, 0.3607112874, -1.768642232, -0.202649463, -0.01653546255,
590  -0.1315184346, -0.1010429096, 0.1637281089, -0.8013080665, -1.915594805, -3.652633751e-14,
591  -2.106812279, -0.3174159872, -2.798953942, -0.2389235612, 4.203062081, 7.690383574e-08,
592  -2.237403378, -1.095805402e-08, 0.05937990357, 0.2648259866, -0.0799285254, 0.196742672,
593  -0.9570615415, -1.670885652e-14, 0.8559739634, 5.551115123e-17, -0.0009312358165, -0.0002137464399,
594  -0.00188481782, -0.00010560363, -0.1690957318, 4.685551946e-11, 0.2135706414, 4.229011585e-11,
595  -0.0003234768609, -0.0008471008737, 0.0004354250735, -0.000629313951, 0.2632686421, -3.330669074e-16,
596  -0.0009312358165, 0.0002137464399, 0.6904253737, -5.551115123e-17, -2.29895663e-05, -6.759071163e-05,
597  8.571154293e-05, -2.029129991e-05, -0.0002755957962, 6.15929991e-05, 0.01997752724, -0.03959823579,
598  -0.03441868427, -0.03703684606, -0.000228978625, 5.163016816e-05, -0.00188481782, 0.00010560363,
599  -2.29895663e-05, 6.759071163e-05, 0.6910137193, 5.551115123e-17, 0.0001789239062, -9.720928471e-06,
600  -0.0005431524679, 3.125021981e-05, -0.0370646689, -0.03445096274, 0.03964615719, -0.02000771091,
601  -0.0004552815242, 2.596515681e-05, -0.1690957318, -4.68551864e-11, 8.571154293e-05, 2.029129991e-05,
602  0.0001789239062, 9.720928471e-06, 0.6510614284, 1.665334537e-16, -0.3264829535, -1.12007903e-11,
603  -1.628936668e-05, -5.104526467e-05, 2.191572411e-05, -3.7934521e-05, -0.1113259098, -1.13588583e-11,
604  -0.1398697498, 1.676108652e-08, 0.009839305325, -0.03110873861, -0.2743146761, 0.001115813703,
605  0.8078232309, 7.208311281e-09, -0.3939381356, 7.358152143e-10, 0.01322489101, 0.01202117753,
606  -0.01780140861, 0.008930678923, -0.1179882345, 4.657475086e-09, 0.1004197767, -0.3872746911,
607  -0.4631185606, -0.1267234616, 0.251940383, 0.3432366261, 0.1357909209, 0.1882467853,
608  -0.04657813737, -0.1481318675, -0.01754424248, -0.008096072777, 0.04046503369, 0.04321394303,
609  0.01594660416, -0.1297942252, -0.1351702458, -0.2877115343, 0.5308175294, -0.1497099916,
610  -0.302802717, 0.3542104112, -0.1827815947, 0.1398510191, 0.06269657586, -0.1100491288,
611  0.04046503056, -0.0383308965, -0.04195025864, 0.008096072041, -0.02146496291, -0.09642585851,
612  1.583358633, 5.551115123e-15, 0.5400010285, 0.09966977815, 0.8788816574, 0.06123897411,
613  -2.364549544, -2.291840495e-08, 1.339153745, 3.484556643e-09, -0.02316768222, -0.0762631529,
614  0.03118495745, -0.05665686084, 0.6842069691, 1.720845688e-15,
615  };
616  }
617  else if constexpr (std::is_floating_point<OneBodyDensityMatrices::Value>::value)
618  {
619  data = {
620  0.4279869817, -0.0005724911282, -0.0009952107251, -0.08454786586, 0.1067853207, -0.0005852888673,
621  -9.694443874e-05, 0.131634321, -0.0005724911282, 0.3452101128, 0.07732509099, 5.299452569e-05,
622  -0.0001685944122, -0.01345751284, -0.03081850011, -0.0001403043966, -0.0009952107251, 0.07732509099,
623  0.3455042856, 9.430886761e-05, -0.0002871986719, -0.03308358791, 0.006219573346, -0.0002406233405,
624  -0.08454786586, 5.299452569e-05, 9.430886761e-05, 0.3255307142, -0.1632414767, -3.366691551e-05,
625  -8.009937089e-06, -0.0556629549, 0.4974537841, 0.04423363063, -0.03284416792, -0.5270959104,
626  0.3214980027, 0.001616059574, -0.005749976824, 0.1941805869, 1.008798472, 0.9822905645,
627  1.251738209, -1.948318603, 1.056353583, -0.1335433475, -0.04654405447, 0.4784832548,
628  1.017610724, 1.176381561, 1.414874433, -1.963489414, 1.064676762, -0.1520816495,
629  -0.05908794817, 0.4825180877, -0.9577974025, -1.212114133, -1.518938752, 2.101531079,
630  -1.118701694, 0.1621029451, 0.05840707331, -0.4785307708, 0.4279869817, -0.0005724911282,
631  -0.0009952107251, -0.08454786586, 0.1067853207, -0.0005852888673, -9.694443874e-05, 0.131634321,
632  -0.0005724911282, 0.3452101128, 0.07732509099, 5.299452569e-05, -0.0001685944122, -0.01345751284,
633  -0.03081850011, -0.0001403043966, -0.0009952107251, 0.07732509099, 0.3455042856, 9.430886761e-05,
634  -0.0002871986719, -0.03308358791, 0.006219573346, -0.0002406233405, -0.08454786586, 5.299452569e-05,
635  9.430886761e-05, 0.3255307142, -0.1632414767, -3.366691551e-05, -8.00993709e-06, -0.0556629549,
636  -0.06993488329, -0.01063469842, -0.1365994325, 0.4039116172, -0.1969690692, 0.01262303497,
637  -0.004435365783, -0.05899411958, 0.2438472339, -0.1415980475, -0.09314308211, -0.0262279326,
638  0.05077686603, 0.01074587779, 0.0101183853, 0.07287041467, 0.07627064426, 0.3044594466,
639  -0.2645758779, -0.1613163063, 0.08637285105, 0.01857458999, -0.04049312511, 0.0374804478,
640  0.7916793164, 0.3198354033, 0.4700603157, -1.182274783, 0.6695768741, -0.04971541756,
641  -0.01273595169, 0.3421034845,
642  };
643  }
644  break;
645  }
648  {
649  data = {
650  0.8261941614, 0,
651  -0.04207477301, -0.01042641052,
652  -0.09193944978, -0.004771495112,
653  -0.1577639198, 2.251655112e-09,
654  0.2529218804, -3.528230208e-10,
655  0.01395852521, -0.01381912806,
656  -0.01878889069, -0.0102664046,
657  0.2669242402, -1.665334537e-16,
658  -0.04207477301, 0.01042641052,
659  0.7237235295, -2.220446049e-16,
660  -0.1362791782, 0.02883929086,
661  -0.0345711543, -0.01578687034,
662  0.07089423603, 0.01585967312,
663  0.06197413691, -0.0673375704,
664  -0.07846982731, -0.0428311471,
665  0.03401975498, 0.01636340115,
666  -0.09193944978, 0.004771495112,
667  -0.1362791782, -0.02883929086,
668  0.4726910113, 3.330669074e-16,
669  0.1392076512, 0.003920553607,
670  -0.1398496516, -0.008039785222,
671  0.02312496202, 0.01413794966,
672  -0.04626371454, 0.02462349008,
673  -0.1442914533, -0.003858018529,
674  -0.1577639198, -2.251654918e-09,
675  -0.0345711543, 0.01578687034,
676  0.1392076512, -0.003920553607,
677  0.6406077933, -2.220446049e-16,
678  -0.3821572706, -6.950154341e-11,
679  -0.03309714608, -0.006729454624,
680  0.04455046481, -0.004999419621,
681  -0.06308963653, -2.521567177e-09,
682  0.9566798997, 3.993940545e-08,
683  0.1287007623, -0.03913606336,
684  -0.3450997639, 0.01459536453,
685  -1.049841949, -8.448363198e-08,
686  0.7807572133, 5.406439579e-08,
687  0.0636014061, -0.0144126824,
688  -0.08561076546, -0.01070736661,
689  0.3453233488, 7.91976143e-09,
690  -0.1652425535, -1.673204981,
691  -0.2864309551, -1.827550661,
692  0.1092013893, -0.4672139335,
693  0.4181616684, 3.856599084,
694  -0.2678120582, -2.450374785,
695  -0.1661115247, -0.5147466032,
696  0.03950582141, 0.7348504126,
697  -0.01780951734, -0.3040424203,
698  0.2224251057, -1.243046364,
699  0.4872309759, -1.334778155,
700  0.2833545206, -0.2677130338,
701  -0.5628674436, 2.865118995,
702  0.3604889296, -1.820416191,
703  0.03950573518, -0.4243869589,
704  -0.1899388751, 0.5147465492,
705  0.02397254984, -0.2258771884,
706  -1.317725416, -2.642330799e-14,
707  -1.83826268, -0.06023003292,
708  -0.5311043355, -0.2084686282,
709  3.662896282, 2.914568187e-08,
710  -2.270465021, 2.617466888e-08,
711  -0.5283222051, 0.08433855146,
712  0.7111490664, 0.06265614476,
713  -0.1773974531, -2.775557562e-15,
714  0.8261941614, -5.551115123e-17,
715  -0.04207477301, -0.01042641052,
716  -0.09193944978, -0.004771495112,
717  -0.1577639198, 2.251654183e-09,
718  0.2529218804, -3.528226322e-10,
719  0.01395852521, -0.01381912806,
720  -0.01878889069, -0.0102664046,
721  0.2669242402, -1.942890293e-16,
722  -0.04207477301, 0.01042641052,
723  0.7237235295, -2.220446049e-16,
724  -0.1362791782, 0.02883929086,
725  -0.0345711543, -0.01578687034,
726  0.07089423603, 0.01585967312,
727  0.06197413691, -0.0673375704,
728  -0.07846982731, -0.0428311471,
729  0.03401975498, 0.01636340115,
730  -0.09193944978, 0.004771495112,
731  -0.1362791782, -0.02883929086,
732  0.4726910113, 5.551115123e-17,
733  0.1392076512, 0.003920553607,
734  -0.1398496516, -0.008039785222,
735  0.02312496202, 0.01413794966,
736  -0.04626371454, 0.02462349008,
737  -0.1442914533, -0.003858018529,
738  -0.1577639198, -2.251654516e-09,
739  -0.0345711543, 0.01578687034,
740  0.1392076512, -0.003920553607,
741  0.6406077933, -1.665334537e-16,
742  -0.3821572706, -6.950134912e-11,
743  -0.03309714608, -0.006729454624,
744  0.04455046481, -0.004999419621,
745  -0.06308963653, -2.521567073e-09,
746  -0.09246990704, 1.714692167e-08,
747  0.01718168811, -0.002306490387,
748  -0.0203382379, 0.001948481938,
749  0.7391743085, 3.440059682e-09,
750  -0.4038538237, 3.37227038e-09,
751  -0.05196998174, -0.01900034084,
752  0.06995427856, -0.01411562006,
753  0.001586062406, 3.942820508e-09,
754  0.09257244551, -0.4150339965,
755  -0.543569044, -0.1978266927,
756  0.2876563907, 0.3016008129,
757  0.2041868448, 0.2689153029,
758  -0.1456849905, -0.2702324047,
759  -0.06303854185, 0.03455118382,
760  0.02035291547, 0.04212529665,
761  -0.04182725494, -0.2125070964,
762  -0.1246073282, -0.308334295,
763  0.6467124191, -0.2090012216,
764  -0.3333005141, 0.3431996583,
765  -0.2748460697, 0.1997807403,
766  0.1960995565, -0.2007592471,
767  0.02035290751, -0.06296446567,
768  -0.07531415043, -0.03455118719,
769  0.05630165983, -0.15787436,
770  1.360996296, 4.440892099e-15,
771  0.460662518, -0.009109012253,
772  -0.08032289324, 0.05224157011,
773  -2.1641054, -3.065196319e-09,
774  1.447670416, -1.103484504e-08,
775  0.2132345387, -0.01793768528,
776  -0.2870247183, -0.01332608765,
777  0.369682822, 9.436895709e-16,
778  };
779  }
780  else if constexpr (std::is_floating_point<OneBodyDensityMatrices::Value>::value)
781  {
782  data = {
783  0.4130970807, -0.02625059177, -0.04835547244, -0.07888195875, 0.12646094, 6.96985738e-05,
784  -0.01452764765, 0.1334621201, -0.02625059177, 0.3466032331, -0.001161163685, -0.009392146568,
785  0.02751728009, 0.0007624905835, -0.06528656996, 0.008828176914, -0.04835547244, -0.001161163685,
786  0.2210869708, 0.06764354093, -0.06590493546, 0.02539089034, -0.01991866217, -0.07021671737,
787  -0.07888195875, -0.009392146568, 0.06764354093, 0.3203038942, -0.1910786333, -0.01991330156,
788  0.01977552422, -0.031544817, 0.4783399299, 0.04478231453, -0.165252207, -0.5249209266,
789  0.3903785766, 0.02459435045, -0.04815905069, 0.1726616705, 0.7539812136, 0.7595813826,
790  0.2502049139, -1.719218703, 1.091281367, 0.1944111372, -0.3327445041, 0.1431164515,
791  0.7327357349, 0.9257821561, 0.3266874516, -1.713993226, 1.090452555, 0.2048993189,
792  -0.372436296, 0.1249248691, -0.6588627081, -0.9492463565, -0.3697864819, 1.831448156,
793  -1.135232497, -0.2219918268, 0.3869026056, -0.08869872655, 0.4130970807, -0.02625059177,
794  -0.04835547244, -0.07888195875, 0.12646094, 6.96985738e-05, -0.01452764765, 0.1334621201,
795  -0.02625059177, 0.3466032331, -0.001161163685, -0.009392146568, 0.02751728009, 0.0007624905835,
796  -0.06528656996, 0.008828176914, -0.04835547244, -0.001161163685, 0.2210869708, 0.06764354093,
797  -0.06590493546, 0.02539089034, -0.01991866217, -0.07021671737, -0.07888195875, -0.009392146568,
798  0.06764354093, 0.3203038942, -0.1910786333, -0.01991330156, 0.01977552422, -0.031544817,
799  -0.0462349621, 0.007437619475, -0.009194881695, 0.3695871528, -0.2019269101, -0.03548515978,
800  0.02791932722, 0.0007930292316, 0.253803221, -0.1431755519, -0.0640034998, -0.03236422918,
801  0.06227370761, -0.006365722925, 0.02063497821, 0.08533992072, 0.0918634834, 0.3878849668,
802  -0.2614830065, -0.2373134048, 0.1984294012, -0.01545314238, -0.06281061419, 0.1070880099,
803  0.6804981478, 0.2257767529, -0.01404066156, -1.082052701, 0.7238352026, 0.09764842672,
804  -0.150175403, 0.184841411,
805  };
806  }
807  break;
808  }
811  {
812  data = {
813  0.9924480859, -1.110223025e-16,
814  -0.1885291696, 0.0319276602,
815  0.2815361269, -0.02138019278,
816  0.1872112957, -3.684201194e-09,
817  -0.5244043425, 1.095947622e-08,
818  -0.03918278828, -0.3023406232,
819  0.05274204952, -0.2246129078,
820  0.01383924214, 9.714451465e-17,
821  -0.1885291696, -0.0319276602,
822  0.9180621234, 2.775557562e-16,
823  0.2402287308, -0.01308718541,
824  0.2057882223, 0.004465819328,
825  0.1024448325, 0.02374868021,
826  -0.05743147279, -0.07637345821,
827  0.03130436103, -0.08626204101,
828  0.1630605924, -0.02521124324,
829  0.2815361269, 0.02138019278,
830  0.2402287308, 0.01308718541,
831  1.031979911, 3.330669074e-16,
832  -0.03937932123, -0.02333745214,
833  -0.2094143097, -0.01161775883,
834  -0.1348871245, -0.190165311,
835  0.1601241654, -0.149665044,
836  0.2223111845, -0.01849191565,
837  0.1872112957, 3.684200389e-09,
838  0.2057882223, -0.004465819328,
839  -0.03937932123, 0.02333745214,
840  1.323421938, 4.440892099e-16,
841  0.2905844345, 1.315278547e-08,
842  -0.1594450613, -0.2611690287,
843  0.2146213636, -0.1940259416,
844  0.2116494581, 6.064937663e-09,
845  0.5632626352, 1.703273084e-08,
846  -0.3488988503, 0.02573820756,
847  0.2269574287, -0.03956679718,
848  -1.637494571, -1.812725284e-07,
849  -0.8314244571, -8.494759768e-08,
850  0.1951723949, 0.1024546704,
851  -0.2627122381, 0.0761148891,
852  -0.282499705, -2.532976356e-09,
853  -0.09350615209, -0.756293056,
854  0.1463790836, -2.234477269,
855  -0.5935525371, -5.517868252,
856  0.829208332, 7.707654298,
857  0.3647990603, 3.119019826,
858  0.1734726648, -0.3225198405,
859  0.1526169414, 0.4999332675,
860  -0.07866503442, -0.05924377889,
861  0.1258641643, -0.5618602341,
862  1.097431538, -1.4619276,
863  1.305834496, -3.987176014,
864  -1.116157862, 5.726119334,
865  -0.4910385854, 2.317161517,
866  0.1526170644, -0.3054089574,
867  0.08142330397, 0.3225199255,
868  0.105887215, -0.04401298315,
869  -0.4339186561, -1.629252289e-14,
870  -2.305548843, -0.6188863722,
871  -5.45730064, -0.2614613878,
872  7.525382195, 1.317113143e-07,
873  2.880837963, 2.432125967e-08,
874  -0.325402542, -0.228642532,
875  0.4380088526, -0.1698614745,
876  -0.1214395281, -3.462508058e-15,
877  0.9924480859, 5.551115123e-17,
878  -0.1885291696, 0.0319276602,
879  0.2815361269, -0.02138019278,
880  0.1872112957, -3.684200667e-09,
881  -0.5244043425, 1.095947622e-08,
882  -0.03918278828, -0.3023406232,
883  0.05274204952, -0.2246129078,
884  0.01383924214, -2.975050761e-16,
885  -0.1885291696, -0.0319276602,
886  0.9180621234, -1.665334537e-16,
887  0.2402287308, -0.01308718541,
888  0.2057882223, 0.004465819328,
889  0.1024448325, 0.02374868021,
890  -0.05743147279, -0.07637345821,
891  0.03130436103, -0.08626204101,
892  0.1630605924, -0.02521124324,
893  0.2815361269, 0.02138019278,
894  0.2402287308, 0.01308718541,
895  1.031979911, 0,
896  -0.03937932123, -0.02333745214,
897  -0.2094143097, -0.01161775883,
898  -0.1348871245, -0.190165311,
899  0.1601241654, -0.149665044,
900  0.2223111845, -0.01849191565,
901  0.1872112957, 3.68420057e-09,
902  0.2057882223, -0.004465819328,
903  -0.03937932123, 0.02333745214,
904  1.323421938, 0,
905  0.2905844345, 1.315278556e-08,
906  -0.1594450613, -0.2611690287,
907  0.2146213636, -0.1940259416,
908  0.2116494581, 6.064937469e-09,
909  0.2051059875, 2.820113419e-08,
910  0.1614723089, -0.04906840403,
911  -0.4326817327, 0.01831176767,
912  1.706267288, 1.643257963e-08,
913  0.4080117321, 4.721028735e-09,
914  -0.1551120856, -0.2788399487,
915  0.2087889614, -0.2071539023,
916  0.1814356421, -1.855115922e-09,
917  0.4257831978, -0.1752088462,
918  -0.4998292759, 0.07462493825,
919  0.2373005349, 0.3540303805,
920  0.2248224859, 0.09288554326,
921  -0.1494463458, 0.1278007236,
922  -0.05568689306, -0.2034175338,
923  0.1004071574, 0.007263799659,
924  0.03200635492, 0.1228505216,
925  -0.5731263129, -0.130165,
926  0.5767337619, -0.002042114292,
927  -0.330578729, 0.3705292719,
928  -0.3026226697, 0.06900584231,
929  0.2011625579, 0.09494483239,
930  0.1004071511, 0.1154251684,
931  -0.1162463317, 0.2034175291,
932  -0.04308221061, 0.09126728585,
933  0.7806523758, 3.275157923e-15,
934  0.1160951465, 0.2239745488,
935  1.974993197, 0.01316587446,
936  -4.058777904, -4.039943313e-08,
937  -1.782306544, -2.124291576e-08,
938  0.2934590981, 0.1825596902,
939  -0.3950112733, 0.1356259601,
940  -0.2695462139, -7.21644966e-16,
941  };
942  }
943  else if constexpr (std::is_floating_point<OneBodyDensityMatrices::Value>::value)
944  {
945  data = {
946  0.9408447786, -0.3935847362, -0.2468586928, 0.1346424361, -0.2790550847, -0.2728613819, 0.118510686,
947  0.373636656, -0.3935847362, 1.465454239, 0.3906355741, 0.3826898337, -0.03726984102, -0.07963536415,
948  -0.4819140243, -0.2130898108, -0.2468586928, 0.3906355741, 0.8424465027, 0.2589273731, -0.03575380357,
949  -0.253719428, 0.1136907332, 0.09985411935, 0.1346424361, 0.3826898337, 0.2589273731, 1.097297774,
950  -0.2458047841, -0.1441440472, -0.1503264601, 0.09909758853, 0.6175832301, -0.6985723605, -0.6243271971,
951  -1.376534064, 0.09219592002, -0.02237084205, 0.2231280489, 0.1476634744, -1.208204987, 1.712369561,
952  1.715275159, -4.409355333, 0.8508272302, -0.494084464, 0.3242451383, -0.2368202448, -1.480471134,
953  2.537211763, 2.19567287, -4.176679345, 0.8267380683, -0.6223330608, 0.1422109575, -0.3020474178,
954  1.853750075, -2.613447154, -2.388944783, 4.650543191, -0.9996739752, 0.5654611738, -0.2023390067,
955  0.4167421543, 0.9408447786, -0.3935847362, -0.2468586928, 0.1346424361, -0.2790550847, -0.2728613819,
956  0.118510686, 0.373636656, -0.3935847362, 1.465454239, 0.3906355741, 0.3826898337, -0.03726984102,
957  -0.07963536415, -0.4819140243, -0.2130898108, -0.2468586928, 0.3906355741, 0.8424465027, 0.2589273731,
958  -0.03575380357, -0.253719428, 0.1136907332, 0.09985411935, 0.1346424361, 0.3826898337, 0.2589273731,
959  1.097297774, -0.2458047841, -0.1441440472, -0.1503264601, 0.09909758853, 0.3273111356, 0.3790594202,
960  -0.01488345199, 1.314728254, -0.321881638, -0.1063667686, -0.2555566944, 0.1032133755, 0.7484376476,
961  -0.8135360105, -0.4255269663, -0.03246449107, -0.1597379358, -0.09287811809, 0.2144758734, 0.2846013811,
962  -0.1713163625, 1.015300497, -0.563966848, -0.3464659878, 0.087960779, 0.2142716444, -0.5686552903,
963  -0.3532239607, 0.1204139568, -0.3049230367, 0.1268687562, -2.916771087, 0.4142601541, -0.209253697,
964  0.4634440052, 0.1000471674,
965  };
966  }
967  break;
968  }
969  }
970  return data;
971 }
972 
973 template<typename T>
975 {
976  Data data;
978  {
979  data = {
980  2.9291234638, 0.0000000000, -0.2132989360, 0.0745674664, 0.6575312763, -0.0241892145, 0.2532299248,
981  -0.0000000109, -1.1240579737, 0.0000000028, -0.1397854178, 0.1645671318, 0.1881583956, 0.1222590747,
982  -0.0019595386, 0.0000000000, -0.2132989360, -0.0745674664, 3.1431572436, -0.0000000000, 0.4556655323,
983  -0.0032110369, 0.0525287831, 0.0018092687, -0.1435885589, 0.0345014350, 0.1530394508, 0.2388837829,
984  -0.1875535298, 0.1030355338, -0.1497527705, -0.0024706633, 0.6575312763, 0.0241892145, 0.4556655323,
985  0.0032110369, 3.1711078068, 0.0000000000, -0.0159540742, -0.0059570373, -0.3042315457, 0.0162837209,
986  -0.2911526549, 0.0615081334, 0.4400568181, 0.0798583623, 0.0217861975, 0.0169827532, 0.2532299248,
987  0.0000000109, 0.0525287831, -0.0018092687, -0.0159540742, 0.0059570373, 3.4070503585, 0.0000000000,
988  0.8103626817, 0.0000000056, -0.0952865678, -0.2711773981, 0.1282606987, -0.2014612680, 0.4847746576,
989  -0.0000000011, -6.0458119707, 0.0000000067, 5.8517166367, -2.3147839762, -20.4116215167, 0.6636152886,
990  -7.3052152124, 0.0000005920, -5.9833739861, 0.0000001586, 0.4080623680, -5.2593437195, -0.5492733308,
991  -3.9072369333, 3.6961289593, 0.0000000740, 0.9493241103, 1.5333691280, -1.1621879823, -3.5401562485,
992  2.3633267850, 2.7990104110, 3.6921767198, -1.5637638719, 0.9327053281, 0.8055850013, -1.3612795276,
993  -0.3390817263, -0.6755117159, 0.5662076117, 0.3702207032, -1.8456735112, -1.2778398461, 1.1391606465,
994  0.8788742485, -3.0888907845, -2.2904629312, 2.2770067320, -4.9698620698, -1.1617407073, -1.2554700698,
995  0.5984800818, -0.6755119158, -0.3616944760, -0.9538525628, 0.3390815677, -0.4983361691, -1.3711756489,
996  -4.0798709203, -0.0000000000, 2.9144861374, -1.2298470587, -10.8447153218, 0.3305180969, -10.0697430530,
997  0.0000001859, -3.4632636632, 0.0000000419, 0.1108403565, -1.6391623703, -0.1491970188, -1.2177557551,
998  -0.5534367488, 0.0000000000, 2.9291234638, -0.0000000000, -0.2132989360, 0.0745674664, 0.6575312763,
999  -0.0241892145, 0.2532299248, -0.0000000109, -1.1240579737, 0.0000000028, -0.1397854178, 0.1645671318,
1000  0.1881583956, 0.1222590747, -0.0019595386, 0.0000000000, -0.2132989360, -0.0745674664, 3.1431572436,
1001  0.0000000000, 0.4556655323, -0.0032110369, 0.0525287831, 0.0018092687, -0.1435885589, 0.0345014350,
1002  0.1530394508, 0.2388837829, -0.1875535298, 0.1030355338, -0.1497527705, -0.0024706633, 0.6575312763,
1003  0.0241892145, 0.4556655323, 0.0032110369, 3.1711078068, 0.0000000000, -0.0159540742, -0.0059570373,
1004  -0.3042315457, 0.0162837209, -0.2911526549, 0.0615081334, 0.4400568181, 0.0798583623, 0.0217861975,
1005  0.0169827532, 0.2532299248, 0.0000000109, 0.0525287831, -0.0018092687, -0.0159540742, 0.0059570373,
1006  3.4070503585, -0.0000000000, 0.8103626817, 0.0000000056, -0.0952865678, -0.2711773981, 0.1282606987,
1007  -0.2014612680, 0.4847746576, -0.0000000011, -0.0535485614, 0.0000000456, -0.2634589697, -0.4273305284,
1008  -3.7681731806, -0.0298776483, 6.3003941549, -0.0000001048, 2.5146140015, -0.0000000790, 0.2019179652,
1009  0.3783014377, -0.2717920962, 0.2810452661, -0.2195164577, -0.0000000335, 0.8248473372, -0.0302759687,
1010  0.0564839927, 0.0627061848, 0.3233263517, -2.6372267985, -0.7480447015, 4.8497448204, -0.3290319572,
1011  1.3407680257, 1.0988018851, -0.0721739419, 0.8942277438, 0.9612852986, -0.2616888674, 0.6620697600,
1012  -1.1102874022, -0.0224924193, 0.5422928928, -0.0307012816, -0.4413043984, -1.9024908391, 1.0069071996,
1013  3.6029399491, 0.4428940955, 0.9960743917, 0.8942278339, -0.9177543595, 0.5594585213, 0.0721740209,
1014  0.3522468123, 0.4918604753, -2.1530368296, 0.0000000000, -0.4899613886, 0.0492231475, 0.4340461159,
1015  -0.0555640559, -5.8227729586, -0.0000000037, -0.9649526524, -0.0000000466, 0.6122771455, 0.4591252615,
1016  -0.8241569084, 0.3410902494, -0.4924760260, 0.0000000000,
1017  };
1018  }
1019  else if constexpr (std::is_floating_point<OneBodyDensityMatrices::Value>::value)
1020  {
1021  data = {
1022  2.9418870316, -0.0730367697, 0.2825942788, -0.0482128086, -0.4076353925, 0.0565848266, 1.1815582427,
1023  0.1397691924, -0.0730367697, 3.4212182000, 0.8908985700, 0.2575937162, -0.1332837838, 0.5354984573,
1024  0.0223089765, -0.4943417103, 0.2825942788, 0.8908985700, 3.2846114781, 0.2100203992, -0.2969335200,
1025  -0.0668884727, 0.7393865326, 0.0744676157, -0.0482128086, 0.2575937162, 0.2100203992, 2.7804644672,
1026  0.1479959718, -0.2108628317, -0.1058624903, 0.1094655214, -10.4814218220, -2.7285795717, -29.9852578068,
1027  0.7835860171, 3.1523292003, -6.8371790968, -5.0309586919, 7.6730591933, -0.5950386079, 2.9929684361,
1028  0.1106432525, 6.0237785519, -0.6494115783, 0.8415771440, 0.5782321647, 0.4242339082, -3.6209466698,
1029  1.7496134966, -6.8014176733, -1.2619245515, 0.8164606643, -0.6836185979, -1.7382138253, 1.3061371959,
1030  -5.8646902639, -2.1311266822, -17.1140445213, -5.8610844090, 2.2307905583, -4.2035012366, -3.0671418156,
1031  3.5223740553, 2.9418870316, -0.0730367697, 0.2825942788, -0.0482128086, -0.4076353925, 0.0565848266,
1032  1.1815582427, 0.1397691924, -0.0730367697, 3.4212182000, 0.8908985700, 0.2575937162, -0.1332837838,
1033  0.5354984573, 0.0223089765, -0.4943417103, 0.2825942788, 0.8908985700, 3.2846114781, 0.2100203992,
1034  -0.2969335200, -0.0668884727, 0.7393865326, 0.0744676157, -0.0482128086, 0.2575937162, 0.2100203992,
1035  2.7804644672, 0.1479959718, -0.2108628317, -0.1058624903, 0.1094655214, 0.8730451012, 1.0231657997,
1036  -2.5208434318, 3.9689595744, 0.5349419029, 0.5691939230, -0.9648631708, -0.9783541221, 0.3220945960,
1037  -0.8712038765, 1.5289446738, -4.6535453897, 1.2120349971, 0.7055562347, 2.1390028670, -0.7203237486,
1038  -1.3012372435, -0.5715963470, 0.5311792927, -1.7615315077, 0.5661913908, 0.2251863434, 0.8706914019,
1039  -0.1236678876, -2.9956864970, -1.7155487648, -1.2654380257, -4.2335858349, 1.4018920220, -0.0193915347,
1040  -0.1188045474, -0.0461947842,
1041  };
1042  }
1043  return data;
1044 }
1045 
1046 template<typename T>
1048 {
1049  Data data;
1051  {
1052  data = {0.80788614, 0.03850586, 0.31091131, 0.38205625, 0.44777905, 0.21091672, 0.10743698, 0.96152914, 0.20791861,
1053  0.6312172, 0.45626258, 0.08119136, 0.90431926, 0.11375654, 0.8067544, 0.54595911, 0.84150056, 0.20621477,
1054  0.03414278, 0.29305023, 0.45018883, 0.20085125, 0.40840351, 0.66825399, 0.84792934, 0.16337382, 0.29096146,
1055  0.39040606, 0.83718991, 0.36255192, 0.75091973, 0.46448582, 0.52401399, 0.90176318, 0.38023557, 0.71028634,
1056  0.84000996, 0.23419118, 0.30553855, 0.72794856, 0.32839433, 0.09785895, 0.4900262, 0.03030494, 0.14355491,
1057  0.57876809, 0.75103524, 0.02254518, 0.66510853, 0.03025319, 0.14589509, 0.01467558, 0.48894704, 0.04131938,
1058  0.78907133, 0.99612016, 0.96833451, 0.67448307, 0.93273006, 0.20283358, 0.73427328, 0.51123229, 0.58574252,
1059  0.75568488, 0.9315693, 0.31038253, 0.34644132, 0.36777262, 0.63775962, 0.05487314, 0.83054504, 0.73965924,
1060  0.37751022, 0.93023354, 0.48824517, 0.03106788, 0.91065432, 0.65471251, 0.00795302, 0.31126833, 0.58215106,
1061  0.89782226, 0.33120933, 0.38434518, 0.94334496, 0.36956736, 0.30646061, 0.07810828, 0.97154079, 0.02514768,
1062  0.84117586, 0.23971438, 0.22623238, 0.35594978, 0.36437615, 0.77427557, 0.53121471, 0.53240224, 0.39277132,
1063  0.69702301, 0.00654816, 0.94311127, 0.99594941, 0.49187803, 0.10536416, 0.97302689, 0.44064925, 0.28800946,
1064  0.29735792, 0.09520762, 0.16488844, 0.53418347, 0.65144419, 0.03973145, 0.53885783, 0.7358722, 0.45309111,
1065  0.79190158, 0.22632916, 0.20165111, 0.42914622, 0.22774585, 0.85961036, 0.56634649, 0.10872035, 0.69251128,
1066  0.2021623, 0.67509793, 0.52598371, 0.46793523, 0.91830822, 0.77020753, 0.09795741, 0.76007691, 0.50678007,
1067  0.4795003, 0.66071835, 0.32442191, 0.23907092, 0.00210992, 0.75524376, 0.94234713, 0.5690928, 0.34679526,
1068  0.30868398, 0.77709926, 0.51133303, 0.20918474, 0.39361117, 0.7962174, 0.48195751, 0.99470092, 0.30140288,
1069  0.54417926, 0.83619258, 0.66607407, 0.63815312, 0.84754301, 0.46217861, 0.40382941, 0.8669163, 0.80892846,
1070  0.78275811, 0.08359574, 0.74127963, 0.83093771, 0.88526743, 0.83504267, 0.97853459, 0.21692576, 0.45164991,
1071  0.9074328, 0.87753109, 0.85276772, 0.53516365, 0.53168515, 0.40269904, 0.89566943, 0.37372088, 0.73002752,
1072  0.35593885, 0.73251118, 0.77375582, 0.12760254, 0.8455303, 0.91220549, 0.97342524, 0.6901983, 0.01956135,
1073  0.98966668, 0.89508114, 0.86107002, 0.85139379, 0.26868744, 0.49119517, 0.79074938, 0.91333433, 0.16790021,
1074  0.14352574, 0.97134471, 0.69411371, 0.28605858, 0.70411151, 0.36656519, 0.70878013, 0.21327726, 0.34290399,
1075  0.84309746, 0.90860334, 0.97362624, 0.05473755, 0.71643348, 0.14711903, 0.38781449, 0.784074, 0.40246134,
1076  0.40066814, 0.93058349, 0.43298608, 0.11167385, 0.27113968, 0.33209627, 0.40601194, 0.81328762, 0.68107437,
1077  0.46367926, 0.13120107, 0.38408714, 0.64249068, 0.68798637, 0.35231959, 0.98679773, 0.12638461, 0.75466016,
1078  0.97212161, 0.15569373, 0.55338423, 0.2814492, 0.88983892, 0.33155614, 0.25340461, 0.02949572, 0.08162776,
1079  0.49678983, 0.59962038, 0.20915831, 0.7750513, 0.6575729, 0.50223288, 0.37927361, 0.09678806, 0.22351711,
1080  0.127808, 0.16958427, 0.26687417, 0.37408405};
1081  }
1082  else if constexpr (std::is_floating_point<OneBodyDensityMatrices::Value>::value)
1083  {
1084  data = {0.19381403, 0.71283387, 0.74185289, 0.01593182, 0.4237967, 0.30713958, 0.04877154, 0.32839586, 0.55966269,
1085  0.86166257, 0.13164395, 0.18967966, 0.52716562, 0.75501117, 0.89099129, 0.91181301, 0.56299847, 0.84040689,
1086  0.42003362, 0.54417536, 0.91786624, 0.24179691, 0.80097938, 0.52312247, 0.49135109, 0.06204171, 0.32918368,
1087  0.19109569, 0.27910843, 0.13869181, 0.39276429, 0.84432737, 0.55296738, 0.58336158, 0.42912991, 0.73847165,
1088  0.70821508, 0.44336715, 0.34394486, 0.81857914, 0.17582283, 0.85206451, 0.00962322, 0.6018576, 0.77709871,
1089  0.9363809, 0.77303195, 0.79399817, 0.64872202, 0.90561921, 0.34909147, 0.5382709, 0.4735065, 0.97592345,
1090  0.64042891, 0.98233348, 0.61072865, 0.99648271, 0.93723708, 0.12341335, 0.87404106, 0.52492966, 0.50025206,
1091  0.32956586, 0.35388674, 0.41701219, 0.35787114, 0.78154075, 0.19389593, 0.56085759, 0.42076409, 0.45505835,
1092  0.13691315, 0.92741853, 0.65416634, 0.01324141, 0.70580805, 0.36063625, 0.20206282, 0.04019175, 0.3161708,
1093  0.8021294, 0.47419179, 0.58339627, 0.94680233, 0.14275504, 0.51723762, 0.88195736, 0.02861162, 0.54720941,
1094  0.47704361, 0.72112318, 0.71249342, 0.57327699, 0.82174918, 0.65460258, 0.58492448, 0.0654615, 0.23514782,
1095  0.56317195, 0.99078012, 0.16018222, 0.98232388, 0.48021303, 0.45997915, 0.22901306, 0.30486665, 0.47519321,
1096  0.11869839, 0.25773838, 0.30733499, 0.03014402, 0.41846284, 0.51370103, 0.14486378, 0.91023931, 0.45369315,
1097  0.18793261, 0.51507439, 0.46019929, 0.67773434, 0.20830221, 0.59268401, 0.48456955, 0.9678142, 0.50709602,
1098  0.85130517, 0.60737725};
1099  }
1100  return data;
1101 }
1102 
1103 } // namespace testing
1104 } // namespace qmcplusplus
class that handles xmlDoc
Definition: Libxml2Doc.h:76
bool open(const std::filesystem::path &fname, unsigned flags=H5F_ACC_RDWR)
open a file
The operator<< functions provide output for data structures that can be used to directly initialize t...
const std::string & string() const
Return std::string for use with HDF5.
Definition: hdf_path.h:67
helper functions for EinsplineSetBuilder
Definition: Configuration.h:43
if(!okay) throw std xmlNodePtr node
Per crowd Estimator for OneBodyDensityMatrices aka 1RDM DensityMatrices1B.
MakeReturn< UnaryNode< FnFabs, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t abs(const Vector< T1, C1 > &l)
class ProjectData
Definition: ProjectData.h:36
static ParticleSetPool make_diamondC_1x1x1(Communicate *c)
TEST_CASE("complex_helper", "[type_traits]")
void checkData(Real *ref_in, Real *test_in, size_t size)
Checking approximate equality for complex valued data as two reals is not consistent with testing pra...
void close()
close all the open groups and file
Definition: hdf_archive.cpp:38
static WaveFunctionPool make_diamondC_1x1x1(const RuntimeOptions &runtime_options, Communicate *comm, ParticleSetPool &particle_pool)
xmlNodePtr getRoot()
Definition: Libxml2Doc.h:88
void testGenerateSamplesForSpinor(ValidOneBodyDensityMatricesInput::valid input, OneBodyDensityMatrices &obdm, ParticleSet &pset_target, StdRandom< T > &rng)
constexpr bool generate_test_data
set to true to generate new random R for particles.
static WaveFunctionPool make_O2_spinor(const RuntimeOptions &runtime_options, Communicate *comm, ParticleSetPool &particle_pool)
Some ParticleSet functions use the global Random so we need some helper functions to avoid interminan...
ProjectData test_project("test", ProjectData::DriverVersion::BATCH)
class to handle hdf file
Definition: hdf_archive.h:51
void readSlabReshaped(T &data, const std::array< IT, RANK > &shape, const std::string &aname)
read file dataset with a specific shape into a container and check status
Definition: hdf_archive.h:321
const char walkers[]
Definition: HDFVersion.h:36
void testGenerateSamples(ValidOneBodyDensityMatricesInput::valid input, OneBodyDensityMatrices &obdm, ParticleSet &pset_target, StdRandom< T > &rng)
Attaches a unit to a Vector for IO.
Communicate * Controller
Global Communicator for a process.
Definition: Communicate.cpp:35
This wrapper is to allow us to leave the user facing operator<< for classes alone.
std::unique_ptr< OperatorEstBase > spawnCrowdClone() const override
standard interface
static RefVector< T > convertUPtrToRefVector(const UPtrVector< T > &ptr_list)
convert a vector of std::unique_ptrs<T> to a refvector<T>
const RuntimeOptions & getRuntimeOptions() const noexcept
#define OHMMS_DIM
Definition: config.h:64
void generateSamples(const Real weight, ParticleSet &pset_target, RandomBase< FullPrecReal > &rng, int steps=0)
Dispatch method to difference methods of generating samples.
Wrapping information on parallelism.
Definition: Communicate.h:68
MakeReturn< BinaryNode< FnPow, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t, typename CreateLeaf< Vector< T2, C2 > >::Leaf_t > >::Expression_t pow(const Vector< T1, C1 > &l, const Vector< T2, C2 > &r)
CrystalLattice< OHMMS_PRECISION, OHMMS_DIM > lattice
void evaluateMatrix(ParticleSet &pset_target, TrialWaveFunction &psi_target, const MCPWalker &walker, RandomBase< FullPrecReal > &rng)
static ParticleSetPool make_O2_spinor(Communicate *c)
Specialized paritlce class for atomistic simulations.
Definition: ParticleSet.h:55
Native representation for DensityMatrices1B Estimator&#39;s inputs.
REQUIRE(std::filesystem::exists(filename))
This a subclass for runtime errors that will occur on all ranks.
A minimally functional wrapper for the since c++11 <random>
OneBodyDensityMatricesInput obdmi(node)
void implAccumulate(const RefVector< MCPWalker > &walkers, const RefVector< ParticleSet > &psets, const RefVector< TrialWaveFunction > &wfns, RandomBase< FullPrecReal > &rng)
Unfortunate design RandomGenerator type aliasing and virtual inheritance requires this for testing...
void testEvaluateMatrix(OneBodyDensityMatrices &obdm, ParticleSet &pset, TrialWaveFunction &trial_wavefunction, MCPWalker &walker, StdRandom< T > &rng)
no change test for evaluateMatrix.
Sampling sampling_
Sampling method, this derived values in input_.
static const std::unordered_map< std::string, std::any > lookup_input_enum_value
mapping for enumerated options of OneBodyDensityMatrices This data object is the basis of input enum ...
Declaration of a TrialWaveFunction.
static std::string reverseLookupInputEnumMap(ENUM_T enum_val, const std::unordered_map< std::string, std::any > &enum_map)
Get string represtation of enum class type value from enum_val.
Definition: InputSection.h:148
std::vector< QMCT::RealType > Data
void init(int iseed_in) override
Definition: StdRandom.h:64
std::vector< std::reference_wrapper< T > > RefVector
bool create(const std::filesystem::path &fname, unsigned flags=H5F_ACC_TRUNC)
create a file
Class to represent a many-body trial wave function.
bool parseFromString(const std::string_view data)
Definition: Libxml2Doc.cpp:204
Walker< QMCTraits, PtclOnLatticeTraits > MCPWalker
Definition: test_walker.cpp:31
CHECK(log_values[0]==ComplexApprox(std::complex< double >{ 5.603777579195571, -6.1586603331188225 }))
OEBAccessor oeba(sdn)
break encapsulation of data_ by OperatorEstBase only for testing!
void testAccumulate(OneBodyDensityMatrices &obdm, RefVector< MCPWalker > &walkers, RefVector< ParticleSet > &psets, RefVector< TrialWaveFunction > &twfcs, StdRandom< T > &rng)
no change test for accumulate
void testCopyConstructor(const OneBodyDensityMatrices &obdm)
void write(hdf_archive &file)
Write to previously registered observable_helper hdf5 wrapper.
void registerOperatorEstimator(hdf_archive &file) override
create and tie OperatorEstimator&#39;s observable_helper hdf5 wrapper to stat.h5 file ...
SpinDensityNew original(std::move(sdi), lattice, species_set)
QMCTraits::FullPrecRealType value_type
A container class to represent a walker.
Definition: Walker.h:49
SpeciesSet makeSpeciesSet(const SpeciesCases species_case)
Walker< QMCTraits, PtclOnLatticeTraits > MCPWalker