14 #include "kernels.hpp"
16 #include "../general/forall.hpp"
18 #if defined(MFEM_USE_SUNDIALS)
20 #if defined(MFEM_USE_MPI)
21 #include <nvector/nvector_parallel.h>
25 #ifdef MFEM_USE_OPENMP
41 const int s = v.
Size();
45 MFEM_ASSERT(!v.
data.
Empty(),
"invalid source vector");
62 for (i = 0; i < np; i++)
70 for (i = 0; i < np; i++)
72 for (j = 0; j < dim[i]; j++)
77 if (!*in[i] && errno == ERANGE)
89 for (
int i = 0; i <
size; i++)
94 if (!in && errno == ERANGE)
114 #ifdef MFEM_USE_LEGACY_OPENMP
115 #pragma omp parallel for reduction(+:dot)
117 for (
int i = 0; i <
size; i++)
119 dot +=
data[i] * v[i];
139 const bool use_dev =
UseDevice() || vuse;
142 if (use_dev) {
Write(); }
151 data = std::move(v.data);
162 auto y =
Write(use_dev);
163 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = value;);
172 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] *= c;);
178 MFEM_ASSERT(
size == v.
size,
"incompatible Vectors!");
183 auto x = v.
Read(use_dev);
184 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] *= x[i];);
192 const double m = 1.0/c;
194 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] *= m;);
200 MFEM_ASSERT(
size == v.
size,
"incompatible Vectors!");
205 auto x = v.
Read(use_dev);
206 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] /= x[i];);
215 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] -= c;);
221 MFEM_ASSERT(
size == v.
size,
"incompatible Vectors!");
226 auto x = v.
Read(use_dev);
227 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] -= x[i];);
236 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] += c;);
242 MFEM_ASSERT(
size == v.
size,
"incompatible Vectors!");
247 auto x = v.
Read(use_dev);
248 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] += x[i];);
254 MFEM_ASSERT(
size == Va.
size,
"incompatible Vectors!");
261 auto x = Va.
Read(use_dev);
262 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] += a * x[i];);
269 MFEM_ASSERT(
size == Va.
size,
"incompatible Vectors!");
273 auto x = Va.
Read(use_dev);
274 auto y =
Write(use_dev);
275 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = a * x[i];);
281 MFEM_ASSERT(v.
Size() + offset <=
size,
"invalid sub-vector");
283 const int vs = v.
Size();
284 const double *vp = v.
data;
285 double *
p =
data + offset;
286 for (
int i = 0; i < vs; i++)
297 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = -y[i];);
303 "incompatible Vectors!");
305 #if !defined(MFEM_USE_LEGACY_OPENMP)
307 const int N = v.
size;
309 auto x1 = v1.
Read(use_dev);
310 auto x2 = v2.
Read(use_dev);
311 auto y = v.
Write(use_dev);
312 MFEM_FORALL_SWITCH(use_dev, i, N, y[i] = x1[i] + x2[i];);
314 #pragma omp parallel for
315 for (
int i = 0; i < v.
size; i++)
325 "incompatible Vectors!");
331 else if (alpha == 1.0)
337 #if !defined(MFEM_USE_LEGACY_OPENMP)
339 const int N = v.
size;
341 auto d_x = v1.
Read(use_dev);
342 auto d_y = v2.
Read(use_dev);
343 auto d_z = v.
Write(use_dev);
344 MFEM_FORALL_SWITCH(use_dev, i, N, d_z[i] = d_x[i] + alpha * d_y[i];);
346 const double *v1p = v1.
data, *v2p = v2.
data;
348 const int s = v.
size;
349 #pragma omp parallel for
350 for (
int i = 0; i <
s; i++)
352 vp[i] = v1p[i] + alpha*v2p[i];
361 "incompatible Vectors!");
373 #if !defined(MFEM_USE_LEGACY_OPENMP)
375 const int N = x.
size;
377 auto xd = x.
Read(use_dev);
378 auto yd = y.
Read(use_dev);
379 auto zd = z.
Write(use_dev);
380 MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = a * (xd[i] + yd[i]););
382 const double *xp = x.
data;
383 const double *yp = y.
data;
385 const int s = x.
size;
386 #pragma omp parallel for
387 for (
int i = 0; i <
s; i++)
389 zp[i] = a * (xp[i] + yp[i]);
399 "incompatible Vectors!");
425 #if !defined(MFEM_USE_LEGACY_OPENMP)
427 const int N = x.
size;
429 auto xd = x.
Read(use_dev);
430 auto yd = y.
Read(use_dev);
431 auto zd = z.
Write(use_dev);
432 MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = a * xd[i] + b * yd[i];);
434 const double *xp = x.
data;
435 const double *yp = y.
data;
437 const int s = x.
size;
438 #pragma omp parallel for
439 for (
int i = 0; i <
s; i++)
441 zp[i] = a * xp[i] + b * yp[i];
450 "incompatible Vectors!");
452 #if !defined(MFEM_USE_LEGACY_OPENMP)
454 const int N = x.
size;
456 auto xd = x.
Read(use_dev);
457 auto yd = y.
Read(use_dev);
458 auto zd = z.
Write(use_dev);
459 MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = xd[i] - yd[i];);
461 const double *xp = x.
data;
462 const double *yp = y.
data;
464 const int s = x.
size;
465 #pragma omp parallel for
466 for (
int i = 0; i <
s; i++)
468 zp[i] = xp[i] - yp[i];
476 "incompatible Vectors!");
488 #if !defined(MFEM_USE_LEGACY_OPENMP)
490 const int N = x.
size;
492 auto xd = x.
Read(use_dev);
493 auto yd = y.
Read(use_dev);
494 auto zd = z.
Write(use_dev);
495 MFEM_FORALL_SWITCH(use_dev, i, N, zd[i] = a * (xd[i] - yd[i]););
497 const double *xp = x.
data;
498 const double *yp = y.
data;
500 const int s = x.
size;
501 #pragma omp parallel for
502 for (
int i = 0; i <
s; i++)
504 zp[i] = a * (xp[i] - yp[i]);
513 "incompatible Vectors!");
518 auto l = lo.
Read(use_dev);
519 auto h = hi.
Read(use_dev);
520 auto m =
Write(use_dev);
521 MFEM_FORALL_SWITCH(use_dev, i, N,
527 else if (m[i] > h[i])
536 const int n = dofs.
Size();
539 auto d_y = elemvect.
Write(use_dev);
540 auto d_X =
Read(use_dev);
541 auto d_dofs = dofs.
Read(use_dev);
542 MFEM_FORALL_SWITCH(use_dev, i, n,
544 const int dof_i = d_dofs[i];
545 d_y[i] = dof_i >= 0 ? d_X[dof_i] : -d_X[-dof_i-1];
552 const int n = dofs.
Size();
553 for (
int i = 0; i < n; i++)
555 const int j = dofs[i];
556 elem_data[i] = (j >= 0) ?
data[j] : -
data[-1-j];
563 const int n = dofs.
Size();
566 auto d_dofs = dofs.
Read(use_dev);
567 MFEM_FORALL_SWITCH(use_dev, i, n,
569 const int j = d_dofs[i];
583 MFEM_ASSERT(dofs.
Size() == elemvect.
Size(),
584 "Size mismatch: length of dofs is " << dofs.
Size()
585 <<
", length of elemvect is " << elemvect.
Size());
588 const int n = dofs.
Size();
591 auto d_y = elemvect.
Read(use_dev);
592 auto d_dofs = dofs.
Read(use_dev);
593 MFEM_FORALL_SWITCH(use_dev, i, n,
595 const int dof_i = d_dofs[i];
602 d_X[-1-dof_i] = -d_y[i];
611 const int n = dofs.
Size();
612 for (
int i = 0; i < n; i++)
614 const int j= dofs[i];
628 MFEM_ASSERT(dofs.
Size() == elemvect.
Size(),
"Size mismatch: "
629 "length of dofs is " << dofs.
Size() <<
630 ", length of elemvect is " << elemvect.
Size());
633 const int n = dofs.
Size();
634 auto d_y = elemvect.
Read(use_dev);
636 auto d_dofs = dofs.
Read(use_dev);
637 MFEM_FORALL_SWITCH(use_dev, i, n,
639 const int j = d_dofs[i];
654 const int n = dofs.
Size();
655 for (
int i = 0; i < n; i++)
657 const int j = dofs[i];
672 MFEM_ASSERT(dofs.
Size() == elemvect.
Size(),
"Size mismatch: "
673 "length of dofs is " << dofs.
Size() <<
674 ", length of elemvect is " << elemvect.
Size());
677 const int n = dofs.
Size();
679 auto d_x = elemvect.
Read(use_dev);
680 auto d_dofs = dofs.
Read(use_dev);
681 MFEM_FORALL_SWITCH(use_dev, i, n,
683 const int j = d_dofs[i];
686 d_y[j] += a * d_x[i];
690 d_y[-1-j] -= a * d_x[i];
698 const int n = dofs.
Size();
700 Vector dofs_vals(n, use_dev ?
704 auto d_dofs_vals = dofs_vals.
Write(use_dev);
705 auto d_dofs = dofs.
Read(use_dev);
706 MFEM_FORALL_SWITCH(use_dev, i, n, d_dofs_vals[i] = d_data[d_dofs[i]];);
707 MFEM_FORALL_SWITCH(use_dev, i, N, d_data[i] = val;);
708 MFEM_FORALL_SWITCH(use_dev, i, n, d_data[d_dofs[i]] = d_dofs_vals[i];);
713 if (!
size) {
return; }
723 if ( i % width == 0 )
735 #ifdef MFEM_USE_ADIOS2
737 const std::string& variable_name)
const
739 if (!
size) {
return; }
741 os.engine.Put(variable_name, &
data[0] );
748 std::ios::fmtflags old_fmt = os.flags();
749 os.setf(std::ios::scientific);
750 std::streamsize old_prec = os.precision(14);
755 for (i = 0; i <
size; i++)
760 os.precision(old_prec);
766 os <<
"size: " <<
size <<
'\n';
769 os <<
"hash: " << hf.
GetHash() <<
'\n';
775 const double max = (double)(RAND_MAX) + 1.;
783 srand((
unsigned)seed);
786 for (
int i = 0; i <
size; i++)
788 data[i] = std::abs(rand()/max);
805 return std::abs(
data[0]);
814 for (
int i = 0; i <
size; i++)
816 max = std::max(std::abs(
data[i]), max);
825 for (
int i = 0; i <
size; i++)
827 sum += std::abs(
data[i]);
834 MFEM_ASSERT(p > 0.0,
"Vector::Normlp");
856 return std::abs(
data[0]);
862 for (
int i = 0; i <
size; i++)
866 const double absdata = std::abs(
data[i]);
867 if (scale <= absdata)
869 sum = 1.0 + sum *
std::pow(scale / absdata, p);
873 sum +=
std::pow(absdata / scale, p);
876 return scale *
std::pow(sum, 1.0/p);
887 double max =
data[0];
889 for (
int i = 1; i <
size; i++)
904 const double *h_data = this->
HostRead();
905 for (
int i = 0; i <
size; i++)
914 static __global__
void cuKernelMin(
const int N,
double *gdsr,
const double *x)
916 __shared__
double s_min[MFEM_CUDA_BLOCKS];
917 const int n = blockDim.x*blockIdx.x + threadIdx.x;
918 if (n>=N) {
return; }
919 const int bid = blockIdx.x;
920 const int tid = threadIdx.x;
921 const int bbd = bid*blockDim.x;
922 const int rid = bbd+tid;
924 for (
int workers=blockDim.x>>1; workers>0; workers>>=1)
927 if (tid >= workers) {
continue; }
928 if (rid >= N) {
continue; }
929 const int dualTid = tid + workers;
930 if (dualTid >= N) {
continue; }
931 const int rdd = bbd+dualTid;
932 if (rdd >= N) {
continue; }
933 if (dualTid >= blockDim.x) {
continue; }
934 s_min[tid] = fmin(s_min[tid], s_min[dualTid]);
936 if (tid==0) { gdsr[bid] = s_min[0]; }
939 static Array<double> cuda_reduce_buf;
941 static double cuVectorMin(
const int N,
const double *X)
943 const int tpb = MFEM_CUDA_BLOCKS;
944 const int blockSize = MFEM_CUDA_BLOCKS;
945 const int gridSize = (N+blockSize-1)/blockSize;
946 const int min_sz = (N%tpb)==0? (N/tpb) : (1+N/tpb);
947 cuda_reduce_buf.
SetSize(min_sz);
948 Memory<double> &buf = cuda_reduce_buf.
GetMemory();
950 cuKernelMin<<<gridSize,blockSize>>>(N, d_min, X);
951 MFEM_GPU_CHECK(cudaGetLastError());
954 for (
int i = 0; i < min_sz; i++) { min = fmin(min, h_min[i]); }
958 static __global__
void cuKernelDot(
const int N,
double *gdsr,
959 const double *x,
const double *y)
961 __shared__
double s_dot[MFEM_CUDA_BLOCKS];
962 const int n = blockDim.x*blockIdx.x + threadIdx.x;
963 if (n>=N) {
return; }
964 const int bid = blockIdx.x;
965 const int tid = threadIdx.x;
966 const int bbd = bid*blockDim.x;
967 const int rid = bbd+tid;
968 s_dot[tid] = x[n] * y[n];
969 for (
int workers=blockDim.x>>1; workers>0; workers>>=1)
972 if (tid >= workers) {
continue; }
973 if (rid >= N) {
continue; }
974 const int dualTid = tid + workers;
975 if (dualTid >= N) {
continue; }
976 const int rdd = bbd+dualTid;
977 if (rdd >= N) {
continue; }
978 if (dualTid >= blockDim.x) {
continue; }
979 s_dot[tid] += s_dot[dualTid];
981 if (tid==0) { gdsr[bid] = s_dot[0]; }
984 static double cuVectorDot(
const int N,
const double *X,
const double *Y)
986 const int tpb = MFEM_CUDA_BLOCKS;
987 const int blockSize = MFEM_CUDA_BLOCKS;
988 const int gridSize = (N+blockSize-1)/blockSize;
989 const int dot_sz = (N%tpb)==0? (N/tpb) : (1+N/tpb);
991 Memory<double> &buf = cuda_reduce_buf.
GetMemory();
993 cuKernelDot<<<gridSize,blockSize>>>(N, d_dot, X, Y);
994 MFEM_GPU_CHECK(cudaGetLastError());
997 for (
int i = 0; i < dot_sz; i++) { dot += h_dot[i]; }
1000 #endif // MFEM_USE_CUDA
1003 static __global__
void hipKernelMin(
const int N,
double *gdsr,
const double *x)
1005 __shared__
double s_min[MFEM_HIP_BLOCKS];
1006 const int n = hipBlockDim_x*hipBlockIdx_x + hipThreadIdx_x;
1007 if (n>=N) {
return; }
1008 const int bid = hipBlockIdx_x;
1009 const int tid = hipThreadIdx_x;
1010 const int bbd = bid*hipBlockDim_x;
1011 const int rid = bbd+tid;
1013 for (
int workers=hipBlockDim_x>>1; workers>0; workers>>=1)
1016 if (tid >= workers) {
continue; }
1017 if (rid >= N) {
continue; }
1018 const int dualTid = tid + workers;
1019 if (dualTid >= N) {
continue; }
1020 const int rdd = bbd+dualTid;
1021 if (rdd >= N) {
continue; }
1022 if (dualTid >= hipBlockDim_x) {
continue; }
1023 s_min[tid] = fmin(s_min[tid], s_min[dualTid]);
1025 if (tid==0) { gdsr[bid] = s_min[0]; }
1028 static Array<double> cuda_reduce_buf;
1030 static double hipVectorMin(
const int N,
const double *X)
1032 const int tpb = MFEM_HIP_BLOCKS;
1033 const int blockSize = MFEM_HIP_BLOCKS;
1034 const int gridSize = (N+blockSize-1)/blockSize;
1035 const int min_sz = (N%tpb)==0 ? (N/tpb) : (1+N/tpb);
1036 cuda_reduce_buf.
SetSize(min_sz);
1037 Memory<double> &buf = cuda_reduce_buf.
GetMemory();
1039 hipLaunchKernelGGL(hipKernelMin,gridSize,blockSize,0,0,N,d_min,X);
1040 MFEM_GPU_CHECK(hipGetLastError());
1043 for (
int i = 0; i < min_sz; i++) { min = fmin(min, h_min[i]); }
1047 static __global__
void hipKernelDot(
const int N,
double *gdsr,
1048 const double *x,
const double *y)
1050 __shared__
double s_dot[MFEM_HIP_BLOCKS];
1051 const int n = hipBlockDim_x*hipBlockIdx_x + hipThreadIdx_x;
1052 if (n>=N) {
return; }
1053 const int bid = hipBlockIdx_x;
1054 const int tid = hipThreadIdx_x;
1055 const int bbd = bid*hipBlockDim_x;
1056 const int rid = bbd+tid;
1057 s_dot[tid] = x[n] * y[n];
1058 for (
int workers=hipBlockDim_x>>1; workers>0; workers>>=1)
1061 if (tid >= workers) {
continue; }
1062 if (rid >= N) {
continue; }
1063 const int dualTid = tid + workers;
1064 if (dualTid >= N) {
continue; }
1065 const int rdd = bbd+dualTid;
1066 if (rdd >= N) {
continue; }
1067 if (dualTid >= hipBlockDim_x) {
continue; }
1068 s_dot[tid] += s_dot[dualTid];
1070 if (tid==0) { gdsr[bid] = s_dot[0]; }
1073 static double hipVectorDot(
const int N,
const double *X,
const double *Y)
1075 const int tpb = MFEM_HIP_BLOCKS;
1076 const int blockSize = MFEM_HIP_BLOCKS;
1077 const int gridSize = (N+blockSize-1)/blockSize;
1078 const int dot_sz = (N%tpb)==0 ? (N/tpb) : (1+N/tpb);
1079 cuda_reduce_buf.
SetSize(dot_sz);
1080 Memory<double> &buf = cuda_reduce_buf.
GetMemory();
1082 hipLaunchKernelGGL(hipKernelDot,gridSize,blockSize,0,0,N,d_dot,X,Y);
1083 MFEM_GPU_CHECK(hipGetLastError());
1086 for (
int i = 0; i < dot_sz; i++) { dot += h_dot[i]; }
1089 #endif // MFEM_USE_HIP
1093 MFEM_ASSERT(
size == v.
size,
"incompatible Vectors!");
1094 if (
size == 0) {
return 0.0; }
1097 #if defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP) || defined(MFEM_USE_OPENMP)
1098 auto m_data =
Read(use_dev);
1102 auto v_data = v.
Read(use_dev);
1104 if (!use_dev) {
goto vector_dot_cpu; }
1106 #ifdef MFEM_USE_OCCA
1109 return occa::linalg::dot<double,double,double>(
1114 #ifdef MFEM_USE_CUDA
1117 return cuVectorDot(
size, m_data, v_data);
1124 return hipVectorDot(
size, m_data, v_data);
1128 #ifdef MFEM_USE_OPENMP
1131 #define MFEM_USE_OPENMP_DETERMINISTIC_DOT
1132 #ifdef MFEM_USE_OPENMP_DETERMINISTIC_DOT
1135 #pragma omp parallel
1137 const int nt = omp_get_num_threads();
1140 const int tid = omp_get_thread_num();
1141 const int stride = (
size + nt - 1)/nt;
1142 const int start = tid*stride;
1143 const int stop = std::min(start + stride,
size);
1144 double my_dot = 0.0;
1145 for (
int i = start; i < stop; i++)
1147 my_dot += m_data[i] * v_data[i];
1150 th_dot(tid) = my_dot;
1152 return th_dot.
Sum();
1156 #pragma omp parallel for reduction(+:prod)
1157 for (
int i = 0; i <
size; i++)
1159 prod += m_data[i] * v_data[i];
1162 #endif // MFEM_USE_OPENMP_DETERMINISTIC_DOT
1164 #endif // MFEM_USE_OPENMP
1168 auto v_data_ = v.
Read();
1169 auto m_data_ =
Read();
1172 auto d_dot = dot.
Write();
1174 MFEM_FORALL(i, N, d_dot[0] += m_data_[i] * v_data_[i];);
1187 auto m_data =
Read(use_dev);
1189 if (!use_dev) {
goto vector_min_cpu; }
1191 #ifdef MFEM_USE_OCCA
1198 #ifdef MFEM_USE_CUDA
1201 return cuVectorMin(
size, m_data);
1208 return hipVectorMin(
size, m_data);
1212 #ifdef MFEM_USE_OPENMP
1215 double minimum = m_data[0];
1216 #pragma omp parallel for reduction(min:minimum)
1217 for (
int i = 0; i <
size; i++)
1219 minimum = std::min(minimum, m_data[i]);
1228 auto m_data_ =
Read();
1233 MFEM_FORALL(i, N, d_min[0] = (d_min[0]<m_data_[i])?d_min[0]:m_data_[i];);
1239 double minimum =
data[0];
1240 for (
int i = 1; i <
size; i++)
1242 if (m_data[i] < minimum)
1244 minimum = m_data[i];
1251 #ifdef MFEM_USE_SUNDIALS
1255 N_Vector_ID nvid = N_VGetVectorID(nv);
1259 case SUNDIALS_NVEC_SERIAL:
1263 case SUNDIALS_NVEC_PARALLEL:
1268 MFEM_ABORT(
"N_Vector type " << nvid <<
" is not supported");
1274 MFEM_ASSERT(nv,
"N_Vector handle is NULL");
1275 N_Vector_ID nvid = N_VGetVectorID(nv);
1279 case SUNDIALS_NVEC_SERIAL:
1280 MFEM_ASSERT(NV_OWN_DATA_S(nv) == SUNFALSE,
"invalid serial N_Vector");
1281 NV_DATA_S(nv) =
data;
1282 NV_LENGTH_S(nv) =
size;
1285 case SUNDIALS_NVEC_PARALLEL:
1286 MFEM_ASSERT(NV_OWN_DATA_P(nv) == SUNFALSE,
"invalid parallel N_Vector");
1287 NV_DATA_P(nv) =
data;
1288 NV_LOCLENGTH_P(nv) =
size;
1289 if (global_length == 0)
1291 global_length = NV_GLOBLENGTH_P(nv);
1293 if (global_length == 0 && global_length !=
size)
1295 MPI_Comm sundials_comm = NV_COMM_P(nv);
1296 long local_size =
size;
1297 MPI_Allreduce(&local_size, &global_length, 1, MPI_LONG,
1298 MPI_SUM,sundials_comm);
1301 NV_GLOBLENGTH_P(nv) = global_length;
1305 MFEM_ABORT(
"N_Vector type " << nvid <<
" is not supported");
1309 #endif // MFEM_USE_SUNDIALS
Hash function for data sequences.
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)
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 GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
virtual double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
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.
std::string GetHash() const
Return the hash string for the current sequence and reset (clear) the sequence.
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.
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
void add(const Vector &v1, const Vector &v2, Vector &v)
HashFunction & AppendDoubles(const double *doubles, size_t num_doubles)
Add a sequence of doubles for hashing, given as a c-array.
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.
double Normlp(double p) const
Returns the l_p norm of 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...
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
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 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...
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
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.
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
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
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
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. ...
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()
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
void PrintHash(std::ostream &out) const
Print the Vector size and hash of its data.
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.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
[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)