MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
transfer.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, 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 #include "transfer.hpp"
13 #include "../general/forall.hpp"
14 
15 namespace mfem
16 {
17 
19  const FiniteElementSpace& hFESpace_)
20  : Operator(hFESpace_.GetVSize(), lFESpace_.GetVSize())
21 {
22  if (lFESpace_.FEColl() == hFESpace_.FEColl())
23  {
25  hFESpace_.GetTransferOperator(lFESpace_, P);
26  P.SetOperatorOwner(false);
27  opr = P.Ptr();
28  }
29  else if (lFESpace_.GetMesh()->GetNE() > 0
30  && hFESpace_.GetMesh()->GetNE() > 0
31  && dynamic_cast<const TensorBasisElement*>(lFESpace_.GetFE(0))
32  && dynamic_cast<const TensorBasisElement*>(hFESpace_.GetFE(0)))
33  {
34  opr = new TensorProductPRefinementTransferOperator(lFESpace_, hFESpace_);
35  }
36  else
37  {
38  opr = new PRefinementTransferOperator(lFESpace_, hFESpace_);
39  }
40 }
41 
43 
44 void TransferOperator::Mult(const Vector& x, Vector& y) const
45 {
46  opr->Mult(x, y);
47 }
48 
50 {
51  opr->MultTranspose(x, y);
52 }
53 
55  const FiniteElementSpace& lFESpace_, const FiniteElementSpace& hFESpace_)
56  : Operator(hFESpace_.GetVSize(), lFESpace_.GetVSize()), lFESpace(lFESpace_),
57  hFESpace(hFESpace_)
58 {
59 }
60 
62 
64 {
65  Mesh* mesh = hFESpace.GetMesh();
66  Array<int> l_dofs, h_dofs, l_vdofs, h_vdofs;
67  DenseMatrix loc_prol;
68  Vector subY, subX;
69 
70  Geometry::Type cached_geom = Geometry::INVALID;
71  const FiniteElement* h_fe = NULL;
72  const FiniteElement* l_fe = NULL;
74 
75  int vdim = lFESpace.GetVDim();
76 
77  for (int i = 0; i < mesh->GetNE(); i++)
78  {
79  hFESpace.GetElementDofs(i, h_dofs);
80  lFESpace.GetElementDofs(i, l_dofs);
81 
83  if (geom != cached_geom)
84  {
85  h_fe = hFESpace.GetFE(i);
86  l_fe = lFESpace.GetFE(i);
88  h_fe->GetTransferMatrix(*l_fe, T, loc_prol);
89  subY.SetSize(loc_prol.Height());
90  cached_geom = geom;
91  }
92 
93  for (int vd = 0; vd < vdim; vd++)
94  {
95  l_dofs.Copy(l_vdofs);
96  lFESpace.DofsToVDofs(vd, l_vdofs);
97  h_dofs.Copy(h_vdofs);
98  hFESpace.DofsToVDofs(vd, h_vdofs);
99  x.GetSubVector(l_vdofs, subX);
100  loc_prol.Mult(subX, subY);
101  y.SetSubVector(h_vdofs, subY);
102  }
103  }
104 }
105 
107  Vector& y) const
108 {
109  y = 0.0;
110 
111  Mesh* mesh = hFESpace.GetMesh();
112  Array<int> l_dofs, h_dofs, l_vdofs, h_vdofs;
113  DenseMatrix loc_prol;
114  Vector subY, subX;
115 
116  Array<char> processed(hFESpace.GetVSize());
117  processed = 0;
118 
119  Geometry::Type cached_geom = Geometry::INVALID;
120  const FiniteElement* h_fe = NULL;
121  const FiniteElement* l_fe = NULL;
123 
124  int vdim = lFESpace.GetVDim();
125 
126  for (int i = 0; i < mesh->GetNE(); i++)
127  {
128  hFESpace.GetElementDofs(i, h_dofs);
129  lFESpace.GetElementDofs(i, l_dofs);
130 
131  const Geometry::Type geom = mesh->GetElementBaseGeometry(i);
132  if (geom != cached_geom)
133  {
134  h_fe = hFESpace.GetFE(i);
135  l_fe = lFESpace.GetFE(i);
137  h_fe->GetTransferMatrix(*l_fe, T, loc_prol);
138  loc_prol.Transpose();
139  subY.SetSize(loc_prol.Height());
140  cached_geom = geom;
141  }
142 
143  for (int vd = 0; vd < vdim; vd++)
144  {
145  l_dofs.Copy(l_vdofs);
146  lFESpace.DofsToVDofs(vd, l_vdofs);
147  h_dofs.Copy(h_vdofs);
148  hFESpace.DofsToVDofs(vd, h_vdofs);
149 
150  x.GetSubVector(h_vdofs, subX);
151  for (int p = 0; p < h_dofs.Size(); ++p)
152  {
153  if (processed[lFESpace.DecodeDof(h_dofs[p])])
154  {
155  subX[p] = 0.0;
156  }
157  }
158 
159  loc_prol.Mult(subX, subY);
160  y.AddElementVector(l_vdofs, subY);
161  }
162 
163  for (int p = 0; p < h_dofs.Size(); ++p)
164  {
165  processed[lFESpace.DecodeDof(h_dofs[p])] = 1;
166  }
167  }
168 }
169 
172  const FiniteElementSpace& lFESpace_,
173  const FiniteElementSpace& hFESpace_)
174  : Operator(hFESpace_.GetVSize(), lFESpace_.GetVSize()), lFESpace(lFESpace_),
175  hFESpace(hFESpace_)
176 {
177  // Assuming the same element type
178  Mesh* mesh = lFESpace.GetMesh();
179  dim = mesh->Dimension();
180  if (mesh->GetNE() == 0)
181  {
182  return;
183  }
184  const FiniteElement& el = *lFESpace.GetFE(0);
185 
186  const TensorBasisElement* ltel =
187  dynamic_cast<const TensorBasisElement*>(&el);
188  MFEM_VERIFY(ltel, "Low order FE space must be tensor product space");
189 
190  const TensorBasisElement* htel =
191  dynamic_cast<const TensorBasisElement*>(hFESpace.GetFE(0));
192  MFEM_VERIFY(htel, "High order FE space must be tensor product space");
193  const Array<int>& hdofmap = htel->GetDofMap();
194 
195  const IntegrationRule& ir = hFESpace.GetFE(0)->GetNodes();
196  IntegrationRule irLex = ir;
197 
198  // The quadrature points, or equivalently, the dofs of the high order space
199  // must be sorted in lexicographical order
200  for (int i = 0; i < ir.GetNPoints(); ++i)
201  {
202  irLex.IntPoint(i) = ir.IntPoint(hdofmap[i]);
203  }
204 
205  NE = lFESpace.GetNE();
206  const DofToQuad& maps = el.GetDofToQuad(irLex, DofToQuad::TENSOR);
207 
208  D1D = maps.ndof;
209  Q1D = maps.nqpt;
210  B = maps.B;
211  Bt = maps.Bt;
212 
213  elem_restrict_lex_l =
215 
216  MFEM_VERIFY(elem_restrict_lex_l,
217  "Low order ElementRestriction not available");
218 
219  elem_restrict_lex_h =
221 
222  MFEM_VERIFY(elem_restrict_lex_h,
223  "High order ElementRestriction not available");
224 
225  localL.SetSize(elem_restrict_lex_l->Height(), Device::GetMemoryType());
226  localH.SetSize(elem_restrict_lex_h->Height(), Device::GetMemoryType());
227  localL.UseDevice(true);
228  localH.UseDevice(true);
229 
230  MFEM_VERIFY(dynamic_cast<const ElementRestriction*>(elem_restrict_lex_h),
231  "High order element restriction is of unsupported type");
232 
233  mask.SetSize(localH.Size(), Device::GetMemoryType());
234  static_cast<const ElementRestriction*>(elem_restrict_lex_h)
235  ->BooleanMask(mask);
236  mask.UseDevice(true);
237 }
238 
239 namespace TransferKernels
240 {
241 void Prolongation2D(const int NE, const int D1D, const int Q1D,
242  const Vector& localL, Vector& localH,
243  const Array<double>& B, const Vector& mask)
244 {
245  auto x_ = Reshape(localL.Read(), D1D, D1D, NE);
246  auto y_ = Reshape(localH.ReadWrite(), Q1D, Q1D, NE);
247  auto B_ = Reshape(B.Read(), Q1D, D1D);
248  auto m_ = Reshape(mask.Read(), Q1D, Q1D, NE);
249 
250  localH = 0.0;
251 
252  MFEM_FORALL(e, NE,
253  {
254  for (int dy = 0; dy < D1D; ++dy)
255  {
256  double sol_x[MAX_Q1D];
257  for (int qy = 0; qy < Q1D; ++qy)
258  {
259  sol_x[qy] = 0.0;
260  }
261  for (int dx = 0; dx < D1D; ++dx)
262  {
263  const double s = x_(dx, dy, e);
264  for (int qx = 0; qx < Q1D; ++qx)
265  {
266  sol_x[qx] += B_(qx, dx) * s;
267  }
268  }
269  for (int qy = 0; qy < Q1D; ++qy)
270  {
271  const double d2q = B_(qy, dy);
272  for (int qx = 0; qx < Q1D; ++qx)
273  {
274  y_(qx, qy, e) += d2q * sol_x[qx];
275  }
276  }
277  }
278  for (int qy = 0; qy < Q1D; ++qy)
279  {
280  for (int qx = 0; qx < Q1D; ++qx)
281  {
282  y_(qx, qy, e) *= m_(qx, qy, e);
283  }
284  }
285  });
286 }
287 
288 void Prolongation3D(const int NE, const int D1D, const int Q1D,
289  const Vector& localL, Vector& localH,
290  const Array<double>& B, const Vector& mask)
291 {
292  auto x_ = Reshape(localL.Read(), D1D, D1D, D1D, NE);
293  auto y_ = Reshape(localH.ReadWrite(), Q1D, Q1D, Q1D, NE);
294  auto B_ = Reshape(B.Read(), Q1D, D1D);
295  auto m_ = Reshape(mask.Read(), Q1D, Q1D, Q1D, NE);
296 
297  localH = 0.0;
298 
299  MFEM_FORALL(e, NE,
300  {
301  for (int dz = 0; dz < D1D; ++dz)
302  {
303  double sol_xy[MAX_Q1D][MAX_Q1D];
304  for (int qy = 0; qy < Q1D; ++qy)
305  {
306  for (int qx = 0; qx < Q1D; ++qx)
307  {
308  sol_xy[qy][qx] = 0.0;
309  }
310  }
311  for (int dy = 0; dy < D1D; ++dy)
312  {
313  double sol_x[MAX_Q1D];
314  for (int qx = 0; qx < Q1D; ++qx)
315  {
316  sol_x[qx] = 0;
317  }
318  for (int dx = 0; dx < D1D; ++dx)
319  {
320  const double s = x_(dx, dy, dz, e);
321  for (int qx = 0; qx < Q1D; ++qx)
322  {
323  sol_x[qx] += B_(qx, dx) * s;
324  }
325  }
326  for (int qy = 0; qy < Q1D; ++qy)
327  {
328  const double wy = B_(qy, dy);
329  for (int qx = 0; qx < Q1D; ++qx)
330  {
331  sol_xy[qy][qx] += wy * sol_x[qx];
332  }
333  }
334  }
335  for (int qz = 0; qz < Q1D; ++qz)
336  {
337  const double wz = B_(qz, dz);
338  for (int qy = 0; qy < Q1D; ++qy)
339  {
340  for (int qx = 0; qx < Q1D; ++qx)
341  {
342  y_(qx, qy, qz, e) += wz * sol_xy[qy][qx];
343  }
344  }
345  }
346  }
347  for (int qz = 0; qz < Q1D; ++qz)
348  {
349  for (int qy = 0; qy < Q1D; ++qy)
350  {
351  for (int qx = 0; qx < Q1D; ++qx)
352  {
353  y_(qx, qy, qz, e) *= m_(qx, qy, qz, e);
354  }
355  }
356  }
357  });
358 }
359 
360 void Restriction2D(const int NE, const int D1D, const int Q1D,
361  const Vector& localH, Vector& localL,
362  const Array<double>& Bt, const Vector& mask)
363 {
364  auto x_ = Reshape(localH.Read(), Q1D, Q1D, NE);
365  auto y_ = Reshape(localL.ReadWrite(), D1D, D1D, NE);
366  auto Bt_ = Reshape(Bt.Read(), D1D, Q1D);
367  auto m_ = Reshape(mask.Read(), Q1D, Q1D, NE);
368 
369  localL = 0.0;
370 
371  MFEM_FORALL(e, NE,
372  {
373  for (int qy = 0; qy < Q1D; ++qy)
374  {
375  double sol_x[MAX_D1D];
376  for (int dx = 0; dx < D1D; ++dx)
377  {
378  sol_x[dx] = 0.0;
379  }
380  for (int qx = 0; qx < Q1D; ++qx)
381  {
382  const double s = m_(qx, qy, e) * x_(qx, qy, e);
383  for (int dx = 0; dx < D1D; ++dx)
384  {
385  sol_x[dx] += Bt_(dx, qx) * s;
386  }
387  }
388  for (int dy = 0; dy < D1D; ++dy)
389  {
390  const double q2d = Bt_(dy, qy);
391  for (int dx = 0; dx < D1D; ++dx)
392  {
393  y_(dx, dy, e) += q2d * sol_x[dx];
394  }
395  }
396  }
397  });
398 }
399 void Restriction3D(const int NE, const int D1D, const int Q1D,
400  const Vector& localH, Vector& localL,
401  const Array<double>& Bt, const Vector& mask)
402 {
403  auto x_ = Reshape(localH.Read(), Q1D, Q1D, Q1D, NE);
404  auto y_ = Reshape(localL.ReadWrite(), D1D, D1D, D1D, NE);
405  auto Bt_ = Reshape(Bt.Read(), D1D, Q1D);
406  auto m_ = Reshape(mask.Read(), Q1D, Q1D, Q1D, NE);
407 
408  localL = 0.0;
409 
410  MFEM_FORALL(e, NE,
411  {
412  for (int qz = 0; qz < Q1D; ++qz)
413  {
414  double sol_xy[MAX_D1D][MAX_D1D];
415  for (int dy = 0; dy < D1D; ++dy)
416  {
417  for (int dx = 0; dx < D1D; ++dx)
418  {
419  sol_xy[dy][dx] = 0;
420  }
421  }
422  for (int qy = 0; qy < Q1D; ++qy)
423  {
424  double sol_x[MAX_D1D];
425  for (int dx = 0; dx < D1D; ++dx)
426  {
427  sol_x[dx] = 0;
428  }
429  for (int qx = 0; qx < Q1D; ++qx)
430  {
431  const double s = m_(qx, qy, qz, e) * x_(qx, qy, qz, e);
432  for (int dx = 0; dx < D1D; ++dx)
433  {
434  sol_x[dx] += Bt_(dx, qx) * s;
435  }
436  }
437  for (int dy = 0; dy < D1D; ++dy)
438  {
439  const double wy = Bt_(dy, qy);
440  for (int dx = 0; dx < D1D; ++dx)
441  {
442  sol_xy[dy][dx] += wy * sol_x[dx];
443  }
444  }
445  }
446  for (int dz = 0; dz < D1D; ++dz)
447  {
448  const double wz = Bt_(dz, qz);
449  for (int dy = 0; dy < D1D; ++dy)
450  {
451  for (int dx = 0; dx < D1D; ++dx)
452  {
453  y_(dx, dy, dz, e) += wz * sol_xy[dy][dx];
454  }
455  }
456  }
457  }
458  });
459 }
460 } // namespace TransferKernels
461 
464 {
465 }
466 
468  Vector& y) const
469 {
470  if (lFESpace.GetMesh()->GetNE() == 0)
471  {
472  return;
473  }
474 
475  elem_restrict_lex_l->Mult(x, localL);
476  if (dim == 2)
477  {
478  TransferKernels::Prolongation2D(NE, D1D, Q1D, localL, localH, B, mask);
479  }
480  else if (dim == 3)
481  {
482  TransferKernels::Prolongation3D(NE, D1D, Q1D, localL, localH, B, mask);
483  }
484  else
485  {
486  MFEM_ABORT("TensorProductPRefinementTransferOperator::Mult not "
487  "implemented for dim = "
488  << dim);
489  }
490  elem_restrict_lex_h->MultTranspose(localH, y);
491 }
492 
494  Vector& y) const
495 {
496  if (lFESpace.GetMesh()->GetNE() == 0)
497  {
498  return;
499  }
500 
501  elem_restrict_lex_h->Mult(x, localH);
502  if (dim == 2)
503  {
504  TransferKernels::Restriction2D(NE, D1D, Q1D, localH, localL, Bt, mask);
505  }
506  else if (dim == 3)
507  {
508  TransferKernels::Restriction3D(NE, D1D, Q1D, localH, localL, Bt, mask);
509  }
510  else
511  {
512  MFEM_ABORT("TensorProductPRefinementTransferOperator::MultTranspose not "
513  "implemented for dim = "
514  << dim);
515  }
516  elem_restrict_lex_l->MultTranspose(localL, y);
517 }
518 
519 #ifdef MFEM_USE_MPI
521  ParFiniteElementSpace& lFESpace_,
522  const ParFiniteElementSpace& hFESpace_)
523  : Operator(hFESpace_.GetTrueVSize(), lFESpace_.GetTrueVSize()),
524  lFESpace(lFESpace_),
525  hFESpace(hFESpace_)
526 {
527  localTransferOperator = new TransferOperator(lFESpace_, hFESpace_);
528 
529  tmpL.SetSize(lFESpace_.GetVSize());
530  tmpH.SetSize(hFESpace_.GetVSize());
531 
532  hFESpace.GetRestrictionMatrix()->BuildTranspose();
533 }
534 
536 {
537  delete localTransferOperator;
538 }
539 
540 void TrueTransferOperator::Mult(const Vector& x, Vector& y) const
541 {
542  lFESpace.GetProlongationMatrix()->Mult(x, tmpL);
543  localTransferOperator->Mult(tmpL, tmpH);
544  hFESpace.GetRestrictionMatrix()->Mult(tmpH, y);
545 }
546 
548 {
549  hFESpace.GetRestrictionMatrix()->MultTranspose(x, tmpH);
550  localTransferOperator->MultTranspose(tmpH, tmpL);
551  lFESpace.GetProlongationMatrix()->MultTranspose(tmpL, y);
552 }
553 #endif
554 
555 } // namespace mfem
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:245
Abstract class for all finite elements.
Definition: fe.hpp:235
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:521
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:400
void BuildTranspose() const
Build and store internally the transpose of this matrix which will be used in the methods AddMultTran...
Definition: sparsemat.cpp:773
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
Definition: transfer.cpp:493
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition: fe.hpp:157
virtual ~TensorProductPRefinementTransferOperator()
Destructor.
Definition: transfer.cpp:463
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:895
Matrix-free transfer operator between finite element spaces on the same mesh.
Definition: transfer.hpp:57
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:459
const Geometry::Type geom
Definition: ex1.cpp:40
static int DecodeDof(int dof)
Helpers to remove encoded sign from a DOF.
Definition: fespace.hpp:173
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:89
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
Definition: transfer.cpp:547
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition: array.hpp:831
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:495
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition: fe.hpp:169
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
virtual ~TransferOperator()
Destructor.
Definition: transfer.cpp:42
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
Definition: fe.hpp:165
~TrueTransferOperator()
Destructor.
Definition: transfer.cpp:535
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a vector x corresponding to the coarse space to the vector y corresp...
Definition: transfer.cpp:44
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
Abstract parallel finite element space.
Definition: pfespace.hpp:28
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:834
virtual ~PRefinementTransferOperator()
Destructor.
Definition: transfer.cpp:61
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: operator.hpp:92
const Operator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:894
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition: dtensor.hpp:134
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe.hpp:312
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:427
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:248
ID for the base class Operator, i.e. any type.
Definition: operator.hpp:239
double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:388
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
Definition: transfer.cpp:49
const int MAX_Q1D
Definition: forall.hpp:27
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:65
TrueTransferOperator(const ParFiniteElementSpace &lFESpace_, const ParFiniteElementSpace &hFESpace_)
Constructs a transfer operator working on true degrees of freedom from from lFESpace to hFESpace...
Definition: transfer.cpp:520
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:314
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
Definition: device.hpp:265
TensorProductPRefinementTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator from lFESpace to hFESpace which have different FE collections.
Definition: transfer.cpp:171
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:1696
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition: fe.hpp:382
void Restriction2D(const int NE, const int D1D, const int Q1D, const Vector &localH, Vector &localL, const Array< double > &Bt, const Vector &mask)
Definition: transfer.cpp:360
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:290
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:389
int Dimension() const
Definition: mesh.hpp:788
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:587
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
void GetTransferOperator(const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
Construct and return an Operator that can be used to transfer GridFunction data from coarse_fes...
Definition: fespace.cpp:2167
Array< double > Bt
Transpose of B.
Definition: fe.hpp:186
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
Definition: transfer.cpp:106
void Transpose()
(*this) = (*this)^t
Definition: densemat.cpp:1378
void SetIdentityTransformation(Geometry::Type GeomType)
Set the FiniteElement Geometry for the reference elements being used.
Definition: eltrans.cpp:371
Matrix-free transfer operator between finite element spaces.
Definition: transfer.hpp:26
Operator that converts FiniteElementSpace L-vectors to E-vectors.
Definition: restriction.hpp:31
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
Matrix-free transfer operator between finite element spaces on the same mesh exploiting the tensor pr...
Definition: transfer.hpp:87
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:594
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a vector x corresponding to the coarse space to the vector y corresp...
Definition: transfer.cpp:467
void Restriction3D(const int NE, const int D1D, const int Q1D, const Vector &localH, Vector &localL, const Array< double > &Bt, const Vector &mask)
Definition: transfer.cpp:399
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Definition: fe.hpp:131
Array< double > B
Basis functions evaluated at quadrature points.
Definition: fe.hpp:180
const Array< int > & GetDofMap() const
Get an Array&lt;int&gt; that maps lexicographically ordered indices to the indices of the respective nodes/...
Definition: fe.hpp:2029
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.cpp:376
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:112
A standard isoparametric element transformation.
Definition: eltrans.hpp:348
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:372
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a vector x corresponding to the coarse space to the vector y corresp...
Definition: transfer.cpp:63
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:339
PRefinementTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator from lFESpace to hFESpace which have different FE collections.
Definition: transfer.cpp:54
const int MAX_D1D
Definition: forall.hpp:26
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:721
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:1798
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a true dof vector x to a true dof vector y.
Definition: transfer.cpp:540
Lexicographic ordering for tensor-product FiniteElements.
Vector data type.
Definition: vector.hpp:51
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:414
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:175
virtual void GetTransferMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this f...
Definition: fe.cpp:123
Abstract operator.
Definition: operator.hpp:24
void Prolongation2D(const int NE, const int D1D, const int Q1D, const Vector &localL, Vector &localH, const Array< double > &B, const Vector &mask)
Definition: transfer.cpp:241
TransferOperator(const FiniteElementSpace &lFESpace, const FiniteElementSpace &hFESpace)
Constructs a transfer operator from lFESpace to hFESpace.
Definition: transfer.cpp:18
void Prolongation3D(const int NE, const int D1D, const int Q1D, const Vector &localL, Vector &localH, const Array< double > &B, const Vector &mask)
Definition: transfer.cpp:288
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Definition: fespace.cpp:108