MFEM  v3.3.2
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 
16 #if defined(MFEM_USE_SUNDIALS) && defined(MFEM_USE_MPI)
17 #include <nvector/nvector_parallel.h>
18 #include <nvector/nvector_parhyp.h>
19 #endif
20 
21 #include <iostream>
22 #include <iomanip>
23 #include <cmath>
24 #include <cstdlib>
25 #include <ctime>
26 #include <limits>
27 
28 namespace mfem
29 {
30 
32 {
33  int s = v.Size();
34  if (s > 0)
35  {
36  allocsize = size = s;
37  data = new double[s];
38  for (int i = 0; i < s; i++)
39  {
40  data[i] = v(i);
41  }
42  }
43  else
44  {
45  allocsize = size = 0;
46  data = NULL;
47  }
48 }
49 
50 void Vector::Load(std::istream **in, int np, int *dim)
51 {
52  int i, j, s;
53 
54  s = 0;
55  for (i = 0; i < np; i++)
56  {
57  s += dim[i];
58  }
59 
60  SetSize(s);
61 
62  int p = 0;
63  for (i = 0; i < np; i++)
64  for (j = 0; j < dim[i]; j++)
65  {
66  *in[i] >> data[p++];
67  }
68 }
69 
70 void Vector::Load(std::istream &in, int Size)
71 {
72  SetSize(Size);
73 
74  for (int i = 0; i < size; i++)
75  {
76  in >> data[i];
77  }
78 }
79 
80 double &Vector::Elem(int i)
81 {
82  return operator()(i);
83 }
84 
85 const double &Vector::Elem(int i) const
86 {
87  return operator()(i);
88 }
89 
90 double Vector::operator*(const double *v) const
91 {
92  int s = size;
93  const double *d = data;
94  double prod = 0.0;
95 #ifdef MFEM_USE_OPENMP
96  #pragma omp parallel for reduction(+:prod)
97 #endif
98  for (int i = 0; i < s; i++)
99  {
100  prod += d[i] * v[i];
101  }
102  return prod;
103 }
104 
105 double Vector::operator*(const Vector &v) const
106 {
107 #ifdef MFEM_DEBUG
108  if (v.size != size)
109  {
110  mfem_error("Vector::operator*(const Vector &) const");
111  }
112 #endif
113 
114  return operator*(v.data);
115 }
116 
117 Vector &Vector::operator=(const double *v)
118 {
119  for (int i = 0; i < size; i++)
120  {
121  data[i] = v[i];
122  }
123  return *this;
124 }
125 
127 {
128  SetSize(v.Size());
129  for (int i = 0; i < size; i++)
130  {
131  data[i] = v.data[i];
132  }
133  return *this;
134 }
135 
136 Vector &Vector::operator=(double value)
137 {
138  int i, s = size;
139  double *p = data, v = value;
140  for (i = 0; i < s; i++)
141  {
142  *(p++) = v;
143  }
144  return *this;
145 }
146 
148 {
149  for (int i = 0; i < size; i++)
150  {
151  data[i] *= c;
152  }
153  return *this;
154 }
155 
157 {
158  double m = 1.0/c;
159  for (int i = 0; i < size; i++)
160  {
161  data[i] *= m;
162  }
163  return *this;
164 }
165 
167 {
168  for (int i = 0; i < size; i++)
169  {
170  data[i] -= c;
171  }
172  return *this;
173 }
174 
176 {
177 #ifdef MFEM_DEBUG
178  if (size != v.size)
179  {
180  mfem_error("Vector::operator-=(const Vector &)");
181  }
182 #endif
183  for (int i = 0; i < size; i++)
184  {
185  data[i] -= v(i);
186  }
187  return *this;
188 }
189 
191 {
192 #ifdef MFEM_DEBUG
193  if (size != v.size)
194  {
195  mfem_error("Vector::operator+=(const Vector &)");
196  }
197 #endif
198  for (int i = 0; i < size; i++)
199  {
200  data[i] += v(i);
201  }
202  return *this;
203 }
204 
205 Vector &Vector::Add(const double a, const Vector &Va)
206 {
207 #ifdef MFEM_DEBUG
208  if (size != Va.size)
209  {
210  mfem_error("Vector::Add(const double, const Vector &)");
211  }
212 #endif
213  if (a != 0.0)
214  {
215  for (int i = 0; i < size; i++)
216  {
217  data[i] += a * Va(i);
218  }
219  }
220  return *this;
221 }
222 
223 Vector &Vector::Set(const double a, const Vector &Va)
224 {
225 #ifdef MFEM_DEBUG
226  if (size != Va.size)
227  {
228  mfem_error("Vector::Set(const double, const Vector &)");
229  }
230 #endif
231  for (int i = 0; i < size; i++)
232  {
233  data[i] = a * Va(i);
234  }
235  return *this;
236 }
237 
238 void Vector::SetVector(const Vector &v, int offset)
239 {
240  int vs = v.Size();
241  double *vp = v.data, *p = data + offset;
242 
243 #ifdef MFEM_DEBUG
244  if (offset+vs > size)
245  {
246  mfem_error("Vector::SetVector(const Vector &, int)");
247  }
248 #endif
249 
250  for (int i = 0; i < vs; i++)
251  {
252  p[i] = vp[i];
253  }
254 }
255 
257 {
258  for (int i = 0; i < size; i++)
259  {
260  data[i] = -data[i];
261  }
262 }
263 
264 void add(const Vector &v1, const Vector &v2, Vector &v)
265 {
266 #ifdef MFEM_DEBUG
267  if (v.size != v1.size || v.size != v2.size)
268  {
269  mfem_error("add(Vector &v1, Vector &v2, Vector &v)");
270  }
271 #endif
272 
273 #ifdef MFEM_USE_OPENMP
274  #pragma omp parallel for
275 #endif
276  for (int i = 0; i < v.size; i++)
277  {
278  v.data[i] = v1.data[i] + v2.data[i];
279  }
280 }
281 
282 void add(const Vector &v1, double alpha, const Vector &v2, Vector &v)
283 {
284 #ifdef MFEM_DEBUG
285  if (v.size != v1.size || v.size != v2.size)
286  {
287  mfem_error ("add(Vector &v1, double alpha, Vector &v2, Vector &v)");
288  }
289 #endif
290  if (alpha == 0.0)
291  {
292  v = v1;
293  }
294  else if (alpha == 1.0)
295  {
296  add(v1, v2, v);
297  }
298  else
299  {
300  const double *v1p = v1.data, *v2p = v2.data;
301  double *vp = v.data;
302  int s = v.size;
303 #ifdef MFEM_USE_OPENMP
304  #pragma omp parallel for
305 #endif
306  for (int i = 0; i < s; i++)
307  {
308  vp[i] = v1p[i] + alpha*v2p[i];
309  }
310  }
311 }
312 
313 void add(const double a, const Vector &x, const Vector &y, Vector &z)
314 {
315 #ifdef MFEM_DEBUG
316  if (x.size != y.size || x.size != z.size)
317  mfem_error ("add(const double a, const Vector &x, const Vector &y,"
318  " Vector &z)");
319 #endif
320  if (a == 0.0)
321  {
322  z = 0.0;
323  }
324  else if (a == 1.0)
325  {
326  add(x, y, z);
327  }
328  else
329  {
330  const double *xp = x.data;
331  const double *yp = y.data;
332  double *zp = z.data;
333  int s = x.size;
334 
335 #ifdef MFEM_USE_OPENMP
336  #pragma omp parallel for
337 #endif
338  for (int i = 0; i < s; i++)
339  {
340  zp[i] = a * (xp[i] + yp[i]);
341  }
342  }
343 }
344 
345 void add(const double a, const Vector &x,
346  const double b, const Vector &y, Vector &z)
347 {
348 #ifdef MFEM_DEBUG
349  if (x.size != y.size || x.size != z.size)
350  mfem_error("add(const double a, const Vector &x,\n"
351  " const double b, const Vector &y, Vector &z)");
352 #endif
353  if (a == 0.0)
354  {
355  z.Set(b, y);
356  }
357  else if (b == 0.0)
358  {
359  z.Set(a, x);
360  }
361  else if (a == 1.0)
362  {
363  add(x, b, y, z);
364  }
365  else if (b == 1.0)
366  {
367  add(y, a, x, z);
368  }
369  else if (a == b)
370  {
371  add(a, x, y, z);
372  }
373  else
374  {
375  const double *xp = x.data;
376  const double *yp = y.data;
377  double *zp = z.data;
378  int s = x.size;
379 
380 #ifdef MFEM_USE_OPENMP
381  #pragma omp parallel for
382 #endif
383  for (int i = 0; i < s; i++)
384  {
385  zp[i] = a * xp[i] + b * yp[i];
386  }
387  }
388 }
389 
390 void subtract(const Vector &x, const Vector &y, Vector &z)
391 {
392 #ifdef MFEM_DEBUG
393  if (x.size != y.size || x.size != z.size)
394  {
395  mfem_error ("subtract(const Vector &, const Vector &, Vector &)");
396  }
397 #endif
398  const double *xp = x.data;
399  const double *yp = y.data;
400  double *zp = z.data;
401  int s = x.size;
402 
403 #ifdef MFEM_USE_OPENMP
404  #pragma omp parallel for
405 #endif
406  for (int i = 0; i < s; i++)
407  {
408  zp[i] = xp[i] - yp[i];
409  }
410 }
411 
412 void subtract(const double a, const Vector &x, const Vector &y, Vector &z)
413 {
414 #ifdef MFEM_DEBUG
415  if (x.size != y.size || x.size != z.size)
416  mfem_error("subtract(const double a, const Vector &x,"
417  " const Vector &y, Vector &z)");
418 #endif
419 
420  if (a == 0.)
421  {
422  z = 0.;
423  }
424  else if (a == 1.)
425  {
426  subtract(x, y, z);
427  }
428  else
429  {
430  const double *xp = x.data;
431  const double *yp = y.data;
432  double *zp = z.data;
433  int s = x.size;
434 
435 #ifdef MFEM_USE_OPENMP
436  #pragma omp parallel for
437 #endif
438  for (int i = 0; i < s; i++)
439  {
440  zp[i] = a * (xp[i] - yp[i]);
441  }
442  }
443 }
444 
445 void Vector::median(const Vector &lo, const Vector &hi)
446 {
447  double *v = data;
448 
449  for (int i = 0; i < size; i++)
450  {
451  if (v[i] < lo[i])
452  {
453  v[i] = lo[i];
454  }
455  else if (v[i] > hi[i])
456  {
457  v[i] = hi[i];
458  }
459  }
460 }
461 
462 void Vector::GetSubVector(const Array<int> &dofs, Vector &elemvect) const
463 {
464  int i, j, n = dofs.Size();
465 
466  elemvect.SetSize (n);
467 
468  for (i = 0; i < n; i++)
469  {
470  if ((j=dofs[i]) >= 0)
471  {
472  elemvect(i) = data[j];
473  }
474  else
475  {
476  elemvect(i) = -data[-1-j];
477  }
478  }
479 }
480 
481 void Vector::GetSubVector(const Array<int> &dofs, double *elem_data) const
482 {
483  int i, j, n = dofs.Size();
484 
485  for (i = 0; i < n; i++)
486  {
487  if ((j=dofs[i]) >= 0)
488  {
489  elem_data[i] = data[j];
490  }
491  else
492  {
493  elem_data[i] = -data[-1-j];
494  }
495  }
496 }
497 
498 void Vector::SetSubVector(const Array<int> &dofs, const double value)
499 {
500  const int n = dofs.Size();
501 
502  for (int i = 0; i < n; i++)
503  {
504  const int j = dofs[i];
505  if (j >= 0)
506  {
507  data[j] = value;
508  }
509  else
510  {
511  data[-1-j] = -value;
512  }
513  }
514 }
515 
516 void Vector::SetSubVector(const Array<int> &dofs, const Vector &elemvect)
517 {
518  int i, j, n = dofs.Size();
519 
520  for (i = 0; i < n; i++)
521  {
522  if ((j=dofs[i]) >= 0)
523  {
524  data[j] = elemvect(i);
525  }
526  else
527  {
528  data[-1-j] = -elemvect(i);
529  }
530  }
531 }
532 
533 void Vector::SetSubVector(const Array<int> &dofs, double *elem_data)
534 {
535  int i, j, n = dofs.Size();
536 
537  for (i = 0; i < n; i++)
538  {
539  if ((j=dofs[i]) >= 0)
540  {
541  data[j] = elem_data[i];
542  }
543  else
544  {
545  data[-1-j] = -elem_data[i];
546  }
547  }
548 }
549 
550 void Vector::AddElementVector(const Array<int> &dofs, const Vector &elemvect)
551 {
552  MFEM_ASSERT(dofs.Size() == elemvect.Size(), "");
553  int i, j, n = dofs.Size();
554 
555  for (i = 0; i < n; i++)
556  {
557  if ((j=dofs[i]) >= 0)
558  {
559  data[j] += elemvect(i);
560  }
561  else
562  {
563  data[-1-j] -= elemvect(i);
564  }
565  }
566 }
567 
568 void Vector::AddElementVector(const Array<int> &dofs, double *elem_data)
569 {
570  int i, j, n = dofs.Size();
571 
572  for (i = 0; i < n; i++)
573  {
574  if ((j = dofs[i]) >= 0)
575  {
576  data[j] += elem_data[i];
577  }
578  else
579  {
580  data[-1-j] -= elem_data[i];
581  }
582  }
583 }
584 
585 void Vector::AddElementVector(const Array<int> &dofs, const double a,
586  const Vector &elemvect)
587 {
588  int i, j, n = dofs.Size();
589 
590  for (i = 0; i < n; i++)
591  if ((j=dofs[i]) >= 0)
592  {
593  data[j] += a * elemvect(i);
594  }
595  else
596  {
597  data[-1-j] -= a * elemvect(i);
598  }
599 }
600 
601 void Vector::SetSubVectorComplement(const Array<int> &dofs, const double val)
602 {
603  Vector dofs_vals;
604  GetSubVector(dofs, dofs_vals);
605  operator=(val);
606  SetSubVector(dofs, dofs_vals);
607 }
608 
609 void Vector::Print(std::ostream &out, int width) const
610 {
611  if (!size) { return; }
612 
613  for (int i = 0; 1; )
614  {
615  out << data[i];
616  i++;
617  if (i == size)
618  {
619  break;
620  }
621  if ( i % width == 0 )
622  {
623  out << '\n';
624  }
625  else
626  {
627  out << ' ';
628  }
629  }
630  out << '\n';
631 }
632 
633 void Vector::Print_HYPRE(std::ostream &out) const
634 {
635  int i;
636  std::ios::fmtflags old_fmt = out.flags();
637  out.setf(std::ios::scientific);
638  std::streamsize old_prec = out.precision(14);
639 
640  out << size << '\n'; // number of rows
641 
642  for (i = 0; i < size; i++)
643  {
644  out << data[i] << '\n';
645  }
646 
647  out.precision(old_prec);
648  out.flags(old_fmt);
649 }
650 
651 void Vector::Randomize(int seed)
652 {
653  // static unsigned int seed = time(0);
654  const double max = (double)(RAND_MAX) + 1.;
655 
656  if (seed == 0)
657  {
658  seed = (int)time(0);
659  }
660 
661  // srand(seed++);
662  srand((unsigned)seed);
663 
664  for (int i = 0; i < size; i++)
665  {
666  data[i] = std::abs(rand()/max);
667  }
668 }
669 
670 double Vector::Norml2() const
671 {
672  // Scale entries of Vector on the fly, using algorithms from
673  // std::hypot() and LAPACK's drm2. This scaling ensures that the
674  // argument of each call to std::pow is <= 1 to avoid overflow.
675  if (0 == size)
676  {
677  return 0.0;
678  } // end if 0 == size
679 
680  if (1 == size)
681  {
682  return std::abs(data[0]);
683  } // end if 1 == size
684 
685  double scale = 0.0;
686  double sum = 0.0;
687 
688  for (int i = 0; i < size; i++)
689  {
690  if (data[i] != 0.0)
691  {
692  const double absdata = std::abs(data[i]);
693  if (scale <= absdata)
694  {
695  const double sqr_arg = scale / absdata;
696  sum = 1.0 + sum * (sqr_arg * sqr_arg);
697  scale = absdata;
698  continue;
699  } // end if scale <= absdata
700  const double sqr_arg = absdata / scale;
701  sum += (sqr_arg * sqr_arg); // else scale > absdata
702  } // end if data[i] != 0
703  }
704  return scale * std::sqrt(sum);
705 }
706 
707 double Vector::Normlinf() const
708 {
709  double max = 0.0;
710  for (int i = 0; i < size; i++)
711  {
712  max = std::max(std::abs(data[i]), max);
713  }
714  return max;
715 }
716 
717 double Vector::Norml1() const
718 {
719  double sum = 0.0;
720  for (int i = 0; i < size; i++)
721  {
722  sum += std::abs(data[i]);
723  }
724  return sum;
725 }
726 
727 double Vector::Normlp(double p) const
728 {
729  MFEM_ASSERT(p > 0.0, "Vector::Normlp");
730  if (p == 1.0)
731  {
732  return Norml1();
733  }
734  if (p == 2.0)
735  {
736  return Norml2();
737  }
738  if (p < std::numeric_limits<double>::infinity())
739  {
740  // Scale entries of Vector on the fly, using algorithms from
741  // std::hypot() and LAPACK's drm2. This scaling ensures that the
742  // argument of each call to std::pow is <= 1 to avoid overflow.
743  if (0 == size)
744  {
745  return 0.0;
746  } // end if 0 == size
747 
748  if (1 == size)
749  {
750  return std::abs(data[0]);
751  } // end if 1 == size
752 
753  double scale = 0.0;
754  double sum = 0.0;
755 
756  for (int i = 0; i < size; i++)
757  {
758  if (data[i] != 0.0)
759  {
760  const double absdata = std::abs(data[i]);
761  if (scale <= absdata)
762  {
763  sum = 1.0 + sum * std::pow(scale / absdata, p);
764  scale = absdata;
765  continue;
766  } // end if scale <= absdata
767  sum += std::pow(absdata / scale, p); // else scale > absdata
768  } // end if data[i] != 0
769  }
770  return scale * std::pow(sum, 1.0/p);
771  } // end if p < std::numeric_limits<double>::infinity()
772 
773  return Normlinf(); // else p >= std::numeric_limits<double>::infinity()
774 }
775 
776 double Vector::Max() const
777 {
778  double max = data[0];
779 
780  for (int i = 1; i < size; i++)
781  if (data[i] > max)
782  {
783  max = data[i];
784  }
785 
786  return max;
787 }
788 
789 double Vector::Min() const
790 {
791  double min = data[0];
792 
793  for (int i = 1; i < size; i++)
794  if (data[i] < min)
795  {
796  min = data[i];
797  }
798 
799  return min;
800 }
801 
802 double Vector::Sum() const
803 {
804  double sum = 0.0;
805 
806  for (int i = 0; i < size; i++)
807  {
808  sum += data[i];
809  }
810 
811  return sum;
812 }
813 
814 #ifdef MFEM_USE_SUNDIALS
815 
816 #ifndef SUNTRUE
817 #define SUNTRUE TRUE
818 #endif
819 #ifndef SUNFALSE
820 #define SUNFALSE FALSE
821 #endif
822 
823 Vector::Vector(N_Vector nv)
824 {
825  N_Vector_ID nvid = N_VGetVectorID(nv);
826  switch (nvid)
827  {
828  case SUNDIALS_NVEC_SERIAL:
829  SetDataAndSize(NV_DATA_S(nv), NV_LENGTH_S(nv));
830  break;
831 #ifdef MFEM_USE_MPI
832  case SUNDIALS_NVEC_PARALLEL:
833  SetDataAndSize(NV_DATA_P(nv), NV_LOCLENGTH_P(nv));
834  break;
835  case SUNDIALS_NVEC_PARHYP:
836  {
837  hypre_Vector *hpv_local = N_VGetVector_ParHyp(nv)->local_vector;
838  SetDataAndSize(hpv_local->data, hpv_local->size);
839  break;
840  }
841 #endif
842  default:
843  MFEM_ABORT("N_Vector type " << nvid << " is not supported");
844  }
845 }
846 
847 void Vector::ToNVector(N_Vector &nv)
848 {
849  MFEM_ASSERT(nv, "N_Vector handle is NULL");
850  N_Vector_ID nvid = N_VGetVectorID(nv);
851  switch (nvid)
852  {
853  case SUNDIALS_NVEC_SERIAL:
854  MFEM_ASSERT(NV_OWN_DATA_S(nv) == SUNFALSE, "invalid serial N_Vector");
855  NV_DATA_S(nv) = data;
856  NV_LENGTH_S(nv) = size;
857  break;
858 #ifdef MFEM_USE_MPI
859  case SUNDIALS_NVEC_PARALLEL:
860  MFEM_ASSERT(NV_OWN_DATA_P(nv) == SUNFALSE, "invalid parallel N_Vector");
861  NV_DATA_P(nv) = data;
862  NV_LOCLENGTH_P(nv) = size;
863  break;
864  case SUNDIALS_NVEC_PARHYP:
865  {
866  hypre_Vector *hpv_local = N_VGetVector_ParHyp(nv)->local_vector;
867  MFEM_ASSERT(hpv_local->owns_data == false, "invalid hypre N_Vector");
868  hpv_local->data = data;
869  hpv_local->size = size;
870  break;
871  }
872 #endif
873  default:
874  MFEM_ABORT("N_Vector type " << nvid << " is not supported");
875  }
876 }
877 
878 #endif // MFEM_USE_SUNDIALS
879 
880 }
int Size() const
Logical size of the array.
Definition: array.hpp:110
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:238
Vector()
Default constructor for Vector. Sets size = 0 and data = NULL.
Definition: vector.hpp:51
double & Elem(int i)
Access Vector entries. Index i = 0 .. size-1.
Definition: vector.cpp:80
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:320
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:670
double & operator()(int i)
Access Vector entries using () for 0-based indexing.
Definition: vector.hpp:349
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Definition: vector.cpp:462
int Size() const
Returns the size of the vector.
Definition: vector.hpp:113
double Normlinf() const
Returns the l_infinity norm of the vector.
Definition: vector.cpp:707
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:651
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:264
double operator*(const double *) const
Dot product with a double * array.
Definition: vector.cpp:90
Vector & operator=(const double *v)
Copy Size() entries from v.
Definition: vector.cpp:117
void Load(std::istream **in, int np, int *dim)
Reads a vector from multiple files.
Definition: vector.cpp:50
double Normlp(double p) const
Returns the l_p norm of the vector.
Definition: vector.cpp:727
int dim
Definition: ex3.cpp:47
void median(const Vector &lo, const Vector &hi)
v = median(v,lo,hi) entrywise. Implementation assumes lo &lt;= hi.
Definition: vector.cpp:445
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add (element) subvector to the vector.
Definition: vector.cpp:550
Vector & operator/=(double c)
Definition: vector.cpp:156
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:601
Vector & operator*=(double c)
Definition: vector.cpp:147
double Min() const
Returns the minimal element of the vector.
Definition: vector.cpp:789
double Norml1() const
Returns the l_1 norm of the vector.
Definition: vector.cpp:717
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
Definition: vector.cpp:609
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:390
void mfem_error(const char *msg)
Definition: error.cpp:107
void Print_HYPRE(std::ostream &out) const
Prints vector to stream out in HYPRE_Vector format.
Definition: vector.cpp:633
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Definition: vector.hpp:94
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:223
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:205
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:776
Vector & operator-=(double c)
Definition: vector.cpp:166
const double alpha
Definition: ex15.cpp:337
Vector data type.
Definition: vector.hpp:41
int allocsize
Definition: vector.hpp:45
Vector & operator+=(const Vector &v)
Definition: vector.cpp:190
double * data
Definition: vector.hpp:46
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
double Sum() const
Return the sum of the vector entries.
Definition: vector.cpp:802
virtual N_Vector ToNVector()
Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_SERIAL.
Definition: vector.hpp:275
void Neg()
(*this) = -(*this)
Definition: vector.cpp:256