QMCPACK
test_Vector.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) 2021 QMCPACK developers.
6 //
7 // File developed by: Mark Dewing, markdewing@gmail.com, University of Illinois at Urbana-Champaign
8 // Peter Doak, doakpw@ornl.gov, Oak Ridge National Laboratory
9 //
10 // File created by: Mark Dewing, markdewing@gmail.com, University of Illinois at Urbana-Champaign
11 //////////////////////////////////////////////////////////////////////////////////////
12 
13 
14 #include "catch.hpp"
15 #include <cstdio>
16 #include <vector>
17 #include <string>
18 #include "OhmmsPETE/OhmmsVector.h"
19 #include "OhmmsPETE/TinyVector.h"
20 
21 using std::string;
22 
23 namespace qmcplusplus
24 {
25 TEST_CASE("vector", "[OhmmsPETE]")
26 {
27  using vec_t = Vector<double>;
28  vec_t A(3);
29  vec_t B(3);
30 
31  // iterator
32  vec_t::iterator ia = A.begin();
33  for (; ia != A.end(); ia++)
34  {
35  *ia = 1.0;
36  }
37 
38  // Assignment. This eventually generates a call to 'evaluate' in OhmmsVector.h
39  // To do: pointer to tutorial on expression template techniques
40  B = A;
41 
42  // *= operator in OhmmVectorOperators.h
43  B *= 3.1;
44 
45  CHECK(B[0] == Approx(3.1));
46  CHECK(B[1] == Approx(3.1));
47  CHECK(B[2] == Approx(3.1));
48  REQUIRE(B == B);
49  REQUIRE(!(B == A));
50 
51  vec_t C(2);
52  REQUIRE(A != B);
53  REQUIRE(A != B);
54 }
55 
56 TEST_CASE("Vector simple intializer list", "[OhmmsPETE]")
57 {
58  //empty list should work
59  Vector<int> vec_int{};
60  Vector<double> vec_double{5.0, 4.0, 3.0, 2.0, 1.0};
61  CHECK(vec_double[0] == Approx(5.0));
62  CHECK(vec_double[4] == Approx(1.0));
63 }
64 
65 TEST_CASE("Vector nested intializer list", "[OhmmsPETE]")
66 {
67  Vector<TinyVector<double, 3>> vec_tinyd3({{1, 2, 3}, {4, 5, 6}, {7, 8, 9}});
68  CHECK(vec_tinyd3[1][1] == 5);
69  CHECK(vec_tinyd3[2][0] == 7);
70  Vector<std::pair<int, int>> vec_pair{{1, 2}, {3, 4}};
71  CHECK(vec_pair[0].first == 1);
72  CHECK(vec_pair[1].second == 4);
73 }
74 
75 TEST_CASE("VectorViewer", "[OhmmsPETE]")
76 {
77  int a[3];
78  a[0] = 2;
79  a[1] = 4;
80  a[2] = -5;
81  Vector<int> view_a(a, 3);
82 
83  REQUIRE(view_a.size() == 3);
84 
85  // operator[]
86  REQUIRE(view_a[0] == 2);
87  REQUIRE(view_a[1] == 4);
88  REQUIRE(view_a[2] == -5);
89 
90  // operator[] returning a reference
91  view_a[1] = 42;
92 
93  REQUIRE(a[0] == 2);
94  REQUIRE(a[1] == 42);
95  REQUIRE(a[2] == -5);
96 
97  // TODO: add optional bounds checking to accesses via operator[]
98 }
99 
100 TEST_CASE("NestedContainers", "[OhmmsPETE]")
101 {
102  Vector<std::vector<int>> vec_of_vecs(2);
103  vec_of_vecs[0].push_back(123);
104  vec_of_vecs.resize(5);
105  vec_of_vecs[0].clear();
106  vec_of_vecs[0].push_back(123);
107  vec_of_vecs.resize(0);
108  vec_of_vecs.resize(3);
109  vec_of_vecs[0].push_back(123);
110  CHECK(vec_of_vecs[0].back() == 123);
111 
112  Vector<std::vector<int>> vec_copy(vec_of_vecs);
113  REQUIRE(vec_copy.size() == 3);
114  REQUIRE(vec_copy[0].size() == 1);
115  CHECK(vec_copy[0].back() == 123);
116 
117  Vector<std::vector<int>> vec_assign;
118  vec_assign = vec_of_vecs;
119  REQUIRE(vec_copy.size() == 3);
120  REQUIRE(vec_copy[0].size() == 1);
121  CHECK(vec_copy[0].back() == 123);
122 }
123 
124 } // namespace qmcplusplus
void resize(size_type n, Type_t val=Type_t())
Resize the container.
Definition: OhmmsVector.h:166
helper functions for EinsplineSetBuilder
Definition: Configuration.h:43
TEST_CASE("complex_helper", "[type_traits]")
size_type size() const
return the current size
Definition: OhmmsVector.h:162
REQUIRE(std::filesystem::exists(filename))
Declaraton of Vector<T,Alloc> Manage memory through Alloc directly and allow referencing an existing ...
CHECK(log_values[0]==ComplexApprox(std::complex< double >{ 5.603777579195571, -6.1586603331188225 }))
void clear()
clear
Definition: OhmmsVector.h:191
double B(double x, int k, int i, const std::vector< double > &t)