MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
vector.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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 "../general/forall.hpp"
16#include "../general/hash.hpp"
17#include "../general/scan.hpp"
18#include "vector.hpp"
19
20#ifdef MFEM_USE_OPENMP
21#include <omp.h>
22#endif
23
24#include <iostream>
25#include <cmath>
26#include <ctime>
27
28namespace mfem
29{
30
31/**
32 * Reducer for helping to compute L2-norms. Given two partial results:
33 * a0 = sum_i (|v_i|/a1)^2
34 * b0 = sum_j (|v_j|/b1)^2 (j disjoint from i for vector v)
35 * computes:
36 * a1 = max(a1, b1)
37 * a0 = (a1 == 0 ? 0 : sum_{k in union(i,j)} (|v_k|/a1)^2)
38 *
39 * This form is resiliant against overflow/underflow, similar to std::hypot
40 */
41struct L2Reducer
42{
43 using value_type = DevicePair<real_t, real_t>;
44 static MFEM_HOST_DEVICE void Join(value_type& a, const value_type &b)
45 {
46 real_t scale = fmax(a.second, b.second);
47 if (scale > 0)
48 {
49 real_t s = a.second / scale;
50 a.first *= s * s;
51 s = b.second / scale;
52 a.first += b.first * s * s;
53 a.second = scale;
54 }
55 }
56
57 static MFEM_HOST_DEVICE void SetInitialValue(value_type &a)
58 {
59 a.first = 0;
60 a.second = 0;
61 }
62};
63
64/**
65 * Reducer for helping to compute Lp-norms. Given two partial results:
66 * a0 = sum_i (|v_i|/a1)^p
67 * b0 = sum_j (|v_j|/b1)^p (j disjoint from i for vector v)
68 * computes:
69 * a1 = max(a1, b1)
70 * a0 = (a1 == 0 ? 0 : sum_{k in union(i,j)} (|v_k|/a1)^p)
71 *
72 * This form is resiliant against overflow/underflow, similar to std::hypot
73 */
74struct LpReducer
75{
76 real_t p;
77 using value_type = DevicePair<real_t, real_t>;
78 MFEM_HOST_DEVICE void Join(value_type& a, const value_type &b) const
79 {
80 real_t scale = fmax(a.second, b.second);
81 if (scale > 0)
82 {
83 a.first = a.first * pow(a.second / scale, p) +
84 b.first * pow(b.second / scale, p);
85 a.second = scale;
86 }
87 }
88
89 static MFEM_HOST_DEVICE void SetInitialValue(value_type &a)
90 {
91 a.first = 0;
92 a.second = 0;
93 }
94};
95
96static Array<real_t>& vector_workspace()
97{
98 static Array<real_t> instance;
99 return instance;
100}
101
102static Array<DevicePair<real_t, real_t>> &Lpvector_workspace()
103{
104 static Array<DevicePair<real_t, real_t>> instance;
105 return instance;
106}
107
109{
110 const int s = v.Size();
111 size = s;
112 if (s > 0)
113 {
114 MFEM_ASSERT(!v.data.Empty(), "invalid source vector");
115 data.New(s, v.data.GetMemoryType());
116 data.CopyFrom(v.data, s);
117 }
118 UseDevice(v.UseDevice());
119}
120
122 : data(std::move(v.data)), size(v.size)
123{
124 v.size = 0;
125}
126
127void Vector::Load(std::istream **in, int np, int *dim)
128{
129 int i, j, s;
130
131 s = 0;
132 for (i = 0; i < np; i++)
133 {
134 s += dim[i];
135 }
136
137 SetSize(s);
138 HostWrite();
139
140 int p = 0;
141 for (i = 0; i < np; i++)
142 {
143 for (j = 0; j < dim[i]; j++)
144 {
145 *in[i] >> data[p++];
146 // Clang's libc++ sets the failbit when (correctly) parsing subnormals,
147 // so we reset the failbit here.
148 if (!*in[i] && errno == ERANGE)
149 {
150 in[i]->clear();
151 }
152 }
153 }
154}
155
156void Vector::Load(std::istream &in, int Size)
157{
158 SetSize(Size);
159 HostWrite();
160
161 for (int i = 0; i < size; i++)
162 {
163 in >> data[i];
164 // Clang's libc++ sets the failbit when (correctly) parsing subnormals,
165 // so we reset the failbit here.
166 if (!in && errno == ERANGE)
167 {
168 in.clear();
169 }
170 }
171}
172
174{
175 return operator()(i);
176}
177
178const real_t &Vector::Elem(int i) const
179{
180 return operator()(i);
181}
182
184{
185 HostRead();
186 real_t dot = 0.0;
187#ifdef MFEM_USE_LEGACY_OPENMP
188 #pragma omp parallel for reduction(+:dot)
189#endif
190 for (int i = 0; i < size; i++)
191 {
192 dot += data[i] * v[i];
193 }
194 return dot;
195}
196
198{
200 return *this;
201}
202
204{
205#if 0
206 SetSize(v.Size(), v.data.GetMemoryType());
207 data.CopyFrom(v.data, v.Size());
208 UseDevice(v.UseDevice());
209#else
210 SetSize(v.Size());
211 const bool vuse = v.UseDevice();
212 const bool use_dev = UseDevice() || vuse;
213 if (use_dev != vuse) { v.UseDevice(use_dev); }
214 // keep 'data' where it is, unless 'use_dev' is true
215 if (use_dev) { Write(); }
216 data.CopyFrom(v.data, v.Size());
217 if (use_dev != vuse) { v.UseDevice(vuse); }
218#endif
219 return *this;
220}
221
223{
224 if (this != &v)
225 {
226 v.Swap(*this);
227 v.Destroy();
228 }
229 return *this;
230}
231
233{
234 const bool use_dev = UseDevice();
235 const int N = size;
236 auto y = Write(use_dev);
237 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] = value; });
238 return *this;
239}
240
242{
243 const bool use_dev = UseDevice();
244 const int N = size;
245 auto y = ReadWrite(use_dev);
246 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] *= c; });
247 return *this;
248}
249
251{
252 MFEM_ASSERT(size == v.size, "incompatible Vectors!");
253
254 const bool use_dev = UseDevice() || v.UseDevice();
255 const int N = size;
256 const auto x = v.Read(use_dev);
257 auto y = ReadWrite(use_dev);
258 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] *= x[i]; });
259 return *this;
260}
261
263{
264 const bool use_dev = UseDevice();
265 const int N = size;
266 const real_t m = 1.0/c;
267 auto y = ReadWrite(use_dev);
268 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] *= m; });
269 return *this;
270}
271
273{
274 MFEM_ASSERT(size == v.size, "incompatible Vectors!");
275
276 const bool use_dev = UseDevice() || v.UseDevice();
277 const int N = size;
278 const auto x = v.Read(use_dev);
279 auto y = ReadWrite(use_dev);
280 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] /= x[i]; });
281 return *this;
282}
283
285{
286 const bool use_dev = UseDevice();
287 const int N = size;
288 auto y = ReadWrite(use_dev);
289 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] -= c; });
290 return *this;
291}
292
294{
295 MFEM_ASSERT(size == v.size, "incompatible Vectors!");
296
297 const bool use_dev = UseDevice() || v.UseDevice();
298 const int N = size;
299 const auto x = v.Read(use_dev);
300 auto y = ReadWrite(use_dev);
301 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] -= x[i]; });
302 return *this;
303}
304
306{
307 const bool use_dev = UseDevice();
308 const int N = size;
309 auto y = ReadWrite(use_dev);
310 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] += c; });
311 return *this;
312}
313
315{
316 MFEM_ASSERT(size == v.size, "incompatible Vectors!");
317
318 const bool use_dev = UseDevice() || v.UseDevice();
319 const int N = size;
320 const auto x = v.Read(use_dev);
321 auto y = ReadWrite(use_dev);
322 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] += x[i]; });
323 return *this;
324}
325
326Vector &Vector::Add(const real_t a, const Vector &Va)
328 MFEM_ASSERT(size == Va.size, "incompatible Vectors!");
329
330 if (a != 0.0)
331 {
332 const int N = size;
333 const bool use_dev = UseDevice() || Va.UseDevice();
334 const auto x = Va.Read(use_dev);
335 auto y = ReadWrite(use_dev);
336 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] += a * x[i]; });
337 }
338 return *this;
339}
340
341Vector &Vector::Set(const real_t a, const Vector &Va)
342{
343 MFEM_ASSERT(size == Va.size, "incompatible Vectors!");
344
345 const bool use_dev = UseDevice() || Va.UseDevice();
346 const int N = size;
347 const auto x = Va.Read(use_dev);
348 auto y = Write(use_dev);
349 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] = a * x[i]; });
350 return *this;
351}
352
353void Vector::SetVector(const Vector &v, int offset)
354{
355 MFEM_ASSERT(v.Size() + offset <= size, "invalid sub-vector");
356
357 const bool use_dev = UseDevice() || v.UseDevice();
358 const int vs = v.Size();
359 const auto vp = v.Read(use_dev);
360 // Use read+write access for *this - we only modify some of its entries
361 auto p = ReadWrite(use_dev) + offset;
362 mfem::forall_switch(use_dev, vs, [=] MFEM_HOST_DEVICE (int i) { p[i] = vp[i]; });
363}
364
365void Vector::AddSubVector(const Vector &v, int offset)
366{
367 MFEM_ASSERT(v.Size() + offset <= size, "invalid sub-vector");
368
369 const bool use_dev = UseDevice() || v.UseDevice();
370 const int vs = v.Size();
371 const auto vp = v.Read(use_dev);
372 auto p = ReadWrite(use_dev) + offset;
373 mfem::forall_switch(use_dev, vs, [=] MFEM_HOST_DEVICE (int i) { p[i] += vp[i]; });
374}
375
377{
378 const bool use_dev = UseDevice();
379 const int N = size;
380 auto y = ReadWrite(use_dev);
381 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] = -y[i]; });
382}
383
385{
386 const bool use_dev = UseDevice();
387 const int N = size;
388 auto y = ReadWrite(use_dev);
389 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] = 1.0/y[i]; });
390}
391
393{
394 const bool use_dev = UseDevice();
395 const int N = size;
396 auto y = ReadWrite(use_dev);
397 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i)
398 {
399 y[i] = std::abs(y[i]);
400 });
401}
402
404{
405 const bool use_dev = UseDevice();
406 const int N = size;
407 auto y = ReadWrite(use_dev);
408 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i)
409 {
410 y[i] = std::pow(y[i], p);
411 });
412}
413
414void add(const Vector &v1, const Vector &v2, Vector &v)
415{
416 MFEM_ASSERT(v.size == v1.size && v.size == v2.size,
417 "incompatible Vectors!");
418
419#if !defined(MFEM_USE_LEGACY_OPENMP)
420 const bool use_dev = v1.UseDevice() || v2.UseDevice() || v.UseDevice();
421 const int N = v.size;
422 // Note: get read access first, in case v is the same as v1/v2.
423 const auto x1 = v1.Read(use_dev);
424 const auto x2 = v2.Read(use_dev);
425 auto y = v.Write(use_dev);
426 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { y[i] = x1[i] + x2[i]; });
427#else
428 #pragma omp parallel for
429 for (int i = 0; i < v.size; i++)
430 {
431 v.data[i] = v1.data[i] + v2.data[i];
432 }
433#endif
434}
435
436void add(const Vector &v1, real_t alpha, const Vector &v2, Vector &v)
437{
438 MFEM_ASSERT(v.size == v1.size && v.size == v2.size,
439 "incompatible Vectors!");
440
441 if (alpha == 0.0)
442 {
443 v = v1;
444 }
445 else if (alpha == 1.0)
446 {
447 add(v1, v2, v);
448 }
449 else
450 {
451#if !defined(MFEM_USE_LEGACY_OPENMP)
452 const bool use_dev = v1.UseDevice() || v2.UseDevice() || v.UseDevice();
453 const int N = v.size;
454 // Note: get read access first, in case v is the same as v1/v2.
455 const auto d_x = v1.Read(use_dev);
456 const auto d_y = v2.Read(use_dev);
457 auto d_z = v.Write(use_dev);
458 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i)
459 {
460 d_z[i] = d_x[i] + alpha * d_y[i];
461 });
462#else
463 const real_t *v1p = v1.data, *v2p = v2.data;
464 real_t *vp = v.data;
465 const int s = v.size;
466 #pragma omp parallel for
467 for (int i = 0; i < s; i++)
468 {
469 vp[i] = v1p[i] + alpha*v2p[i];
470 }
471#endif
472 }
473}
474
475void add(const real_t a, const Vector &x, const Vector &y, Vector &z)
476{
477 MFEM_ASSERT(x.size == y.size && x.size == z.size,
478 "incompatible Vectors!");
479
480 if (a == 0.0)
481 {
482 z = 0.0;
483 }
484 else if (a == 1.0)
485 {
486 add(x, y, z);
487 }
488 else
489 {
490#if !defined(MFEM_USE_LEGACY_OPENMP)
491 const bool use_dev = x.UseDevice() || y.UseDevice() || z.UseDevice();
492 const int N = x.size;
493 // Note: get read access first, in case z is the same as x/y.
494 const auto xd = x.Read(use_dev);
495 const auto yd = y.Read(use_dev);
496 auto zd = z.Write(use_dev);
497 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i)
498 {
499 zd[i] = a * (xd[i] + yd[i]);
500 });
501#else
502 const real_t *xp = x.data;
503 const real_t *yp = y.data;
504 real_t *zp = z.data;
505 const int s = x.size;
506 #pragma omp parallel for
507 for (int i = 0; i < s; i++)
508 {
509 zp[i] = a * (xp[i] + yp[i]);
510 }
511#endif
512 }
513}
514
515void add(const real_t a, const Vector &x,
516 const real_t b, const Vector &y, Vector &z)
517{
518 MFEM_ASSERT(x.size == y.size && x.size == z.size,
519 "incompatible Vectors!");
520
521 if (a == 0.0)
522 {
523 z.Set(b, y);
524 }
525 else if (b == 0.0)
526 {
527 z.Set(a, x);
528 }
529#if 0
530 else if (a == 1.0)
531 {
532 add(x, b, y, z);
533 }
534 else if (b == 1.0)
535 {
536 add(y, a, x, z);
537 }
538 else if (a == b)
539 {
540 add(a, x, y, z);
541 }
542#endif
543 else
544 {
545#if !defined(MFEM_USE_LEGACY_OPENMP)
546 const bool use_dev = x.UseDevice() || y.UseDevice() || z.UseDevice();
547 const int N = x.size;
548 // Note: get read access first, in case z is the same as x/y.
549 const auto xd = x.Read(use_dev);
550 const auto yd = y.Read(use_dev);
551 auto zd = z.Write(use_dev);
552 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i)
553 {
554 zd[i] = a * xd[i] + b * yd[i];
555 });
556#else
557 const real_t *xp = x.data;
558 const real_t *yp = y.data;
559 real_t *zp = z.data;
560 const int s = x.size;
561 #pragma omp parallel for
562 for (int i = 0; i < s; i++)
563 {
564 zp[i] = a * xp[i] + b * yp[i];
565 }
566#endif
567 }
568}
569
570void subtract(const Vector &x, const Vector &y, Vector &z)
571{
572 MFEM_ASSERT(x.size == y.size && x.size == z.size,
573 "incompatible Vectors!");
574
575#if !defined(MFEM_USE_LEGACY_OPENMP)
576 const bool use_dev = x.UseDevice() || y.UseDevice() || z.UseDevice();
577 const int N = x.size;
578 // Note: get read access first, in case z is the same as x/y.
579 const auto xd = x.Read(use_dev);
580 const auto yd = y.Read(use_dev);
581 auto zd = z.Write(use_dev);
582 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i)
583 {
584 zd[i] = xd[i] - yd[i];
585 });
586#else
587 const real_t *xp = x.data;
588 const real_t *yp = y.data;
589 real_t *zp = z.data;
590 const int s = x.size;
591 #pragma omp parallel for
592 for (int i = 0; i < s; i++)
593 {
594 zp[i] = xp[i] - yp[i];
595 }
596#endif
597}
598
599void subtract(const real_t a, const Vector &x, const Vector &y, Vector &z)
600{
601 MFEM_ASSERT(x.size == y.size && x.size == z.size,
602 "incompatible Vectors!");
603
604 if (a == 0.)
605 {
606 z = 0.;
607 }
608 else if (a == 1.)
609 {
610 subtract(x, y, z);
611 }
612 else
613 {
614#if !defined(MFEM_USE_LEGACY_OPENMP)
615 const bool use_dev = x.UseDevice() || y.UseDevice() || z.UseDevice();
616 const int N = x.size;
617 // Note: get read access first, in case z is the same as x/y.
618 const auto xd = x.Read(use_dev);
619 const auto yd = y.Read(use_dev);
620 auto zd = z.Write(use_dev);
621 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i)
622 {
623 zd[i] = a * (xd[i] - yd[i]);
624 });
625#else
626 const real_t *xp = x.data;
627 const real_t *yp = y.data;
628 real_t *zp = z.data;
629 const int s = x.size;
630 #pragma omp parallel for
631 for (int i = 0; i < s; i++)
632 {
633 zp[i] = a * (xp[i] - yp[i]);
634 }
635#endif
636 }
637}
638
639void Vector::cross3D(const Vector &vin, Vector &vout) const
640{
641 HostRead();
642 vin.HostRead();
643 vout.HostWrite();
644 MFEM_VERIFY(size == 3, "Only 3D vectors supported in cross.");
645 MFEM_VERIFY(vin.Size() == 3, "Only 3D vectors supported in cross.");
646 vout.SetSize(3);
647 vout(0) = data[1]*vin(2)-data[2]*vin(1);
648 vout(1) = data[2]*vin(0)-data[0]*vin(2);
649 vout(2) = data[0]*vin(1)-data[1]*vin(0);
650}
651
652void Vector::median(const Vector &lo, const Vector &hi)
653{
654 MFEM_ASSERT(size == lo.size && size == hi.size,
655 "incompatible Vectors!");
656
657 const bool use_dev = UseDevice() || lo.UseDevice() || hi.UseDevice();
658 const int N = size;
659 // Note: get read access first, in case *this is the same as lo/hi.
660 const auto l = lo.Read(use_dev);
661 const auto h = hi.Read(use_dev);
662 auto m = Write(use_dev);
663 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i)
664 {
665 if (m[i] < l[i])
666 {
667 m[i] = l[i];
668 }
669 else if (m[i] > h[i])
670 {
671 m[i] = h[i];
672 }
673 });
674}
675
676void Vector::GetSubVector(const Array<int> &dofs, Vector &elemvect) const
677{
678 const int n = dofs.Size();
679 elemvect.SetSize(n);
680 const bool use_dev = dofs.UseDevice() || elemvect.UseDevice();
681 const auto d_X = Read(use_dev);
682 const auto d_dofs = dofs.Read(use_dev);
683 auto d_y = elemvect.Write(use_dev);
684 mfem::forall_switch(use_dev, n, [=] MFEM_HOST_DEVICE (int i)
685 {
686 const int dof_i = d_dofs[i];
687 d_y[i] = dof_i >= 0 ? d_X[dof_i] : -d_X[-dof_i-1];
688 });
689}
690
691void Vector::GetSubVector(const Array<int> &dofs, real_t *elem_data) const
692{
693 HostRead();
694 const int n = dofs.Size();
695 for (int i = 0; i < n; i++)
696 {
697 const int j = dofs[i];
698 elem_data[i] = (j >= 0) ? data[j] : -data[-1-j];
699 }
700}
701
702void Vector::SetSubVector(const Array<int> &dofs, const real_t value)
703{
704 const bool use_dev = UseDevice() || dofs.UseDevice();
705 const int n = dofs.Size();
706 // Use read+write access for *this - we only modify some of its entries
707 auto d_X = ReadWrite(use_dev);
708 const auto d_dofs = dofs.Read(use_dev);
709 mfem::forall_switch(use_dev, n, [=] MFEM_HOST_DEVICE (int i)
710 {
711 const int j = d_dofs[i];
712 if (j >= 0)
713 {
714 d_X[j] = value;
715 }
716 else
717 {
718 d_X[-1-j] = -value;
719 }
720 });
721}
722
723void Vector::SetSubVectorHost(const Array<int> &dofs, const real_t value)
724{
726 for (int i = 0; i < dofs.Size(); ++i)
727 {
728 const int j = dofs[i];
729 if (j >= 0)
730 {
731 (*this)[j] = value;
732 }
733 else
734 {
735 (*this)[-1-j] = -value;
736 }
737 }
738}
739
740void Vector::SetSubVector(const Array<int> &dofs, const Vector &elemvect)
741{
742 MFEM_ASSERT(dofs.Size() <= elemvect.Size(),
743 "Size mismatch: length of dofs is " << dofs.Size()
744 << ", length of elemvect is " << elemvect.Size());
745
746 const bool use_dev = dofs.UseDevice() || elemvect.UseDevice();
747 const int n = dofs.Size();
748 // Use read+write access for X - we only modify some of its entries
749 auto d_X = ReadWrite(use_dev);
750 const auto d_y = elemvect.Read(use_dev);
751 const auto d_dofs = dofs.Read(use_dev);
752 mfem::forall_switch(use_dev, n, [=] MFEM_HOST_DEVICE (int i)
753 {
754 const int dof_i = d_dofs[i];
755 if (dof_i >= 0)
756 {
757 d_X[dof_i] = d_y[i];
758 }
759 else
760 {
761 d_X[-1-dof_i] = -d_y[i];
762 }
763 });
764}
765
766void Vector::SetSubVector(const Array<int> &dofs, real_t *elem_data)
767{
768 // Use read+write access because we overwrite only part of the data.
770 const int n = dofs.Size();
771 for (int i = 0; i < n; i++)
772 {
773 const int j= dofs[i];
774 if (j >= 0)
775 {
776 operator()(j) = elem_data[i];
777 }
778 else
779 {
780 operator()(-1-j) = -elem_data[i];
781 }
782 }
783}
784
785void Vector::AddElementVector(const Array<int> &dofs, const Vector &elemvect)
786{
787 MFEM_ASSERT(dofs.Size() <= elemvect.Size(), "Size mismatch: "
788 "length of dofs is " << dofs.Size() <<
789 ", length of elemvect is " << elemvect.Size());
790
791 const bool use_dev = dofs.UseDevice() || elemvect.UseDevice();
792 const int n = dofs.Size();
793 const auto d_y = elemvect.Read(use_dev);
794 const auto d_dofs = dofs.Read(use_dev);
795 auto d_X = ReadWrite(use_dev);
796 mfem::forall_switch(use_dev, n, [=] MFEM_HOST_DEVICE (int i)
797 {
798 const int j = d_dofs[i];
799 if (j >= 0)
800 {
801 d_X[j] += d_y[i];
802 }
803 else
804 {
805 d_X[-1-j] -= d_y[i];
806 }
807 });
808}
809
810void Vector::AddElementVector(const Array<int> &dofs, real_t *elem_data)
811{
813 const int n = dofs.Size();
814 for (int i = 0; i < n; i++)
815 {
816 const int j = dofs[i];
817 if (j >= 0)
818 {
819 operator()(j) += elem_data[i];
820 }
821 else
822 {
823 operator()(-1-j) -= elem_data[i];
824 }
825 }
826}
827
829 const Vector &elemvect)
830{
831 MFEM_ASSERT(dofs.Size() <= elemvect.Size(), "Size mismatch: "
832 "length of dofs is " << dofs.Size() <<
833 ", length of elemvect is " << elemvect.Size());
834
835 const bool use_dev = dofs.UseDevice() || elemvect.UseDevice();
836 const int n = dofs.Size();
837 const auto d_x = elemvect.Read(use_dev);
838 const auto d_dofs = dofs.Read(use_dev);
839 auto d_y = ReadWrite(use_dev);
840 mfem::forall_switch(use_dev, n, [=] MFEM_HOST_DEVICE (int i)
841 {
842 const int j = d_dofs[i];
843 if (j >= 0)
844 {
845 d_y[j] += a * d_x[i];
846 }
847 else
848 {
849 d_y[-1-j] -= a * d_x[i];
850 }
851 });
852}
853
855{
856 const bool use_dev = UseDevice() || dofs.UseDevice();
857 const int n = dofs.Size();
858 const int N = size;
859 Vector dofs_vals(n, use_dev ?
862 auto d_data = ReadWrite(use_dev);
863 auto d_dofs_vals = dofs_vals.Write(use_dev);
864 const auto d_dofs = dofs.Read(use_dev);
865 mfem::forall_switch(use_dev, n, [=] MFEM_HOST_DEVICE (int i) { d_dofs_vals[i] = d_data[d_dofs[i]]; });
866 mfem::forall_switch(use_dev, N, [=] MFEM_HOST_DEVICE (int i) { d_data[i] = val; });
867 mfem::forall_switch(use_dev, n, [=] MFEM_HOST_DEVICE (int i) { d_data[d_dofs[i]] = d_dofs_vals[i]; });
868}
869
870void Vector::Print(std::ostream &os, int width) const
871{
872 if (!size) { return; }
873 HostRead();
874 for (int i = 0; 1; )
875 {
876 os << ZeroSubnormal(data[i]);
877 i++;
878 if (i == size)
879 {
880 break;
881 }
882 if ( i % width == 0 )
883 {
884 os << '\n';
885 }
886 else
887 {
888 os << ' ';
889 }
890 }
891 os << '\n';
892}
893
894#ifdef MFEM_USE_ADIOS2
896 const std::string& variable_name) const
897{
898 if (!size) { return; }
899 HostRead();
900 os.engine.Put(variable_name, &data[0] );
901}
902#endif
903
904void Vector::Print_HYPRE(std::ostream &os) const
905{
906 int i;
907 std::ios::fmtflags old_fmt = os.flags();
908 os.setf(std::ios::scientific);
909 std::streamsize old_prec = os.precision(14);
910
911 os << size << '\n'; // number of rows
912
914 for (i = 0; i < size; i++)
915 {
916 os << ZeroSubnormal(data[i]) << '\n';
917 }
918
919 os.precision(old_prec);
920 os.flags(old_fmt);
921}
922
923void Vector::PrintMathematica(std::ostream & os) const
924{
925 std::ios::fmtflags old_fmt = os.flags();
926 os.setf(std::ios::scientific);
927 std::streamsize old_prec = os.precision(14);
928
929 os << "(* Read file into Mathematica using: "
930 << "myVec = Get[\"this_file_name\"] *)\n";
931 os << "{\n";
932
934 for (int i = 0; i < size; i++)
935 {
936 os << "Internal`StringToMReal[\"" << ZeroSubnormal(data[i]) << "\"]";
937 if (i < size - 1) { os << ','; }
938 os << '\n';
939 }
940
941 os << "}\n";
942
943 os.precision(old_prec);
944 os.flags(old_fmt);
945}
946
947void Vector::PrintHash(std::ostream &os) const
948{
949 os << "size: " << size << '\n';
950 HashFunction hf;
952 os << "hash: " << hf.GetHash() << '\n';
953}
954
955void Vector::Randomize(int seed)
956{
957 if (seed == 0) { seed = (int)time(0); }
958
959 srand((unsigned)seed);
960
961 HostWrite();
962 for (int i = 0; i < size; i++)
963 {
964 data[i] = rand_real();
965 }
966}
967
969{
970 // Scale entries of Vector on the fly, using algorithms from
971 // std::hypot() and LAPACK's drm2. This scaling ensures that the
972 // argument of each call to std::pow is <= 1 to avoid overflow.
973 if (size == 0) { return 0.0; }
974
975 const auto m_data = Read(UseDevice());
976 using value_type = DevicePair<real_t, real_t>;
977 value_type res;
978 res.first = 0;
979 res.second = 0;
980 // first compute sum (|m_data|/scale)^2
981 reduce(size, res, [=] MFEM_HOST_DEVICE(int i, value_type &r)
982 {
983 real_t n = fabs(m_data[i]);
984 if (n > 0)
985 {
986 if (r.second <= n)
987 {
988 real_t arg = r.second / n;
989 r.first = r.first * (arg * arg) + 1;
990 r.second = n;
991 }
992 else
993 {
994 real_t arg = n / r.second;
995 r.first += arg * arg;
996 }
997 }
998 },
999 L2Reducer{}, UseDevice(), Lpvector_workspace());
1000 // final answer
1001 return res.second * sqrt(res.first);
1002}
1003
1005{
1006 if (size == 0) { return 0; }
1007
1008 real_t res = 0;
1009 const auto m_data = Read(UseDevice());
1010 reduce(size, res, [=] MFEM_HOST_DEVICE(int i, real_t &r)
1011 {
1012 r = fmax(r, fabs(m_data[i]));
1013 },
1014 MaxReducer<real_t> {}, UseDevice(), vector_workspace());
1015 return res;
1016}
1017
1019{
1020 if (size == 0) { return 0.0; }
1021
1022 real_t res = 0;
1023 const auto m_data = Read(UseDevice());
1024 reduce(size, res, [=] MFEM_HOST_DEVICE(int i, real_t &r)
1025 {
1026 r += fabs(m_data[i]);
1027 },
1028 SumReducer<real_t> {}, UseDevice(), vector_workspace());
1029 return res;
1030}
1031
1033{
1034 MFEM_ASSERT(p > 0.0, "Vector::Normlp");
1035
1036 if (p == 1.0) { return Norml1(); }
1037
1038 if (p == 2.0) { return Norml2(); }
1039
1040 if (p < infinity())
1041 {
1042 // Scale entries of Vector on the fly, using algorithms from
1043 // std::hypot() and LAPACK's drm2. This scaling ensures that the
1044 // argument of each call to std::pow is <= 1 to avoid overflow.
1045 if (size == 0) { return 0.0; }
1046
1047 using value_type = DevicePair<real_t, real_t>;
1048 value_type res;
1049 res.first = 0;
1050 res.second = 0;
1051 const auto m_data = Read(UseDevice());
1052 // first compute sum (|m_data|/scale)^p
1053 reduce(size, res, [=] MFEM_HOST_DEVICE(int i, value_type &r)
1054 {
1055 real_t n = fabs(m_data[i]);
1056 if (n > 0)
1057 {
1058 if (r.second <= n)
1059 {
1060 real_t arg = r.second / n;
1061 r.first = r.first * pow(arg, p) + 1;
1062 r.second = n;
1063 }
1064 else
1065 {
1066 real_t arg = n / r.second;
1067 r.first += pow(arg, p);
1068 }
1069 }
1070 },
1071 LpReducer{p}, UseDevice(), Lpvector_workspace());
1072 // final answer
1073 return res.second * pow(res.first, 1.0 / p);
1074 } // end if p < infinity()
1075
1076 return Normlinf(); // else p >= infinity()
1077}
1078
1080{
1081 MFEM_ASSERT(size == v.size, "incompatible Vectors!");
1082
1083 if (size == 0) { return 0.0; }
1084
1085 const bool use_dev = UseDevice() || v.UseDevice();
1086 const auto m_data = Read(use_dev), v_data = v.Read(use_dev);
1087
1088 // If OCCA is enabled, it handles all selected backends
1089#ifdef MFEM_USE_OCCA
1090 if (use_dev && DeviceCanUseOcca())
1091 {
1092 return occa::linalg::dot<real_t, real_t, real_t>(
1094 }
1095#endif
1096
1097 const auto compute_dot = [&]()
1098 {
1099 real_t res = 0;
1100 reduce(size, res, [=] MFEM_HOST_DEVICE (int i, real_t &r)
1101 {
1102 r += m_data[i] * v_data[i];
1103 },
1104 SumReducer<real_t> {}, use_dev, vector_workspace());
1105 return res;
1106 };
1107
1108 // Device backends have top priority
1109 if (Device::Allows(Backend::DEVICE_MASK)) { return compute_dot(); }
1110
1111 // Special path for OpenMP
1112#ifdef MFEM_USE_OPENMP
1113 if (use_dev && Device::Allows(Backend::OMP_MASK))
1114 {
1115 // By default, use a deterministic way of computing the dot product
1116#define MFEM_USE_OPENMP_DETERMINISTIC_DOT
1117#ifdef MFEM_USE_OPENMP_DETERMINISTIC_DOT
1118 static Vector th_dot;
1119 #pragma omp parallel
1120 {
1121 const int nt = omp_get_num_threads();
1122 #pragma omp master
1123 th_dot.SetSize(nt);
1124 const int tid = omp_get_thread_num();
1125 const int stride = (size + nt - 1) / nt;
1126 const int start = tid * stride;
1127 const int stop = std::min(start + stride, size);
1128 real_t my_dot = 0.0;
1129 for (int i = start; i < stop; i++)
1130 {
1131 my_dot += m_data[i] * v_data[i];
1132 }
1133 #pragma omp barrier
1134 th_dot(tid) = my_dot;
1135 }
1136 return th_dot.Sum();
1137#else
1138 // The standard way of computing the dot product is non-deterministic
1139 real_t prod = 0.0;
1140 #pragma omp parallel for reduction(+ : prod)
1141 for (int i = 0; i < size; i++)
1142 {
1143 prod += m_data[i] * v_data[i];
1144 }
1145 return prod;
1146#endif // MFEM_USE_OPENMP_DETERMINISTIC_DOT
1147 }
1148#endif // MFEM_USE_OPENMP
1149
1150 // All other CPU backends
1151 return compute_dot();
1152}
1153
1155{
1156 if (size == 0) { return infinity(); }
1157
1158 const auto use_dev = UseDevice();
1159 const auto m_data = Read(use_dev);
1160
1161#ifdef MFEM_USE_OCCA
1162 if (use_dev && DeviceCanUseOcca())
1163 {
1164 return occa::linalg::min<real_t,real_t>(OccaMemoryRead(data, size));
1165 }
1166#endif
1167
1168 const auto compute_min = [&]()
1169 {
1170 real_t res = infinity();
1171 reduce(size, res, [=] MFEM_HOST_DEVICE(int i, real_t &r)
1172 {
1173 r = fmin(r, m_data[i]);
1174 },
1175 MinReducer<real_t> {}, use_dev, vector_workspace());
1176 return res;
1177 };
1178
1179 // Device backends have top priority
1180 if (Device::Allows(Backend::DEVICE_MASK)) { return compute_min(); }
1181
1182 // Special path for OpenMP
1183#ifdef MFEM_USE_OPENMP
1184 if (use_dev && Device::Allows(Backend::OMP_MASK))
1185 {
1186 real_t minimum = m_data[0];
1187 #pragma omp parallel for reduction(min:minimum)
1188 for (int i = 0; i < size; i++)
1189 {
1190 minimum = std::min(minimum, m_data[i]);
1191 }
1192 return minimum;
1193 }
1194#endif
1195
1196 // All other CPU backends
1197 return compute_min();
1198}
1199
1201{
1202 if (size == 0) { return -infinity(); }
1203
1204 const auto use_dev = UseDevice();
1205 const auto m_data = Read(use_dev);
1206
1207#ifdef MFEM_USE_OCCA
1208 if (use_dev && DeviceCanUseOcca())
1209 {
1210 return occa::linalg::max<real_t, real_t>(OccaMemoryRead(data, size));
1211 }
1212#endif
1213
1214 const auto compute_max = [&]()
1215 {
1216 real_t res = -infinity();
1217 reduce(size, res, [=] MFEM_HOST_DEVICE(int i, real_t &r)
1218 {
1219 r = fmax(r, m_data[i]);
1220 },
1221 MaxReducer<real_t> {}, use_dev, vector_workspace());
1222 return res;
1223 };
1224
1225 // Device backends have top priority
1226 if (Device::Allows(Backend::DEVICE_MASK)) { return compute_max(); }
1227
1228 // Special path for OpenMP
1229#ifdef MFEM_USE_OPENMP
1230 if (use_dev && Device::Allows(Backend::OMP_MASK))
1231 {
1232 real_t maximum = m_data[0];
1233 #pragma omp parallel for reduction(max : maximum)
1234 for (int i = 0; i < size; i++)
1235 {
1236 maximum = fmax(maximum, m_data[i]);
1237 }
1238 return maximum;
1239 }
1240#endif
1241
1242 // All other CPU backends
1243 return compute_max();
1244}
1245
1247{
1248 if (size == 0) { return 0.0; }
1249
1250 real_t res = 0;
1251 const auto m_data = Read(UseDevice());
1252 reduce(size, res, [=] MFEM_HOST_DEVICE(int i, real_t &r)
1253 {
1254 r += m_data[i];
1255 },
1256 SumReducer<real_t> {}, UseDevice(), vector_workspace());
1257 return res;
1258}
1259
1260void Vector::DeleteAt(const Array<int> &indices)
1261{
1262 if (indices.Size())
1263 {
1264 const bool use_dev = UseDevice();
1265
1266 // extra entry for number of selected out
1267 Array<int> workspace(size + 1);
1268 const auto d_flag = workspace.Write(use_dev);
1269 mfem::forall_switch(use_dev, size,
1270 [=] MFEM_HOST_DEVICE(int i) { d_flag[i] = true; });
1271 const auto d_indices = indices.Read(use_dev);
1272 mfem::forall_switch(use_dev, indices.Size(), [=] MFEM_HOST_DEVICE(int i)
1273 {
1274 // fine as long as indices are unique; to support non-unique indices
1275 // assignment to d_flag must be atomic
1276 d_flag[d_indices[i]] = false;
1277 });
1278
1279 Vector copy(*this);
1280 auto d_in = copy.Read(use_dev);
1281 auto d_out = Write(use_dev);
1282 CopyFlagged(use_dev, d_in, d_flag, d_out, d_flag + size, size);
1283
1284 // assumes indices are unique
1285 size -= indices.Size();
1286 }
1287}
1288
1289} // namespace mfem
int Size() const
Return the logical size of the array.
Definition array.hpp:166
bool UseDevice() const
Return the device flag of the Memory object used by the Array.
Definition array.hpp:151
T * Write(bool on_dev=true)
Shortcut for mfem::Write(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:389
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:381
static MemoryType GetHostMemoryType()
Get the current Host MemoryType. This is the MemoryType used by most MFEM classes when allocating mem...
Definition device.hpp:268
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:262
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition device.hpp:277
Hash function for data sequences.
Definition hash.hpp:532
std::string GetHash() const
Return the hash string for the current sequence and reset (clear) the sequence.
Definition hash.cpp:60
HashFunction & AppendDoubles(const real_t *doubles, size_t num_doubles)
Add a sequence of doubles for hashing, given as a c-array.
Definition hash.hpp:584
void CopyFromHost(const T *src, int size)
Copy size entries from the host pointer src to *this.
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...
bool Empty() const
Return true if the Memory object is empty, see Reset().
const T * Read(MemoryClass mc, int size) const
Get read-only access to the memory with the given MemoryClass.
void CopyFrom(const Memory &src, int size)
Copy size entries from src to *this.
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
Vector data type.
Definition vector.hpp:82
void Randomize(int seed=0)
Set random values in the vector.
Definition vector.cpp:955
void PrintHash(std::ostream &out) const
Print the Vector size and hash of its data.
Definition vector.cpp:947
void Pow(const real_t p)
(*this)(i) = pow((*this)(i), p)
Definition vector.cpp:403
real_t & Elem(int i)
Access Vector entries. Index i = 0 .. size-1.
Definition vector.cpp:173
void SetVector(const Vector &v, int offset)
(*this)[i + offset] = v[i]
Definition vector.cpp:353
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:524
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:520
void median(const Vector &lo, const Vector &hi)
v = median(v,lo,hi) entrywise. Implementation assumes lo <= hi.
Definition vector.cpp:652
void Neg()
(*this) = -(*this)
Definition vector.cpp:376
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
Definition vector.cpp:870
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:536
real_t operator*(const real_t *v) const
Definition vector.cpp:183
real_t Normlinf() const
Returns the l_infinity norm of the vector.
Definition vector.cpp:1004
real_t Norml1() const
Returns the l_1 norm of the vector.
Definition vector.cpp:1018
void SetSubVectorHost(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value (always on host).
Definition vector.cpp:723
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
Definition vector.cpp:702
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:785
Memory< real_t > data
Definition vector.hpp:85
void AddSubVector(const Vector &v, int offset)
(*this)[i + offset] += v[i]
Definition vector.cpp:365
void Swap(Vector &other)
Swap the contents of two Vectors.
Definition vector.hpp:697
void DeleteAt(const Array< int > &indices)
Definition vector.cpp:1260
Vector & operator*=(real_t c)
Definition vector.cpp:241
real_t Norml2() const
Returns the l2 norm of the vector.
Definition vector.cpp:968
void Load(std::istream **in, int np, int *dim)
Reads a vector from multiple files.
Definition vector.cpp:127
Vector & Set(const real_t a, const Vector &x)
(*this) = a * x
Definition vector.cpp:341
real_t Max() const
Returns the maximal element of the vector.
Definition vector.cpp:1200
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition vector.hpp:148
int Size() const
Returns the size of the vector.
Definition vector.hpp:234
void Abs()
(*this)(i) = abs((*this)(i))
Definition vector.cpp:392
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition vector.hpp:145
real_t Sum() const
Return the sum of the vector entries.
Definition vector.cpp:1246
void PrintMathematica(std::ostream &out=mfem::out) const
Prints vector as a List for importing into Mathematica.
Definition vector.cpp:923
void Print_HYPRE(std::ostream &out) const
Prints vector to stream out in HYPRE_Vector format.
Definition vector.cpp:904
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:584
void Reciprocal()
(*this)(i) = 1.0 / (*this)(i)
Definition vector.cpp:384
virtual real_t * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:532
real_t Normlp(real_t p) const
Returns the l_p norm of the vector.
Definition vector.cpp:1032
Vector & operator-=(real_t c)
Definition vector.cpp:284
Vector & operator=(const real_t *v)
Copy Size() entries from v.
Definition vector.cpp:197
void SetSubVectorComplement(const Array< int > &dofs, const real_t val)
Set all vector entries NOT in the dofs Array to the given val.
Definition vector.cpp:854
virtual real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:540
real_t Min() const
Returns the minimal element of the vector.
Definition vector.cpp:1154
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:676
Vector & operator+=(real_t c)
Definition vector.cpp:305
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:528
Vector & Add(const real_t a, const Vector &Va)
(*this) += a * Va
Definition vector.cpp:326
real_t & operator()(int i)
Access Vector entries using () for 0-based indexing.
Definition vector.hpp:681
void cross3D(const Vector &vin, Vector &vout) const
Definition vector.cpp:639
Vector & operator/=(real_t c)
Definition vector.cpp:262
const real_t alpha
Definition ex15.cpp:369
int dim
Definition ex24.cpp:53
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
mfem::real_t real_t
MFEM_HOST_DEVICE dual< value_type, gradient_type > pow(dual< value_type, gradient_type > a, dual< value_type, gradient_type > b)
implementation of a (dual) raised to the b (dual) power
Definition dual.hpp:374
real_t infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition vector.hpp:47
void reduce(int N, T &res, B &&body, const R &reducer, bool use_dev, Array< T > &workspace)
Performs a 1D reduction on the range [0,N). res initial value and where the result will be written....
Definition reducers.hpp:532
void CopyFlagged(bool use_dev, InputIt d_in, FlagIt d_flags, OutputIt d_out, NumSelectedIt d_num_selected_out, size_t num_items)
Equivalent to *d_num_selected_out = std::copy_if(d_in, d_in+num_items, d_out, [=](auto iter){ return ...
Definition scan.hpp:238
real_t rand_real()
Generate a random real_t number in the interval [0,1) using rand().
Definition vector.hpp:61
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition vector.cpp:414
T ZeroSubnormal(T val)
Definition vector.hpp:548
bool DeviceCanUseOcca()
Function that determines if an OCCA kernel should be used, based on the current mfem::Device configur...
Definition occa.hpp:69
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition vector.cpp:570
float real_t
Definition config.hpp:46
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....
Definition occa.hpp:37
void forall_switch(bool use_dev, int N, lambda &&body)
Definition forall.hpp:919
real_t p(const Vector &x, real_t t)
@ DEVICE_MASK
Biwise-OR of all device backends.
Definition device.hpp:99
@ OMP_MASK
Biwise-OR of all OpenMP backends.
Definition device.hpp:95
Pair of values which can be used in device code.
Definition reducers.hpp:27
a = max(a,b)
Definition reducers.hpp:145
a = min(a,b)
Definition reducers.hpp:98