MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
quadinterpolator_face.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, 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 
13 #include "../general/forall.hpp"
14 #include "../linalg/dtensor.hpp"
15 #include "../linalg/kernels.hpp"
16 
17 namespace mfem
18 {
19 
20 /// Return the sign to apply to the normals on each face to point from e1 to e2.
21 static void GetSigns(const FiniteElementSpace &fes, const FaceType type,
22  Array<bool> &signs)
23 {
24  const Mesh &mesh = *fes.GetMesh();
25  const int dim = mesh.SpaceDimension();
26  int face_id;
27  int f_ind = 0;
28  for (int f = 0; f < mesh.GetNumFacesWithGhost(); ++f)
29  {
30  Mesh::FaceInformation face = mesh.GetFaceInformation(f);
31  face_id = face.element[0].local_face_id;
32  if (face.IsNonconformingCoarse())
33  {
34  // We skip nonconforming coarse-fine faces as they are treated
35  // by the corresponding nonconforming fine-coarse faces.
36  continue;
37  }
38  else if ( face.IsOfFaceType(type) )
39  {
40  if (dim==2)
41  {
42  if (face_id==2 || face_id==3)
43  {
44  signs[f_ind] = true;
45  }
46  else
47  {
48  signs[f_ind] = false;
49  }
50  }
51  else if (dim==3)
52  {
53  if (face_id==0 || face_id==3 || face_id==4)
54  {
55  signs[f_ind] = true;
56  }
57  else
58  {
59  signs[f_ind] = false;
60  }
61  }
62  f_ind++;
63  }
64  }
65 }
66 
68  const FiniteElementSpace &fes,
69  const IntegrationRule &ir, FaceType type_)
70  : type(type_), nf(fes.GetNFbyType(type)), signs(nf)
71 {
72  fespace = &fes;
73  IntRule = &ir;
74  use_tensor_products = true; // not implemented yet (not used)
75 
76  if (fespace->GetNE() == 0) { return; }
77  GetSigns(*fespace, type, signs);
78  const FiniteElement *fe = fespace->GetFE(0);
79  const ScalarFiniteElement *sfe =
80  dynamic_cast<const ScalarFiniteElement*>(fe);
81  const TensorBasisElement *tfe =
82  dynamic_cast<const TensorBasisElement*>(fe);
83  MFEM_VERIFY(sfe != NULL, "Only scalar finite elements are supported");
84  MFEM_VERIFY(tfe != NULL &&
87  "Only Gauss-Lobatto and Bernstein basis are supported in "
88  "FaceQuadratureInterpolator.");
89 }
90 
91 template<const int T_VDIM, const int T_ND1D, const int T_NQ1D>
93  const int NF,
94  const int vdim,
95  const DofToQuad &maps,
96  const Array<bool> &signs,
97  const Vector &f_vec,
98  Vector &q_val,
99  Vector &q_der,
100  Vector &q_det,
101  Vector &q_nor,
102  const int eval_flags)
103 {
104  const int nd1d = maps.ndof;
105  const int nq1d = maps.nqpt;
106  const int ND1D = T_ND1D ? T_ND1D : nd1d;
107  const int NQ1D = T_NQ1D ? T_NQ1D : nq1d;
108  const int VDIM = T_VDIM ? T_VDIM : vdim;
109  MFEM_VERIFY(ND1D <= MAX_ND1D, "");
110  MFEM_VERIFY(NQ1D <= MAX_NQ1D, "");
111  MFEM_VERIFY(VDIM == 2 || !(eval_flags & DETERMINANTS), "");
112  auto B = Reshape(maps.B.Read(), NQ1D, ND1D);
113  auto G = Reshape(maps.G.Read(), NQ1D, ND1D);
114  auto F = Reshape(f_vec.Read(), ND1D, VDIM, NF);
115  auto sign = signs.Read();
116  auto val = Reshape(q_val.Write(), NQ1D, VDIM, NF);
117  // auto der = Reshape(q_der.Write(), NQ1D, VDIM, NF); // only tangential der
118  auto det = Reshape(q_det.Write(), NQ1D, NF);
119  auto n = Reshape(q_nor.Write(), NQ1D, VDIM, NF);
120  MFEM_VERIFY(eval_flags | DERIVATIVES,
121  "Derivatives on the faces are not yet supported.");
122  // If Gauss-Lobatto
123  MFEM_FORALL(f, NF,
124  {
125  const int ND1D = T_ND1D ? T_ND1D : nd1d;
126  const int NQ1D = T_NQ1D ? T_NQ1D : nq1d;
127  const int VDIM = T_VDIM ? T_VDIM : vdim;
128  constexpr int max_ND1D = T_ND1D ? T_ND1D : MAX_ND1D;
129  constexpr int max_VDIM = T_VDIM ? T_VDIM : MAX_VDIM2D;
130  double r_F[max_ND1D][max_VDIM];
131  for (int d = 0; d < ND1D; d++)
132  {
133  for (int c = 0; c < VDIM; c++)
134  {
135  r_F[d][c] = F(d,c,f);
136  }
137  }
138  for (int q = 0; q < NQ1D; ++q)
139  {
140  if (eval_flags & VALUES)
141  {
142  double ed[max_VDIM];
143  for (int c = 0; c < VDIM; c++) { ed[c] = 0.0; }
144  for (int d = 0; d < ND1D; ++d)
145  {
146  const double b = B(q,d);
147  for (int c = 0; c < VDIM; c++) { ed[c] += b*r_F[d][c]; }
148  }
149  for (int c = 0; c < VDIM; c++) { val(q,c,f) = ed[c]; }
150  }
151  if ((eval_flags & DERIVATIVES)
152  || (eval_flags & DETERMINANTS)
153  || (eval_flags & NORMALS))
154  {
155  double D[max_VDIM];
156  for (int i = 0; i < VDIM; i++) { D[i] = 0.0; }
157  for (int d = 0; d < ND1D; ++d)
158  {
159  const double w = G(q,d);
160  for (int c = 0; c < VDIM; c++)
161  {
162  double s_e = r_F[d][c];
163  D[c] += s_e * w;
164  }
165  }
166  if (VDIM == 2 &&
167  ((eval_flags & NORMALS)
168  || (eval_flags & DETERMINANTS)))
169  {
170  const double norm = sqrt(D[0]*D[0]+D[1]*D[1]);
171  if (eval_flags & DETERMINANTS)
172  {
173  det(q,f) = norm;
174  }
175  if (eval_flags & NORMALS)
176  {
177  const double s = sign[f] ? -1.0 : 1.0;
178  n(q,0,f) = s*D[1]/norm;
179  n(q,1,f) = -s*D[0]/norm;
180  }
181  }
182  }
183  }
184  });
185 }
186 
187 template<const int T_VDIM, const int T_ND1D, const int T_NQ1D>
189  const int NF,
190  const int vdim,
191  const DofToQuad &maps,
192  const Array<bool> &signs,
193  const Vector &e_vec,
194  Vector &q_val,
195  Vector &q_der,
196  Vector &q_det,
197  Vector &q_nor,
198  const int eval_flags)
199 {
200  const int nd1d = maps.ndof;
201  const int nq1d = maps.nqpt;
202  const int ND1D = T_ND1D ? T_ND1D : nd1d;
203  const int NQ1D = T_NQ1D ? T_NQ1D : nq1d;
204  const int VDIM = T_VDIM ? T_VDIM : vdim;
205  MFEM_VERIFY(ND1D <= MAX_ND1D, "");
206  MFEM_VERIFY(NQ1D <= MAX_NQ1D, "");
207  MFEM_VERIFY(VDIM == 3 || !(eval_flags & DETERMINANTS), "");
208  auto B = Reshape(maps.B.Read(), NQ1D, ND1D);
209  auto G = Reshape(maps.G.Read(), NQ1D, ND1D);
210  auto F = Reshape(e_vec.Read(), ND1D, ND1D, VDIM, NF);
211  auto sign = signs.Read();
212  auto val = Reshape(q_val.Write(), NQ1D, NQ1D, VDIM, NF);
213  // auto der = Reshape(q_der.Write(), NQ1D, VDIM, 3, NF);
214  auto det = Reshape(q_det.Write(), NQ1D, NQ1D, NF);
215  auto nor = Reshape(q_nor.Write(), NQ1D, NQ1D, 3, NF);
216  MFEM_VERIFY(eval_flags | DERIVATIVES,
217  "Derivatives on the faces are not yet supported.");
218  MFEM_FORALL(f, NF,
219  {
220  const int ND1D = T_ND1D ? T_ND1D : nd1d;
221  const int NQ1D = T_NQ1D ? T_NQ1D : nq1d;
222  const int VDIM = T_VDIM ? T_VDIM : vdim;
223  constexpr int max_ND1D = T_ND1D ? T_ND1D : MAX_ND1D;
224  constexpr int max_NQ1D = T_NQ1D ? T_NQ1D : MAX_NQ1D;
225  constexpr int max_VDIM = T_VDIM ? T_VDIM : MAX_VDIM3D;
226  double r_F[max_ND1D][max_ND1D][max_VDIM];
227  for (int d1 = 0; d1 < ND1D; d1++)
228  {
229  for (int d2 = 0; d2 < ND1D; d2++)
230  {
231  for (int c = 0; c < VDIM; c++)
232  {
233  r_F[d1][d2][c] = F(d1,d2,c,f);
234  }
235  }
236  }
237  if (eval_flags & VALUES)
238  {
239  double Bu[max_NQ1D][max_ND1D][max_VDIM];
240  for (int d2 = 0; d2 < ND1D; ++d2)
241  {
242  for (int q = 0; q < NQ1D; ++q)
243  {
244  for (int c = 0; c < VDIM; c++) { Bu[q][d2][c] = 0.0; }
245  for (int d1 = 0; d1 < ND1D; ++d1)
246  {
247  const double b = B(q,d1);
248  for (int c = 0; c < VDIM; c++)
249  {
250  Bu[q][d2][c] += b*r_F[d1][d2][c];
251  }
252  }
253  }
254  }
255  double BBu[max_NQ1D][max_NQ1D][max_VDIM];
256  for (int q2 = 0; q2 < NQ1D; ++q2)
257  {
258  for (int q1 = 0; q1 < NQ1D; ++q1)
259  {
260  for (int c = 0; c < VDIM; c++) { BBu[q2][q1][c] = 0.0; }
261  for (int d2 = 0; d2 < ND1D; ++d2)
262  {
263  const double b = B(q2,d2);
264  for (int c = 0; c < VDIM; c++)
265  {
266  BBu[q2][q1][c] += b*Bu[q1][d2][c];
267  }
268  }
269  for (int c = 0; c < VDIM; c++)
270  {
271  val(q1,q2,c,f) = BBu[q2][q1][c];
272  }
273  }
274  }
275  }
276  if ((eval_flags & DERIVATIVES)
277  || (eval_flags & DETERMINANTS)
278  || (eval_flags & NORMALS))
279  {
280  // We only compute the tangential derivatives
281  double Gu[max_NQ1D][max_ND1D][max_VDIM];
282  double Bu[max_NQ1D][max_ND1D][max_VDIM];
283  for (int d2 = 0; d2 < ND1D; ++d2)
284  {
285  for (int q = 0; q < NQ1D; ++q)
286  {
287  for (int c = 0; c < VDIM; c++)
288  {
289  Gu[q][d2][c] = 0.0;
290  Bu[q][d2][c] = 0.0;
291  }
292  for (int d1 = 0; d1 < ND1D; ++d1)
293  {
294  const double b = B(q,d1);
295  const double g = G(q,d1);
296  for (int c = 0; c < VDIM; c++)
297  {
298  const double u = r_F[d1][d2][c];
299  Gu[q][d2][c] += g*u;
300  Bu[q][d2][c] += b*u;
301  }
302  }
303  }
304  }
305  double BGu[max_NQ1D][max_NQ1D][max_VDIM];
306  double GBu[max_NQ1D][max_NQ1D][max_VDIM];
307  for (int q2 = 0; q2 < NQ1D; ++q2)
308  {
309  for (int q1 = 0; q1 < NQ1D; ++q1)
310  {
311  for (int c = 0; c < VDIM; c++)
312  {
313  BGu[q2][q1][c] = 0.0;
314  GBu[q2][q1][c] = 0.0;
315  }
316  for (int d2 = 0; d2 < ND1D; ++d2)
317  {
318  const double b = B(q2,d2);
319  const double g = G(q2,d2);
320  for (int c = 0; c < VDIM; c++)
321  {
322  BGu[q2][q1][c] += b*Gu[q1][d2][c];
323  GBu[q2][q1][c] += g*Bu[q1][d2][c];
324  }
325  }
326  }
327  }
328  if (VDIM == 3 && ((eval_flags & NORMALS) ||
329  (eval_flags & DETERMINANTS)))
330  {
331  double n[3];
332  for (int q2 = 0; q2 < NQ1D; ++q2)
333  {
334  for (int q1 = 0; q1 < NQ1D; ++q1)
335  {
336  const double s = sign[f] ? -1.0 : 1.0;
337  n[0] = s*( BGu[q2][q1][1]*GBu[q2][q1][2]-GBu[q2][q1][1]*
338  BGu[q2][q1][2] );
339  n[1] = s*(-BGu[q2][q1][0]*GBu[q2][q1][2]+GBu[q2][q1][0]*
340  BGu[q2][q1][2] );
341  n[2] = s*( BGu[q2][q1][0]*GBu[q2][q1][1]-GBu[q2][q1][0]*
342  BGu[q2][q1][1] );
343  const double norm = sqrt(n[0]*n[0]+n[1]*n[1]+n[2]*n[2]);
344  if (eval_flags & DETERMINANTS) { det(q1,q2,f) = norm; }
345  if (eval_flags & NORMALS)
346  {
347  nor(q1,q2,0,f) = n[0]/norm;
348  nor(q1,q2,1,f) = n[1]/norm;
349  nor(q1,q2,2,f) = n[2]/norm;
350  }
351  }
352  }
353  }
354  }
355  });
356 }
357 
359  const Vector &e_vec, unsigned eval_flags,
360  Vector &q_val, Vector &q_der, Vector &q_det, Vector &q_nor) const
361 {
362  if (nf == 0) { return; }
363  const int vdim = fespace->GetVDim();
364  const int dim = fespace->GetMesh()->Dimension();
365  const FiniteElement *fe =
367  const IntegrationRule *ir = IntRule;
368  const DofToQuad &maps = fe->GetDofToQuad(*ir, DofToQuad::TENSOR);
369  const int nd1d = maps.ndof;
370  const int nq1d = maps.nqpt;
371  void (*eval_func)(
372  const int NF,
373  const int vdim,
374  const DofToQuad &maps,
375  const Array<bool> &signs,
376  const Vector &e_vec,
377  Vector &q_val,
378  Vector &q_der,
379  Vector &q_det,
380  Vector &q_nor,
381  const int eval_flags) = NULL;
382  if (vdim == 1)
383  {
384  if (dim == 2)
385  {
386  switch (10*nd1d + nq1d)
387  {
388  // Q0
389  case 11: eval_func = &Eval2D<1,1,1>; break;
390  case 12: eval_func = &Eval2D<1,1,2>; break;
391  // Q1
392  case 22: eval_func = &Eval2D<1,2,2>; break;
393  case 23: eval_func = &Eval2D<1,2,3>; break;
394  // Q2
395  case 33: eval_func = &Eval2D<1,3,3>; break;
396  case 34: eval_func = &Eval2D<1,3,4>; break;
397  // Q3
398  case 44: eval_func = &Eval2D<1,4,4>; break;
399  case 45: eval_func = &Eval2D<1,4,5>; break;
400  case 46: eval_func = &Eval2D<1,4,6>; break;
401  // Q4
402  case 55: eval_func = &Eval2D<1,5,5>; break;
403  case 56: eval_func = &Eval2D<1,5,6>; break;
404  case 57: eval_func = &Eval2D<1,5,7>; break;
405  case 58: eval_func = &Eval2D<1,5,8>; break;
406  }
407  if (nq1d >= 10 || !eval_func)
408  {
409  eval_func = &Eval2D<1>;
410  }
411  }
412  else if (dim == 3)
413  {
414  switch (10*nd1d + nq1d)
415  {
416  // Q0
417  case 11: eval_func = &Eval3D<1,1,1>; break;
418  case 12: eval_func = &Eval3D<1,1,2>; break;
419  // Q1
420  case 22: eval_func = &Eval3D<1,2,2>; break;
421  case 23: eval_func = &Eval3D<1,2,3>; break;
422  // Q2
423  case 33: eval_func = &Eval3D<1,3,3>; break;
424  case 34: eval_func = &Eval3D<1,3,4>; break;
425  // Q3
426  case 44: eval_func = &Eval3D<1,4,4>; break;
427  case 45: eval_func = &Eval3D<1,4,5>; break;
428  case 46: eval_func = &Eval3D<1,4,6>; break;
429  // Q4
430  case 55: eval_func = &Eval3D<1,5,5>; break;
431  case 56: eval_func = &Eval3D<1,5,6>; break;
432  }
433  if (nq1d >= 10 || !eval_func)
434  {
435  eval_func = &Eval3D<1>;
436  }
437  }
438  }
439  else if (vdim == dim)
440  {
441  if (dim == 2)
442  {
443  switch (10*nd1d + nq1d)
444  {
445  // Q1
446  case 22: eval_func = &Eval2D<2,2,2>; break;
447  case 23: eval_func = &Eval2D<2,2,3>; break;
448  // Q2
449  case 33: eval_func = &Eval2D<2,3,3>; break;
450  case 34: eval_func = &Eval2D<2,3,4>; break;
451  // Q3
452  case 44: eval_func = &Eval2D<2,4,4>; break;
453  case 45: eval_func = &Eval2D<2,4,5>; break;
454  case 46: eval_func = &Eval2D<2,4,6>; break;
455  // Q4
456  case 55: eval_func = &Eval2D<2,5,5>; break;
457  case 56: eval_func = &Eval2D<2,5,6>; break;
458  case 57: eval_func = &Eval2D<2,5,7>; break;
459  case 58: eval_func = &Eval2D<2,5,8>; break;
460  }
461  if (nq1d >= 10 || !eval_func)
462  {
463  eval_func = &Eval2D<2>;
464  }
465  }
466  else if (dim == 3)
467  {
468  switch (10*nd1d + nq1d)
469  {
470  // Q1
471  case 22: eval_func = &Eval3D<3,2,2>; break;
472  case 23: eval_func = &Eval3D<3,2,3>; break;
473  // Q2
474  case 33: eval_func = &Eval3D<3,3,3>; break;
475  case 34: eval_func = &Eval3D<3,3,4>; break;
476  // Q3
477  case 44: eval_func = &Eval3D<3,4,4>; break;
478  case 45: eval_func = &Eval3D<3,4,5>; break;
479  case 46: eval_func = &Eval3D<3,4,6>; break;
480  // Q4
481  case 55: eval_func = &Eval3D<3,5,5>; break;
482  case 56: eval_func = &Eval3D<3,5,6>; break;
483  }
484  if (nq1d >= 10 || !eval_func)
485  {
486  eval_func = &Eval3D<3>;
487  }
488  }
489  }
490  if (eval_func)
491  {
492  eval_func(nf, vdim, maps, signs, e_vec,
493  q_val, q_der, q_det, q_nor, eval_flags);
494  }
495  else
496  {
497  MFEM_ABORT("case not supported yet");
498  }
499 }
500 
502  const Vector &e_vec, Vector &q_val) const
503 {
504  Vector q_der, q_det, q_nor;
505  Mult(e_vec, VALUES, q_val, q_der, q_det, q_nor);
506 }
507 
508 } // namespace mfem
Abstract class for all finite elements.
Definition: fe_base.hpp:235
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
static void Eval3D(const int NF, const int vdim, const DofToQuad &maps, const Array< bool > &signs, const Vector &e_vec, Vector &q_val, Vector &q_der, Vector &q_det, Vector &q_nor, const int eval_flags)
Template compute kernel for 3D.
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition: fe_base.hpp:162
FaceQuadratureInterpolator(const FiniteElementSpace &fes, const IntegrationRule &ir, FaceType type)
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element &#39;i&#39; to the given &#39;geom_type&#39;.
Definition: fespace.cpp:3175
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition: fe_base.hpp:174
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
Definition: fe_base.hpp:170
Geometry::Type GetFaceBaseGeometry(int i) const
Definition: mesh.hpp:1059
const IntegrationRule * IntRule
Not owned.
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition: dtensor.hpp:131
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:590
Class for finite elements with basis functions that return scalar values.
Definition: fe_base.hpp:629
int GetBasisType() const
Definition: fe_base.hpp:1173
FaceType
Definition: mesh.hpp:45
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:433
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:446
double b
Definition: lissajous.cpp:42
void Mult(const Vector &e_vec, unsigned eval_flags, Vector &q_val, Vector &q_der, Vector &q_det, Vector &q_nor) const
Interpolate the E-vector e_vec to quadrature points.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:304
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:557
int Dimension() const
Definition: mesh.hpp:999
void Values(const Vector &e_vec, Vector &q_val) const
Interpolate the values of the E-vector e_vec at quadrature points.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Definition: fe_base.hpp:136
Evaluate the values at quadrature points.
Array< double > B
Basis functions evaluated at quadrature points.
Definition: fe_base.hpp:185
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
Definition: fe_base.cpp:366
FDualNumber< tbase > sqrt(const FDualNumber< tbase > &f)
sqrt([dual number])
Definition: fdual.hpp:600
int dim
Definition: ex24.cpp:53
Array< double > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Definition: fe_base.hpp:206
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2783
static void Eval2D(const int NF, const int vdim, const DofToQuad &maps, const Array< bool > &signs, const Vector &e_vec, Vector &q_val, Vector &q_der, Vector &q_det, Vector &q_nor, const int eval_flags)
Template compute kernel for 2D.
Bernstein polynomials.
Definition: fe_base.hpp:32
Evaluate the derivatives at quadrature points.
Vector data type.
Definition: vector.hpp:60
void f_vec(const Vector &xvec, Vector &f)
Definition: lor_mms.hpp:68
RefCoord s[3]
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
Assuming the derivative at quadrature points form a matrix, this flag can be used to compute and stor...
const FiniteElementSpace * fespace
Not owned.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:438