MFEM  v4.6.0
Finite element discretization library
bilinearform_ext.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 // Implementations of classes FABilinearFormExtension, EABilinearFormExtension,
13 // PABilinearFormExtension and MFBilinearFormExtension.
14 
15 #include "../general/forall.hpp"
16 #include "bilinearform.hpp"
17 #include "pbilinearform.hpp"
18 #include "pgridfunc.hpp"
19 #include "ceed/interface/util.hpp"
20 
21 namespace mfem
22 {
23 
25  : Operator(form->Size()), a(form)
26 {
27  // empty
28 }
29 
31 {
32  return a->GetProlongation();
33 }
34 
36 {
37  return a->GetRestriction();
38 }
39 
40 // Data and methods for partially-assembled bilinear forms
42  : BilinearFormExtension(form),
43  trial_fes(a->FESpace()),
44  test_fes(a->FESpace())
45 {
46  elem_restrict = NULL;
47  int_face_restrict_lex = NULL;
48  bdr_face_restrict_lex = NULL;
49 }
50 
52 {
53  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
54  const int integratorCount = integrators.Size();
55  for (int i = 0; i < integratorCount; ++i)
56  {
57  integrators[i]->AssembleMF(*a->FESpace());
58  }
59 
60  MFEM_VERIFY(a->GetBBFI()->Size() == 0, "AddBoundaryIntegrator is not "
61  "currently supported in MFBilinearFormExtension");
62 }
63 
65 {
66  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
67 
68  const int iSz = integrators.Size();
70  {
71  localY = 0.0;
72  for (int i = 0; i < iSz; ++i)
73  {
74  integrators[i]->AssembleDiagonalMF(localY);
75  }
76  const ElementRestriction* H1elem_restrict =
77  dynamic_cast<const ElementRestriction*>(elem_restrict);
78  if (H1elem_restrict)
79  {
80  H1elem_restrict->MultTransposeUnsigned(localY, y);
81  }
82  else
83  {
85  }
86  }
87  else
88  {
89  y.UseDevice(true); // typically this is a large vector, so store on device
90  y = 0.0;
91  for (int i = 0; i < iSz; ++i)
92  {
93  integrators[i]->AssembleDiagonalMF(y);
94  }
95  }
96 }
97 
99 {
100  FiniteElementSpace *fes = a->FESpace();
101  height = width = fes->GetVSize();
102  trial_fes = fes;
103  test_fes = fes;
104 
105  elem_restrict = nullptr;
106  int_face_restrict_lex = nullptr;
107  bdr_face_restrict_lex = nullptr;
108 }
109 
111  OperatorHandle &A)
112 {
113  Operator *oper;
114  Operator::FormSystemOperator(ess_tdof_list, oper);
115  A.Reset(oper); // A will own oper
116 }
117 
119  Vector &x, Vector &b,
120  OperatorHandle &A,
121  Vector &X, Vector &B,
122  int copy_interior)
123 {
124  Operator *oper;
125  Operator::FormLinearSystem(ess_tdof_list, x, b, oper, X, B, copy_interior);
126  A.Reset(oper); // A will own oper
127 }
128 
130 {
131  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
132 
133  const int iSz = integrators.Size();
135  {
136  y.UseDevice(true); // typically this is a large vector, so store on device
137  y = 0.0;
138  for (int i = 0; i < iSz; ++i)
139  {
140  integrators[i]->AddMultMF(x, y);
141  }
142  }
143  else
144  {
146  localY = 0.0;
147  for (int i = 0; i < iSz; ++i)
148  {
149  integrators[i]->AddMultMF(localX, localY);
150  }
152  }
153 
154  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
155  const int iFISz = intFaceIntegrators.Size();
156  if (int_face_restrict_lex && iFISz>0)
157  {
159  if (int_face_X.Size()>0)
160  {
161  int_face_Y = 0.0;
162  for (int i = 0; i < iFISz; ++i)
163  {
164  intFaceIntegrators[i]->AddMultMF(int_face_X, int_face_Y);
165  }
167  }
168  }
169 
170  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
171  const int bFISz = bdrFaceIntegrators.Size();
172  if (bdr_face_restrict_lex && bFISz>0)
173  {
175  if (bdr_face_X.Size()>0)
176  {
177  bdr_face_Y = 0.0;
178  for (int i = 0; i < bFISz; ++i)
179  {
180  bdrFaceIntegrators[i]->AddMultMF(bdr_face_X, bdr_face_Y);
181  }
183  }
184  }
185 }
186 
188 {
189  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
190  const int iSz = integrators.Size();
191  if (elem_restrict)
192  {
194  localY = 0.0;
195  for (int i = 0; i < iSz; ++i)
196  {
197  integrators[i]->AddMultTransposeMF(localX, localY);
198  }
200  }
201  else
202  {
203  y.UseDevice(true);
204  y = 0.0;
205  for (int i = 0; i < iSz; ++i)
206  {
207  integrators[i]->AddMultTransposeMF(x, y);
208  }
209  }
210 
211  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
212  const int iFISz = intFaceIntegrators.Size();
213  if (int_face_restrict_lex && iFISz>0)
214  {
216  if (int_face_X.Size()>0)
217  {
218  int_face_Y = 0.0;
219  for (int i = 0; i < iFISz; ++i)
220  {
221  intFaceIntegrators[i]->AddMultTransposeMF(int_face_X, int_face_Y);
222  }
224  }
225  }
226 
227  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
228  const int bFISz = bdrFaceIntegrators.Size();
229  if (bdr_face_restrict_lex && bFISz>0)
230  {
232  if (bdr_face_X.Size()>0)
233  {
234  bdr_face_Y = 0.0;
235  for (int i = 0; i < bFISz; ++i)
236  {
237  bdrFaceIntegrators[i]->AddMultTransposeMF(bdr_face_X, bdr_face_Y);
238  }
240  }
241  }
242 }
243 
244 // Data and methods for partially-assembled bilinear forms
246  : BilinearFormExtension(form),
247  trial_fes(a->FESpace()),
248  test_fes(a->FESpace())
249 {
250  elem_restrict = NULL;
251  int_face_restrict_lex = NULL;
252  bdr_face_restrict_lex = NULL;
253 }
254 
256 {
257  if ( Device::Allows(Backend::CEED_MASK) ) { return; }
262  if (elem_restrict)
263  {
266  localY.UseDevice(true); // ensure 'localY = 0.0' is done on device
267 
268  // Gather the attributes on the host from all the elements
269  const Mesh &mesh = *trial_fes->GetMesh();
270  elem_attributes.SetSize(mesh.GetNE());
271  for (int i = 0; i < mesh.GetNE(); ++i)
272  {
273  elem_attributes[i] = mesh.GetAttribute(i);
274  }
275  }
276 
277  // Construct face restriction operators only if the bilinear form has
278  // interior or boundary face integrators
279  if (int_face_restrict_lex == NULL && a->GetFBFI()->Size() > 0)
280  {
286  int_face_Y.UseDevice(true); // ensure 'int_face_Y = 0.0' is done on device
287  }
288 
289  const bool has_bdr_integs = (a->GetBFBFI()->Size() > 0 ||
290  a->GetBBFI()->Size() > 0);
291  if (bdr_face_restrict_lex == NULL && has_bdr_integs)
292  {
296  m);
299  bdr_face_Y.UseDevice(true); // ensure 'faceBoundY = 0.0' is done on device
300 
301  const Mesh &mesh = *trial_fes->GetMesh();
302  // See LinearFormExtension::Update for explanation of f_to_be logic.
303  std::unordered_map<int,int> f_to_be;
304  for (int i = 0; i < mesh.GetNBE(); ++i)
305  {
306  const int f = mesh.GetBdrElementEdgeIndex(i);
307  f_to_be[f] = i;
308  }
309  const int nf_bdr = trial_fes->GetNFbyType(FaceType::Boundary);
310  bdr_attributes.SetSize(nf_bdr);
311  int f_ind = 0;
312  int missing_bdr_elems = 0;
313  for (int f = 0; f < mesh.GetNumFaces(); ++f)
314  {
316  {
317  continue;
318  }
319  int attribute = 1; // default value
320  if (f_to_be.find(f) != f_to_be.end())
321  {
322  const int be = f_to_be[f];
323  attribute = mesh.GetBdrAttribute(be);
324  }
325  else
326  {
327  // If a boundary face does not correspond to the a boundary element,
328  // we assign it the default attribute of 1. We also generate a
329  // warning at runtime with the number of such missing elements.
330  ++missing_bdr_elems;
331  }
332  bdr_attributes[f_ind] = attribute;
333  ++f_ind;
334  }
335  if (missing_bdr_elems)
336  {
337  MFEM_WARNING("Missing " << missing_bdr_elems << " boundary elements "
338  "for boundary faces.");
339  }
340  }
341 }
342 
344 {
346 
347  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
348  for (BilinearFormIntegrator *integ : integrators)
349  {
350  if (integ->Patchwise())
351  {
352  MFEM_VERIFY(a->FESpace()->GetNURBSext(),
353  "Patchwise integration requires a NURBS FE space");
354  integ->AssembleNURBSPA(*a->FESpace());
355  }
356  else
357  {
358  integ->AssemblePA(*a->FESpace());
359  }
360  }
361 
362  Array<BilinearFormIntegrator*> &bdr_integrators = *a->GetBBFI();
363  for (BilinearFormIntegrator *integ : bdr_integrators)
364  {
365  integ->AssemblePABoundary(*a->FESpace());
366  }
367 
368  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
369  for (BilinearFormIntegrator *integ : intFaceIntegrators)
370  {
371  integ->AssemblePAInteriorFaces(*a->FESpace());
372  }
373 
374  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
375  for (BilinearFormIntegrator *integ : bdrFaceIntegrators)
376  {
377  integ->AssemblePABoundaryFaces(*a->FESpace());
378  }
379 }
380 
382 {
383  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
384 
385  const int iSz = integrators.Size();
387  {
388  if (iSz > 0)
389  {
390  localY = 0.0;
391  for (int i = 0; i < iSz; ++i)
392  {
393  integrators[i]->AssembleDiagonalPA(localY);
394  }
395  const ElementRestriction* H1elem_restrict =
396  dynamic_cast<const ElementRestriction*>(elem_restrict);
397  if (H1elem_restrict)
398  {
399  H1elem_restrict->MultTransposeUnsigned(localY, y);
400  }
401  else
402  {
404  }
405  }
406  else
407  {
408  y = 0.0;
409  }
410  }
411  else
412  {
413  y.UseDevice(true); // typically this is a large vector, so store on device
414  y = 0.0;
415  for (int i = 0; i < iSz; ++i)
416  {
417  integrators[i]->AssembleDiagonalPA(y);
418  }
419  }
420 
421  Array<BilinearFormIntegrator*> &bdr_integs = *a->GetBBFI();
422  const int n_bdr_integs = bdr_integs.Size();
423  if (bdr_face_restrict_lex && n_bdr_integs > 0)
424  {
425  bdr_face_Y = 0.0;
426  for (int i = 0; i < n_bdr_integs; ++i)
427  {
428  bdr_integs[i]->AssembleDiagonalPA(bdr_face_Y);
429  }
431  }
432 }
433 
435 {
436  FiniteElementSpace *fes = a->FESpace();
437  height = width = fes->GetVSize();
438  trial_fes = fes;
439  test_fes = fes;
440 
441  elem_restrict = nullptr;
442  int_face_restrict_lex = nullptr;
443  bdr_face_restrict_lex = nullptr;
444 }
445 
447  OperatorHandle &A)
448 {
449  Operator *oper;
450  Operator::FormSystemOperator(ess_tdof_list, oper);
451  A.Reset(oper); // A will own oper
452 }
453 
455  Vector &x, Vector &b,
456  OperatorHandle &A,
457  Vector &X, Vector &B,
458  int copy_interior)
459 {
460  Operator *oper;
461  Operator::FormLinearSystem(ess_tdof_list, x, b, oper, X, B, copy_interior);
462  A.Reset(oper); // A will own oper
463 }
464 
466 {
467  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
468 
469  const int iSz = integrators.Size();
470 
471  bool allPatchwise = true;
472  bool somePatchwise = false;
473 
474  for (int i = 0; i < iSz; ++i)
475  {
476  if (integrators[i]->Patchwise())
477  {
478  somePatchwise = true;
479  }
480  else
481  {
482  allPatchwise = false;
483  }
484  }
485 
486  MFEM_VERIFY(!(somePatchwise && !allPatchwise),
487  "All or none of the integrators should be patchwise");
488 
489  if (DeviceCanUseCeed() || !elem_restrict || allPatchwise)
490  {
491  y.UseDevice(true); // typically this is a large vector, so store on device
492  y = 0.0;
493  for (int i = 0; i < iSz; ++i)
494  {
495  if (integrators[i]->Patchwise())
496  {
497  integrators[i]->AddMultNURBSPA(x, y);
498  }
499  else
500  {
501  integrators[i]->AddMultPA(x, y);
502  }
503  }
504  }
505  else
506  {
507  if (iSz)
508  {
509  Array<Array<int>*> &elem_markers = *a->GetDBFI_Marker();
511  localY = 0.0;
512  for (int i = 0; i < iSz; ++i)
513  {
514  AddMultWithMarkers(*integrators[i], localX, elem_markers[i], elem_attributes,
515  false, localY);
516  }
518  }
519  else
520  {
521  y = 0.0;
522  }
523  }
524 
525  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
526  const int iFISz = intFaceIntegrators.Size();
527  if (int_face_restrict_lex && iFISz>0)
528  {
530  if (int_face_X.Size()>0)
531  {
532  int_face_Y = 0.0;
533  for (int i = 0; i < iFISz; ++i)
534  {
535  intFaceIntegrators[i]->AddMultPA(int_face_X, int_face_Y);
536  }
538  }
539  }
540 
541  Array<BilinearFormIntegrator*> &bdr_integs = *a->GetBBFI();
542  Array<BilinearFormIntegrator*> &bdr_face_integs = *a->GetBFBFI();
543  const int n_bdr_integs = bdr_integs.Size();
544  const int n_bdr_face_integs = bdr_face_integs.Size();
545  const bool has_bdr_integs = (n_bdr_face_integs > 0 || n_bdr_integs > 0);
546  if (bdr_face_restrict_lex && has_bdr_integs)
547  {
548  Array<Array<int>*> &bdr_markers = *a->GetBBFI_Marker();
549  Array<Array<int>*> &bdr_face_markers = *a->GetBFBFI_Marker();
551  if (bdr_face_X.Size()>0)
552  {
553  bdr_face_Y = 0.0;
554  for (int i = 0; i < n_bdr_integs; ++i)
555  {
556  AddMultWithMarkers(*bdr_integs[i], bdr_face_X, bdr_markers[i], bdr_attributes,
557  false, bdr_face_Y);
558  }
559  for (int i = 0; i < n_bdr_face_integs; ++i)
560  {
561  AddMultWithMarkers(*bdr_face_integs[i], bdr_face_X, bdr_face_markers[i],
562  bdr_attributes, false, bdr_face_Y);
563  }
565  }
566  }
567 }
568 
570 {
571  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
572  const int iSz = integrators.Size();
573  if (elem_restrict)
574  {
575  Array<Array<int>*> &elem_markers = *a->GetDBFI_Marker();
577  localY = 0.0;
578  for (int i = 0; i < iSz; ++i)
579  {
580  AddMultWithMarkers(*integrators[i], localX, elem_markers[i], elem_attributes,
581  true, localY);
582  }
584  }
585  else
586  {
587  y.UseDevice(true);
588  y = 0.0;
589  for (int i = 0; i < iSz; ++i)
590  {
591  integrators[i]->AddMultTransposePA(x, y);
592  }
593  }
594 
595  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
596  const int iFISz = intFaceIntegrators.Size();
597  if (int_face_restrict_lex && iFISz>0)
598  {
600  if (int_face_X.Size()>0)
601  {
602  int_face_Y = 0.0;
603  for (int i = 0; i < iFISz; ++i)
604  {
605  intFaceIntegrators[i]->AddMultTransposePA(int_face_X, int_face_Y);
606  }
608  }
609  }
610 
611  Array<BilinearFormIntegrator*> &bdr_integs = *a->GetBBFI();
612  Array<BilinearFormIntegrator*> &bdr_face_integs = *a->GetBFBFI();
613  const int n_bdr_integs = bdr_integs.Size();
614  const int n_bdr_face_integs = bdr_face_integs.Size();
615  const bool has_bdr_integs = (n_bdr_face_integs > 0 || n_bdr_integs > 0);
616  if (bdr_face_restrict_lex && has_bdr_integs)
617  {
618  Array<Array<int>*> &bdr_markers = *a->GetBBFI_Marker();
619  Array<Array<int>*> &bdr_face_markers = *a->GetBFBFI_Marker();
620 
622  if (bdr_face_X.Size() > 0)
623  {
624  bdr_face_Y = 0.0;
625  for (int i = 0; i < n_bdr_integs; ++i)
626  {
627  AddMultWithMarkers(*bdr_integs[i], bdr_face_X, bdr_markers[i], bdr_attributes,
628  true, bdr_face_Y);
629  }
630  for (int i = 0; i < n_bdr_face_integs; ++i)
631  {
632  AddMultWithMarkers(*bdr_face_integs[i], bdr_face_X, bdr_face_markers[i],
633  bdr_attributes, true, bdr_face_Y);
634  }
636  }
637  }
638 }
639 
640 // Compute kernels for PABilinearFormExtension::AddMultWithMarkers.
641 // Cannot be in member function with non-public visibility.
642 static void AddWithMarkers_(
643  const int ne,
644  const int nd,
645  const Vector &x,
646  const Array<int> &markers,
647  const Array<int> &attributes,
648  Vector &y)
649 {
650  const auto d_x = Reshape(x.Read(), nd, ne);
651  const auto d_m = Reshape(markers.Read(), markers.Size());
652  const auto d_attr = Reshape(attributes.Read(), ne);
653  auto d_y = Reshape(y.ReadWrite(), nd, ne);
654  mfem::forall(ne, [=] MFEM_HOST_DEVICE (int e)
655  {
656  const int attr = d_attr[e];
657  if (d_m[attr - 1] == 0) { return; }
658  for (int i = 0; i < nd; ++i)
659  {
660  d_y(i, e) += d_x(i, e);
661  }
662  });
663 }
664 
666  const BilinearFormIntegrator &integ,
667  const Vector &x,
668  const Array<int> *markers,
669  const Array<int> &attributes,
670  const bool transpose,
671  Vector &y) const
672 {
673  if (markers)
674  {
675  tmp_evec.SetSize(y.Size());
676  tmp_evec = 0.0;
677  if (transpose) { integ.AddMultTransposePA(x, tmp_evec); }
678  else { integ.AddMultPA(x, tmp_evec); }
679  const int ne = attributes.Size();
680  const int nd = x.Size() / ne;
681  AddWithMarkers_(ne, nd, tmp_evec, *markers, attributes, y);
682  }
683  else
684  {
685  if (transpose) { integ.AddMultTransposePA(x, y); }
686  else { integ.AddMultPA(x, y); }
687  }
688 }
689 
690 // Data and methods for element-assembled bilinear forms
692  : PABilinearFormExtension(form),
693  factorize_face_terms(false)
694 {
695  if (form->FESpace()->IsDGSpace() && form->FESpace()->Conforming())
696  {
697  factorize_face_terms = true;
698  }
699 }
700 
702 {
704 
705  ne = trial_fes->GetMesh()->GetNE();
706  elemDofs = trial_fes->GetFE(0)->GetDof();
707 
709  ea_data.UseDevice(true);
710 
711  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
712  const int integratorCount = integrators.Size();
713  if ( integratorCount == 0 )
714  {
715  ea_data = 0.0;
716  }
717  for (int i = 0; i < integratorCount; ++i)
718  {
719  integrators[i]->AssembleEA(*a->FESpace(), ea_data, i);
720  }
721 
722  faceDofs = trial_fes ->
723  GetTraceElement(0, trial_fes->GetMesh()->GetFaceGeometry(0)) ->
724  GetDof();
725 
726  MFEM_VERIFY(a->GetBBFI()->Size() == 0,
727  "Element assembly does not support AddBoundaryIntegrator yet.");
728 
729  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
730  const int intFaceIntegratorCount = intFaceIntegrators.Size();
731  if (intFaceIntegratorCount>0)
732  {
736  }
737  for (int i = 0; i < intFaceIntegratorCount; ++i)
738  {
739  intFaceIntegrators[i]->AssembleEAInteriorFaces(*a->FESpace(),
740  ea_data_int,
741  ea_data_ext,
742  i);
743  }
744 
745  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
746  const int boundFaceIntegratorCount = bdrFaceIntegrators.Size();
747  if (boundFaceIntegratorCount>0)
748  {
751  ea_data_bdr = 0.0;
752  }
753  for (int i = 0; i < boundFaceIntegratorCount; ++i)
754  {
755  bdrFaceIntegrators[i]->AssembleEABoundaryFaces(*a->FESpace(),ea_data_bdr,i);
756  }
757 
759  {
760  auto restFint = dynamic_cast<const L2FaceRestriction*>(int_face_restrict_lex);
762  }
764  {
765  auto restFbdr = dynamic_cast<const L2FaceRestriction*>(bdr_face_restrict_lex);
767  }
768 }
769 
771 {
772  // Apply the Element Restriction
773  const bool useRestrict = !DeviceCanUseCeed() && elem_restrict;
774  if (!useRestrict)
775  {
776  y.UseDevice(true); // typically this is a large vector, so store on device
777  y = 0.0;
778  }
779  else
780  {
782  localY = 0.0;
783  }
784  // Apply the Element Matrices
785  {
786  const int NDOFS = elemDofs;
787  auto X = Reshape(useRestrict?localX.Read():x.Read(), NDOFS, ne);
788  auto Y = Reshape(useRestrict?localY.ReadWrite():y.ReadWrite(), NDOFS, ne);
789  auto A = Reshape(ea_data.Read(), NDOFS, NDOFS, ne);
790  mfem::forall(ne*NDOFS, [=] MFEM_HOST_DEVICE (int glob_j)
791  {
792  const int e = glob_j/NDOFS;
793  const int j = glob_j%NDOFS;
794  double res = 0.0;
795  for (int i = 0; i < NDOFS; i++)
796  {
797  res += A(i, j, e)*X(i, e);
798  }
799  Y(j, e) += res;
800  });
801  // Apply the Element Restriction transposed
802  if (useRestrict)
803  {
805  }
806  }
807 
808  // Treatment of interior faces
809  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
810  const int iFISz = intFaceIntegrators.Size();
811  if (int_face_restrict_lex && iFISz>0)
812  {
813  // Apply the Interior Face Restriction
815  if (int_face_X.Size()>0)
816  {
817  int_face_Y = 0.0;
818  // Apply the interior face matrices
819  const int NDOFS = faceDofs;
820  auto X = Reshape(int_face_X.Read(), NDOFS, 2, nf_int);
821  auto Y = Reshape(int_face_Y.ReadWrite(), NDOFS, 2, nf_int);
823  {
824  auto A_int = Reshape(ea_data_int.Read(), NDOFS, NDOFS, 2, nf_int);
825  mfem::forall(nf_int*NDOFS, [=] MFEM_HOST_DEVICE (int glob_j)
826  {
827  const int f = glob_j/NDOFS;
828  const int j = glob_j%NDOFS;
829  double res = 0.0;
830  for (int i = 0; i < NDOFS; i++)
831  {
832  res += A_int(i, j, 0, f)*X(i, 0, f);
833  }
834  Y(j, 0, f) += res;
835  res = 0.0;
836  for (int i = 0; i < NDOFS; i++)
837  {
838  res += A_int(i, j, 1, f)*X(i, 1, f);
839  }
840  Y(j, 1, f) += res;
841  });
842  }
843  auto A_ext = Reshape(ea_data_ext.Read(), NDOFS, NDOFS, 2, nf_int);
844  mfem::forall(nf_int*NDOFS, [=] MFEM_HOST_DEVICE (int glob_j)
845  {
846  const int f = glob_j/NDOFS;
847  const int j = glob_j%NDOFS;
848  double res = 0.0;
849  for (int i = 0; i < NDOFS; i++)
850  {
851  res += A_ext(i, j, 0, f)*X(i, 0, f);
852  }
853  Y(j, 1, f) += res;
854  res = 0.0;
855  for (int i = 0; i < NDOFS; i++)
856  {
857  res += A_ext(i, j, 1, f)*X(i, 1, f);
858  }
859  Y(j, 0, f) += res;
860  });
861  // Apply the Interior Face Restriction transposed
863  }
864  }
865 
866  // Treatment of boundary faces
867  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
868  const int bFISz = bdrFaceIntegrators.Size();
869  if (!factorize_face_terms && bdr_face_restrict_lex && bFISz>0)
870  {
871  // Apply the Boundary Face Restriction
873  if (bdr_face_X.Size()>0)
874  {
875  bdr_face_Y = 0.0;
876  // Apply the boundary face matrices
877  const int NDOFS = faceDofs;
878  auto X = Reshape(bdr_face_X.Read(), NDOFS, nf_bdr);
879  auto Y = Reshape(bdr_face_Y.ReadWrite(), NDOFS, nf_bdr);
880  auto A = Reshape(ea_data_bdr.Read(), NDOFS, NDOFS, nf_bdr);
881  mfem::forall(nf_bdr*NDOFS, [=] MFEM_HOST_DEVICE (int glob_j)
882  {
883  const int f = glob_j/NDOFS;
884  const int j = glob_j%NDOFS;
885  double res = 0.0;
886  for (int i = 0; i < NDOFS; i++)
887  {
888  res += A(i, j, f)*X(i, f);
889  }
890  Y(j, f) += res;
891  });
892  // Apply the Boundary Face Restriction transposed
894  }
895  }
896 }
897 
899 {
900  // Apply the Element Restriction
901  const bool useRestrict = !DeviceCanUseCeed() && elem_restrict;
902  if (!useRestrict)
903  {
904  y.UseDevice(true); // typically this is a large vector, so store on device
905  y = 0.0;
906  }
907  else
908  {
910  localY = 0.0;
911  }
912  // Apply the Element Matrices transposed
913  {
914  const int NDOFS = elemDofs;
915  auto X = Reshape(useRestrict?localX.Read():x.Read(), NDOFS, ne);
916  auto Y = Reshape(useRestrict?localY.ReadWrite():y.ReadWrite(), NDOFS, ne);
917  auto A = Reshape(ea_data.Read(), NDOFS, NDOFS, ne);
918  mfem::forall(ne*NDOFS, [=] MFEM_HOST_DEVICE (int glob_j)
919  {
920  const int e = glob_j/NDOFS;
921  const int j = glob_j%NDOFS;
922  double res = 0.0;
923  for (int i = 0; i < NDOFS; i++)
924  {
925  res += A(j, i, e)*X(i, e);
926  }
927  Y(j, e) += res;
928  });
929  // Apply the Element Restriction transposed
930  if (useRestrict)
931  {
933  }
934  }
935 
936  // Treatment of interior faces
937  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
938  const int iFISz = intFaceIntegrators.Size();
939  if (int_face_restrict_lex && iFISz>0)
940  {
941  // Apply the Interior Face Restriction
943  if (int_face_X.Size()>0)
944  {
945  int_face_Y = 0.0;
946  // Apply the interior face matrices transposed
947  const int NDOFS = faceDofs;
948  auto X = Reshape(int_face_X.Read(), NDOFS, 2, nf_int);
949  auto Y = Reshape(int_face_Y.ReadWrite(), NDOFS, 2, nf_int);
951  {
952  auto A_int = Reshape(ea_data_int.Read(), NDOFS, NDOFS, 2, nf_int);
953  mfem::forall(nf_int*NDOFS, [=] MFEM_HOST_DEVICE (int glob_j)
954  {
955  const int f = glob_j/NDOFS;
956  const int j = glob_j%NDOFS;
957  double res = 0.0;
958  for (int i = 0; i < NDOFS; i++)
959  {
960  res += A_int(j, i, 0, f)*X(i, 0, f);
961  }
962  Y(j, 0, f) += res;
963  res = 0.0;
964  for (int i = 0; i < NDOFS; i++)
965  {
966  res += A_int(j, i, 1, f)*X(i, 1, f);
967  }
968  Y(j, 1, f) += res;
969  });
970  }
971  auto A_ext = Reshape(ea_data_ext.Read(), NDOFS, NDOFS, 2, nf_int);
972  mfem::forall(nf_int*NDOFS, [=] MFEM_HOST_DEVICE (int glob_j)
973  {
974  const int f = glob_j/NDOFS;
975  const int j = glob_j%NDOFS;
976  double res = 0.0;
977  for (int i = 0; i < NDOFS; i++)
978  {
979  res += A_ext(j, i, 1, f)*X(i, 0, f);
980  }
981  Y(j, 1, f) += res;
982  res = 0.0;
983  for (int i = 0; i < NDOFS; i++)
984  {
985  res += A_ext(j, i, 0, f)*X(i, 1, f);
986  }
987  Y(j, 0, f) += res;
988  });
989  // Apply the Interior Face Restriction transposed
991  }
992  }
993 
994  // Treatment of boundary faces
995  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
996  const int bFISz = bdrFaceIntegrators.Size();
997  if (!factorize_face_terms && bdr_face_restrict_lex && bFISz>0)
998  {
999  // Apply the Boundary Face Restriction
1001  if (bdr_face_X.Size()>0)
1002  {
1003  bdr_face_Y = 0.0;
1004  // Apply the boundary face matrices transposed
1005  const int NDOFS = faceDofs;
1006  auto X = Reshape(bdr_face_X.Read(), NDOFS, nf_bdr);
1007  auto Y = Reshape(bdr_face_Y.ReadWrite(), NDOFS, nf_bdr);
1008  auto A = Reshape(ea_data_bdr.Read(), NDOFS, NDOFS, nf_bdr);
1009  mfem::forall(nf_bdr*NDOFS, [=] MFEM_HOST_DEVICE (int glob_j)
1010  {
1011  const int f = glob_j/NDOFS;
1012  const int j = glob_j%NDOFS;
1013  double res = 0.0;
1014  for (int i = 0; i < NDOFS; i++)
1015  {
1016  res += A(j, i, f)*X(i, f);
1017  }
1018  Y(j, f) += res;
1019  });
1020  // Apply the Boundary Face Restriction transposed
1022  }
1023  }
1024 }
1025 
1026 // Data and methods for fully-assembled bilinear forms
1028  : EABilinearFormExtension(form),
1029  mat(a->mat)
1030 {
1031 #ifdef MFEM_USE_MPI
1032  ParFiniteElementSpace *pfes = nullptr;
1033  if ( a->GetFBFI()->Size()>0 &&
1034  (pfes = dynamic_cast<ParFiniteElementSpace*>(form->FESpace())) )
1035  {
1036  pfes->ExchangeFaceNbrData();
1037  }
1038 #endif
1039 }
1040 
1042 {
1044  FiniteElementSpace &fes = *a->FESpace();
1045  int width = fes.GetVSize();
1046  int height = fes.GetVSize();
1047  bool keep_nbr_block = false;
1048 #ifdef MFEM_USE_MPI
1049  ParFiniteElementSpace *pfes = nullptr;
1050  if ( a->GetFBFI()->Size()>0 &&
1051  (pfes = dynamic_cast<ParFiniteElementSpace*>(&fes)) )
1052  {
1053  pfes->ExchangeFaceNbrData();
1054  width += pfes->GetFaceNbrVSize();
1055  dg_x.SetSize(width);
1056  ParBilinearForm *pb = nullptr;
1057  if ((pb = dynamic_cast<ParBilinearForm*>(a)) && (pb->keep_nbr_block))
1058  {
1059  height += pfes->GetFaceNbrVSize();
1060  dg_y.SetSize(height);
1061  keep_nbr_block = true;
1062  }
1063  }
1064 #endif
1065  if (a->mat) // We reuse the sparse matrix memory
1066  {
1067  if (fes.IsDGSpace())
1068  {
1069  const L2ElementRestriction *restE =
1070  static_cast<const L2ElementRestriction*>(elem_restrict);
1071  const L2FaceRestriction *restF =
1072  static_cast<const L2FaceRestriction*>(int_face_restrict_lex);
1073  MFEM_VERIFY(
1074  fes.Conforming(),
1075  "Full Assembly not yet supported on NCMesh.");
1076  // 1. Fill J and Data
1077  // 1.1 Fill J and Data with Elem ea_data
1078  restE->FillJAndData(ea_data, *mat);
1079  // 1.2 Fill J and Data with Face ea_data_ext
1080  if (restF) { restF->FillJAndData(ea_data_ext, *mat, keep_nbr_block); }
1081  // 1.3 Shift indirections in I back to original
1082  auto I = mat->HostReadWriteI();
1083  for (int i = height; i > 0; i--)
1084  {
1085  I[i] = I[i-1];
1086  }
1087  I[0] = 0;
1088  }
1089  else
1090  {
1091  const ElementRestriction &rest =
1092  static_cast<const ElementRestriction&>(*elem_restrict);
1093  rest.FillJAndData(ea_data, *mat);
1094  }
1095  }
1096  else // We create, compute the sparsity, and fill the sparse matrix
1097  {
1098  mat = new SparseMatrix;
1099  mat->OverrideSize(height, width);
1100  if (fes.IsDGSpace())
1101  {
1102  const L2ElementRestriction *restE =
1103  static_cast<const L2ElementRestriction*>(elem_restrict);
1104  const L2FaceRestriction *restF =
1105  static_cast<const L2FaceRestriction*>(int_face_restrict_lex);
1106  MFEM_VERIFY(
1107  fes.Conforming(),
1108  "Full Assembly not yet supported on NCMesh.");
1109  // 1. Fill I
1110  mat->GetMemoryI().New(height+1, mat->GetMemoryI().GetMemoryType());
1111  // 1.1 Increment with restE
1112  restE->FillI(*mat);
1113  // 1.2 Increment with restF
1114  if (restF) { restF->FillI(*mat, keep_nbr_block); }
1115  // 1.3 Sum the non-zeros in I
1116  auto h_I = mat->HostReadWriteI();
1117  int cpt = 0;
1118  for (int i = 0; i < height; i++)
1119  {
1120  const int nnz = h_I[i];
1121  h_I[i] = cpt;
1122  cpt += nnz;
1123  }
1124  const int nnz = cpt;
1125  h_I[height] = nnz;
1126  mat->GetMemoryJ().New(nnz, mat->GetMemoryJ().GetMemoryType());
1127  mat->GetMemoryData().New(nnz, mat->GetMemoryData().GetMemoryType());
1128  // 2. Fill J and Data
1129  // 2.1 Fill J and Data with Elem ea_data
1130  restE->FillJAndData(ea_data, *mat);
1131  // 2.2 Fill J and Data with Face ea_data_ext
1132  if (restF) { restF->FillJAndData(ea_data_ext, *mat, keep_nbr_block); }
1133  // 2.3 Shift indirections in I back to original
1134  auto I = mat->HostReadWriteI();
1135  for (int i = height; i > 0; i--)
1136  {
1137  I[i] = I[i-1];
1138  }
1139  I[0] = 0;
1140  }
1141  else // continuous Galerkin case
1142  {
1143  const ElementRestriction &rest =
1144  static_cast<const ElementRestriction&>(*elem_restrict);
1145  rest.FillSparseMatrix(ea_data, *mat);
1146  }
1147  a->mat = mat;
1148  }
1149  if ( a->sort_sparse_matrix )
1150  {
1151  a->mat->SortColumnIndices();
1152  }
1153 }
1154 
1155 
1157 {
1158 #ifdef MFEM_USE_MPI
1159  if ( auto pa = dynamic_cast<ParBilinearForm*>(a) )
1160  {
1161  pa->ParallelRAP(*pa->mat, A);
1162  }
1163  else
1164 #endif
1165  {
1166  a->SerialRAP(A);
1167  }
1168 }
1169 
1171  OperatorHandle &A)
1172 {
1173  MFEM_VERIFY(a->diag_policy == DiagonalPolicy::DIAG_ONE,
1174  "Only DiagonalPolicy::DIAG_ONE supported with"
1175  " FABilinearFormExtension.");
1176 #ifdef MFEM_USE_MPI
1177  if ( dynamic_cast<ParBilinearForm*>(a) )
1178  {
1179  A.As<HypreParMatrix>()->EliminateBC(ess_dofs,
1180  DiagonalPolicy::DIAG_ONE);
1181  }
1182  else
1183 #endif
1184  {
1185  A.As<SparseMatrix>()->EliminateBC(ess_dofs,
1186  DiagonalPolicy::DIAG_ONE);
1187  }
1188 }
1189 
1191  OperatorHandle &A)
1192 {
1193  RAP(A);
1194  EliminateBC(ess_dofs, A);
1195 }
1196 
1198  Vector &x, Vector &b,
1199  OperatorHandle &A,
1200  Vector &X, Vector &B,
1201  int copy_interior)
1202 {
1203  Operator *A_out;
1204  Operator::FormLinearSystem(ess_tdof_list, x, b, A_out, X, B, copy_interior);
1205  delete A_out;
1206  FormSystemMatrix(ess_tdof_list, A);
1207 }
1208 
1210 {
1211 #ifdef MFEM_USE_MPI
1212  const ParFiniteElementSpace *pfes;
1213  if ( (pfes = dynamic_cast<const ParFiniteElementSpace*>(test_fes)) )
1214  {
1215  // DG Prolongation
1216  ParGridFunction x_gf;
1217  x_gf.MakeRef(const_cast<ParFiniteElementSpace*>(pfes),
1218  const_cast<Vector&>(x),0);
1219  x_gf.ExchangeFaceNbrData();
1220  Vector &shared_x = x_gf.FaceNbrData();
1221  const int local_size = a->FESpace()->GetVSize();
1222  auto dg_x_ptr = dg_x.Write();
1223  auto x_ptr = x.Read();
1224  mfem::forall(local_size, [=] MFEM_HOST_DEVICE (int i)
1225  {
1226  dg_x_ptr[i] = x_ptr[i];
1227  });
1228  const int shared_size = shared_x.Size();
1229  auto shared_x_ptr = shared_x.Read();
1230  mfem::forall(shared_size, [=] MFEM_HOST_DEVICE (int i)
1231  {
1232  dg_x_ptr[local_size+i] = shared_x_ptr[i];
1233  });
1234  ParBilinearForm *pform = nullptr;
1235  if ((pform = dynamic_cast<ParBilinearForm*>(a)) && (pform->keep_nbr_block))
1236  {
1237  mat->Mult(dg_x, dg_y);
1238  // DG Restriction
1239  auto dg_y_ptr = dg_y.Read();
1240  auto y_ptr = y.ReadWrite();
1241  mfem::forall(local_size, [=] MFEM_HOST_DEVICE (int i)
1242  {
1243  y_ptr[i] += dg_y_ptr[i];
1244  });
1245  }
1246  else
1247  {
1248  mat->Mult(dg_x, y);
1249  }
1250  }
1251  else
1252 #endif
1253  {
1254  mat->Mult(x, y);
1255  }
1256 }
1257 
1259 {
1260  if ( a->GetFBFI()->Size()>0 )
1261  {
1262  DGMult(x, y);
1263  }
1264  else
1265  {
1266  mat->Mult(x, y);
1267  }
1268 }
1269 
1271 {
1272 #ifdef MFEM_USE_MPI
1273  const ParFiniteElementSpace *pfes;
1274  if ( (pfes = dynamic_cast<const ParFiniteElementSpace*>(test_fes)) )
1275  {
1276  // DG Prolongation
1277  ParGridFunction x_gf;
1278  x_gf.MakeRef(const_cast<ParFiniteElementSpace*>(pfes),
1279  const_cast<Vector&>(x),0);
1280  x_gf.ExchangeFaceNbrData();
1281  Vector &shared_x = x_gf.FaceNbrData();
1282  const int local_size = a->FESpace()->GetVSize();
1283  auto dg_x_ptr = dg_x.Write();
1284  auto x_ptr = x.Read();
1285  mfem::forall(local_size, [=] MFEM_HOST_DEVICE (int i)
1286  {
1287  dg_x_ptr[i] = x_ptr[i];
1288  });
1289  const int shared_size = shared_x.Size();
1290  auto shared_x_ptr = shared_x.Read();
1291  mfem::forall(shared_size, [=] MFEM_HOST_DEVICE (int i)
1292  {
1293  dg_x_ptr[local_size+i] = shared_x_ptr[i];
1294  });
1295  ParBilinearForm *pb = nullptr;
1296  if ((pb = dynamic_cast<ParBilinearForm*>(a)) && (pb->keep_nbr_block))
1297  {
1298  mat->MultTranspose(dg_x, dg_y);
1299  // DG Restriction
1300  auto dg_y_ptr = dg_y.Read();
1301  auto y_ptr = y.ReadWrite();
1302  mfem::forall(local_size, [=] MFEM_HOST_DEVICE (int i)
1303  {
1304  y_ptr[i] += dg_y_ptr[i];
1305  });
1306  }
1307  else
1308  {
1309  mat->MultTranspose(dg_x, y);
1310  }
1311  }
1312  else
1313 #endif
1314  {
1315  mat->MultTranspose(x, y);
1316  }
1317 }
1318 
1320 {
1321  if ( a->GetFBFI()->Size()>0 )
1322  {
1323  DGMultTranspose(x, y);
1324  }
1325  else
1326  {
1327  mat->MultTranspose(x, y);
1328  }
1329 }
1330 
1331 
1333  : Operator(form->Height(), form->Width()), a(form)
1334 {
1335  // empty
1336 }
1337 
1339 {
1340  return a->GetProlongation();
1341 }
1342 
1344 {
1345  return a->GetRestriction();
1346 }
1347 
1349 {
1350  return a->GetOutputProlongation();
1351 }
1352 
1354 {
1355  return a->GetOutputRestriction();
1356 }
1357 
1358 // Data and methods for partially-assembled bilinear forms
1359 
1361  MixedBilinearForm *form)
1363  trial_fes(form->TrialFESpace()),
1364  test_fes(form->TestFESpace()),
1365  elem_restrict_trial(NULL),
1366  elem_restrict_test(NULL)
1367 {
1368  Update();
1369 }
1370 
1372 {
1373  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1374  const int integratorCount = integrators.Size();
1375  for (int i = 0; i < integratorCount; ++i)
1376  {
1377  integrators[i]->AssemblePA(*trial_fes, *test_fes);
1378  }
1379  MFEM_VERIFY(a->GetBBFI()->Size() == 0,
1380  "Partial assembly does not support AddBoundaryIntegrator yet.");
1381  MFEM_VERIFY(a->GetTFBFI()->Size() == 0,
1382  "Partial assembly does not support AddTraceFaceIntegrator yet.");
1383  MFEM_VERIFY(a->GetBTFBFI()->Size() == 0,
1384  "Partial assembly does not support AddBdrTraceFaceIntegrator yet.");
1385 }
1386 
1388 {
1389  trial_fes = a->TrialFESpace();
1390  test_fes = a->TestFESpace();
1391  height = test_fes->GetVSize();
1392  width = trial_fes->GetVSize();
1397  if (elem_restrict_trial)
1398  {
1399  localTrial.UseDevice(true);
1402  }
1403  if (elem_restrict_test)
1404  {
1405  localTest.UseDevice(true); // ensure 'localY = 0.0' is done on device
1407  }
1408 }
1409 
1411  const Array<int> &trial_tdof_list,
1412  const Array<int> &test_tdof_list,
1413  OperatorHandle &A)
1414 {
1415  Operator * oper;
1416  Operator::FormRectangularSystemOperator(trial_tdof_list, test_tdof_list,
1417  oper);
1418  A.Reset(oper); // A will own oper
1419 }
1420 
1422  const Array<int> &trial_tdof_list,
1423  const Array<int> &test_tdof_list,
1424  Vector &x, Vector &b,
1425  OperatorHandle &A,
1426  Vector &X, Vector &B)
1427 {
1428  Operator *oper;
1429  Operator::FormRectangularLinearSystem(trial_tdof_list, test_tdof_list, x, b,
1430  oper, X, B);
1431  A.Reset(oper); // A will own oper
1432 }
1433 
1435  const Operator *elem_restrict_x,
1436  const Vector &x,
1437  Vector &localX,
1438  const Operator *elem_restrict_y,
1439  Vector &y,
1440  Vector &localY,
1441  const double c) const
1442 {
1443  // * G operation: localX = c*local(x)
1444  if (elem_restrict_x)
1445  {
1446  elem_restrict_x->Mult(x, localX);
1447  if (c != 1.0)
1448  {
1449  localX *= c;
1450  }
1451  }
1452  else
1453  {
1454  if (c == 1.0)
1455  {
1456  localX.SyncAliasMemory(x);
1457  }
1458  else
1459  {
1460  localX.Set(c, x);
1461  }
1462  }
1463  if (elem_restrict_y)
1464  {
1465  localY = 0.0;
1466  }
1467  else
1468  {
1469  y.UseDevice(true);
1470  localY.SyncAliasMemory(y);
1471  }
1472 }
1473 
1475 {
1476  y = 0.0;
1477  AddMult(x, y);
1478 }
1479 
1481  const double c) const
1482 {
1483  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1484  const int iSz = integrators.Size();
1485 
1486  // * G operation
1488  elem_restrict_test, y, localTest, c);
1489 
1490  // * B^TDB operation
1491  for (int i = 0; i < iSz; ++i)
1492  {
1493  integrators[i]->AddMultPA(localTrial, localTest);
1494  }
1495 
1496  // * G^T operation
1497  if (elem_restrict_test)
1498  {
1499  tempY.SetSize(y.Size());
1501  y += tempY;
1502  }
1503 }
1504 
1506  Vector &y) const
1507 {
1508  y = 0.0;
1509  AddMultTranspose(x, y);
1510 }
1511 
1513  const double c) const
1514 {
1515  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1516  const int iSz = integrators.Size();
1517 
1518  // * G operation
1521 
1522  // * B^TD^TB operation
1523  for (int i = 0; i < iSz; ++i)
1524  {
1525  integrators[i]->AddMultTransposePA(localTest, localTrial);
1526  }
1527 
1528  // * G^T operation
1529  if (elem_restrict_trial)
1530  {
1531  tempY.SetSize(y.Size());
1533  y += tempY;
1534  }
1535 }
1536 
1538  Vector &diag) const
1539 {
1540  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1541 
1542  const int iSz = integrators.Size();
1543 
1544  if (elem_restrict_trial)
1545  {
1546  const ElementRestriction* H1elem_restrict_trial =
1547  dynamic_cast<const ElementRestriction*>(elem_restrict_trial);
1548  if (H1elem_restrict_trial)
1549  {
1550  H1elem_restrict_trial->MultUnsigned(D, localTrial);
1551  }
1552  else
1553  {
1555  }
1556  }
1557 
1558  if (elem_restrict_test)
1559  {
1560  localTest = 0.0;
1561  for (int i = 0; i < iSz; ++i)
1562  {
1563  if (elem_restrict_trial)
1564  {
1565  integrators[i]->AssembleDiagonalPA_ADAt(localTrial, localTest);
1566  }
1567  else
1568  {
1569  integrators[i]->AssembleDiagonalPA_ADAt(D, localTest);
1570  }
1571  }
1572  const ElementRestriction* H1elem_restrict_test =
1573  dynamic_cast<const ElementRestriction*>(elem_restrict_test);
1574  if (H1elem_restrict_test)
1575  {
1576  H1elem_restrict_test->MultTransposeUnsigned(localTest, diag);
1577  }
1578  else
1579  {
1581  }
1582  }
1583  else
1584  {
1585  diag.UseDevice(true); // typically this is a large vector, so store on device
1586  diag = 0.0;
1587  for (int i = 0; i < iSz; ++i)
1588  {
1589  if (elem_restrict_trial)
1590  {
1591  integrators[i]->AssembleDiagonalPA_ADAt(localTrial, diag);
1592  }
1593  else
1594  {
1595  integrators[i]->AssembleDiagonalPA_ADAt(D, diag);
1596  }
1597  }
1598  }
1599 }
1600 
1602  DiscreteLinearOperator *linop) :
1604 {
1605 }
1606 
1607 const
1609 const
1610 {
1611  return a->GetOutputRestrictionTranspose();
1612 }
1613 
1615 {
1616  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1617  const int integratorCount = integrators.Size();
1618  for (int i = 0; i < integratorCount; ++i)
1619  {
1620  integrators[i]->AssemblePA(*trial_fes, *test_fes);
1621  }
1622 
1623  test_multiplicity.UseDevice(true);
1624  test_multiplicity.SetSize(elem_restrict_test->Width()); // l-vector
1625  Vector ones(elem_restrict_test->Height()); // e-vector
1626  ones = 1.0;
1627 
1628  const ElementRestriction* elem_restrict =
1629  dynamic_cast<const ElementRestriction*>(elem_restrict_test);
1630  if (elem_restrict)
1631  {
1632  elem_restrict->MultTransposeUnsigned(ones, test_multiplicity);
1633  }
1634  else
1635  {
1636  mfem_error("A real ElementRestriction is required in this setting!");
1637  }
1638 
1639  auto tm = test_multiplicity.ReadWrite();
1640  mfem::forall(test_multiplicity.Size(), [=] MFEM_HOST_DEVICE (int i)
1641  {
1642  tm[i] = 1.0 / tm[i];
1643  });
1644 }
1645 
1647  const Vector &x, Vector &y, const double c) const
1648 {
1649  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1650  const int iSz = integrators.Size();
1651 
1652  // * G operation
1654  elem_restrict_test, y, localTest, c);
1655 
1656  // * B^TDB operation
1657  for (int i = 0; i < iSz; ++i)
1658  {
1659  integrators[i]->AddMultPA(localTrial, localTest);
1660  }
1661 
1662  // do a kind of "set" rather than "add" in the below
1663  // operation as compared to the BilinearForm case
1664  // * G^T operation (kind of...)
1665  const ElementRestriction* elem_restrict =
1666  dynamic_cast<const ElementRestriction*>(elem_restrict_test);
1667  if (elem_restrict)
1668  {
1669  tempY.SetSize(y.Size());
1670  elem_restrict->MultLeftInverse(localTest, tempY);
1671  y += tempY;
1672  }
1673  else
1674  {
1675  mfem_error("In this setting you need a real ElementRestriction!");
1676  }
1677 }
1678 
1680  const Vector &x, Vector &y, const double c) const
1681 {
1682  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1683  const int iSz = integrators.Size();
1684 
1685  // do a kind of "set" rather than "add" in the below
1686  // operation as compared to the BilinearForm case
1687  // * G operation (kinda)
1688  Vector xscaled(x);
1689  MFEM_VERIFY(x.Size() == test_multiplicity.Size(), "Input vector of wrong size");
1690  auto xs = xscaled.ReadWrite();
1691  auto tm = test_multiplicity.Read();
1692  mfem::forall(x.Size(), [=] MFEM_HOST_DEVICE (int i)
1693  {
1694  xs[i] *= tm[i];
1695  });
1698 
1699  // * B^TD^TB operation
1700  for (int i = 0; i < iSz; ++i)
1701  {
1702  integrators[i]->AddMultTransposePA(localTest, localTrial);
1703  }
1704 
1705  // * G^T operation
1706  if (elem_restrict_trial)
1707  {
1708  tempY.SetSize(y.Size());
1710  y += tempY;
1711  }
1712  else
1713  {
1714  mfem_error("Trial ElementRestriction not defined");
1715  }
1716 }
1717 
1719  const Array<int>& ess1, const Array<int>& ess2, OperatorHandle &A)
1720 {
1721  const Operator *Pi = this->GetProlongation();
1722  const Operator *RoT = this->GetOutputRestrictionTranspose();
1723  Operator *rap = SetupRAP(Pi, RoT);
1724 
1726  = new RectangularConstrainedOperator(rap, ess1, ess2, rap != this);
1727 
1728  A.Reset(Arco);
1729 }
1730 
1731 } // namespace mfem
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:307
void FillI(SparseMatrix &mat) const
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:93
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:733
void FormRectangularSystemOperator(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Operator *&A)
Return in A a parallel (on truedofs) version of this rectangular operator (including constraints)...
Definition: operator.cpp:235
void AddMultTranspose(const Vector &x, Vector &y, const double c=1.0) const
y += c*A^T*x
Array< BilinearFormIntegrator * > * GetBBFI()
Access all the integrators added with AddBoundaryIntegrator().
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
virtual const Operator * GetRestriction() const
Get the finite element space restriction matrix.
void FormRectangularSystemOperator(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)
Setup OperatorHandle A to contain constrained linear operator.
void Mult(const Vector &x, Vector &y) const override=0
Extract the face degrees of freedom from x into y.
void SetupMultInputs(const Operator *elem_restrict_x, const Vector &x, Vector &localX, const Operator *elem_restrict_y, Vector &y, Vector &localY, const double c) const
Helper function to set up inputs/outputs for Mult or MultTranspose.
const FiniteElementSpace * test_fes
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
virtual const Operator * GetProlongation() const
Get the finite element space prolongation operator.
int GetBdrElementEdgeIndex(int i) const
Definition: mesh.cpp:6588
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
Array< BilinearFormIntegrator * > * GetDBFI()
Access all the integrators added with AddDomainIntegrator().
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:5668
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:517
MFBilinearFormExtension(BilinearForm *form)
Array< BilinearFormIntegrator * > * GetBFBFI()
Access all integrators added with AddBdrFaceIntegrator().
void Assemble()
Assemble at the level given for the BilinearFormExtension subclass.
void Assemble()
Assemble at the level given for the BilinearFormExtension subclass.
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:115
virtual void AddMultTransposeUnsigned(const Vector &x, Vector &y, const double a=1.0) const
Add the face degrees of freedom x to the element degrees of freedom y ignoring the signs from DOF ori...
*Array< Array< int > * > * GetDBFI_Marker()
Access all boundary markers added with AddDomainIntegrator().
Geometry::Type GetFaceGeometry(int i) const
Return the Geometry::Type associated with face i.
Definition: mesh.cpp:1421
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
PAMixedBilinearFormExtension(MixedBilinearForm *form)
virtual const FaceRestriction * GetFaceRestriction(ElementDofOrdering f_ordering, FaceType, L2FaceValues mul=L2FaceValues::DoubleValued) const
Return an Operator that converts L-vectors to E-vectors on each face.
Definition: fespace.cpp:1335
int Size() const
Returns the size of the vector.
Definition: vector.hpp:197
const NURBSExtension * GetNURBSext() const
Definition: fespace.hpp:557
virtual const Operator * GetRestriction() const
Get the finite element space restriction operator.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:453
Array< Array< int > * > * GetBFBFI_Marker()
Access all boundary markers added with AddBdrFaceIntegrator(). If no marker was specified when the in...
void DGMultTranspose(const Vector &x, Vector &y) const
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
Abstract parallel finite element space.
Definition: pfespace.hpp:28
Array< int > elem_attributes
Attributes of all mesh elements.
virtual void FillJAndData(const Vector &fea_data, SparseMatrix &mat, const bool keep_nbr_block=false) const
Fill the J and Data arrays of the SparseMatrix corresponding to the sparsity pattern given by this L2...
BilinearFormExtension(BilinearForm *form)
virtual const Operator * GetOutputProlongation() const
Get the test finite element space prolongation matrix.
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Definition: fespace.hpp:1306
void AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
Array< BilinearFormIntegrator * > * GetFBFI()
Access all integrators added with AddInteriorFaceIntegrator().
virtual const Operator * GetProlongation() const
Get the input finite element space prolongation matrix.
void SortColumnIndices()
Sort the column indices corresponding to each row.
Definition: sparsemat.cpp:457
Memory< double > & GetMemoryData()
Definition: sparsemat.hpp:269
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 ...
virtual const Operator * GetOutputProlongation() const
Get the output finite element space restriction matrix.
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:1089
Memory< int > & GetMemoryI()
Definition: sparsemat.hpp:237
Array< Array< int > * > * GetBBFI_Marker()
Access all boundary markers added with AddBoundaryIntegrator(). If no marker was specified when the i...
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:2841
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 ...
Data and methods for partially-assembled bilinear forms.
Class extending the MixedBilinearForm class to support different AssemblyLevels.
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:1302
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
Definition: bilininteg.cpp:99
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition: fespace.hpp:753
Memory< int > & GetMemoryJ()
Definition: sparsemat.hpp:253
void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const
Assemble the diagonal of ADA^T for a diagonal vector D.
std::function< double(const Vector &)> f(double mass_coeff)
Definition: lor_mms.hpp:30
void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B)
Form a column-constrained linear system using a matrix-free approach.
Definition: operator.cpp:131
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:111
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1190
virtual void AddFaceMatricesToElementMatrices(const Vector &fea_data, Vector &ea_data) const
This methods adds the DG face matrices to the element matrices.
virtual const Operator * GetRestriction() const
Get the finite element space restriction matrix.
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 ...
PADiscreteLinearOperatorExtension(DiscreteLinearOperator *linop)
virtual void AddMultTransposeInPlace(Vector &x, Vector &y) const
Add the face degrees of freedom x to the element degrees of freedom y. Perform the same computation a...
void Assemble()
Partial assembly of all internal integrators.
Data type sparse matrix.
Definition: sparsemat.hpp:50
Native ordering as defined by the FiniteElement.
Operator that extracts Face degrees of freedom for L2 spaces.
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:461
SparseMatrix * mat
Sparse matrix to be associated with the form. Owned.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:154
double b
Definition: lissajous.cpp:42
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition: device.hpp:273
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
Definition: device.hpp:277
BilinearForm * a
Not owned.
void AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
FiniteElementSpace * TestFESpace()
Return the test FE space associated with the BilinearForm.
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void MultLeftInverse(const Vector &x, Vector &y) const
void EliminateBC(const Array< int > &ess_dofs, OperatorHandle &A)
void MultTransposeUnsigned(const Vector &x, Vector &y) const
Compute MultTranspose without applying signs based on DOF orientations.
Array< BilinearFormIntegrator * > * GetBTFBFI()
Access all integrators added with AddBdrTraceFaceIntegrator().
Class extending the BilinearForm class to support different AssemblyLevels.
virtual const Operator * GetOutputRestrictionTranspose() const
Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type...
Definition: operator.hpp:156
void FormRectangularSystemOperator(const Array< int > &, const Array< int > &, OperatorHandle &A)
Setup OperatorHandle A to contain constrained linear operator.
virtual void FillI(SparseMatrix &mat, const bool keep_nbr_block=false) const
Fill the I array of SparseMatrix corresponding to the sparsity pattern given by this L2FaceRestrictio...
FABilinearFormExtension(BilinearForm *form)
const FaceRestriction * bdr_face_restrict_lex
void Assemble()
Partial assembly of all internal integrators.
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
bool IsOfFaceType(FaceType type) const
Return true if the face is of the same type as type.
Definition: mesh.hpp:1712
void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
Array< BilinearFormIntegrator * > * GetTFBFI()
Access all integrators added with AddTraceFaceIntegrator().
Array< BilinearFormIntegrator * > * GetBBFI()
Access all integrators added with AddBoundaryIntegrator().
const Operator * GetOutputRestrictionTranspose() const
Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type...
void Assemble()
Assemble at the level given for the BilinearFormExtension subclass.
Vector & FaceNbrData()
Definition: pgridfunc.hpp:208
void SetupRestrictionOperators(const L2FaceValues m)
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:1196
virtual void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:908
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:555
Data and methods for element-assembled bilinear forms.
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:24
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
void forall(int N, lambda &&body)
Definition: forall.hpp:742
Operator that converts FiniteElementSpace L-vectors to E-vectors.
Definition: restriction.hpp:37
const FiniteElementSpace * test_fes
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
void SerialRAP(OperatorHandle &A)
Compute serial RAP operator and store it in A as a SparseMatrix.
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 ...
virtual const Operator * GetOutputRestriction() const
Get the output finite element space restriction matrix.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:687
void Assemble()
Assemble at the level given for the BilinearFormExtension subclass.
bool DeviceCanUseCeed()
Function that determines if a CEED kernel should be used, based on the current mfem::Device configura...
Definition: util.cpp:33
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
void MultUnsigned(const Vector &x, Vector &y) const
Compute Mult without applying signs based on DOF orientations.
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
Form a constrained linear system using a matrix-free approach.
Definition: operator.cpp:114
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
void Mult(const Vector &x, Vector &y) const
y = A*x
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
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:258
void FillSparseMatrix(const Vector &mat_ea, SparseMatrix &mat) const
Fill a Sparse Matrix with Element Matrices.
DiagonalPolicy diag_policy
const FaceRestriction * int_face_restrict_lex
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:263
EABilinearFormExtension(BilinearForm *form)
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:1086
double a
Definition: lissajous.cpp:41
MixedBilinearForm * a
Not owned.
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
const FaceRestriction * bdr_face_restrict_lex
const FaceRestriction * int_face_restrict_lex
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:469
Rectangular Operator for imposing essential boundary conditions on the input space using only the act...
Definition: operator.hpp:938
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
void SyncAliasMemory(const Vector &v) const
Update the alias memory location of the vector to match v.
Definition: vector.hpp:238
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
const FiniteElementSpace * test_fes
void FormSystemOperator(const Array< int > &ess_tdof_list, Operator *&A)
Return in A a parallel (on truedofs) version of this square operator.
Definition: operator.cpp:227
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:74
void Update()
Update internals for when a new MixedBilinearForm is given to this class.
void AddMult(const Vector &x, Vector &y, const double c=1.0) const
y += c*A*x
void AddMult(const Vector &x, Vector &y, const double c=1.0) const
y += c*A*x
void OverrideSize(int height_, int width_)
Sets the height and width of the matrix.
Definition: sparsemat.cpp:297
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:709
Lexicographic ordering for tensor-product FiniteElements.
Class for parallel bilinear form.
Operator * SetupRAP(const Operator *Pi, const Operator *Po)
Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P"...
Definition: operator.cpp:168
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
FaceInformation GetFaceInformation(int f) const
Definition: mesh.cpp:1138
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
Data and methods for partially-assembled mixed bilinear forms.
Vector data type.
Definition: vector.hpp:58
bool IsDGSpace() const
Return whether or not the space is discontinuous (L2)
Definition: fespace.hpp:1274
void DGMult(const Vector &x, Vector &y) const
void MultTranspose(const Vector &x, Vector &y) const
y = A^T*x
Array< BilinearFormIntegrator * > * GetDBFI()
Access all integrators added with AddDomainIntegrator().
int * HostReadWriteI()
Definition: sparsemat.hpp:249
virtual const Operator * GetOutputRestriction() const
Get the test finite element space restriction matrix.
void AddMultWithMarkers(const BilinearFormIntegrator &integ, const Vector &x, const Array< int > *markers, const Array< int > &attributes, const bool transpose, Vector &y) const
Accumulate the action (or transpose) of the integrator on x into y, taking into account the (possibly...
Class for parallel grid function.
Definition: pgridfunc.hpp:32
virtual const Operator * GetRestriction() const
Get the input finite element space restriction matrix.
const FiniteElementSpace * trial_fes
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
Bitwise-OR of all CEED backends.
Definition: device.hpp:94
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition: dtensor.hpp:131
FiniteElementSpace * TrialFESpace()
Return the trial FE space associated with the BilinearForm.
const FiniteElementSpace * trial_fes
virtual void AddMultTransposePA(const Vector &x, Vector &y) const
Method for partially assembled transposed action.
Definition: bilininteg.cpp:111
bool Conforming() const
Definition: fespace.hpp:561
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
MixedBilinearFormExtension(MixedBilinearForm *form)
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Operator that converts L2 FiniteElementSpace L-vectors to E-vectors.
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
const FiniteElementSpace * trial_fes
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...
void AddMultTranspose(const Vector &x, Vector &y, const double c=1.0) const
y += c*A^T*x
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:145