16 #include "../general/forall.hpp"
18 #if defined(MFEM_USE_SUNDIALS)
20 #if defined(MFEM_USE_MPI)
21 #include <nvector/nvector_parallel.h>
37 const int s = v.
Size();
40 MFEM_ASSERT(!v.
data.
Empty(),
"invalid source vector");
58 for (i = 0; i < np; i++)
66 for (i = 0; i < np; i++)
68 for (j = 0; j < dim[i]; j++)
73 if (!*in[i] && errno == ERANGE)
85 for (
int i = 0; i <
size; i++)
90 if (!in && errno == ERANGE)
110 #ifdef MFEM_USE_LEGACY_OPENMP
111 #pragma omp parallel for reduction(+:dot)
113 for (
int i = 0; i <
size; i++)
115 dot +=
data[i] * v[i];
137 if (use_dev) {
Write(); }
147 auto y =
Write(use_dev);
148 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = value;);
157 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] *= c;);
165 const double m = 1.0/c;
167 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] *= m;);
176 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] -= c;);
182 MFEM_ASSERT(
size == v.
size,
"incompatible Vectors!");
187 auto x = v.
Read(use_dev);
188 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] -= x[i];);
197 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] += c;);
203 MFEM_ASSERT(
size == v.
size,
"incompatible Vectors!");
208 auto x = v.
Read(use_dev);
209 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] += x[i];);
215 MFEM_ASSERT(
size == Va.
size,
"incompatible Vectors!");
222 auto x = Va.
Read(use_dev);
223 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] += a * x[i];);
230 MFEM_ASSERT(
size == Va.
size,
"incompatible Vectors!");
234 auto x = Va.
Read(use_dev);
235 auto y =
Write(use_dev);
236 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = a * x[i];);
242 MFEM_ASSERT(v.
Size() + offset <=
size,
"invalid sub-vector");
244 const int vs = v.
Size();
245 const double *vp = v.
data;
246 double *
p =
data + offset;
247 for (
int i = 0; i < vs; i++)
258 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = -y[i];);
264 "incompatible Vectors!");
266 #if !defined(MFEM_USE_LEGACY_OPENMP)
268 const int N = v.
size;
270 auto x1 = v1.
Read(use_dev);
271 auto x2 = v2.
Read(use_dev);
272 auto y = v.
Write(use_dev);
273 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = x1[i] + x2[i];);
275 #pragma omp parallel for
276 for (
int i = 0; i < v.
size; i++)
286 "incompatible Vectors!");
292 else if (alpha == 1.0)
298 #if !defined(MFEM_USE_LEGACY_OPENMP)
300 const int N = v.
size;
302 auto d_x = v1.
Read(use_dev);
303 auto d_y = v2.
Read(use_dev);
304 auto d_z = v.
Write(use_dev);
305 MFEM_FORALL_SWITCH(use_dev, i, N, d_z[i] = d_x[i] + alpha * d_y[i];);
307 const double *v1p = v1.
data, *v2p = v2.
data;
309 const int s = v.
size;
310 #pragma omp parallel for
311 for (
int i = 0; i < s; i++)
313 vp[i] = v1p[i] + alpha*v2p[i];
322 "incompatible Vectors!");
334 #if !defined(MFEM_USE_LEGACY_OPENMP)
336 const int N = x.
size;
338 auto xd = x.
Read(use_dev);
339 auto yd = y.
Read(use_dev);
340 auto zd = z.
Write(use_dev);
341 MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = a * (xd[i] + yd[i]););
343 const double *xp = x.
data;
344 const double *yp = y.
data;
346 const int s = x.
size;
347 #pragma omp parallel for
348 for (
int i = 0; i < s; i++)
350 zp[i] = a * (xp[i] + yp[i]);
360 "incompatible Vectors!");
386 #if !defined(MFEM_USE_LEGACY_OPENMP)
388 const int N = x.
size;
390 auto xd = x.
Read(use_dev);
391 auto yd = y.
Read(use_dev);
392 auto zd = z.
Write(use_dev);
393 MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = a * xd[i] + b * yd[i];);
395 const double *xp = x.
data;
396 const double *yp = y.
data;
398 const int s = x.
size;
399 #pragma omp parallel for
400 for (
int i = 0; i < s; i++)
402 zp[i] = a * xp[i] + b * yp[i];
411 "incompatible Vectors!");
413 #if !defined(MFEM_USE_LEGACY_OPENMP)
415 const int N = x.
size;
417 auto xd = x.
Read(use_dev);
418 auto yd = y.
Read(use_dev);
419 auto zd = z.
Write(use_dev);
420 MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = xd[i] - yd[i];);
422 const double *xp = x.
data;
423 const double *yp = y.
data;
425 const int s = x.
size;
426 #pragma omp parallel for
427 for (
int i = 0; i < s; i++)
429 zp[i] = xp[i] - yp[i];
437 "incompatible Vectors!");
449 #if !defined(MFEM_USE_LEGACY_OPENMP)
451 const int N = x.
size;
453 auto xd = x.
Read(use_dev);
454 auto yd = y.
Read(use_dev);
455 auto zd = z.
Write(use_dev);
456 MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = a * (xd[i] - yd[i]););
458 const double *xp = x.
data;
459 const double *yp = y.
data;
461 const int s = x.
size;
462 #pragma omp parallel for
463 for (
int i = 0; i < s; i++)
465 zp[i] = a * (xp[i] - yp[i]);
474 "incompatible Vectors!");
479 auto l = lo.
Read(use_dev);
480 auto h = hi.
Read(use_dev);
481 auto m =
Write(use_dev);
482 MFEM_FORALL_SWITCH(use_dev, i, N,
488 else if (m[i] > h[i])
497 const int n = dofs.
Size();
500 auto d_y = elemvect.
Write(use_dev);
501 auto d_X =
Read(use_dev);
502 auto d_dofs = dofs.
Read(use_dev);
503 MFEM_FORALL_SWITCH(use_dev, i, n,
505 const int dof_i = d_dofs[i];
506 d_y[i] = dof_i >= 0 ? d_X[dof_i] : -d_X[-dof_i-1];
513 const int n = dofs.
Size();
514 for (
int i = 0; i < n; i++)
516 const int j = dofs[i];
517 elem_data[i] = (j >= 0) ?
data[j] : -
data[-1-j];
524 const int n = dofs.
Size();
527 auto d_dofs = dofs.
Read(use_dev);
528 MFEM_FORALL_SWITCH(use_dev, i, n,
530 const int j = d_dofs[i];
544 MFEM_ASSERT(dofs.
Size() == elemvect.
Size(),
545 "Size mismatch: length of dofs is " << dofs.
Size()
546 <<
", length of elemvect is " << elemvect.
Size());
549 const int n = dofs.
Size();
552 auto d_y = elemvect.
Read(use_dev);
553 auto d_dofs = dofs.
Read(use_dev);
554 MFEM_FORALL_SWITCH(use_dev, i, n,
556 const int dof_i = d_dofs[i];
563 d_X[-1-dof_i] = -d_y[i];
572 const int n = dofs.
Size();
573 for (
int i = 0; i < n; i++)
575 const int j= dofs[i];
589 MFEM_ASSERT(dofs.
Size() == elemvect.
Size(),
"Size mismatch: "
590 "length of dofs is " << dofs.
Size() <<
591 ", length of elemvect is " << elemvect.
Size());
594 const int n = dofs.
Size();
595 auto d_y = elemvect.
Read(use_dev);
597 auto d_dofs = dofs.
Read(use_dev);
598 MFEM_FORALL_SWITCH(use_dev, i, n,
600 const int j = d_dofs[i];
615 const int n = dofs.
Size();
616 for (
int i = 0; i < n; i++)
618 const int j = dofs[i];
633 MFEM_ASSERT(dofs.
Size() == elemvect.
Size(),
"Size mismatch: "
634 "length of dofs is " << dofs.
Size() <<
635 ", length of elemvect is " << elemvect.
Size());
638 const int n = dofs.
Size();
640 auto d_x = elemvect.
Read(use_dev);
641 auto d_dofs = dofs.
Read(use_dev);
642 MFEM_FORALL_SWITCH(use_dev, i, n,
644 const int j = d_dofs[i];
647 d_y[j] += a * d_x[i];
651 d_y[-1-j] -= a * d_x[i];
659 const int n = dofs.
Size();
661 Vector dofs_vals(n, use_dev ?
665 auto d_dofs_vals = dofs_vals.
Write(use_dev);
666 auto d_dofs = dofs.
Read(use_dev);
667 MFEM_FORALL_SWITCH(use_dev, i, n, d_dofs_vals[i] = d_data[d_dofs[i]];);
668 MFEM_FORALL_SWITCH(use_dev, i, N, d_data[i] = val;);
669 MFEM_FORALL_SWITCH(use_dev, i, n, d_data[d_dofs[i]] = d_dofs_vals[i];);
674 if (!
size) {
return; }
684 if ( i % width == 0 )
696 #ifdef MFEM_USE_ADIOS2
698 const std::string& variable_name)
const
700 if (!
size) {
return; }
702 out.engine.Put(variable_name, &
data[0] );
709 std::ios::fmtflags old_fmt = out.flags();
710 out.setf(std::ios::scientific);
711 std::streamsize old_prec = out.precision(14);
716 for (i = 0; i <
size; i++)
721 out.precision(old_prec);
728 const double max = (double)(RAND_MAX) + 1.;
736 srand((
unsigned)seed);
738 for (
int i = 0; i <
size; i++)
740 data[i] = std::abs(rand()/max);
757 return std::abs(
data[0]);
765 for (
int i = 0; i <
size; i++)
767 max = std::max(std::abs(
data[i]), max);
775 for (
int i = 0; i <
size; i++)
777 sum += std::abs(
data[i]);
784 MFEM_ASSERT(p > 0.0,
"Vector::Normlp");
806 return std::abs(
data[0]);
812 for (
int i = 0; i <
size; i++)
816 const double absdata = std::abs(
data[i]);
817 if (scale <= absdata)
819 sum = 1.0 + sum * std::pow(scale / absdata, p);
823 sum += std::pow(absdata / scale, p);
826 return scale * std::pow(sum, 1.0/p);
836 double max =
data[0];
838 for (
int i = 1; i <
size; i++)
853 const double *h_data = this->
HostRead();
854 for (
int i = 0; i <
size; i++)
863 static __global__
void cuKernelMin(
const int N,
double *gdsr,
const double *x)
865 __shared__
double s_min[MFEM_CUDA_BLOCKS];
866 const int n = blockDim.x*blockIdx.x + threadIdx.x;
867 if (n>=N) {
return; }
868 const int bid = blockIdx.x;
869 const int tid = threadIdx.x;
870 const int bbd = bid*blockDim.x;
871 const int rid = bbd+tid;
873 for (
int workers=blockDim.x>>1; workers>0; workers>>=1)
876 if (tid >= workers) {
continue; }
877 if (rid >= N) {
continue; }
878 const int dualTid = tid + workers;
879 if (dualTid >= N) {
continue; }
880 const int rdd = bbd+dualTid;
881 if (rdd >= N) {
continue; }
882 if (dualTid >= blockDim.x) {
continue; }
883 s_min[tid] = fmin(s_min[tid], s_min[dualTid]);
885 if (tid==0) { gdsr[bid] = s_min[0]; }
888 static Array<double> cuda_reduce_buf;
890 static double cuVectorMin(
const int N,
const double *X)
892 const int tpb = MFEM_CUDA_BLOCKS;
893 const int blockSize = MFEM_CUDA_BLOCKS;
894 const int gridSize = (N+blockSize-1)/blockSize;
895 const int min_sz = (N%tpb)==0? (N/tpb) : (1+N/tpb);
896 cuda_reduce_buf.
SetSize(min_sz);
897 Memory<double> &buf = cuda_reduce_buf.
GetMemory();
899 cuKernelMin<<<gridSize,blockSize>>>(N, d_min, X);
900 MFEM_GPU_CHECK(cudaGetLastError());
903 for (
int i = 0; i < min_sz; i++) { min = fmin(min, h_min[i]); }
907 static __global__
void cuKernelDot(
const int N,
double *gdsr,
908 const double *x,
const double *y)
910 __shared__
double s_dot[MFEM_CUDA_BLOCKS];
911 const int n = blockDim.x*blockIdx.x + threadIdx.x;
912 if (n>=N) {
return; }
913 const int bid = blockIdx.x;
914 const int tid = threadIdx.x;
915 const int bbd = bid*blockDim.x;
916 const int rid = bbd+tid;
917 s_dot[tid] = x[n] * y[n];
918 for (
int workers=blockDim.x>>1; workers>0; workers>>=1)
921 if (tid >= workers) {
continue; }
922 if (rid >= N) {
continue; }
923 const int dualTid = tid + workers;
924 if (dualTid >= N) {
continue; }
925 const int rdd = bbd+dualTid;
926 if (rdd >= N) {
continue; }
927 if (dualTid >= blockDim.x) {
continue; }
928 s_dot[tid] += s_dot[dualTid];
930 if (tid==0) { gdsr[bid] = s_dot[0]; }
933 static double cuVectorDot(
const int N,
const double *X,
const double *Y)
935 const int tpb = MFEM_CUDA_BLOCKS;
936 const int blockSize = MFEM_CUDA_BLOCKS;
937 const int gridSize = (N+blockSize-1)/blockSize;
938 const int dot_sz = (N%tpb)==0? (N/tpb) : (1+N/tpb);
940 Memory<double> &buf = cuda_reduce_buf.
GetMemory();
942 cuKernelDot<<<gridSize,blockSize>>>(N, d_dot, X, Y);
943 MFEM_GPU_CHECK(cudaGetLastError());
946 for (
int i = 0; i < dot_sz; i++) { dot += h_dot[i]; }
949 #endif // MFEM_USE_CUDA
952 static __global__
void hipKernelMin(
const int N,
double *gdsr,
const double *x)
954 __shared__
double s_min[MFEM_HIP_BLOCKS];
955 const int n = hipBlockDim_x*hipBlockIdx_x + hipThreadIdx_x;
956 if (n>=N) {
return; }
957 const int bid = hipBlockIdx_x;
958 const int tid = hipThreadIdx_x;
959 const int bbd = bid*hipBlockDim_x;
960 const int rid = bbd+tid;
962 for (
int workers=hipBlockDim_x>>1; workers>0; workers>>=1)
965 if (tid >= workers) {
continue; }
966 if (rid >= N) {
continue; }
967 const int dualTid = tid + workers;
968 if (dualTid >= N) {
continue; }
969 const int rdd = bbd+dualTid;
970 if (rdd >= N) {
continue; }
971 if (dualTid >= hipBlockDim_x) {
continue; }
972 s_min[tid] = fmin(s_min[tid], s_min[dualTid]);
974 if (tid==0) { gdsr[bid] = s_min[0]; }
977 static Array<double> cuda_reduce_buf;
979 static double hipVectorMin(
const int N,
const double *X)
981 const int tpb = MFEM_HIP_BLOCKS;
982 const int blockSize = MFEM_HIP_BLOCKS;
983 const int gridSize = (N+blockSize-1)/blockSize;
984 const int min_sz = (N%tpb)==0 ? (N/tpb) : (1+N/tpb);
985 cuda_reduce_buf.
SetSize(min_sz);
986 Memory<double> &buf = cuda_reduce_buf.
GetMemory();
988 hipLaunchKernelGGL(hipKernelMin,gridSize,blockSize,0,0,N,d_min,X);
989 MFEM_GPU_CHECK(hipGetLastError());
992 for (
int i = 0; i < min_sz; i++) { min = fmin(min, h_min[i]); }
996 static __global__
void hipKernelDot(
const int N,
double *gdsr,
997 const double *x,
const double *y)
999 __shared__
double s_dot[MFEM_HIP_BLOCKS];
1000 const int n = hipBlockDim_x*hipBlockIdx_x + hipThreadIdx_x;
1001 if (n>=N) {
return; }
1002 const int bid = hipBlockIdx_x;
1003 const int tid = hipThreadIdx_x;
1004 const int bbd = bid*hipBlockDim_x;
1005 const int rid = bbd+tid;
1006 s_dot[tid] = x[n] * y[n];
1007 for (
int workers=hipBlockDim_x>>1; workers>0; workers>>=1)
1010 if (tid >= workers) {
continue; }
1011 if (rid >= N) {
continue; }
1012 const int dualTid = tid + workers;
1013 if (dualTid >= N) {
continue; }
1014 const int rdd = bbd+dualTid;
1015 if (rdd >= N) {
continue; }
1016 if (dualTid >= hipBlockDim_x) {
continue; }
1017 s_dot[tid] += s_dot[dualTid];
1019 if (tid==0) { gdsr[bid] = s_dot[0]; }
1022 static double hipVectorDot(
const int N,
const double *X,
const double *Y)
1024 const int tpb = MFEM_HIP_BLOCKS;
1025 const int blockSize = MFEM_HIP_BLOCKS;
1026 const int gridSize = (N+blockSize-1)/blockSize;
1027 const int dot_sz = (N%tpb)==0 ? (N/tpb) : (1+N/tpb);
1028 cuda_reduce_buf.
SetSize(dot_sz);
1029 Memory<double> &buf = cuda_reduce_buf.
GetMemory();
1031 hipLaunchKernelGGL(hipKernelDot,gridSize,blockSize,0,0,N,d_dot,X,Y);
1032 MFEM_GPU_CHECK(hipGetLastError());
1035 for (
int i = 0; i < dot_sz; i++) { dot += h_dot[i]; }
1038 #endif // MFEM_USE_HIP
1042 MFEM_ASSERT(
size == v.
size,
"incompatible Vectors!");
1045 #if defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP) || defined(MFEM_USE_OPENMP)
1046 auto m_data =
Read(use_dev);
1050 auto v_data = v.
Read(use_dev);
1052 if (!use_dev) {
goto vector_dot_cpu; }
1054 #ifdef MFEM_USE_OCCA
1057 return occa::linalg::dot<double,double,double>(
1062 #ifdef MFEM_USE_CUDA
1065 return cuVectorDot(
size, m_data, v_data);
1072 return hipVectorDot(
size, m_data, v_data);
1076 #ifdef MFEM_USE_OPENMP
1080 #pragma omp parallel for reduction(+:prod)
1081 for (
int i = 0; i <
size; i++)
1083 prod += m_data[i] * v_data[i];
1091 auto v_data = v.
Read();
1092 auto m_data =
Read();
1095 auto d_dot = dot.
Write();
1097 MFEM_FORALL(i, N, d_dot[0] += m_data[i] * v_data[i];);
1110 auto m_data =
Read(use_dev);
1112 if (!use_dev) {
goto vector_min_cpu; }
1114 #ifdef MFEM_USE_OCCA
1121 #ifdef MFEM_USE_CUDA
1124 return cuVectorMin(
size, m_data);
1131 return hipVectorMin(
size, m_data);
1135 #ifdef MFEM_USE_OPENMP
1138 double minimum = m_data[0];
1139 #pragma omp parallel for reduction(min:minimum)
1140 for (
int i = 0; i <
size; i++)
1142 minimum = std::min(minimum, m_data[i]);
1151 auto m_data =
Read();
1156 MFEM_FORALL(i, N, d_min[0] = (d_min[0]<m_data[i])?d_min[0]:m_data[i];);
1162 double minimum =
data[0];
1163 for (
int i = 1; i <
size; i++)
1165 if (m_data[i] < minimum)
1167 minimum = m_data[i];
1174 #ifdef MFEM_USE_SUNDIALS
1178 N_Vector_ID nvid = N_VGetVectorID(nv);
1182 case SUNDIALS_NVEC_SERIAL:
1186 case SUNDIALS_NVEC_PARALLEL:
1191 MFEM_ABORT(
"N_Vector type " << nvid <<
" is not supported");
1197 MFEM_ASSERT(nv,
"N_Vector handle is NULL");
1198 N_Vector_ID nvid = N_VGetVectorID(nv);
1202 case SUNDIALS_NVEC_SERIAL:
1203 MFEM_ASSERT(NV_OWN_DATA_S(nv) == SUNFALSE,
"invalid serial N_Vector");
1204 NV_DATA_S(nv) =
data;
1205 NV_LENGTH_S(nv) =
size;
1208 case SUNDIALS_NVEC_PARALLEL:
1209 MFEM_ASSERT(NV_OWN_DATA_P(nv) == SUNFALSE,
"invalid parallel N_Vector");
1210 NV_DATA_P(nv) =
data;
1211 NV_LOCLENGTH_P(nv) =
size;
1212 if (global_length == 0)
1214 global_length = NV_GLOBLENGTH_P(nv);
1216 if (global_length == 0 && global_length !=
size)
1218 MPI_Comm sundials_comm = NV_COMM_P(nv);
1219 long local_size =
size;
1220 MPI_Allreduce(&local_size, &global_length, 1, MPI_LONG,
1221 MPI_SUM,sundials_comm);
1224 NV_GLOBLENGTH_P(nv) = global_length;
1228 MFEM_ABORT(
"N_Vector type " << nvid <<
" is not supported");
1232 #endif // MFEM_USE_SUNDIALS
int Size() const
Return the logical size of the array.
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
void SetVector(const Vector &v, int offset)
double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Vector()
Default constructor for Vector. Sets size = 0 and data = NULL.
Memory< T > & GetMemory()
Return a reference to the Memory object used by the Array.
double & Elem(int i)
Access Vector entries. Index i = 0 .. size-1.
Device memory; using CUDA or HIP *Malloc and *Free.
void SetSize(int s)
Resize the vector to size s.
double Norml2() const
Returns the l2 norm of the vector.
double & operator()(int i)
Access Vector entries using () for 0-based indexing.
Biwise-OR of all HIP backends.
void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
int Size() const
Returns the size of the vector.
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.
double Normlinf() const
Returns the l_infinity norm of the vector.
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.
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
void add(const Vector &v1, const Vector &v2, Vector &v)
double operator*(const double *) const
Dot product with a double * array.
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...
Vector & operator=(const double *v)
Copy Size() entries from v.
void Load(std::istream **in, int np, int *dim)
Reads a vector from multiple files.
const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
double Normlp(double p) const
Returns the l_p norm of the vector.
double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
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...
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
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.
void median(const Vector &lo, const Vector &hi)
v = median(v,lo,hi) entrywise. Implementation assumes lo <= hi.
Biwise-OR of all OpenMP backends.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
MFEM_HOST_DEVICE double Norml2(const int size, const T *data)
Returns the l2 norm of the Vector with given size and data.
void Reset()
Reset the memory to be empty, ensuring that Delete() will be a no-op.
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof...
Vector & operator/=(double c)
Biwise-OR of all CUDA backends.
void SetSubVectorComplement(const Array< int > &dofs, const double val)
Set all vector entries NOT in the dofs Array to the given val.
Vector & operator*=(double c)
double Min() const
Returns the minimal element of the vector.
Vector & operator+=(double c)
double Norml1() const
Returns the l_1 norm of the vector.
double p(const Vector &x, double t)
static MemoryType GetHostMemoryType()
Get the current Host MemoryType. This is the MemoryType used by most MFEM classes when allocating mem...
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
void subtract(const Vector &x, const Vector &y, Vector &z)
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void Print_HYPRE(std::ostream &out) const
Prints vector to stream out in HYPRE_Vector format.
static bool Allows(unsigned long b_mask)
Return true if any of the backends in the backend mask, b_mask, are allowed.
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Host memory; using new[] and delete[].
virtual MFEM_DEPRECATED N_Vector ToNVector()
(DEPRECATED) Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_SERIAL. ...
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
double Max() const
Returns the maximal element of the vector.
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
bool Empty() const
Return true if the Memory object is empty, see Reset().
T * ReadWrite(MemoryClass mc, int size)
Get read-write access to the memory with the given MemoryClass.
Vector & operator-=(double c)
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
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.
[device] Debug backend: host memory is READ/WRITE protected while a device is in use. It allows to test the "device" code-path (using separate host/device memory pools and host <-> device transfers) without any GPU hardware. As 'DEBUG' is sometimes used as a macro, _DEVICE has been added to avoid conflicts.
void Neg()
(*this) = -(*this)