MFEM  v4.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, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 // Implementation of data type vector
13 
14 #include "vector.hpp"
15 #include "dtensor.hpp"
16 #include "../general/forall.hpp"
17 
18 #if defined(MFEM_USE_SUNDIALS) && defined(MFEM_USE_MPI)
19 #include <nvector/nvector_parallel.h>
20 #include <nvector/nvector_parhyp.h>
21 #endif
22 
23 #include <iostream>
24 #include <iomanip>
25 #include <cmath>
26 #include <cstdlib>
27 #include <ctime>
28 #include <limits>
29 
30 namespace mfem
31 {
32 
34 {
35  const int s = v.Size();
36  if (s > 0)
37  {
38  MFEM_ASSERT(!v.data.Empty(), "invalid source vector");
39  size = s;
40  data.New(s, v.data.GetMemoryType());
41  data.CopyFrom(v.data, s);
42  }
43  else
44  {
45  size = 0;
46  data.Reset();
47  }
48  UseDevice(v.UseDevice());
49 }
50 
51 void Vector::Load(std::istream **in, int np, int *dim)
52 {
53  int i, j, s;
54 
55  s = 0;
56  for (i = 0; i < np; i++)
57  {
58  s += dim[i];
59  }
60 
61  SetSize(s);
62 
63  int p = 0;
64  for (i = 0; i < np; i++)
65  {
66  for (j = 0; j < dim[i]; j++)
67  {
68  *in[i] >> data[p++];
69  }
70  }
71 }
72 
73 void Vector::Load(std::istream &in, int Size)
74 {
75  SetSize(Size);
76 
77  for (int i = 0; i < size; i++)
78  {
79  in >> data[i];
80  }
81 }
82 
83 double &Vector::Elem(int i)
84 {
85  return operator()(i);
86 }
87 
88 const double &Vector::Elem(int i) const
89 {
90  return operator()(i);
91 }
92 
93 double Vector::operator*(const double *v) const
94 {
95  double dot = 0.0;
96 #ifdef MFEM_USE_LEGACY_OPENMP
97  #pragma omp parallel for reduction(+:dot)
98 #endif
99  for (int i = 0; i < size; i++)
100  {
101  dot += data[i] * v[i];
102  }
103  return dot;
104 }
105 
106 Vector &Vector::operator=(const double *v)
107 {
108  data.CopyFromHost(v, size);
109  return *this;
110 }
111 
113 {
114 #if 0
115  SetSize(v.Size(), v.data.GetMemoryType());
116  data.CopyFrom(v.data, v.Size());
117  UseDevice(v.UseDevice());
118 #else
119  SetSize(v.Size());
120  const bool use_dev = UseDevice() || v.UseDevice();
121  v.UseDevice(use_dev);
122  // keep 'data' where it is, unless 'use_dev' is true
123  if (use_dev) { Write(); }
124  data.CopyFrom(v.data, v.Size());
125 #endif
126  return *this;
127 }
128 
129 Vector &Vector::operator=(double value)
130 {
131  const bool use_dev = UseDevice();
132  const int N = size;
133  auto y = Write(use_dev);
134  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = value;);
135  return *this;
136 }
137 
139 {
140  const bool use_dev = UseDevice();
141  const int N = size;
142  auto y = ReadWrite(use_dev);
143  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] *= c;);
144  return *this;
145 }
146 
148 {
149  const bool use_dev = UseDevice();
150  const int N = size;
151  const double m = 1.0/c;
152  auto y = ReadWrite(use_dev);
153  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] *= m;);
154  return *this;
155 }
156 
158 {
159  const bool use_dev = UseDevice();
160  const int N = size;
161  auto y = ReadWrite(use_dev);
162  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] -= c;);
163  return *this;
164 }
165 
167 {
168  MFEM_ASSERT(size == v.size, "incompatible Vectors!");
169 
170  const bool use_dev = UseDevice() || v.UseDevice();
171  const int N = size;
172  auto y = ReadWrite(use_dev);
173  auto x = v.Read(use_dev);
174  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] -= x[i];);
175  return *this;
176 }
177 
179 {
180  MFEM_ASSERT(size == v.size, "incompatible Vectors!");
181 
182  const bool use_dev = UseDevice() || v.UseDevice();
183  const int N = size;
184  auto y = ReadWrite(use_dev);
185  auto x = v.Read(use_dev);
186  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] += x[i];);
187  return *this;
188 }
189 
190 Vector &Vector::Add(const double a, const Vector &Va)
191 {
192  MFEM_ASSERT(size == Va.size, "incompatible Vectors!");
193 
194  if (a != 0.0)
195  {
196  const int N = size;
197  const bool use_dev = UseDevice() || Va.UseDevice();
198  auto y = ReadWrite(use_dev);
199  auto x = Va.Read(use_dev);
200  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] += a * x[i];);
201  }
202  return *this;
203 }
204 
205 Vector &Vector::Set(const double a, const Vector &Va)
206 {
207  MFEM_ASSERT(size == Va.size, "incompatible Vectors!");
208 
209  const bool use_dev = UseDevice() || Va.UseDevice();
210  const int N = size;
211  auto x = Va.Read(use_dev);
212  auto y = Write(use_dev);
213  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = a * x[i];);
214  return *this;
215 }
216 
217 void Vector::SetVector(const Vector &v, int offset)
218 {
219  MFEM_ASSERT(v.Size() + offset <= size, "invalid sub-vector");
220 
221  const int vs = v.Size();
222  const double *vp = v.data;
223  double *p = data + offset;
224  for (int i = 0; i < vs; i++)
225  {
226  p[i] = vp[i];
227  }
228 }
229 
231 {
232  const bool use_dev = UseDevice();
233  const int N = size;
234  auto y = ReadWrite(use_dev);
235  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = -y[i];);
236 }
237 
238 void add(const Vector &v1, const Vector &v2, Vector &v)
239 {
240  MFEM_ASSERT(v.size == v1.size && v.size == v2.size,
241  "incompatible Vectors!");
242 
243 #if !defined(MFEM_USE_LEGACY_OPENMP)
244  const bool use_dev = v1.UseDevice() || v2.UseDevice() || v.UseDevice();
245  const int N = v.size;
246  // Note: get read access first, in case v is the same as v1/v2.
247  auto x1 = v1.Read(use_dev);
248  auto x2 = v2.Read(use_dev);
249  auto y = v.Write(use_dev);
250  MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = x1[i] + x2[i];);
251 #else
252  #pragma omp parallel for
253  for (int i = 0; i < v.size; i++)
254  {
255  v.data[i] = v1.data[i] + v2.data[i];
256  }
257 #endif
258 }
259 
260 void add(const Vector &v1, double alpha, const Vector &v2, Vector &v)
261 {
262  MFEM_ASSERT(v.size == v1.size && v.size == v2.size,
263  "incompatible Vectors!");
264 
265  if (alpha == 0.0)
266  {
267  v = v1;
268  }
269  else if (alpha == 1.0)
270  {
271  add(v1, v2, v);
272  }
273  else
274  {
275 #if !defined(MFEM_USE_LEGACY_OPENMP)
276  const bool use_dev = v1.UseDevice() || v2.UseDevice() || v.UseDevice();
277  const int N = v.size;
278  // Note: get read access first, in case v is the same as v1/v2.
279  auto d_x = v1.Read(use_dev);
280  auto d_y = v2.Read(use_dev);
281  auto d_z = v.Write(use_dev);
282  MFEM_FORALL_SWITCH(use_dev, i, N, d_z[i] = d_x[i] + alpha * d_y[i];);
283 #else
284  const double *v1p = v1.data, *v2p = v2.data;
285  double *vp = v.data;
286  const int s = v.size;
287  #pragma omp parallel for
288  for (int i = 0; i < s; i++)
289  {
290  vp[i] = v1p[i] + alpha*v2p[i];
291  }
292 #endif
293  }
294 }
295 
296 void add(const double a, const Vector &x, const Vector &y, Vector &z)
297 {
298  MFEM_ASSERT(x.size == y.size && x.size == z.size,
299  "incompatible Vectors!");
300 
301  if (a == 0.0)
302  {
303  z = 0.0;
304  }
305  else if (a == 1.0)
306  {
307  add(x, y, z);
308  }
309  else
310  {
311 #if !defined(MFEM_USE_LEGACY_OPENMP)
312  const bool use_dev = x.UseDevice() || y.UseDevice() || z.UseDevice();
313  const int N = x.size;
314  // Note: get read access first, in case z is the same as x/y.
315  auto xd = x.Read(use_dev);
316  auto yd = y.Read(use_dev);
317  auto zd = z.Write(use_dev);
318  MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = a * (xd[i] + yd[i]););
319 #else
320  const double *xp = x.data;
321  const double *yp = y.data;
322  double *zp = z.data;
323  const int s = x.size;
324  #pragma omp parallel for
325  for (int i = 0; i < s; i++)
326  {
327  zp[i] = a * (xp[i] + yp[i]);
328  }
329 #endif
330  }
331 }
332 
333 void add(const double a, const Vector &x,
334  const double b, const Vector &y, Vector &z)
335 {
336  MFEM_ASSERT(x.size == y.size && x.size == z.size,
337  "incompatible Vectors!");
338 
339  if (a == 0.0)
340  {
341  z.Set(b, y);
342  }
343  else if (b == 0.0)
344  {
345  z.Set(a, x);
346  }
347 #if 0
348  else if (a == 1.0)
349  {
350  add(x, b, y, z);
351  }
352  else if (b == 1.0)
353  {
354  add(y, a, x, z);
355  }
356  else if (a == b)
357  {
358  add(a, x, y, z);
359  }
360 #endif
361  else
362  {
363 #if !defined(MFEM_USE_LEGACY_OPENMP)
364  const bool use_dev = x.UseDevice() || y.UseDevice() || z.UseDevice();
365  const int N = x.size;
366  // Note: get read access first, in case z is the same as x/y.
367  auto xd = x.Read(use_dev);
368  auto yd = y.Read(use_dev);
369  auto zd = z.Write(use_dev);
370  MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = a * xd[i] + b * yd[i];);
371 #else
372  const double *xp = x.data;
373  const double *yp = y.data;
374  double *zp = z.data;
375  const int s = x.size;
376  #pragma omp parallel for
377  for (int i = 0; i < s; i++)
378  {
379  zp[i] = a * xp[i] + b * yp[i];
380  }
381 #endif
382  }
383 }
384 
385 void subtract(const Vector &x, const Vector &y, Vector &z)
386 {
387  MFEM_ASSERT(x.size == y.size && x.size == z.size,
388  "incompatible Vectors!");
389 
390 #if !defined(MFEM_USE_LEGACY_OPENMP)
391  const bool use_dev = x.UseDevice() || y.UseDevice() || z.UseDevice();
392  const int N = x.size;
393  // Note: get read access first, in case z is the same as x/y.
394  auto xd = x.Read(use_dev);
395  auto yd = y.Read(use_dev);
396  auto zd = z.Write(use_dev);
397  MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = xd[i] - yd[i];);
398 #else
399  const double *xp = x.data;
400  const double *yp = y.data;
401  double *zp = z.data;
402  const int s = x.size;
403  #pragma omp parallel for
404  for (int i = 0; i < s; i++)
405  {
406  zp[i] = xp[i] - yp[i];
407  }
408 #endif
409 }
410 
411 void subtract(const double a, const Vector &x, const Vector &y, Vector &z)
412 {
413  MFEM_ASSERT(x.size == y.size && x.size == z.size,
414  "incompatible Vectors!");
415 
416  if (a == 0.)
417  {
418  z = 0.;
419  }
420  else if (a == 1.)
421  {
422  subtract(x, y, z);
423  }
424  else
425  {
426 #if !defined(MFEM_USE_LEGACY_OPENMP)
427  const bool use_dev = x.UseDevice() || y.UseDevice() || z.UseDevice();
428  const int N = x.size;
429  // Note: get read access first, in case z is the same as x/y.
430  auto xd = x.Read(use_dev);
431  auto yd = y.Read(use_dev);
432  auto zd = z.Write(use_dev);
433  MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = a * (xd[i] - yd[i]););
434 #else
435  const double *xp = x.data;
436  const double *yp = y.data;
437  double *zp = z.data;
438  const int s = x.size;
439  #pragma omp parallel for
440  for (int i = 0; i < s; i++)
441  {
442  zp[i] = a * (xp[i] - yp[i]);
443  }
444 #endif
445  }
446 }
447 
448 void Vector::median(const Vector &lo, const Vector &hi)
449 {
450  MFEM_ASSERT(size == lo.size && size == hi.size,
451  "incompatible Vectors!");
452 
453  const bool use_dev = UseDevice() || lo.UseDevice() || hi.UseDevice();
454  const int N = size;
455  // Note: get read access first, in case *this is the same as lo/hi.
456  auto l = lo.Read(use_dev);
457  auto h = hi.Read(use_dev);
458  auto m = Write(use_dev);
459  MFEM_FORALL_SWITCH(use_dev, i, N,
460  {
461  if (m[i] < l[i])
462  {
463  m[i] = l[i];
464  }
465  else if (m[i] > h[i])
466  {
467  m[i] = h[i];
468  }
469  });
470 }
471 
472 void Vector::GetSubVector(const Array<int> &dofs, Vector &elemvect) const
473 {
474  const int n = dofs.Size();
475  elemvect.SetSize(n);
476  const bool use_dev = dofs.UseDevice() || elemvect.UseDevice();
477  auto d_y = elemvect.Write(use_dev);
478  auto d_X = Read(use_dev);
479  auto d_dofs = dofs.Read(use_dev);
480  MFEM_FORALL_SWITCH(use_dev, i, n,
481  {
482  const int dof_i = d_dofs[i];
483  d_y[i] = dof_i >= 0 ? d_X[dof_i] : -d_X[-dof_i-1];
484  });
485 }
486 
487 void Vector::GetSubVector(const Array<int> &dofs, double *elem_data) const
488 {
490  const int n = dofs.Size();
491  for (int i = 0; i < n; i++)
492  {
493  const int j = dofs[i];
494  elem_data[i] = (j >= 0) ? data[j] : -data[-1-j];
495  }
496 }
497 
498 void Vector::SetSubVector(const Array<int> &dofs, const double value)
499 {
500  const bool use_dev = dofs.UseDevice();
501  const int n = dofs.Size();
502  // Use read+write access for *this - we only modify some of its entries
503  auto d_X = ReadWrite(use_dev);
504  auto d_dofs = dofs.Read(use_dev);
505  MFEM_FORALL_SWITCH(use_dev, i, n,
506  {
507  const int j = d_dofs[i];
508  if (j >= 0)
509  {
510  d_X[j] = value;
511  }
512  else
513  {
514  d_X[-1-j] = -value;
515  }
516  });
517 }
518 
519 void Vector::SetSubVector(const Array<int> &dofs, const Vector &elemvect)
520 {
521  MFEM_ASSERT(dofs.Size() == elemvect.Size(),
522  "Size mismatch: length of dofs is " << dofs.Size()
523  << ", length of elemvect is " << elemvect.Size());
524 
525  const bool use_dev = dofs.UseDevice() || elemvect.UseDevice();
526  const int n = dofs.Size();
527  // Use read+write access for X - we only modify some of its entries
528  auto d_X = ReadWrite(use_dev);
529  auto d_y = elemvect.Read(use_dev);
530  auto d_dofs = dofs.Read(use_dev);
531  MFEM_FORALL_SWITCH(use_dev, i, n,
532  {
533  const int dof_i = d_dofs[i];
534  if (dof_i >= 0)
535  {
536  d_X[dof_i] = d_y[i];
537  }
538  else
539  {
540  d_X[-1-dof_i] = -d_y[i];
541  }
542  });
543 }
544 
545 void Vector::SetSubVector(const Array<int> &dofs, double *elem_data)
546 {
547  // Use read+write access because we overwrite only part of the data.
549  const int n = dofs.Size();
550  for (int i = 0; i < n; i++)
551  {
552  const int j= dofs[i];
553  if (j >= 0)
554  {
555  operator()(j) = elem_data[i];
556  }
557  else
558  {
559  operator()(-1-j) = -elem_data[i];
560  }
561  }
562 }
563 
564 void Vector::AddElementVector(const Array<int> &dofs, const Vector &elemvect)
565 {
566  MFEM_ASSERT(dofs.Size() == elemvect.Size(), "Size mismatch: "
567  "length of dofs is " << dofs.Size() <<
568  ", length of elemvect is " << elemvect.Size());
569 
570  const bool use_dev = dofs.UseDevice() || elemvect.UseDevice();
571  const int n = dofs.Size();
572  auto d_y = elemvect.Read(use_dev);
573  auto d_X = ReadWrite(use_dev);
574  auto d_dofs = dofs.Read(use_dev);
575  MFEM_FORALL_SWITCH(use_dev, i, n,
576  {
577  const int j = d_dofs[i];
578  if (j >= 0)
579  {
580  d_X[j] += d_y[i];
581  }
582  else
583  {
584  d_X[-1-j] -= d_y[i];
585  }
586  });
587 }
588 
589 void Vector::AddElementVector(const Array<int> &dofs, double *elem_data)
590 {
592  const int n = dofs.Size();
593  for (int i = 0; i < n; i++)
594  {
595  const int j = dofs[i];
596  if (j >= 0)
597  {
598  operator()(j) += elem_data[i];
599  }
600  else
601  {
602  operator()(-1-j) -= elem_data[i];
603  }
604  }
605 }
606 
607 void Vector::AddElementVector(const Array<int> &dofs, const double a,
608  const Vector &elemvect)
609 {
610  MFEM_ASSERT(dofs.Size() == elemvect.Size(), "Size mismatch: "
611  "length of dofs is " << dofs.Size() <<
612  ", length of elemvect is " << elemvect.Size());
613 
614  const bool use_dev = dofs.UseDevice() || elemvect.UseDevice();
615  const int n = dofs.Size();
616  auto d_y = ReadWrite(use_dev);
617  auto d_x = elemvect.Read(use_dev);
618  auto d_dofs = dofs.Read(use_dev);
619  MFEM_FORALL_SWITCH(use_dev, i, n,
620  {
621  const int j = d_dofs[i];
622  if (j >= 0)
623  {
624  d_y[j] += a * d_x[i];
625  }
626  else
627  {
628  d_y[-1-j] -= a * d_x[i];
629  }
630  });
631 }
632 
633 void Vector::SetSubVectorComplement(const Array<int> &dofs, const double val)
634 {
635  const bool use_dev = UseDevice() || dofs.UseDevice();
636  const int n = dofs.Size();
637  const int N = size;
638  Vector dofs_vals(n, use_dev ? Device::GetMemoryType() : MemoryType::HOST);
639  auto d_data = ReadWrite(use_dev);
640  auto d_dofs_vals = dofs_vals.Write(use_dev);
641  auto d_dofs = dofs.Read(use_dev);
642  MFEM_FORALL_SWITCH(use_dev, i, n, d_dofs_vals[i] = d_data[d_dofs[i]];);
643  MFEM_FORALL_SWITCH(use_dev, i, N, d_data[i] = val;);
644  MFEM_FORALL_SWITCH(use_dev, i, n, d_data[d_dofs[i]] = d_dofs_vals[i];);
645 }
646 
647 void Vector::Print(std::ostream &out, int width) const
648 {
649  if (!size) { return; }
651  for (int i = 0; 1; )
652  {
653  out << data[i];
654  i++;
655  if (i == size)
656  {
657  break;
658  }
659  if ( i % width == 0 )
660  {
661  out << '\n';
662  }
663  else
664  {
665  out << ' ';
666  }
667  }
668  out << '\n';
669 }
670 
671 void Vector::Print_HYPRE(std::ostream &out) const
672 {
673  int i;
674  std::ios::fmtflags old_fmt = out.flags();
675  out.setf(std::ios::scientific);
676  std::streamsize old_prec = out.precision(14);
677 
678  out << size << '\n'; // number of rows
679 
680  data.Read(MemoryClass::HOST, size);
681  for (i = 0; i < size; i++)
682  {
683  out << data[i] << '\n';
684  }
685 
686  out.precision(old_prec);
687  out.flags(old_fmt);
688 }
689 
690 void Vector::Randomize(int seed)
691 {
692  // static unsigned int seed = time(0);
693  const double max = (double)(RAND_MAX) + 1.;
694 
695  if (seed == 0)
696  {
697  seed = (int)time(0);
698  }
699 
700  // srand(seed++);
701  srand((unsigned)seed);
702 
703  for (int i = 0; i < size; i++)
704  {
705  data[i] = std::abs(rand()/max);
706  }
707 }
708 
709 double Vector::Norml2() const
710 {
711  // Scale entries of Vector on the fly, using algorithms from
712  // std::hypot() and LAPACK's drm2. This scaling ensures that the
713  // argument of each call to std::pow is <= 1 to avoid overflow.
714  if (0 == size)
715  {
716  return 0.0;
717  } // end if 0 == size
718 
719  if (1 == size)
720  {
721  return std::abs(data[0]);
722  } // end if 1 == size
723 
724  double scale = 0.0;
725  double sum = 0.0;
726 
727  for (int i = 0; i < size; i++)
728  {
729  if (data[i] != 0.0)
730  {
731  const double absdata = std::abs(data[i]);
732  if (scale <= absdata)
733  {
734  const double sqr_arg = scale / absdata;
735  sum = 1.0 + sum * (sqr_arg * sqr_arg);
736  scale = absdata;
737  continue;
738  } // end if scale <= absdata
739  const double sqr_arg = absdata / scale;
740  sum += (sqr_arg * sqr_arg); // else scale > absdata
741  } // end if data[i] != 0
742  }
743  return scale * std::sqrt(sum);
744 }
745 
746 double Vector::Normlinf() const
747 {
748  double max = 0.0;
749  for (int i = 0; i < size; i++)
750  {
751  max = std::max(std::abs(data[i]), max);
752  }
753  return max;
754 }
755 
756 double Vector::Norml1() const
757 {
758  double sum = 0.0;
759  for (int i = 0; i < size; i++)
760  {
761  sum += std::abs(data[i]);
762  }
763  return sum;
764 }
765 
766 double Vector::Normlp(double p) const
767 {
768  MFEM_ASSERT(p > 0.0, "Vector::Normlp");
769 
770  if (p == 1.0)
771  {
772  return Norml1();
773  }
774  if (p == 2.0)
775  {
776  return Norml2();
777  }
778  if (p < infinity())
779  {
780  // Scale entries of Vector on the fly, using algorithms from
781  // std::hypot() and LAPACK's drm2. This scaling ensures that the
782  // argument of each call to std::pow is <= 1 to avoid overflow.
783  if (0 == size)
784  {
785  return 0.0;
786  } // end if 0 == size
787 
788  if (1 == size)
789  {
790  return std::abs(data[0]);
791  } // end if 1 == size
792 
793  double scale = 0.0;
794  double sum = 0.0;
795 
796  for (int i = 0; i < size; i++)
797  {
798  if (data[i] != 0.0)
799  {
800  const double absdata = std::abs(data[i]);
801  if (scale <= absdata)
802  {
803  sum = 1.0 + sum * std::pow(scale / absdata, p);
804  scale = absdata;
805  continue;
806  } // end if scale <= absdata
807  sum += std::pow(absdata / scale, p); // else scale > absdata
808  } // end if data[i] != 0
809  }
810  return scale * std::pow(sum, 1.0/p);
811  } // end if p < infinity()
812 
813  return Normlinf(); // else p >= infinity()
814 }
815 
816 double Vector::Max() const
817 {
818  if (size == 0) { return -infinity(); }
819 
820  double max = data[0];
821 
822  for (int i = 1; i < size; i++)
823  {
824  if (data[i] > max)
825  {
826  max = data[i];
827  }
828  }
829 
830  return max;
831 }
832 
833 double Vector::Sum() const
834 {
835  double sum = 0.0;
836 
837  for (int i = 0; i < size; i++)
838  {
839  sum += data[i];
840  }
841 
842  return sum;
843 }
844 
845 #ifdef MFEM_USE_CUDA
846 static __global__ void cuKernelMin(const int N, double *gdsr, const double *x)
847 {
848  __shared__ double s_min[MFEM_CUDA_BLOCKS];
849  const int n = blockDim.x*blockIdx.x + threadIdx.x;
850  if (n>=N) { return; }
851  const int bid = blockIdx.x;
852  const int tid = threadIdx.x;
853  const int bbd = bid*blockDim.x;
854  const int rid = bbd+tid;
855  s_min[tid] = x[n];
856  for (int workers=blockDim.x>>1; workers>0; workers>>=1)
857  {
858  __syncthreads();
859  if (tid >= workers) { continue; }
860  if (rid >= N) { continue; }
861  const int dualTid = tid + workers;
862  if (dualTid >= N) { continue; }
863  const int rdd = bbd+dualTid;
864  if (rdd >= N) { continue; }
865  if (dualTid >= blockDim.x) { continue; }
866  s_min[tid] = fmin(s_min[tid], s_min[dualTid]);
867  }
868  if (tid==0) { gdsr[bid] = s_min[0]; }
869 }
870 
871 static Array<double> cuda_reduce_buf;
872 
873 static double cuVectorMin(const int N, const double *X)
874 {
875  const int tpb = MFEM_CUDA_BLOCKS;
876  const int blockSize = MFEM_CUDA_BLOCKS;
877  const int gridSize = (N+blockSize-1)/blockSize;
878  const int min_sz = (N%tpb)==0? (N/tpb) : (1+N/tpb);
879  cuda_reduce_buf.SetSize(min_sz);
880  Memory<double> &buf = cuda_reduce_buf.GetMemory();
881  double *d_min = buf.Write(MemoryClass::CUDA, min_sz);
882  cuKernelMin<<<gridSize,blockSize>>>(N, d_min, X);
883  MFEM_CUDA_CHECK(cudaGetLastError());
884  const double *h_min = buf.Read(MemoryClass::HOST, min_sz);
886  for (int i = 0; i < min_sz; i++) { min = fmin(min, h_min[i]); }
887  return min;
888 }
889 
890 static __global__ void cuKernelDot(const int N, double *gdsr,
891  const double *x, const double *y)
892 {
893  __shared__ double s_dot[MFEM_CUDA_BLOCKS];
894  const int n = blockDim.x*blockIdx.x + threadIdx.x;
895  if (n>=N) { return; }
896  const int bid = blockIdx.x;
897  const int tid = threadIdx.x;
898  const int bbd = bid*blockDim.x;
899  const int rid = bbd+tid;
900  s_dot[tid] = x[n] * y[n];
901  for (int workers=blockDim.x>>1; workers>0; workers>>=1)
902  {
903  __syncthreads();
904  if (tid >= workers) { continue; }
905  if (rid >= N) { continue; }
906  const int dualTid = tid + workers;
907  if (dualTid >= N) { continue; }
908  const int rdd = bbd+dualTid;
909  if (rdd >= N) { continue; }
910  if (dualTid >= blockDim.x) { continue; }
911  s_dot[tid] += s_dot[dualTid];
912  }
913  if (tid==0) { gdsr[bid] = s_dot[0]; }
914 }
915 
916 static double cuVectorDot(const int N, const double *X, const double *Y)
917 {
918  const int tpb = MFEM_CUDA_BLOCKS;
919  const int blockSize = MFEM_CUDA_BLOCKS;
920  const int gridSize = (N+blockSize-1)/blockSize;
921  const int dot_sz = (N%tpb)==0? (N/tpb) : (1+N/tpb);
922  cuda_reduce_buf.SetSize(dot_sz);
923  Memory<double> &buf = cuda_reduce_buf.GetMemory();
924  double *d_dot = buf.Write(MemoryClass::CUDA, dot_sz);
925  cuKernelDot<<<gridSize,blockSize>>>(N, d_dot, X, Y);
926  MFEM_CUDA_CHECK(cudaGetLastError());
927  const double *h_dot = buf.Read(MemoryClass::HOST, dot_sz);
928  double dot = 0.0;
929  for (int i = 0; i < dot_sz; i++) { dot += h_dot[i]; }
930  return dot;
931 }
932 #endif // MFEM_USE_CUDA
933 
934 double Vector::operator*(const Vector &v) const
935 {
936  MFEM_ASSERT(size == v.size, "incompatible Vectors!");
937 
938  const bool use_dev = UseDevice() || v.UseDevice();
939 #if defined(MFEM_USE_CUDA) || defined(MFEM_USE_OPENMP)
940  auto m_data = Read(use_dev);
941 #else
942  Read(use_dev);
943 #endif
944  auto v_data = v.Read(use_dev);
945 
946  if (!use_dev) { goto vector_dot_cpu; }
947 
948 #ifdef MFEM_USE_OCCA
949  if (DeviceCanUseOcca())
950  {
951  return occa::linalg::dot<double,double,double>(
953  }
954 #endif
955 
956 #ifdef MFEM_USE_CUDA
958  {
959  return cuVectorDot(size, m_data, v_data);
960  }
961 #endif
962 
963 #ifdef MFEM_USE_OPENMP
965  {
966  double prod = 0.0;
967  #pragma omp parallel for reduction(+:prod)
968  for (int i = 0; i < size; i++)
969  {
970  prod += m_data[i] * v_data[i];
971  }
972  return prod;
973  }
974 #endif
975 
976 vector_dot_cpu:
977  return operator*(v_data);
978 }
979 
980 double Vector::Min() const
981 {
982  if (size == 0) { return infinity(); }
983 
984  const bool use_dev = UseDevice();
985  auto m_data = Read(use_dev);
986 
987  if (!use_dev) { goto vector_min_cpu; }
988 
989 #ifdef MFEM_USE_OCCA
990  if (DeviceCanUseOcca())
991  {
992  return occa::linalg::min<double,double>(OccaMemoryRead(data, size));
993  }
994 #endif
995 
996 #ifdef MFEM_USE_CUDA
998  {
999  return cuVectorMin(size, m_data);
1000  }
1001 #endif
1002 
1003 #ifdef MFEM_USE_OPENMP
1005  {
1006  double minimum = m_data[0];
1007  #pragma omp parallel for reduction(min:minimum)
1008  for (int i = 0; i < size; i++)
1009  {
1010  minimum = std::min(minimum, m_data[i]);
1011  }
1012  return minimum;
1013  }
1014 #endif
1015 
1016 vector_min_cpu:
1017  double minimum = data[0];
1018  for (int i = 1; i < size; i++)
1019  {
1020  if (m_data[i] < minimum)
1021  {
1022  minimum = m_data[i];
1023  }
1024  }
1025  return minimum;
1026 }
1027 
1028 
1029 #ifdef MFEM_USE_SUNDIALS
1030 
1031 #ifndef SUNTRUE
1032 #define SUNTRUE TRUE
1033 #endif
1034 #ifndef SUNFALSE
1035 #define SUNFALSE FALSE
1036 #endif
1037 
1038 Vector::Vector(N_Vector nv)
1039 {
1040  N_Vector_ID nvid = N_VGetVectorID(nv);
1041  switch (nvid)
1042  {
1043  case SUNDIALS_NVEC_SERIAL:
1044  SetDataAndSize(NV_DATA_S(nv), NV_LENGTH_S(nv));
1045  break;
1046 #ifdef MFEM_USE_MPI
1047  case SUNDIALS_NVEC_PARALLEL:
1048  SetDataAndSize(NV_DATA_P(nv), NV_LOCLENGTH_P(nv));
1049  break;
1050  case SUNDIALS_NVEC_PARHYP:
1051  {
1052  hypre_Vector *hpv_local = N_VGetVector_ParHyp(nv)->local_vector;
1053  SetDataAndSize(hpv_local->data, hpv_local->size);
1054  break;
1055  }
1056 #endif
1057  default:
1058  MFEM_ABORT("N_Vector type " << nvid << " is not supported");
1059  }
1060 }
1061 
1062 void Vector::ToNVector(N_Vector &nv)
1063 {
1064  MFEM_ASSERT(nv, "N_Vector handle is NULL");
1065  N_Vector_ID nvid = N_VGetVectorID(nv);
1066  switch (nvid)
1067  {
1068  case SUNDIALS_NVEC_SERIAL:
1069  MFEM_ASSERT(NV_OWN_DATA_S(nv) == SUNFALSE, "invalid serial N_Vector");
1070  NV_DATA_S(nv) = data;
1071  NV_LENGTH_S(nv) = size;
1072  break;
1073 #ifdef MFEM_USE_MPI
1074  case SUNDIALS_NVEC_PARALLEL:
1075  MFEM_ASSERT(NV_OWN_DATA_P(nv) == SUNFALSE, "invalid parallel N_Vector");
1076  NV_DATA_P(nv) = data;
1077  NV_LOCLENGTH_P(nv) = size;
1078  break;
1079  case SUNDIALS_NVEC_PARHYP:
1080  {
1081  hypre_Vector *hpv_local = N_VGetVector_ParHyp(nv)->local_vector;
1082  MFEM_ASSERT(hpv_local->owns_data == false, "invalid hypre N_Vector");
1083  hpv_local->data = data;
1084  hpv_local->size = size;
1085  break;
1086  }
1087 #endif
1088  default:
1089  MFEM_ABORT("N_Vector type " << nvid << " is not supported");
1090  }
1091 }
1092 
1093 #endif // MFEM_USE_SUNDIALS
1094 
1095 }
int Size() const
Logical size of the array.
Definition: array.hpp:118
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:498
void SetVector(const Vector &v, int offset)
Definition: vector.cpp:217
Vector()
Default constructor for Vector. Sets size = 0 and data = NULL.
Definition: vector.hpp:58
Memory< T > & GetMemory()
Return a reference to the Memory object used by the Array.
Definition: array.hpp:97
Memory< double > data
Definition: vector.hpp:52
double & Elem(int i)
Access Vector entries. Index i = 0 .. size-1.
Definition: vector.cpp:83
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:400
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:709
double & operator()(int i)
Access Vector entries using () for 0-based indexing.
Definition: vector.hpp:467
void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:86
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Definition: vector.cpp:472
int Size() const
Returns the size of the vector.
Definition: vector.hpp:150
Memory types: { CUDA, CUDA_UVM }.
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:103
double Normlinf() const
Returns the l_infinity norm of the vector.
Definition: vector.cpp:746
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:690
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:334
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:238
double operator*(const double *) const
Dot product with a double * array.
Definition: vector.cpp:93
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...
unsigned flags
Vector & operator=(const double *v)
Copy Size() entries from v.
Definition: vector.cpp:106
void Load(std::istream **in, int np, int *dim)
Reads a vector from multiple files.
Definition: vector.cpp:51
double Normlp(double p) const
Returns the l_p norm of the vector.
Definition: vector.cpp:766
int dim
Definition: ex3.cpp:48
double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:342
bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition: vector.hpp:89
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
static MemoryType GetMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition: device.hpp:210
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:448
Biwise-OR of all OpenMP backends.
Definition: device.hpp:69
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:261
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 (element) subvector to the vector.
Definition: vector.cpp:564
Vector & operator/=(double c)
Definition: vector.cpp:147
Biwise-OR of all CUDA backends.
Definition: device.hpp:67
void SetSubVectorComplement(const Array< int > &dofs, const double val)
Set all vector entries NOT in the &#39;dofs&#39; array to the given &#39;val&#39;.
Definition: vector.cpp:633
Vector & operator*=(double c)
Definition: vector.cpp:138
double Min() const
Returns the minimal element of the vector.
Definition: vector.cpp:980
double Norml1() const
Returns the l_1 norm of the vector.
Definition: vector.cpp:756
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
Definition: vector.cpp:647
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:385
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:618
void Print_HYPRE(std::ostream &out) const
Prints vector to stream out in HYPRE_Vector format.
Definition: vector.cpp:671
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:204
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Definition: vector.hpp:125
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:205
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:190
Host memory; using new[] and delete[].
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:326
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:816
void New(int size)
Allocate host memory for size entries with type MemoryType::HOST.
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.
Vector & operator-=(double c)
Definition: vector.cpp:157
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:42
const double alpha
Definition: ex15.cpp:339
Vector data type.
Definition: vector.hpp:48
Vector & operator+=(const Vector &v)
Definition: vector.cpp:178
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:64
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:833
virtual N_Vector ToNVector()
Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_SERIAL.
Definition: vector.hpp:355
void Neg()
(*this) = -(*this)
Definition: vector.cpp:230