MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
vector.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 // Implementation of data type vector
13 
14 #include "kernels.hpp"
15 #include "vector.hpp"
16 #include "../general/forall.hpp"
17 
18 #if defined(MFEM_USE_SUNDIALS)
19 #include "sundials.hpp"
20 #if defined(MFEM_USE_MPI)
21 #include <nvector/nvector_parallel.h>
22 #endif
23 #endif
24 
25 #ifdef MFEM_USE_OPENMP
26 #include <omp.h>
27 #endif
28 
29 #include <iostream>
30 #include <iomanip>
31 #include <cmath>
32 #include <cstdlib>
33 #include <ctime>
34 #include <limits>
35 
36 namespace mfem
37 {
38 
40 {
41  const int s = v.Size();
42  if (s > 0)
43  {
44  MFEM_ASSERT(!v.data.Empty(), "invalid source vector");
45  size = s;
46  data.New(s, v.data.GetMemoryType());
47  data.CopyFrom(v.data, s);
48  }
49  else
50  {
51  size = 0;
52  data.Reset();
53  }
54  UseDevice(v.UseDevice());
55 }
56 
57 void Vector::Load(std::istream **in, int np, int *dim)
58 {
59  int i, j, s;
60 
61  s = 0;
62  for (i = 0; i < np; i++)
63  {
64  s += dim[i];
65  }
66 
67  SetSize(s);
68 
69  int p = 0;
70  for (i = 0; i < np; i++)
71  {
72  for (j = 0; j < dim[i]; j++)
73  {
74  *in[i] >> data[p++];
75  // Clang's libc++ sets the failbit when (correctly) parsing subnormals,
76  // so we reset the failbit here.
77  if (!*in[i] && errno == ERANGE)
78  {
79  in[i]->clear();
80  }
81  }
82  }
83 }
84 
85 void Vector::Load(std::istream &in, int Size)
86 {
87  SetSize(Size);
88 
89  for (int i = 0; i < size; i++)
90  {
91  in >> data[i];
92  // Clang's libc++ sets the failbit when (correctly) parsing subnormals,
93  // so we reset the failbit here.
94  if (!in && errno == ERANGE)
95  {
96  in.clear();
97  }
98  }
99 }
100 
101 double &Vector::Elem(int i)
102 {
103  return operator()(i);
104 }
105 
106 const double &Vector::Elem(int i) const
107 {
108  return operator()(i);
109 }
110 
111 double Vector::operator*(const double *v) const
112 {
113  double dot = 0.0;
114 #ifdef MFEM_USE_LEGACY_OPENMP
115  #pragma omp parallel for reduction(+:dot)
116 #endif
117  for (int i = 0; i < size; i++)
118  {
119  dot += data[i] * v[i];
120  }
121  return dot;
122 }
123 
124 Vector &Vector::operator=(const double *v)
125 {
126  data.CopyFromHost(v, size);
127  return *this;
128 }
129 
131 {
132 #if 0
133  SetSize(v.Size(), v.data.GetMemoryType());
134  data.CopyFrom(v.data, v.Size());
135  UseDevice(v.UseDevice());
136 #else
137  SetSize(v.Size());
138  bool vuse = v.UseDevice();
139  const bool use_dev = UseDevice() || vuse;
140  v.UseDevice(use_dev);
141  // keep 'data' where it is, unless 'use_dev' is true
142  if (use_dev) { Write(); }
143  data.CopyFrom(v.data, v.Size());
144  v.UseDevice(vuse);
145 #endif
146  return *this;
147 }
148 
149 Vector &Vector::operator=(double value)
150 {
151  const bool use_dev = UseDevice();
152  const int N = size;
153  auto y = Write(use_dev);
154  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = value;);
155  return *this;
156 }
157 
159 {
160  const bool use_dev = UseDevice();
161  const int N = size;
162  auto y = ReadWrite(use_dev);
163  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] *= c;);
164  return *this;
165 }
166 
168 {
169  MFEM_ASSERT(size == v.size, "incompatible Vectors!");
170 
171  const bool use_dev = UseDevice() || v.UseDevice();
172  const int N = size;
173  auto y = ReadWrite(use_dev);
174  auto x = v.Read(use_dev);
175  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] *= x[i];);
176  return *this;
177 }
178 
180 {
181  const bool use_dev = UseDevice();
182  const int N = size;
183  const double m = 1.0/c;
184  auto y = ReadWrite(use_dev);
185  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] *= m;);
186  return *this;
187 }
188 
190 {
191  MFEM_ASSERT(size == v.size, "incompatible Vectors!");
192 
193  const bool use_dev = UseDevice() || v.UseDevice();
194  const int N = size;
195  auto y = ReadWrite(use_dev);
196  auto x = v.Read(use_dev);
197  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] /= x[i];);
198  return *this;
199 }
200 
202 {
203  const bool use_dev = UseDevice();
204  const int N = size;
205  auto y = ReadWrite(use_dev);
206  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] -= c;);
207  return *this;
208 }
209 
211 {
212  MFEM_ASSERT(size == v.size, "incompatible Vectors!");
213 
214  const bool use_dev = UseDevice() || v.UseDevice();
215  const int N = size;
216  auto y = ReadWrite(use_dev);
217  auto x = v.Read(use_dev);
218  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] -= x[i];);
219  return *this;
220 }
221 
223 {
224  const bool use_dev = UseDevice();
225  const int N = size;
226  auto y = ReadWrite(use_dev);
227  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] += c;);
228  return *this;
229 }
230 
232 {
233  MFEM_ASSERT(size == v.size, "incompatible Vectors!");
234 
235  const bool use_dev = UseDevice() || v.UseDevice();
236  const int N = size;
237  auto y = ReadWrite(use_dev);
238  auto x = v.Read(use_dev);
239  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] += x[i];);
240  return *this;
241 }
242 
243 Vector &Vector::Add(const double a, const Vector &Va)
244 {
245  MFEM_ASSERT(size == Va.size, "incompatible Vectors!");
246 
247  if (a != 0.0)
248  {
249  const int N = size;
250  const bool use_dev = UseDevice() || Va.UseDevice();
251  auto y = ReadWrite(use_dev);
252  auto x = Va.Read(use_dev);
253  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] += a * x[i];);
254  }
255  return *this;
256 }
257 
258 Vector &Vector::Set(const double a, const Vector &Va)
259 {
260  MFEM_ASSERT(size == Va.size, "incompatible Vectors!");
261 
262  const bool use_dev = UseDevice() || Va.UseDevice();
263  const int N = size;
264  auto x = Va.Read(use_dev);
265  auto y = Write(use_dev);
266  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = a * x[i];);
267  return *this;
268 }
269 
270 void Vector::SetVector(const Vector &v, int offset)
271 {
272  MFEM_ASSERT(v.Size() + offset <= size, "invalid sub-vector");
273 
274  const int vs = v.Size();
275  const double *vp = v.data;
276  double *p = data + offset;
277  for (int i = 0; i < vs; i++)
278  {
279  p[i] = vp[i];
280  }
281 }
282 
284 {
285  const bool use_dev = UseDevice();
286  const int N = size;
287  auto y = ReadWrite(use_dev);
288  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = -y[i];);
289 }
290 
291 void add(const Vector &v1, const Vector &v2, Vector &v)
292 {
293  MFEM_ASSERT(v.size == v1.size && v.size == v2.size,
294  "incompatible Vectors!");
295 
296 #if !defined(MFEM_USE_LEGACY_OPENMP)
297  const bool use_dev = v1.UseDevice() || v2.UseDevice() || v.UseDevice();
298  const int N = v.size;
299  // Note: get read access first, in case v is the same as v1/v2.
300  auto x1 = v1.Read(use_dev);
301  auto x2 = v2.Read(use_dev);
302  auto y = v.Write(use_dev);
303  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = x1[i] + x2[i];);
304 #else
305  #pragma omp parallel for
306  for (int i = 0; i < v.size; i++)
307  {
308  v.data[i] = v1.data[i] + v2.data[i];
309  }
310 #endif
311 }
312 
313 void add(const Vector &v1, double alpha, const Vector &v2, Vector &v)
314 {
315  MFEM_ASSERT(v.size == v1.size && v.size == v2.size,
316  "incompatible Vectors!");
317 
318  if (alpha == 0.0)
319  {
320  v = v1;
321  }
322  else if (alpha == 1.0)
323  {
324  add(v1, v2, v);
325  }
326  else
327  {
328 #if !defined(MFEM_USE_LEGACY_OPENMP)
329  const bool use_dev = v1.UseDevice() || v2.UseDevice() || v.UseDevice();
330  const int N = v.size;
331  // Note: get read access first, in case v is the same as v1/v2.
332  auto d_x = v1.Read(use_dev);
333  auto d_y = v2.Read(use_dev);
334  auto d_z = v.Write(use_dev);
335  MFEM_FORALL_SWITCH(use_dev, i, N, d_z[i] = d_x[i] + alpha * d_y[i];);
336 #else
337  const double *v1p = v1.data, *v2p = v2.data;
338  double *vp = v.data;
339  const int s = v.size;
340  #pragma omp parallel for
341  for (int i = 0; i < s; i++)
342  {
343  vp[i] = v1p[i] + alpha*v2p[i];
344  }
345 #endif
346  }
347 }
348 
349 void add(const double a, const Vector &x, const Vector &y, Vector &z)
350 {
351  MFEM_ASSERT(x.size == y.size && x.size == z.size,
352  "incompatible Vectors!");
353 
354  if (a == 0.0)
355  {
356  z = 0.0;
357  }
358  else if (a == 1.0)
359  {
360  add(x, y, z);
361  }
362  else
363  {
364 #if !defined(MFEM_USE_LEGACY_OPENMP)
365  const bool use_dev = x.UseDevice() || y.UseDevice() || z.UseDevice();
366  const int N = x.size;
367  // Note: get read access first, in case z is the same as x/y.
368  auto xd = x.Read(use_dev);
369  auto yd = y.Read(use_dev);
370  auto zd = z.Write(use_dev);
371  MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = a * (xd[i] + yd[i]););
372 #else
373  const double *xp = x.data;
374  const double *yp = y.data;
375  double *zp = z.data;
376  const int s = x.size;
377  #pragma omp parallel for
378  for (int i = 0; i < s; i++)
379  {
380  zp[i] = a * (xp[i] + yp[i]);
381  }
382 #endif
383  }
384 }
385 
386 void add(const double a, const Vector &x,
387  const double b, const Vector &y, Vector &z)
388 {
389  MFEM_ASSERT(x.size == y.size && x.size == z.size,
390  "incompatible Vectors!");
391 
392  if (a == 0.0)
393  {
394  z.Set(b, y);
395  }
396  else if (b == 0.0)
397  {
398  z.Set(a, x);
399  }
400 #if 0
401  else if (a == 1.0)
402  {
403  add(x, b, y, z);
404  }
405  else if (b == 1.0)
406  {
407  add(y, a, x, z);
408  }
409  else if (a == b)
410  {
411  add(a, x, y, z);
412  }
413 #endif
414  else
415  {
416 #if !defined(MFEM_USE_LEGACY_OPENMP)
417  const bool use_dev = x.UseDevice() || y.UseDevice() || z.UseDevice();
418  const int N = x.size;
419  // Note: get read access first, in case z is the same as x/y.
420  auto xd = x.Read(use_dev);
421  auto yd = y.Read(use_dev);
422  auto zd = z.Write(use_dev);
423  MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = a * xd[i] + b * yd[i];);
424 #else
425  const double *xp = x.data;
426  const double *yp = y.data;
427  double *zp = z.data;
428  const int s = x.size;
429  #pragma omp parallel for
430  for (int i = 0; i < s; i++)
431  {
432  zp[i] = a * xp[i] + b * yp[i];
433  }
434 #endif
435  }
436 }
437 
438 void subtract(const Vector &x, const Vector &y, Vector &z)
439 {
440  MFEM_ASSERT(x.size == y.size && x.size == z.size,
441  "incompatible Vectors!");
442 
443 #if !defined(MFEM_USE_LEGACY_OPENMP)
444  const bool use_dev = x.UseDevice() || y.UseDevice() || z.UseDevice();
445  const int N = x.size;
446  // Note: get read access first, in case z is the same as x/y.
447  auto xd = x.Read(use_dev);
448  auto yd = y.Read(use_dev);
449  auto zd = z.Write(use_dev);
450  MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = xd[i] - yd[i];);
451 #else
452  const double *xp = x.data;
453  const double *yp = y.data;
454  double *zp = z.data;
455  const int s = x.size;
456  #pragma omp parallel for
457  for (int i = 0; i < s; i++)
458  {
459  zp[i] = xp[i] - yp[i];
460  }
461 #endif
462 }
463 
464 void subtract(const double a, const Vector &x, const Vector &y, Vector &z)
465 {
466  MFEM_ASSERT(x.size == y.size && x.size == z.size,
467  "incompatible Vectors!");
468 
469  if (a == 0.)
470  {
471  z = 0.;
472  }
473  else if (a == 1.)
474  {
475  subtract(x, y, z);
476  }
477  else
478  {
479 #if !defined(MFEM_USE_LEGACY_OPENMP)
480  const bool use_dev = x.UseDevice() || y.UseDevice() || z.UseDevice();
481  const int N = x.size;
482  // Note: get read access first, in case z is the same as x/y.
483  auto xd = x.Read(use_dev);
484  auto yd = y.Read(use_dev);
485  auto zd = z.Write(use_dev);
486  MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = a * (xd[i] - yd[i]););
487 #else
488  const double *xp = x.data;
489  const double *yp = y.data;
490  double *zp = z.data;
491  const int s = x.size;
492  #pragma omp parallel for
493  for (int i = 0; i < s; i++)
494  {
495  zp[i] = a * (xp[i] - yp[i]);
496  }
497 #endif
498  }
499 }
500 
501 void Vector::median(const Vector &lo, const Vector &hi)
502 {
503  MFEM_ASSERT(size == lo.size && size == hi.size,
504  "incompatible Vectors!");
505 
506  const bool use_dev = UseDevice() || lo.UseDevice() || hi.UseDevice();
507  const int N = size;
508  // Note: get read access first, in case *this is the same as lo/hi.
509  auto l = lo.Read(use_dev);
510  auto h = hi.Read(use_dev);
511  auto m = Write(use_dev);
512  MFEM_FORALL_SWITCH(use_dev, i, N,
513  {
514  if (m[i] < l[i])
515  {
516  m[i] = l[i];
517  }
518  else if (m[i] > h[i])
519  {
520  m[i] = h[i];
521  }
522  });
523 }
524 
525 void Vector::GetSubVector(const Array<int> &dofs, Vector &elemvect) const
526 {
527  const int n = dofs.Size();
528  elemvect.SetSize(n);
529  const bool use_dev = dofs.UseDevice() || elemvect.UseDevice();
530  auto d_y = elemvect.Write(use_dev);
531  auto d_X = Read(use_dev);
532  auto d_dofs = dofs.Read(use_dev);
533  MFEM_FORALL_SWITCH(use_dev, i, n,
534  {
535  const int dof_i = d_dofs[i];
536  d_y[i] = dof_i >= 0 ? d_X[dof_i] : -d_X[-dof_i-1];
537  });
538 }
539 
540 void Vector::GetSubVector(const Array<int> &dofs, double *elem_data) const
541 {
543  const int n = dofs.Size();
544  for (int i = 0; i < n; i++)
545  {
546  const int j = dofs[i];
547  elem_data[i] = (j >= 0) ? data[j] : -data[-1-j];
548  }
549 }
550 
551 void Vector::SetSubVector(const Array<int> &dofs, const double value)
552 {
553  const bool use_dev = dofs.UseDevice();
554  const int n = dofs.Size();
555  // Use read+write access for *this - we only modify some of its entries
556  auto d_X = ReadWrite(use_dev);
557  auto d_dofs = dofs.Read(use_dev);
558  MFEM_FORALL_SWITCH(use_dev, i, n,
559  {
560  const int j = d_dofs[i];
561  if (j >= 0)
562  {
563  d_X[j] = value;
564  }
565  else
566  {
567  d_X[-1-j] = -value;
568  }
569  });
570 }
571 
572 void Vector::SetSubVector(const Array<int> &dofs, const Vector &elemvect)
573 {
574  MFEM_ASSERT(dofs.Size() == elemvect.Size(),
575  "Size mismatch: length of dofs is " << dofs.Size()
576  << ", length of elemvect is " << elemvect.Size());
577 
578  const bool use_dev = dofs.UseDevice() || elemvect.UseDevice();
579  const int n = dofs.Size();
580  // Use read+write access for X - we only modify some of its entries
581  auto d_X = ReadWrite(use_dev);
582  auto d_y = elemvect.Read(use_dev);
583  auto d_dofs = dofs.Read(use_dev);
584  MFEM_FORALL_SWITCH(use_dev, i, n,
585  {
586  const int dof_i = d_dofs[i];
587  if (dof_i >= 0)
588  {
589  d_X[dof_i] = d_y[i];
590  }
591  else
592  {
593  d_X[-1-dof_i] = -d_y[i];
594  }
595  });
596 }
597 
598 void Vector::SetSubVector(const Array<int> &dofs, double *elem_data)
599 {
600  // Use read+write access because we overwrite only part of the data.
602  const int n = dofs.Size();
603  for (int i = 0; i < n; i++)
604  {
605  const int j= dofs[i];
606  if (j >= 0)
607  {
608  operator()(j) = elem_data[i];
609  }
610  else
611  {
612  operator()(-1-j) = -elem_data[i];
613  }
614  }
615 }
616 
617 void Vector::AddElementVector(const Array<int> &dofs, const Vector &elemvect)
618 {
619  MFEM_ASSERT(dofs.Size() == elemvect.Size(), "Size mismatch: "
620  "length of dofs is " << dofs.Size() <<
621  ", length of elemvect is " << elemvect.Size());
622 
623  const bool use_dev = dofs.UseDevice() || elemvect.UseDevice();
624  const int n = dofs.Size();
625  auto d_y = elemvect.Read(use_dev);
626  auto d_X = ReadWrite(use_dev);
627  auto d_dofs = dofs.Read(use_dev);
628  MFEM_FORALL_SWITCH(use_dev, i, n,
629  {
630  const int j = d_dofs[i];
631  if (j >= 0)
632  {
633  d_X[j] += d_y[i];
634  }
635  else
636  {
637  d_X[-1-j] -= d_y[i];
638  }
639  });
640 }
641 
642 void Vector::AddElementVector(const Array<int> &dofs, double *elem_data)
643 {
645  const int n = dofs.Size();
646  for (int i = 0; i < n; i++)
647  {
648  const int j = dofs[i];
649  if (j >= 0)
650  {
651  operator()(j) += elem_data[i];
652  }
653  else
654  {
655  operator()(-1-j) -= elem_data[i];
656  }
657  }
658 }
659 
660 void Vector::AddElementVector(const Array<int> &dofs, const double a,
661  const Vector &elemvect)
662 {
663  MFEM_ASSERT(dofs.Size() == elemvect.Size(), "Size mismatch: "
664  "length of dofs is " << dofs.Size() <<
665  ", length of elemvect is " << elemvect.Size());
666 
667  const bool use_dev = dofs.UseDevice() || elemvect.UseDevice();
668  const int n = dofs.Size();
669  auto d_y = ReadWrite(use_dev);
670  auto d_x = elemvect.Read(use_dev);
671  auto d_dofs = dofs.Read(use_dev);
672  MFEM_FORALL_SWITCH(use_dev, i, n,
673  {
674  const int j = d_dofs[i];
675  if (j >= 0)
676  {
677  d_y[j] += a * d_x[i];
678  }
679  else
680  {
681  d_y[-1-j] -= a * d_x[i];
682  }
683  });
684 }
685 
686 void Vector::SetSubVectorComplement(const Array<int> &dofs, const double val)
687 {
688  const bool use_dev = UseDevice() || dofs.UseDevice();
689  const int n = dofs.Size();
690  const int N = size;
691  Vector dofs_vals(n, use_dev ?
694  auto d_data = ReadWrite(use_dev);
695  auto d_dofs_vals = dofs_vals.Write(use_dev);
696  auto d_dofs = dofs.Read(use_dev);
697  MFEM_FORALL_SWITCH(use_dev, i, n, d_dofs_vals[i] = d_data[d_dofs[i]];);
698  MFEM_FORALL_SWITCH(use_dev, i, N, d_data[i] = val;);
699  MFEM_FORALL_SWITCH(use_dev, i, n, d_data[d_dofs[i]] = d_dofs_vals[i];);
700 }
701 
702 void Vector::Print(std::ostream &out, int width) const
703 {
704  if (!size) { return; }
706  for (int i = 0; 1; )
707  {
708  out << ZeroSubnormal(data[i]);
709  i++;
710  if (i == size)
711  {
712  break;
713  }
714  if ( i % width == 0 )
715  {
716  out << '\n';
717  }
718  else
719  {
720  out << ' ';
721  }
722  }
723  out << '\n';
724 }
725 
726 #ifdef MFEM_USE_ADIOS2
728  const std::string& variable_name) const
729 {
730  if (!size) { return; }
732  out.engine.Put(variable_name, &data[0] );
733 }
734 #endif
735 
736 void Vector::Print_HYPRE(std::ostream &out) const
737 {
738  int i;
739  std::ios::fmtflags old_fmt = out.flags();
740  out.setf(std::ios::scientific);
741  std::streamsize old_prec = out.precision(14);
742 
743  out << size << '\n'; // number of rows
744 
745  data.Read(MemoryClass::HOST, size);
746  for (i = 0; i < size; i++)
747  {
748  out << ZeroSubnormal(data[i]) << '\n';
749  }
750 
751  out.precision(old_prec);
752  out.flags(old_fmt);
753 }
754 
755 void Vector::PrintHash(std::ostream &out) const
756 {
757  out << "size: " << size << '\n';
758  HashFunction hf;
759  hf.AppendDoubles(HostRead(), size);
760  out << "hash: " << hf.GetHash() << '\n';
761 }
762 
763 void Vector::Randomize(int seed)
764 {
765  // static unsigned int seed = time(0);
766  const double max = (double)(RAND_MAX) + 1.;
767 
768  if (seed == 0)
769  {
770  seed = (int)time(0);
771  }
772 
773  // srand(seed++);
774  srand((unsigned)seed);
775 
776  HostWrite();
777  for (int i = 0; i < size; i++)
778  {
779  data[i] = std::abs(rand()/max);
780  }
781 }
782 
783 double Vector::Norml2() const
784 {
785  // Scale entries of Vector on the fly, using algorithms from
786  // std::hypot() and LAPACK's drm2. This scaling ensures that the
787  // argument of each call to std::pow is <= 1 to avoid overflow.
788  if (0 == size)
789  {
790  return 0.0;
791  } // end if 0 == size
792 
794  if (1 == size)
795  {
796  return std::abs(data[0]);
797  } // end if 1 == size
798  return kernels::Norml2(size, (const double*) data);
799 }
800 
801 double Vector::Normlinf() const
802 {
803  HostRead();
804  double max = 0.0;
805  for (int i = 0; i < size; i++)
806  {
807  max = std::max(std::abs(data[i]), max);
808  }
809  return max;
810 }
811 
812 double Vector::Norml1() const
813 {
814  HostRead();
815  double sum = 0.0;
816  for (int i = 0; i < size; i++)
817  {
818  sum += std::abs(data[i]);
819  }
820  return sum;
821 }
822 
823 double Vector::Normlp(double p) const
824 {
825  MFEM_ASSERT(p > 0.0, "Vector::Normlp");
826 
827  if (p == 1.0)
828  {
829  return Norml1();
830  }
831  if (p == 2.0)
832  {
833  return Norml2();
834  }
835  if (p < infinity())
836  {
837  // Scale entries of Vector on the fly, using algorithms from
838  // std::hypot() and LAPACK's drm2. This scaling ensures that the
839  // argument of each call to std::pow is <= 1 to avoid overflow.
840  if (0 == size)
841  {
842  return 0.0;
843  } // end if 0 == size
844 
845  if (1 == size)
846  {
847  return std::abs(data[0]);
848  } // end if 1 == size
849 
850  double scale = 0.0;
851  double sum = 0.0;
852 
853  for (int i = 0; i < size; i++)
854  {
855  if (data[i] != 0.0)
856  {
857  const double absdata = std::abs(data[i]);
858  if (scale <= absdata)
859  {
860  sum = 1.0 + sum * std::pow(scale / absdata, p);
861  scale = absdata;
862  continue;
863  } // end if scale <= absdata
864  sum += std::pow(absdata / scale, p); // else scale > absdata
865  } // end if data[i] != 0
866  }
867  return scale * std::pow(sum, 1.0/p);
868  } // end if p < infinity()
869 
870  return Normlinf(); // else p >= infinity()
871 }
872 
873 double Vector::Max() const
874 {
875  if (size == 0) { return -infinity(); }
876 
877  HostRead();
878  double max = data[0];
879 
880  for (int i = 1; i < size; i++)
881  {
882  if (data[i] > max)
883  {
884  max = data[i];
885  }
886  }
887 
888  return max;
889 }
890 
891 double Vector::Sum() const
892 {
893  double sum = 0.0;
894 
895  const double *h_data = this->HostRead();
896  for (int i = 0; i < size; i++)
897  {
898  sum += h_data[i];
899  }
900 
901  return sum;
902 }
903 
904 #ifdef MFEM_USE_CUDA
905 static __global__ void cuKernelMin(const int N, double *gdsr, const double *x)
906 {
907  __shared__ double s_min[MFEM_CUDA_BLOCKS];
908  const int n = blockDim.x*blockIdx.x + threadIdx.x;
909  if (n>=N) { return; }
910  const int bid = blockIdx.x;
911  const int tid = threadIdx.x;
912  const int bbd = bid*blockDim.x;
913  const int rid = bbd+tid;
914  s_min[tid] = x[n];
915  for (int workers=blockDim.x>>1; workers>0; workers>>=1)
916  {
917  __syncthreads();
918  if (tid >= workers) { continue; }
919  if (rid >= N) { continue; }
920  const int dualTid = tid + workers;
921  if (dualTid >= N) { continue; }
922  const int rdd = bbd+dualTid;
923  if (rdd >= N) { continue; }
924  if (dualTid >= blockDim.x) { continue; }
925  s_min[tid] = fmin(s_min[tid], s_min[dualTid]);
926  }
927  if (tid==0) { gdsr[bid] = s_min[0]; }
928 }
929 
930 static Array<double> cuda_reduce_buf;
931 
932 static double cuVectorMin(const int N, const double *X)
933 {
934  const int tpb = MFEM_CUDA_BLOCKS;
935  const int blockSize = MFEM_CUDA_BLOCKS;
936  const int gridSize = (N+blockSize-1)/blockSize;
937  const int min_sz = (N%tpb)==0? (N/tpb) : (1+N/tpb);
938  cuda_reduce_buf.SetSize(min_sz);
939  Memory<double> &buf = cuda_reduce_buf.GetMemory();
940  double *d_min = buf.Write(MemoryClass::DEVICE, min_sz);
941  cuKernelMin<<<gridSize,blockSize>>>(N, d_min, X);
942  MFEM_GPU_CHECK(cudaGetLastError());
943  const double *h_min = buf.Read(MemoryClass::HOST, min_sz);
945  for (int i = 0; i < min_sz; i++) { min = fmin(min, h_min[i]); }
946  return min;
947 }
948 
949 static __global__ void cuKernelDot(const int N, double *gdsr,
950  const double *x, const double *y)
951 {
952  __shared__ double s_dot[MFEM_CUDA_BLOCKS];
953  const int n = blockDim.x*blockIdx.x + threadIdx.x;
954  if (n>=N) { return; }
955  const int bid = blockIdx.x;
956  const int tid = threadIdx.x;
957  const int bbd = bid*blockDim.x;
958  const int rid = bbd+tid;
959  s_dot[tid] = x[n] * y[n];
960  for (int workers=blockDim.x>>1; workers>0; workers>>=1)
961  {
962  __syncthreads();
963  if (tid >= workers) { continue; }
964  if (rid >= N) { continue; }
965  const int dualTid = tid + workers;
966  if (dualTid >= N) { continue; }
967  const int rdd = bbd+dualTid;
968  if (rdd >= N) { continue; }
969  if (dualTid >= blockDim.x) { continue; }
970  s_dot[tid] += s_dot[dualTid];
971  }
972  if (tid==0) { gdsr[bid] = s_dot[0]; }
973 }
974 
975 static double cuVectorDot(const int N, const double *X, const double *Y)
976 {
977  const int tpb = MFEM_CUDA_BLOCKS;
978  const int blockSize = MFEM_CUDA_BLOCKS;
979  const int gridSize = (N+blockSize-1)/blockSize;
980  const int dot_sz = (N%tpb)==0? (N/tpb) : (1+N/tpb);
981  cuda_reduce_buf.SetSize(dot_sz, Device::GetDeviceMemoryType());
982  Memory<double> &buf = cuda_reduce_buf.GetMemory();
983  double *d_dot = buf.Write(MemoryClass::DEVICE, dot_sz);
984  cuKernelDot<<<gridSize,blockSize>>>(N, d_dot, X, Y);
985  MFEM_GPU_CHECK(cudaGetLastError());
986  const double *h_dot = buf.Read(MemoryClass::HOST, dot_sz);
987  double dot = 0.0;
988  for (int i = 0; i < dot_sz; i++) { dot += h_dot[i]; }
989  return dot;
990 }
991 #endif // MFEM_USE_CUDA
992 
993 #ifdef MFEM_USE_HIP
994 static __global__ void hipKernelMin(const int N, double *gdsr, const double *x)
995 {
996  __shared__ double s_min[MFEM_HIP_BLOCKS];
997  const int n = hipBlockDim_x*hipBlockIdx_x + hipThreadIdx_x;
998  if (n>=N) { return; }
999  const int bid = hipBlockIdx_x;
1000  const int tid = hipThreadIdx_x;
1001  const int bbd = bid*hipBlockDim_x;
1002  const int rid = bbd+tid;
1003  s_min[tid] = x[n];
1004  for (int workers=hipBlockDim_x>>1; workers>0; workers>>=1)
1005  {
1006  __syncthreads();
1007  if (tid >= workers) { continue; }
1008  if (rid >= N) { continue; }
1009  const int dualTid = tid + workers;
1010  if (dualTid >= N) { continue; }
1011  const int rdd = bbd+dualTid;
1012  if (rdd >= N) { continue; }
1013  if (dualTid >= hipBlockDim_x) { continue; }
1014  s_min[tid] = fmin(s_min[tid], s_min[dualTid]);
1015  }
1016  if (tid==0) { gdsr[bid] = s_min[0]; }
1017 }
1018 
1019 static Array<double> cuda_reduce_buf;
1020 
1021 static double hipVectorMin(const int N, const double *X)
1022 {
1023  const int tpb = MFEM_HIP_BLOCKS;
1024  const int blockSize = MFEM_HIP_BLOCKS;
1025  const int gridSize = (N+blockSize-1)/blockSize;
1026  const int min_sz = (N%tpb)==0 ? (N/tpb) : (1+N/tpb);
1027  cuda_reduce_buf.SetSize(min_sz);
1028  Memory<double> &buf = cuda_reduce_buf.GetMemory();
1029  double *d_min = buf.Write(MemoryClass::DEVICE, min_sz);
1030  hipLaunchKernelGGL(hipKernelMin,gridSize,blockSize,0,0,N,d_min,X);
1031  MFEM_GPU_CHECK(hipGetLastError());
1032  const double *h_min = buf.Read(MemoryClass::HOST, min_sz);
1034  for (int i = 0; i < min_sz; i++) { min = fmin(min, h_min[i]); }
1035  return min;
1036 }
1037 
1038 static __global__ void hipKernelDot(const int N, double *gdsr,
1039  const double *x, const double *y)
1040 {
1041  __shared__ double s_dot[MFEM_HIP_BLOCKS];
1042  const int n = hipBlockDim_x*hipBlockIdx_x + hipThreadIdx_x;
1043  if (n>=N) { return; }
1044  const int bid = hipBlockIdx_x;
1045  const int tid = hipThreadIdx_x;
1046  const int bbd = bid*hipBlockDim_x;
1047  const int rid = bbd+tid;
1048  s_dot[tid] = x[n] * y[n];
1049  for (int workers=hipBlockDim_x>>1; workers>0; workers>>=1)
1050  {
1051  __syncthreads();
1052  if (tid >= workers) { continue; }
1053  if (rid >= N) { continue; }
1054  const int dualTid = tid + workers;
1055  if (dualTid >= N) { continue; }
1056  const int rdd = bbd+dualTid;
1057  if (rdd >= N) { continue; }
1058  if (dualTid >= hipBlockDim_x) { continue; }
1059  s_dot[tid] += s_dot[dualTid];
1060  }
1061  if (tid==0) { gdsr[bid] = s_dot[0]; }
1062 }
1063 
1064 static double hipVectorDot(const int N, const double *X, const double *Y)
1065 {
1066  const int tpb = MFEM_HIP_BLOCKS;
1067  const int blockSize = MFEM_HIP_BLOCKS;
1068  const int gridSize = (N+blockSize-1)/blockSize;
1069  const int dot_sz = (N%tpb)==0 ? (N/tpb) : (1+N/tpb);
1070  cuda_reduce_buf.SetSize(dot_sz);
1071  Memory<double> &buf = cuda_reduce_buf.GetMemory();
1072  double *d_dot = buf.Write(MemoryClass::DEVICE, dot_sz);
1073  hipLaunchKernelGGL(hipKernelDot,gridSize,blockSize,0,0,N,d_dot,X,Y);
1074  MFEM_GPU_CHECK(hipGetLastError());
1075  const double *h_dot = buf.Read(MemoryClass::HOST, dot_sz);
1076  double dot = 0.0;
1077  for (int i = 0; i < dot_sz; i++) { dot += h_dot[i]; }
1078  return dot;
1079 }
1080 #endif // MFEM_USE_HIP
1081 
1082 double Vector::operator*(const Vector &v) const
1083 {
1084  MFEM_ASSERT(size == v.size, "incompatible Vectors!");
1085  if (size == 0) { return 0.0; }
1086 
1087  const bool use_dev = UseDevice() || v.UseDevice();
1088 #if defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP) || defined(MFEM_USE_OPENMP)
1089  auto m_data = Read(use_dev);
1090 #else
1091  Read(use_dev);
1092 #endif
1093  auto v_data = v.Read(use_dev);
1094 
1095  if (!use_dev) { goto vector_dot_cpu; }
1096 
1097 #ifdef MFEM_USE_OCCA
1098  if (DeviceCanUseOcca())
1099  {
1100  return occa::linalg::dot<double,double,double>(
1102  }
1103 #endif
1104 
1105 #ifdef MFEM_USE_CUDA
1107  {
1108  return cuVectorDot(size, m_data, v_data);
1109  }
1110 #endif
1111 
1112 #ifdef MFEM_USE_HIP
1114  {
1115  return hipVectorDot(size, m_data, v_data);
1116  }
1117 #endif
1118 
1119 #ifdef MFEM_USE_OPENMP
1121  {
1122 #define MFEM_USE_OPENMP_DETERMINISTIC_DOT
1123 #ifdef MFEM_USE_OPENMP_DETERMINISTIC_DOT
1124  // By default, use a deterministic way of computing the dot product
1125  static Vector th_dot;
1126  #pragma omp parallel
1127  {
1128  const int nt = omp_get_num_threads();
1129  #pragma omp master
1130  th_dot.SetSize(nt);
1131  const int tid = omp_get_thread_num();
1132  const int stride = (size + nt - 1)/nt;
1133  const int start = tid*stride;
1134  const int stop = std::min(start + stride, size);
1135  double my_dot = 0.0;
1136  for (int i = start; i < stop; i++)
1137  {
1138  my_dot += m_data[i] * v_data[i];
1139  }
1140  #pragma omp barrier
1141  th_dot(tid) = my_dot;
1142  }
1143  return th_dot.Sum();
1144 #else
1145  // The standard way of computing the dot product is non-deterministic
1146  double prod = 0.0;
1147  #pragma omp parallel for reduction(+:prod)
1148  for (int i = 0; i < size; i++)
1149  {
1150  prod += m_data[i] * v_data[i];
1151  }
1152  return prod;
1153 #endif // MFEM_USE_OPENMP_DETERMINISTIC_DOT
1154  }
1155 #endif // MFEM_USE_OPENMP
1157  {
1158  const int N = size;
1159  auto v_data = v.Read();
1160  auto m_data = Read();
1161  Vector dot(1);
1162  dot.UseDevice(true);
1163  auto d_dot = dot.Write();
1164  dot = 0.0;
1165  MFEM_FORALL(i, N, d_dot[0] += m_data[i] * v_data[i];);
1166  dot.HostReadWrite();
1167  return dot[0];
1168  }
1169 vector_dot_cpu:
1170  return operator*(v_data);
1171 }
1172 
1173 double Vector::Min() const
1174 {
1175  if (size == 0) { return infinity(); }
1176 
1177  const bool use_dev = UseDevice();
1178  auto m_data = Read(use_dev);
1179 
1180  if (!use_dev) { goto vector_min_cpu; }
1181 
1182 #ifdef MFEM_USE_OCCA
1183  if (DeviceCanUseOcca())
1184  {
1185  return occa::linalg::min<double,double>(OccaMemoryRead(data, size));
1186  }
1187 #endif
1188 
1189 #ifdef MFEM_USE_CUDA
1191  {
1192  return cuVectorMin(size, m_data);
1193  }
1194 #endif
1195 
1196 #ifdef MFEM_USE_HIP
1198  {
1199  return hipVectorMin(size, m_data);
1200  }
1201 #endif
1202 
1203 #ifdef MFEM_USE_OPENMP
1205  {
1206  double minimum = m_data[0];
1207  #pragma omp parallel for reduction(min:minimum)
1208  for (int i = 0; i < size; i++)
1209  {
1210  minimum = std::min(minimum, m_data[i]);
1211  }
1212  return minimum;
1213  }
1214 #endif
1215 
1217  {
1218  const int N = size;
1219  auto m_data = Read();
1220  Vector min(1);
1221  min = infinity();
1222  min.UseDevice(true);
1223  auto d_min = min.ReadWrite();
1224  MFEM_FORALL(i, N, d_min[0] = (d_min[0]<m_data[i])?d_min[0]:m_data[i];);
1225  min.HostReadWrite();
1226  return min[0];
1227  }
1228 
1229 vector_min_cpu:
1230  double minimum = data[0];
1231  for (int i = 1; i < size; i++)
1232  {
1233  if (m_data[i] < minimum)
1234  {
1235  minimum = m_data[i];
1236  }
1237  }
1238  return minimum;
1239 }
1240 
1241 
1242 #ifdef MFEM_USE_SUNDIALS
1243 
1244 Vector::Vector(N_Vector nv)
1245 {
1246  N_Vector_ID nvid = N_VGetVectorID(nv);
1247 
1248  switch (nvid)
1249  {
1250  case SUNDIALS_NVEC_SERIAL:
1251  SetDataAndSize(NV_DATA_S(nv), NV_LENGTH_S(nv));
1252  break;
1253 #ifdef MFEM_USE_MPI
1254  case SUNDIALS_NVEC_PARALLEL:
1255  SetDataAndSize(NV_DATA_P(nv), NV_LOCLENGTH_P(nv));
1256  break;
1257 #endif
1258  default:
1259  MFEM_ABORT("N_Vector type " << nvid << " is not supported");
1260  }
1261 }
1262 
1263 void Vector::ToNVector(N_Vector &nv, long global_length)
1264 {
1265  MFEM_ASSERT(nv, "N_Vector handle is NULL");
1266  N_Vector_ID nvid = N_VGetVectorID(nv);
1267 
1268  switch (nvid)
1269  {
1270  case SUNDIALS_NVEC_SERIAL:
1271  MFEM_ASSERT(NV_OWN_DATA_S(nv) == SUNFALSE, "invalid serial N_Vector");
1272  NV_DATA_S(nv) = data;
1273  NV_LENGTH_S(nv) = size;
1274  break;
1275 #ifdef MFEM_USE_MPI
1276  case SUNDIALS_NVEC_PARALLEL:
1277  MFEM_ASSERT(NV_OWN_DATA_P(nv) == SUNFALSE, "invalid parallel N_Vector");
1278  NV_DATA_P(nv) = data;
1279  NV_LOCLENGTH_P(nv) = size;
1280  if (global_length == 0)
1281  {
1282  global_length = NV_GLOBLENGTH_P(nv);
1283 
1284  if (global_length == 0 && global_length != size)
1285  {
1286  MPI_Comm sundials_comm = NV_COMM_P(nv);
1287  long local_size = size;
1288  MPI_Allreduce(&local_size, &global_length, 1, MPI_LONG,
1289  MPI_SUM,sundials_comm);
1290  }
1291  }
1292  NV_GLOBLENGTH_P(nv) = global_length;
1293  break;
1294 #endif
1295  default:
1296  MFEM_ABORT("N_Vector type " << nvid << " is not supported");
1297  }
1298 }
1299 
1300 #endif // MFEM_USE_SUNDIALS
1301 
1302 }
Hash function for data sequences.
Definition: hash.hpp:218
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:551
void SetVector(const Vector &v, int offset)
Definition: vector.cpp:270
Vector()
Default constructor for Vector. Sets size = 0 and data = NULL.
Definition: vector.hpp:70
Memory< T > & GetMemory()
Return a reference to the Memory object used by the Array.
Definition: array.hpp:113
Memory< double > data
Definition: vector.hpp:64
double & Elem(int i)
Access Vector entries. Index i = 0 .. size-1.
Definition: vector.cpp:101
Device memory; using CUDA or HIP *Malloc and *Free.
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:783
double & operator()(int i)
Access Vector entries using () for 0-based indexing.
Definition: vector.hpp:599
Biwise-OR of all HIP backends.
Definition: device.hpp:90
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:525
virtual double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:438
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
T * Write(MemoryClass mc, int size)
Get write-only access to the memory with the given MemoryClass.
bool UseDevice() const
Return the device flag of the Memory object used by the Array.
Definition: array.hpp:119
std::string GetHash() const
Return the hash string for the current sequence and reset (clear) the sequence.
Definition: hash.cpp:60
double Normlinf() const
Returns the l_infinity norm of the vector.
Definition: vector.cpp:801
void CopyFrom(const Memory &src, int size)
Copy size entries from src to *this.
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:763
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:108
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:291
HashFunction & AppendDoubles(const double *doubles, size_t num_doubles)
Add a sequence of doubles for hashing, given as a c-array.
Definition: hash.hpp:271
double operator*(const double *) const
Dot product with a double * array.
Definition: vector.cpp:111
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...
T ZeroSubnormal(T val)
Definition: vector.hpp:473
Vector & operator=(const double *v)
Copy Size() entries from v.
Definition: vector.cpp:124
void Load(std::istream **in, int np, int *dim)
Reads a vector from multiple files.
Definition: vector.cpp:57
double Normlp(double p) const
Returns the l_p norm of the vector.
Definition: vector.cpp:823
void CopyFromHost(const T *src, int size)
Copy size entries from the host pointer src to *this.
bool DeviceCanUseOcca()
Function that determines if an OCCA kernel should be used, based on the current mfem::Device configur...
Definition: occa.hpp:69
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:434
double b
Definition: lissajous.cpp:42
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition: device.hpp:273
const occa::memory OccaMemoryRead(const Memory< T > &mem, size_t size)
Wrap a Memory object as occa::memory for read only access with the mfem::Device MemoryClass. The returned occa::memory is associated with the default occa::device used by MFEM.
Definition: occa.hpp:37
void median(const Vector &lo, const Vector &hi)
v = median(v,lo,hi) entrywise. Implementation assumes lo &lt;= hi.
Definition: vector.cpp:501
Biwise-OR of all OpenMP backends.
Definition: device.hpp:92
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:300
MFEM_HOST_DEVICE double Norml2(const int size, const T *data)
Returns the l2 norm of the Vector with given size and data.
Definition: kernels.hpp:133
void Reset()
Reset the memory to be empty, ensuring that Delete() will be a no-op.
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof...
Definition: vector.cpp:617
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:430
Vector & operator/=(double c)
Definition: vector.cpp:179
Biwise-OR of all CUDA backends.
Definition: device.hpp:88
void SetSubVectorComplement(const Array< int > &dofs, const double val)
Set all vector entries NOT in the dofs Array to the given val.
Definition: vector.cpp:686
Vector & operator*=(double c)
Definition: vector.cpp:158
double Min() const
Returns the minimal element of the vector.
Definition: vector.cpp:1173
Vector & operator+=(double c)
Definition: vector.cpp:222
double Norml1() const
Returns the l_1 norm of the vector.
Definition: vector.cpp:812
static MemoryType GetHostMemoryType()
Get the current Host MemoryType. This is the MemoryType used by most MFEM classes when allocating mem...
Definition: device.hpp:264
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
Definition: vector.cpp:702
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:438
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:674
void Print_HYPRE(std::ostream &out) const
Prints vector to stream out in HYPRE_Vector format.
Definition: vector.cpp:736
static bool Allows(unsigned long b_mask)
Return true if any of the backends in the backend mask, b_mask, are allowed.
Definition: device.hpp:258
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Definition: vector.hpp:147
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:258
double a
Definition: lissajous.cpp:41
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:442
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:243
Host memory; using new[] and delete[].
virtual MFEM_DEPRECATED N_Vector ToNVector()
(DEPRECATED) Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_SERIAL. ...
Definition: vector.hpp:455
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:873
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
bool Empty() const
Return true if the Memory object is empty, see Reset().
T * ReadWrite(MemoryClass mc, int size)
Get read-write access to the memory with the given MemoryClass.
int dim
Definition: ex24.cpp:53
Vector & operator-=(double c)
Definition: vector.cpp:201
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:46
const double alpha
Definition: ex15.cpp:369
Vector data type.
Definition: vector.hpp:60
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition: vector.hpp:111
void PrintHash(std::ostream &out) const
Print the Vector size and hash of its data.
Definition: vector.cpp:755
RefCoord s[3]
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
const T * Read(MemoryClass mc, int size) const
Get read-only access to the memory with the given MemoryClass.
double Sum() const
Return the sum of the vector entries.
Definition: vector.cpp:891
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:426
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:446
[device] Debug backend: host memory is READ/WRITE protected while a device is in use. It allows to test the &quot;device&quot; code-path (using separate host/device memory pools and host &lt;-&gt; device transfers) without any GPU hardware. As &#39;DEBUG&#39; is sometimes used as a macro, _DEVICE has been added to avoid conflicts.
Definition: device.hpp:75
void Neg()
(*this) = -(*this)
Definition: vector.cpp:283