Peano
Loading...
Searching...
No Matches
VectorOperations.cpph
Go to the documentation of this file.
1#include <limits>
2
3
4template <int Size, typename Scalar>
5std::vector<Scalar> tarch::la::toSTLVector(const Vector<Size, Scalar>& vector) {
6 std::vector<Scalar> result(Size);
7 for (int i=0; i<Size; i++) result[i] = vector(i);
8 return result;
9}
10
11
12template <int Size, typename Scalar>
13int tarch::la::count(const Vector<Size, Scalar>& vector, const Scalar& value) {
14 int result = 0;
15 for (int d = 0; d < Size; d++) {
16 if (tarch::la::equals(vector(d), value)) {
17 result++;
18 }
19 }
20 return result;
21}
22
23template <int Size, typename Scalar>
25 bool result = true;
26 for (int d = 0; d < Size; d++) {
27 result &= tarch::la::equals(vector(d), vector(0));
28 }
29 return result;
30}
31
32template <int Size, typename Scalar>
34 bool result = false;
35 for (int i = 0; i < Size; i++) {
36 result |= std::isnan(vector(i));
37 }
38 return result;
39}
40
41template <int Size, typename Scalar>
43 bool result = true;
44 for (int i = 0; i < Size; i++) {
45 result &= std::isfinite(vector(i));
46 }
47 return result;
48}
49
50template <int Size, typename Scalar>
51bool tarch::la::contains(const Vector<Size, Scalar>& vector, const Scalar& value) {
52 bool result = false;
53 for (int d = 0; d < Size; d++) {
54 result |= vector(d) == value;
55 }
56 return result;
57}
58
59template <int Size, typename Scalar>
62#ifdef CompilerICC
63#pragma ivdep
64#endif
65 for (int i = 0; i < Size; i++) {
66 assertion1(vector(i) != 0.0, vector.toString());
67 result(i) = 1.0 / vector(i);
68 }
69 return result;
70}
71
72template <int Size, typename Scalar>
74 Scalar result(tarch::la::abs(vector(0)));
75 for (int i = 1; i < Size; i++) {
76 result += tarch::la::abs(vector(i));
77 }
78 return result;
79}
80
81template <int Size>
82double tarch::la::norm1(const Vector<Size, std::complex<double>>& vector) {
83 double result(0.0);
84 for (int i = 0; i < Size; i++) {
85 result += tarch::la::abs(vector(i));
86 }
87 return result;
88}
89
90template <int Size, typename Scalar>
92 Scalar result(vector(0) * vector(0));
93 for (int i = 1; i < Size; i++) {
94 result += vector(i) * vector(i);
95 }
96 return result;
97}
98
99template <int Size, typename Scalar>
101 return std::sqrt(norm2Squared(vector));
102}
103
104template <int Size>
105double tarch::la::norm2Squared(const Vector<Size, std::complex<double>>& vector) {
106 double result(0.0);
107 for (int i = 0; i < Size; i++) {
108 result += (vector(i).real() * vector(i).real());
109 result += (vector(i).img() * vector(i).img());
110 }
111 return result;
112}
113
114template <int Size>
115double tarch::la::norm2(const Vector<Size, std::complex<double>>& vector) {
116 return std::sqrt(norm2Squared(vector));
117}
118
119template <int Size, typename Scalar>
121 Scalar result(0);
122
123 for (int i = 0; i < Size; i++) {
124 result = std::abs(vector(i)) > result ? std::abs(vector(i)) : result;
125 }
126
127 return result;
128}
129
130template <int Size>
131double tarch::la::normMax(const Vector<Size, std::complex<double>>& vector) {
132 double result = 0.0;
133
134 for (int i = 0; i < Size; i++) {
135 result = std::abs(vector(i)) > result ? std::abs(vector(i)) : result;
136 }
137
138 return result;
139}
140
141template <int Size, typename Scalar>
144 for (int i = 0; i < Size; i++) {
145 result(i) = tarch::la::abs(vector(i));
146 }
147 return result;
148}
149
150template <int Size>
151tarch::la::Vector<Size, double> tarch::la::abs(const Vector<Size, std::complex<double>>& vector) {
153 for (int i = 0; i < Size; i++) {
154 result(i) = tarch::la::abs(vector(i));
155 }
156 return result;
157}
158
159template <int Size, typename Scalar>
161 return sum(vector) / static_cast<Scalar>(Size);
162}
163
164template <int Size, typename Scalar>
166 Scalar result = vector(0);
167 for (int i = 1; i < Size; i++) {
168 result += vector(i);
169 }
170 return result;
171}
172
173template <int Size, typename Scalar>
175 Scalar result = vector(0);
176 for (int i = 1; i < Size; i++) {
177 result *= vector(i);
178 }
179 return result;
180}
181
182template <int Size, typename Scalar>
184 int indexMax = 0;
185 for (int i = 1; i < Size; i++) {
186 indexMax = vector(i) > vector(indexMax) ? i : indexMax;
187 }
188 return indexMax;
189}
190
191template <int Size, typename Scalar>
193 int indexMin = 0;
194 for (int i = 1; i < Size; i++) {
195 indexMin = vector(i) < vector(indexMin) ? i : indexMin;
196 }
197 return indexMin;
198}
199
200template <int Size, typename Scalar>
202 Scalar largest = vector(0);
203 for (int i = 1; i < Size; i++) {
204 if (largest < vector(i)) {
205 largest = vector(i);
206 }
207 }
208 return largest;
209}
210
211template <int Size, typename Scalar>
213 Scalar largest = std::abs(vector(0));
214 for (int i = 1; i < Size; i++) {
215 if (largest < std::abs(vector(i))) {
216 largest = std::abs(vector(i));
217 }
218 }
219 return largest;
220}
221
222template <int Size>
223tarch::la::Vector<Size, double> tarch::la::real(const Vector<Size, std::complex<double>>& vector) {
225
226 for (int d = 0; d < Size; d++) {
227 result(d) = vector(d).real();
228 }
229
230 return result;
231}
232
233template <int Size>
234tarch::la::Vector<Size, double> tarch::la::imag(const Vector<Size, std::complex<double>>& vector) {
236
237 for (int d = 0; d < Size; d++) {
238 result(d) = vector(d).imag();
239 }
240
241 return result;
242}
243
244template <int Size>
245double tarch::la::maxReal(const Vector<Size, std::complex<double>>& vector) {
246 double result = std::numeric_limits<double>::min();
247 for (int i = 0; i < Size; i++) {
248 result = result > vector(i).real() ? result : vector(i).real();
249 }
250 return result;
251}
252
253template <int Size>
254double tarch::la::maxImag(const Vector<Size, std::complex<double>>& vector) {
255 double result = std::numeric_limits<double>::min();
256 for (int i = 0; i < Size; i++) {
257 result = result > vector(i).imag() ? result : vector(i).imag();
258 }
259 return result;
260}
261
262template <int Size, typename Scalar>
264 Scalar smallest = vector(0);
265 for (int i = 1; i < Size; i++) {
266 if (smallest > vector(i)) {
267 smallest = vector(i);
268 }
269 }
270 return smallest;
271}
272
273template <int Size, typename Scalar>
274std::ostream& operator<<(std::ostream& os, const tarch::la::Vector<Size, Scalar>& vector) {
275 os << ::toString(vector);
276 return os;
277}
#define assertion1(expr, param)
std::ostream & operator<<(std::ostream &os, const tarch::la::Vector< Size, Scalar > &vector)
Streams the component values into a comma separated representation.
Definition vec.h:7
std::string toString(Filter filter)
Definition convert.cpp:170
CF abs(const CF &cf)
Scalar average(const Vector< Size, Scalar > &vector)
Computes the volume of the tetrahedron spanned by the Cartesian unit vectors scaled by the correspond...
int count(const Vector< Size, Scalar > &vector, const Scalar &value)
Scalar volume(const Vector< Size, Scalar > &vector)
Computes the volume of the tetrahedron spanned by the Cartesian unit vectors scaled by the correspond...
Vector< Size, double > real(const Vector< Size, std::complex< double > > &vector)
Scalar sum(const Matrix< Rows, Cols, Scalar > &matrix)
Computes the sum of all entries of the matrix.
double maxImag(const Vector< Size, std::complex< double > > &vector)
int indexMin(const Vector< Size, Scalar > &vector)
Returns the index of the element with minimal value (NOT absolute value).
bool allEntriesAreTheSame(const Vector< Size, Scalar > &vector)
std::vector< Scalar > toSTLVector(const Vector< Size, Scalar > &vector)
double abs(double value)
Returns the absolute value of a type by redirecting to std::abs.
int indexMax(const Vector< Size, Scalar > &vector)
Returns the index of the element with maximal value (NOT absolute value).
double max(double a, double b, double c)
I need the maximum of three values all the time, to I decided to write a function for this.
Definition Scalar.cpp:8
Matrix< Rows, Cols, Scalar > invertEntries(const Matrix< Rows, Cols, Scalar > &matrix)
int isEntryFinite(const Vector< Size, Scalar > &vector)
Scalar norm2Squared(const Vector< Size, Scalar > &vector)
int isEntryNan(const Vector< Size, Scalar > &vector)
Scalar maxAbs(const Vector< Size, Scalar > &vector)
Returns the element with maximal absolute value.
Vector< Size, double > imag(const Vector< Size, std::complex< double > > &vector)
bool equals(const Matrix< Rows, Cols, Scalar > &lhs, const Matrix< Rows, Cols, Scalar > &rhs, const Scalar &tolerance=NUMERICAL_ZERO_DIFFERENCE)
Compares to matrices on equality by means of a numerical accuracy.
Scalar norm1(const Vector< Size, Scalar > &vector)
Computes the 1-norm of the vector, i.e.
bool contains(const Vector< Size, Scalar > &vector, const Scalar &value)
Scalar normMax(const Vector< Size, Scalar > &vector)
Computes the max-norm of the vector.
Scalar min(const Vector< Size, Scalar > &vector)
Returns the element with minimal value (NOT absolute value).
Scalar norm2(const Vector< Size, Scalar > &vector)
Computes the 2-norm of the vector, i.e.
double maxReal(const Vector< Size, std::complex< double > > &vector)
Simple vector class.
Definition Vector.h:134