MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
bilinearform_ext.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 
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 
62 {
63  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
64 
65  const int iSz = integrators.Size();
67  {
68  localY = 0.0;
69  for (int i = 0; i < iSz; ++i)
70  {
71  integrators[i]->AssembleDiagonalMF(localY);
72  }
73  const ElementRestriction* H1elem_restrict =
74  dynamic_cast<const ElementRestriction*>(elem_restrict);
75  if (H1elem_restrict)
76  {
77  H1elem_restrict->MultTransposeUnsigned(localY, y);
78  }
79  else
80  {
82  }
83  }
84  else
85  {
86  y.UseDevice(true); // typically this is a large vector, so store on device
87  y = 0.0;
88  for (int i = 0; i < iSz; ++i)
89  {
90  integrators[i]->AssembleDiagonalMF(y);
91  }
92  }
93 }
94 
96 {
97  FiniteElementSpace *fes = a->FESpace();
98  height = width = fes->GetVSize();
99  trial_fes = fes;
100  test_fes = fes;
101 
102  elem_restrict = nullptr;
103  int_face_restrict_lex = nullptr;
104  bdr_face_restrict_lex = nullptr;
105 }
106 
108  OperatorHandle &A)
109 {
110  Operator *oper;
111  Operator::FormSystemOperator(ess_tdof_list, oper);
112  A.Reset(oper); // A will own oper
113 }
114 
116  Vector &x, Vector &b,
117  OperatorHandle &A,
118  Vector &X, Vector &B,
119  int copy_interior)
120 {
121  Operator *oper;
122  Operator::FormLinearSystem(ess_tdof_list, x, b, oper, X, B, copy_interior);
123  A.Reset(oper); // A will own oper
124 }
125 
127 {
128  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
129 
130  const int iSz = integrators.Size();
132  {
133  y.UseDevice(true); // typically this is a large vector, so store on device
134  y = 0.0;
135  for (int i = 0; i < iSz; ++i)
136  {
137  integrators[i]->AddMultMF(x, y);
138  }
139  }
140  else
141  {
143  localY = 0.0;
144  for (int i = 0; i < iSz; ++i)
145  {
146  integrators[i]->AddMultMF(localX, localY);
147  }
149  }
150 
151  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
152  const int iFISz = intFaceIntegrators.Size();
153  if (int_face_restrict_lex && iFISz>0)
154  {
156  if (int_face_X.Size()>0)
157  {
158  int_face_Y = 0.0;
159  for (int i = 0; i < iFISz; ++i)
160  {
161  intFaceIntegrators[i]->AddMultMF(int_face_X, int_face_Y);
162  }
164  }
165  }
166 
167  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
168  const int bFISz = bdrFaceIntegrators.Size();
169  if (bdr_face_restrict_lex && bFISz>0)
170  {
172  if (bdr_face_X.Size()>0)
173  {
174  bdr_face_Y = 0.0;
175  for (int i = 0; i < bFISz; ++i)
176  {
177  bdrFaceIntegrators[i]->AddMultMF(bdr_face_X, bdr_face_Y);
178  }
180  }
181  }
182 }
183 
185 {
186  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
187  const int iSz = integrators.Size();
188  if (elem_restrict)
189  {
191  localY = 0.0;
192  for (int i = 0; i < iSz; ++i)
193  {
194  integrators[i]->AddMultTransposeMF(localX, localY);
195  }
197  }
198  else
199  {
200  y.UseDevice(true);
201  y = 0.0;
202  for (int i = 0; i < iSz; ++i)
203  {
204  integrators[i]->AddMultTransposeMF(x, y);
205  }
206  }
207 
208  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
209  const int iFISz = intFaceIntegrators.Size();
210  if (int_face_restrict_lex && iFISz>0)
211  {
213  if (int_face_X.Size()>0)
214  {
215  int_face_Y = 0.0;
216  for (int i = 0; i < iFISz; ++i)
217  {
218  intFaceIntegrators[i]->AddMultTransposeMF(int_face_X, int_face_Y);
219  }
221  }
222  }
223 
224  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
225  const int bFISz = bdrFaceIntegrators.Size();
226  if (bdr_face_restrict_lex && bFISz>0)
227  {
229  if (bdr_face_X.Size()>0)
230  {
231  bdr_face_Y = 0.0;
232  for (int i = 0; i < bFISz; ++i)
233  {
234  bdrFaceIntegrators[i]->AddMultTransposeMF(bdr_face_X, bdr_face_Y);
235  }
237  }
238  }
239 }
240 
241 // Data and methods for partially-assembled bilinear forms
243  : BilinearFormExtension(form),
244  trial_fes(a->FESpace()),
245  test_fes(a->FESpace())
246 {
247  elem_restrict = NULL;
248  int_face_restrict_lex = NULL;
249  bdr_face_restrict_lex = NULL;
250 }
251 
253 {
254  if ( Device::Allows(Backend::CEED_MASK) ) { return; }
259  if (elem_restrict)
260  {
261  localX.SetSize(elem_restrict->Height(), Device::GetDeviceMemoryType());
262  localY.SetSize(elem_restrict->Height(), Device::GetDeviceMemoryType());
263  localY.UseDevice(true); // ensure 'localY = 0.0' is done on device
264  }
265 
266  // Construct face restriction operators only if the bilinear form has
267  // interior or boundary face integrators
268  if (int_face_restrict_lex == NULL && a->GetFBFI()->Size() > 0)
269  {
275  int_face_Y.UseDevice(true); // ensure 'int_face_Y = 0.0' is done on device
276  }
277 
278  if (bdr_face_restrict_lex == NULL && a->GetBFBFI()->Size() > 0)
279  {
283  m);
286  bdr_face_Y.UseDevice(true); // ensure 'faceBoundY = 0.0' is done on device
287  }
288 }
289 
291 {
293 
294  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
295  const int integratorCount = integrators.Size();
296  for (int i = 0; i < integratorCount; ++i)
297  {
298  integrators[i]->AssemblePA(*a->FESpace());
299  }
300 
301  MFEM_VERIFY(a->GetBBFI()->Size() == 0,
302  "Partial assembly does not support AddBoundaryIntegrator yet.");
303 
304  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
305  const int intFaceIntegratorCount = intFaceIntegrators.Size();
306  for (int i = 0; i < intFaceIntegratorCount; ++i)
307  {
308  intFaceIntegrators[i]->AssemblePAInteriorFaces(*a->FESpace());
309  }
310 
311  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
312  const int boundFaceIntegratorCount = bdrFaceIntegrators.Size();
313  for (int i = 0; i < boundFaceIntegratorCount; ++i)
314  {
315  bdrFaceIntegrators[i]->AssemblePABoundaryFaces(*a->FESpace());
316  }
317 }
318 
320 {
321  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
322 
323  const int iSz = integrators.Size();
325  {
326  localY = 0.0;
327  for (int i = 0; i < iSz; ++i)
328  {
329  integrators[i]->AssembleDiagonalPA(localY);
330  }
331  const ElementRestriction* H1elem_restrict =
332  dynamic_cast<const ElementRestriction*>(elem_restrict);
333  if (H1elem_restrict)
334  {
335  H1elem_restrict->MultTransposeUnsigned(localY, y);
336  }
337  else
338  {
340  }
341  }
342  else
343  {
344  y.UseDevice(true); // typically this is a large vector, so store on device
345  y = 0.0;
346  for (int i = 0; i < iSz; ++i)
347  {
348  integrators[i]->AssembleDiagonalPA(y);
349  }
350  }
351 }
352 
354 {
355  FiniteElementSpace *fes = a->FESpace();
356  height = width = fes->GetVSize();
357  trial_fes = fes;
358  test_fes = fes;
359 
360  elem_restrict = nullptr;
361  int_face_restrict_lex = nullptr;
362  bdr_face_restrict_lex = nullptr;
363 }
364 
366  OperatorHandle &A)
367 {
368  Operator *oper;
369  Operator::FormSystemOperator(ess_tdof_list, oper);
370  A.Reset(oper); // A will own oper
371 }
372 
374  Vector &x, Vector &b,
375  OperatorHandle &A,
376  Vector &X, Vector &B,
377  int copy_interior)
378 {
379  Operator *oper;
380  Operator::FormLinearSystem(ess_tdof_list, x, b, oper, X, B, copy_interior);
381  A.Reset(oper); // A will own oper
382 }
383 
385 {
386  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
387 
388  const int iSz = integrators.Size();
390  {
391  y.UseDevice(true); // typically this is a large vector, so store on device
392  y = 0.0;
393  for (int i = 0; i < iSz; ++i)
394  {
395  integrators[i]->AddMultPA(x, y);
396  }
397  }
398  else
399  {
401  localY = 0.0;
402  for (int i = 0; i < iSz; ++i)
403  {
404  integrators[i]->AddMultPA(localX, localY);
405  }
407  }
408 
409  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
410  const int iFISz = intFaceIntegrators.Size();
411  if (int_face_restrict_lex && iFISz>0)
412  {
414  if (int_face_X.Size()>0)
415  {
416  int_face_Y = 0.0;
417  for (int i = 0; i < iFISz; ++i)
418  {
419  intFaceIntegrators[i]->AddMultPA(int_face_X, int_face_Y);
420  }
422  }
423  }
424 
425  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
426  const int bFISz = bdrFaceIntegrators.Size();
427  if (bdr_face_restrict_lex && bFISz>0)
428  {
430  if (bdr_face_X.Size()>0)
431  {
432  bdr_face_Y = 0.0;
433  for (int i = 0; i < bFISz; ++i)
434  {
435  bdrFaceIntegrators[i]->AddMultPA(bdr_face_X, bdr_face_Y);
436  }
438  }
439  }
440 }
441 
443 {
444  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
445  const int iSz = integrators.Size();
446  if (elem_restrict)
447  {
449  localY = 0.0;
450  for (int i = 0; i < iSz; ++i)
451  {
452  integrators[i]->AddMultTransposePA(localX, localY);
453  }
455  }
456  else
457  {
458  y.UseDevice(true);
459  y = 0.0;
460  for (int i = 0; i < iSz; ++i)
461  {
462  integrators[i]->AddMultTransposePA(x, y);
463  }
464  }
465 
466  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
467  const int iFISz = intFaceIntegrators.Size();
468  if (int_face_restrict_lex && iFISz>0)
469  {
471  if (int_face_X.Size()>0)
472  {
473  int_face_Y = 0.0;
474  for (int i = 0; i < iFISz; ++i)
475  {
476  intFaceIntegrators[i]->AddMultTransposePA(int_face_X, int_face_Y);
477  }
479  }
480  }
481 
482  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
483  const int bFISz = bdrFaceIntegrators.Size();
484  if (bdr_face_restrict_lex && bFISz>0)
485  {
487  if (bdr_face_X.Size()>0)
488  {
489  bdr_face_Y = 0.0;
490  for (int i = 0; i < bFISz; ++i)
491  {
492  bdrFaceIntegrators[i]->AddMultTransposePA(bdr_face_X, bdr_face_Y);
493  }
495  }
496  }
497 }
498 
499 // Data and methods for element-assembled bilinear forms
501  : PABilinearFormExtension(form),
502  factorize_face_terms(false)
503 {
504  if (form->FESpace()->IsDGSpace() && form->FESpace()->Conforming())
505  {
506  factorize_face_terms = true;
507  }
508 }
509 
511 {
513 
514  ne = trial_fes->GetMesh()->GetNE();
515  elemDofs = trial_fes->GetFE(0)->GetDof();
516 
518  ea_data.UseDevice(true);
519 
520  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
521  const int integratorCount = integrators.Size();
522  if ( integratorCount == 0 )
523  {
524  ea_data = 0.0;
525  }
526  for (int i = 0; i < integratorCount; ++i)
527  {
528  integrators[i]->AssembleEA(*a->FESpace(), ea_data, i);
529  }
530 
531  faceDofs = trial_fes ->
532  GetTraceElement(0, trial_fes->GetMesh()->GetFaceBaseGeometry(0)) ->
533  GetDof();
534 
535  MFEM_VERIFY(a->GetBBFI()->Size() == 0,
536  "Element assembly does not support AddBoundaryIntegrator yet.");
537 
538  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
539  const int intFaceIntegratorCount = intFaceIntegrators.Size();
540  if (intFaceIntegratorCount>0)
541  {
544  ea_data_ext.SetSize(2*nf_int*faceDofs*faceDofs, Device::GetMemoryType());
545  }
546  for (int i = 0; i < intFaceIntegratorCount; ++i)
547  {
548  intFaceIntegrators[i]->AssembleEAInteriorFaces(*a->FESpace(),
549  ea_data_int,
550  ea_data_ext,
551  i);
552  }
553 
554  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
555  const int boundFaceIntegratorCount = bdrFaceIntegrators.Size();
556  if (boundFaceIntegratorCount>0)
557  {
560  ea_data_bdr = 0.0;
561  }
562  for (int i = 0; i < boundFaceIntegratorCount; ++i)
563  {
564  bdrFaceIntegrators[i]->AssembleEABoundaryFaces(*a->FESpace(),ea_data_bdr,i);
565  }
566 
568  {
569  auto restFint = dynamic_cast<const L2FaceRestriction*>(int_face_restrict_lex);
571  }
573  {
574  auto restFbdr = dynamic_cast<const L2FaceRestriction*>(bdr_face_restrict_lex);
576  }
577 }
578 
580 {
581  // Apply the Element Restriction
582  const bool useRestrict = !DeviceCanUseCeed() && elem_restrict;
583  if (!useRestrict)
584  {
585  y.UseDevice(true); // typically this is a large vector, so store on device
586  y = 0.0;
587  }
588  else
589  {
591  localY = 0.0;
592  }
593  // Apply the Element Matrices
594  {
595  const int NDOFS = elemDofs;
596  auto X = Reshape(useRestrict?localX.Read():x.Read(), NDOFS, ne);
597  auto Y = Reshape(useRestrict?localY.ReadWrite():y.ReadWrite(), NDOFS, ne);
598  auto A = Reshape(ea_data.Read(), NDOFS, NDOFS, ne);
599  MFEM_FORALL(glob_j, ne*NDOFS,
600  {
601  const int e = glob_j/NDOFS;
602  const int j = glob_j%NDOFS;
603  double res = 0.0;
604  for (int i = 0; i < NDOFS; i++)
605  {
606  res += A(i, j, e)*X(i, e);
607  }
608  Y(j, e) += res;
609  });
610  // Apply the Element Restriction transposed
611  if (useRestrict)
612  {
614  }
615  }
616 
617  // Treatment of interior faces
618  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
619  const int iFISz = intFaceIntegrators.Size();
620  if (int_face_restrict_lex && iFISz>0)
621  {
622  // Apply the Interior Face Restriction
624  if (int_face_X.Size()>0)
625  {
626  int_face_Y = 0.0;
627  // Apply the interior face matrices
628  const int NDOFS = faceDofs;
629  auto X = Reshape(int_face_X.Read(), NDOFS, 2, nf_int);
630  auto Y = Reshape(int_face_Y.ReadWrite(), NDOFS, 2, nf_int);
632  {
633  auto A_int = Reshape(ea_data_int.Read(), NDOFS, NDOFS, 2, nf_int);
634  MFEM_FORALL(glob_j, nf_int*NDOFS,
635  {
636  const int f = glob_j/NDOFS;
637  const int j = glob_j%NDOFS;
638  double res = 0.0;
639  for (int i = 0; i < NDOFS; i++)
640  {
641  res += A_int(i, j, 0, f)*X(i, 0, f);
642  }
643  Y(j, 0, f) += res;
644  res = 0.0;
645  for (int i = 0; i < NDOFS; i++)
646  {
647  res += A_int(i, j, 1, f)*X(i, 1, f);
648  }
649  Y(j, 1, f) += res;
650  });
651  }
652  auto A_ext = Reshape(ea_data_ext.Read(), NDOFS, NDOFS, 2, nf_int);
653  MFEM_FORALL(glob_j, nf_int*NDOFS,
654  {
655  const int f = glob_j/NDOFS;
656  const int j = glob_j%NDOFS;
657  double res = 0.0;
658  for (int i = 0; i < NDOFS; i++)
659  {
660  res += A_ext(i, j, 0, f)*X(i, 0, f);
661  }
662  Y(j, 1, f) += res;
663  res = 0.0;
664  for (int i = 0; i < NDOFS; i++)
665  {
666  res += A_ext(i, j, 1, f)*X(i, 1, f);
667  }
668  Y(j, 0, f) += res;
669  });
670  // Apply the Interior Face Restriction transposed
672  }
673  }
674 
675  // Treatment of boundary faces
676  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
677  const int bFISz = bdrFaceIntegrators.Size();
678  if (!factorize_face_terms && bdr_face_restrict_lex && bFISz>0)
679  {
680  // Apply the Boundary Face Restriction
682  if (bdr_face_X.Size()>0)
683  {
684  bdr_face_Y = 0.0;
685  // Apply the boundary face matrices
686  const int NDOFS = faceDofs;
687  auto X = Reshape(bdr_face_X.Read(), NDOFS, nf_bdr);
688  auto Y = Reshape(bdr_face_Y.ReadWrite(), NDOFS, nf_bdr);
689  auto A = Reshape(ea_data_bdr.Read(), NDOFS, NDOFS, nf_bdr);
690  MFEM_FORALL(glob_j, nf_bdr*NDOFS,
691  {
692  const int f = glob_j/NDOFS;
693  const int j = glob_j%NDOFS;
694  double res = 0.0;
695  for (int i = 0; i < NDOFS; i++)
696  {
697  res += A(i, j, f)*X(i, f);
698  }
699  Y(j, f) += res;
700  });
701  // Apply the Boundary Face Restriction transposed
703  }
704  }
705 }
706 
708 {
709  // Apply the Element Restriction
710  const bool useRestrict = !DeviceCanUseCeed() && elem_restrict;
711  if (!useRestrict)
712  {
713  y.UseDevice(true); // typically this is a large vector, so store on device
714  y = 0.0;
715  }
716  else
717  {
719  localY = 0.0;
720  }
721  // Apply the Element Matrices transposed
722  {
723  const int NDOFS = elemDofs;
724  auto X = Reshape(useRestrict?localX.Read():x.Read(), NDOFS, ne);
725  auto Y = Reshape(useRestrict?localY.ReadWrite():y.ReadWrite(), NDOFS, ne);
726  auto A = Reshape(ea_data.Read(), NDOFS, NDOFS, ne);
727  MFEM_FORALL(glob_j, ne*NDOFS,
728  {
729  const int e = glob_j/NDOFS;
730  const int j = glob_j%NDOFS;
731  double res = 0.0;
732  for (int i = 0; i < NDOFS; i++)
733  {
734  res += A(j, i, e)*X(i, e);
735  }
736  Y(j, e) += res;
737  });
738  // Apply the Element Restriction transposed
739  if (useRestrict)
740  {
742  }
743  }
744 
745  // Treatment of interior faces
746  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
747  const int iFISz = intFaceIntegrators.Size();
748  if (int_face_restrict_lex && iFISz>0)
749  {
750  // Apply the Interior Face Restriction
752  if (int_face_X.Size()>0)
753  {
754  int_face_Y = 0.0;
755  // Apply the interior face matrices transposed
756  const int NDOFS = faceDofs;
757  auto X = Reshape(int_face_X.Read(), NDOFS, 2, nf_int);
758  auto Y = Reshape(int_face_Y.ReadWrite(), NDOFS, 2, nf_int);
760  {
761  auto A_int = Reshape(ea_data_int.Read(), NDOFS, NDOFS, 2, nf_int);
762  MFEM_FORALL(glob_j, nf_int*NDOFS,
763  {
764  const int f = glob_j/NDOFS;
765  const int j = glob_j%NDOFS;
766  double res = 0.0;
767  for (int i = 0; i < NDOFS; i++)
768  {
769  res += A_int(j, i, 0, f)*X(i, 0, f);
770  }
771  Y(j, 0, f) += res;
772  res = 0.0;
773  for (int i = 0; i < NDOFS; i++)
774  {
775  res += A_int(j, i, 1, f)*X(i, 1, f);
776  }
777  Y(j, 1, f) += res;
778  });
779  }
780  auto A_ext = Reshape(ea_data_ext.Read(), NDOFS, NDOFS, 2, nf_int);
781  MFEM_FORALL(glob_j, nf_int*NDOFS,
782  {
783  const int f = glob_j/NDOFS;
784  const int j = glob_j%NDOFS;
785  double res = 0.0;
786  for (int i = 0; i < NDOFS; i++)
787  {
788  res += A_ext(j, i, 1, f)*X(i, 0, f);
789  }
790  Y(j, 1, f) += res;
791  res = 0.0;
792  for (int i = 0; i < NDOFS; i++)
793  {
794  res += A_ext(j, i, 0, f)*X(i, 1, f);
795  }
796  Y(j, 0, f) += res;
797  });
798  // Apply the Interior Face Restriction transposed
800  }
801  }
802 
803  // Treatment of boundary faces
804  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
805  const int bFISz = bdrFaceIntegrators.Size();
806  if (!factorize_face_terms && bdr_face_restrict_lex && bFISz>0)
807  {
808  // Apply the Boundary Face Restriction
810  if (bdr_face_X.Size()>0)
811  {
812  bdr_face_Y = 0.0;
813  // Apply the boundary face matrices transposed
814  const int NDOFS = faceDofs;
815  auto X = Reshape(bdr_face_X.Read(), NDOFS, nf_bdr);
816  auto Y = Reshape(bdr_face_Y.ReadWrite(), NDOFS, nf_bdr);
817  auto A = Reshape(ea_data_bdr.Read(), NDOFS, NDOFS, nf_bdr);
818  MFEM_FORALL(glob_j, nf_bdr*NDOFS,
819  {
820  const int f = glob_j/NDOFS;
821  const int j = glob_j%NDOFS;
822  double res = 0.0;
823  for (int i = 0; i < NDOFS; i++)
824  {
825  res += A(j, i, f)*X(i, f);
826  }
827  Y(j, f) += res;
828  });
829  // Apply the Boundary Face Restriction transposed
831  }
832  }
833 }
834 
835 // Data and methods for fully-assembled bilinear forms
837  : EABilinearFormExtension(form),
838  mat(a->mat)
839 {
840 #ifdef MFEM_USE_MPI
841  ParFiniteElementSpace *pfes = nullptr;
842  if ( a->GetFBFI()->Size()>0 &&
843  (pfes = dynamic_cast<ParFiniteElementSpace*>(form->FESpace())) )
844  {
845  pfes->ExchangeFaceNbrData();
846  }
847 #endif
848 }
849 
851 {
853  FiniteElementSpace &fes = *a->FESpace();
854  int width = fes.GetVSize();
855  int height = fes.GetVSize();
856  bool keep_nbr_block = false;
857 #ifdef MFEM_USE_MPI
858  ParFiniteElementSpace *pfes = nullptr;
859  if ( a->GetFBFI()->Size()>0 &&
860  (pfes = dynamic_cast<ParFiniteElementSpace*>(&fes)) )
861  {
862  pfes->ExchangeFaceNbrData();
863  width += pfes->GetFaceNbrVSize();
864  dg_x.SetSize(width);
865  ParBilinearForm *pb = nullptr;
866  if ((pb = dynamic_cast<ParBilinearForm*>(a)) && (pb->keep_nbr_block))
867  {
868  height += pfes->GetFaceNbrVSize();
869  dg_y.SetSize(height);
870  keep_nbr_block = true;
871  }
872  }
873 #endif
874  if (a->mat) // We reuse the sparse matrix memory
875  {
876  if (fes.IsDGSpace())
877  {
878  const L2ElementRestriction *restE =
879  static_cast<const L2ElementRestriction*>(elem_restrict);
880  const L2FaceRestriction *restF =
881  static_cast<const L2FaceRestriction*>(int_face_restrict_lex);
882  MFEM_VERIFY(
883  fes.Conforming(),
884  "Full Assembly not yet supported on NCMesh.");
885  // 1. Fill J and Data
886  // 1.1 Fill J and Data with Elem ea_data
887  restE->FillJAndData(ea_data, *mat);
888  // 1.2 Fill J and Data with Face ea_data_ext
889  if (restF) { restF->FillJAndData(ea_data_ext, *mat, keep_nbr_block); }
890  // 1.3 Shift indirections in I back to original
891  auto I = mat->HostReadWriteI();
892  for (int i = height; i > 0; i--)
893  {
894  I[i] = I[i-1];
895  }
896  I[0] = 0;
897  }
898  else
899  {
900  const ElementRestriction &rest =
901  static_cast<const ElementRestriction&>(*elem_restrict);
902  rest.FillJAndData(ea_data, *mat);
903  }
904  }
905  else // We create, compute the sparsity, and fill the sparse matrix
906  {
907  mat = new SparseMatrix;
908  mat->OverrideSize(height, width);
909  if (fes.IsDGSpace())
910  {
911  const L2ElementRestriction *restE =
912  static_cast<const L2ElementRestriction*>(elem_restrict);
913  const L2FaceRestriction *restF =
914  static_cast<const L2FaceRestriction*>(int_face_restrict_lex);
915  MFEM_VERIFY(
916  fes.Conforming(),
917  "Full Assembly not yet supported on NCMesh.");
918  // 1. Fill I
919  mat->GetMemoryI().New(height+1, mat->GetMemoryI().GetMemoryType());
920  // 1.1 Increment with restE
921  restE->FillI(*mat);
922  // 1.2 Increment with restF
923  if (restF) { restF->FillI(*mat, keep_nbr_block); }
924  // 1.3 Sum the non-zeros in I
925  auto h_I = mat->HostReadWriteI();
926  int cpt = 0;
927  for (int i = 0; i < height; i++)
928  {
929  const int nnz = h_I[i];
930  h_I[i] = cpt;
931  cpt += nnz;
932  }
933  const int nnz = cpt;
934  h_I[height] = nnz;
935  mat->GetMemoryJ().New(nnz, mat->GetMemoryJ().GetMemoryType());
936  mat->GetMemoryData().New(nnz, mat->GetMemoryData().GetMemoryType());
937  // 2. Fill J and Data
938  // 2.1 Fill J and Data with Elem ea_data
939  restE->FillJAndData(ea_data, *mat);
940  // 2.2 Fill J and Data with Face ea_data_ext
941  if (restF) { restF->FillJAndData(ea_data_ext, *mat, keep_nbr_block); }
942  // 2.3 Shift indirections in I back to original
943  auto I = mat->HostReadWriteI();
944  for (int i = height; i > 0; i--)
945  {
946  I[i] = I[i-1];
947  }
948  I[0] = 0;
949  }
950  else // continuous Galerkin case
951  {
952  const ElementRestriction &rest =
953  static_cast<const ElementRestriction&>(*elem_restrict);
954  rest.FillSparseMatrix(ea_data, *mat);
955  }
956  a->mat = mat;
957  }
958  if ( a->sort_sparse_matrix )
959  {
960  a->mat->SortColumnIndices();
961  }
962 }
963 
964 
966 {
967 #ifdef MFEM_USE_MPI
968  if ( auto pa = dynamic_cast<ParBilinearForm*>(a) )
969  {
970  pa->ParallelRAP(*pa->mat, A);
971  }
972  else
973 #endif
974  {
975  a->SerialRAP(A);
976  }
977 }
978 
980  OperatorHandle &A)
981 {
982  MFEM_VERIFY(a->diag_policy == DiagonalPolicy::DIAG_ONE,
983  "Only DiagonalPolicy::DIAG_ONE supported with"
984  " FABilinearFormExtension.");
985 #ifdef MFEM_USE_MPI
986  if ( dynamic_cast<ParBilinearForm*>(a) )
987  {
988  A.As<HypreParMatrix>()->EliminateBC(ess_dofs,
989  DiagonalPolicy::DIAG_ONE);
990  }
991  else
992 #endif
993  {
994  A.As<SparseMatrix>()->EliminateBC(ess_dofs,
995  DiagonalPolicy::DIAG_ONE);
996  }
997 }
998 
1000  OperatorHandle &A)
1001 {
1002  RAP(A);
1003  EliminateBC(ess_dofs, A);
1004 }
1005 
1007  Vector &x, Vector &b,
1008  OperatorHandle &A,
1009  Vector &X, Vector &B,
1010  int copy_interior)
1011 {
1012  Operator *A_out;
1013  Operator::FormLinearSystem(ess_tdof_list, x, b, A_out, X, B, copy_interior);
1014  delete A_out;
1015  FormSystemMatrix(ess_tdof_list, A);
1016 }
1017 
1019 {
1020 #ifdef MFEM_USE_MPI
1021  const ParFiniteElementSpace *pfes;
1022  if ( (pfes = dynamic_cast<const ParFiniteElementSpace*>(test_fes)) )
1023  {
1024  // DG Prolongation
1025  ParGridFunction x_gf;
1026  x_gf.MakeRef(const_cast<ParFiniteElementSpace*>(pfes),
1027  const_cast<Vector&>(x),0);
1028  x_gf.ExchangeFaceNbrData();
1029  Vector &shared_x = x_gf.FaceNbrData();
1030  const int local_size = a->FESpace()->GetVSize();
1031  auto dg_x_ptr = dg_x.Write();
1032  auto x_ptr = x.Read();
1033  MFEM_FORALL(i,local_size,
1034  {
1035  dg_x_ptr[i] = x_ptr[i];
1036  });
1037  const int shared_size = shared_x.Size();
1038  auto shared_x_ptr = shared_x.Read();
1039  MFEM_FORALL(i,shared_size,
1040  {
1041  dg_x_ptr[local_size+i] = shared_x_ptr[i];
1042  });
1043  ParBilinearForm *pform = nullptr;
1044  if ((pform = dynamic_cast<ParBilinearForm*>(a)) && (pform->keep_nbr_block))
1045  {
1046  mat->Mult(dg_x, dg_y);
1047  // DG Restriction
1048  auto dg_y_ptr = dg_y.Read();
1049  auto y_ptr = y.ReadWrite();
1050  MFEM_FORALL(i,local_size,
1051  {
1052  y_ptr[i] += dg_y_ptr[i];
1053  });
1054  }
1055  else
1056  {
1057  mat->Mult(dg_x, y);
1058  }
1059  }
1060  else
1061 #endif
1062  {
1063  mat->Mult(x, y);
1064  }
1065 }
1066 
1068 {
1069  if ( a->GetFBFI()->Size()>0 )
1070  {
1071  DGMult(x, y);
1072  }
1073  else
1074  {
1075  mat->Mult(x, y);
1076  }
1077 }
1078 
1080 {
1081 #ifdef MFEM_USE_MPI
1082  const ParFiniteElementSpace *pfes;
1083  if ( (pfes = dynamic_cast<const ParFiniteElementSpace*>(test_fes)) )
1084  {
1085  // DG Prolongation
1086  ParGridFunction x_gf;
1087  x_gf.MakeRef(const_cast<ParFiniteElementSpace*>(pfes),
1088  const_cast<Vector&>(x),0);
1089  x_gf.ExchangeFaceNbrData();
1090  Vector &shared_x = x_gf.FaceNbrData();
1091  const int local_size = a->FESpace()->GetVSize();
1092  auto dg_x_ptr = dg_x.Write();
1093  auto x_ptr = x.Read();
1094  MFEM_FORALL(i,local_size,
1095  {
1096  dg_x_ptr[i] = x_ptr[i];
1097  });
1098  const int shared_size = shared_x.Size();
1099  auto shared_x_ptr = shared_x.Read();
1100  MFEM_FORALL(i,shared_size,
1101  {
1102  dg_x_ptr[local_size+i] = shared_x_ptr[i];
1103  });
1104  ParBilinearForm *pb = nullptr;
1105  if ((pb = dynamic_cast<ParBilinearForm*>(a)) && (pb->keep_nbr_block))
1106  {
1107  mat->MultTranspose(dg_x, dg_y);
1108  // DG Restriction
1109  auto dg_y_ptr = dg_y.Read();
1110  auto y_ptr = y.ReadWrite();
1111  MFEM_FORALL(i,local_size,
1112  {
1113  y_ptr[i] += dg_y_ptr[i];
1114  });
1115  }
1116  else
1117  {
1118  mat->MultTranspose(dg_x, y);
1119  }
1120  }
1121  else
1122 #endif
1123  {
1124  mat->MultTranspose(x, y);
1125  }
1126 }
1127 
1129 {
1130  if ( a->GetFBFI()->Size()>0 )
1131  {
1132  DGMultTranspose(x, y);
1133  }
1134  else
1135  {
1136  mat->MultTranspose(x, y);
1137  }
1138 }
1139 
1140 
1142  : Operator(form->Height(), form->Width()), a(form)
1143 {
1144  // empty
1145 }
1146 
1148 {
1149  return a->GetProlongation();
1150 }
1151 
1153 {
1154  return a->GetRestriction();
1155 }
1156 
1158 {
1159  return a->GetOutputProlongation();
1160 }
1161 
1163 {
1164  return a->GetOutputRestriction();
1165 }
1166 
1167 // Data and methods for partially-assembled bilinear forms
1168 
1170  MixedBilinearForm *form)
1172  trial_fes(form->TrialFESpace()),
1173  test_fes(form->TestFESpace()),
1174  elem_restrict_trial(NULL),
1175  elem_restrict_test(NULL)
1176 {
1177  Update();
1178 }
1179 
1181 {
1182  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1183  const int integratorCount = integrators.Size();
1184  for (int i = 0; i < integratorCount; ++i)
1185  {
1186  integrators[i]->AssemblePA(*trial_fes, *test_fes);
1187  }
1188  MFEM_VERIFY(a->GetBBFI()->Size() == 0,
1189  "Partial assembly does not support AddBoundaryIntegrator yet.");
1190  MFEM_VERIFY(a->GetTFBFI()->Size() == 0,
1191  "Partial assembly does not support AddTraceFaceIntegrator yet.");
1192  MFEM_VERIFY(a->GetBTFBFI()->Size() == 0,
1193  "Partial assembly does not support AddBdrTraceFaceIntegrator yet.");
1194 }
1195 
1197 {
1198  trial_fes = a->TrialFESpace();
1199  test_fes = a->TestFESpace();
1200  height = test_fes->GetVSize();
1201  width = trial_fes->GetVSize();
1206  if (elem_restrict_trial)
1207  {
1208  localTrial.UseDevice(true);
1211  }
1212  if (elem_restrict_test)
1213  {
1214  localTest.UseDevice(true); // ensure 'localY = 0.0' is done on device
1216  }
1217 }
1218 
1220  const Array<int> &trial_tdof_list,
1221  const Array<int> &test_tdof_list,
1222  OperatorHandle &A)
1223 {
1224  Operator * oper;
1225  Operator::FormRectangularSystemOperator(trial_tdof_list, test_tdof_list,
1226  oper);
1227  A.Reset(oper); // A will own oper
1228 }
1229 
1231  const Array<int> &trial_tdof_list,
1232  const Array<int> &test_tdof_list,
1233  Vector &x, Vector &b,
1234  OperatorHandle &A,
1235  Vector &X, Vector &B)
1236 {
1237  Operator *oper;
1238  Operator::FormRectangularLinearSystem(trial_tdof_list, test_tdof_list, x, b,
1239  oper, X, B);
1240  A.Reset(oper); // A will own oper
1241 }
1242 
1244  const Operator *elem_restrict_x,
1245  const Vector &x,
1246  Vector &localX,
1247  const Operator *elem_restrict_y,
1248  Vector &y,
1249  Vector &localY,
1250  const double c) const
1251 {
1252  // * G operation: localX = c*local(x)
1253  if (elem_restrict_x)
1254  {
1255  elem_restrict_x->Mult(x, localX);
1256  if (c != 1.0)
1257  {
1258  localX *= c;
1259  }
1260  }
1261  else
1262  {
1263  if (c == 1.0)
1264  {
1265  localX.SyncAliasMemory(x);
1266  }
1267  else
1268  {
1269  localX.Set(c, x);
1270  }
1271  }
1272  if (elem_restrict_y)
1273  {
1274  localY = 0.0;
1275  }
1276  else
1277  {
1278  y.UseDevice(true);
1279  localY.SyncAliasMemory(y);
1280  }
1281 }
1282 
1284 {
1285  y = 0.0;
1286  AddMult(x, y);
1287 }
1288 
1290  const double c) const
1291 {
1292  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1293  const int iSz = integrators.Size();
1294 
1295  // * G operation
1297  elem_restrict_test, y, localTest, c);
1298 
1299  // * B^TDB operation
1300  for (int i = 0; i < iSz; ++i)
1301  {
1302  integrators[i]->AddMultPA(localTrial, localTest);
1303  }
1304 
1305  // * G^T operation
1306  if (elem_restrict_test)
1307  {
1308  tempY.SetSize(y.Size());
1310  y += tempY;
1311  }
1312 }
1313 
1315  Vector &y) const
1316 {
1317  y = 0.0;
1318  AddMultTranspose(x, y);
1319 }
1320 
1322  const double c) const
1323 {
1324  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1325  const int iSz = integrators.Size();
1326 
1327  // * G operation
1330 
1331  // * B^TD^TB operation
1332  for (int i = 0; i < iSz; ++i)
1333  {
1334  integrators[i]->AddMultTransposePA(localTest, localTrial);
1335  }
1336 
1337  // * G^T operation
1338  if (elem_restrict_trial)
1339  {
1340  tempY.SetSize(y.Size());
1342  y += tempY;
1343  }
1344 }
1345 
1347  Vector &diag) const
1348 {
1349  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1350 
1351  const int iSz = integrators.Size();
1352 
1353  if (elem_restrict_trial)
1354  {
1355  const ElementRestriction* H1elem_restrict_trial =
1356  dynamic_cast<const ElementRestriction*>(elem_restrict_trial);
1357  if (H1elem_restrict_trial)
1358  {
1359  H1elem_restrict_trial->MultUnsigned(D, localTrial);
1360  }
1361  else
1362  {
1364  }
1365  }
1366 
1367  if (elem_restrict_test)
1368  {
1369  localTest = 0.0;
1370  for (int i = 0; i < iSz; ++i)
1371  {
1372  if (elem_restrict_trial)
1373  {
1374  integrators[i]->AssembleDiagonalPA_ADAt(localTrial, localTest);
1375  }
1376  else
1377  {
1378  integrators[i]->AssembleDiagonalPA_ADAt(D, localTest);
1379  }
1380  }
1381  const ElementRestriction* H1elem_restrict_test =
1382  dynamic_cast<const ElementRestriction*>(elem_restrict_test);
1383  if (H1elem_restrict_test)
1384  {
1385  H1elem_restrict_test->MultTransposeUnsigned(localTest, diag);
1386  }
1387  else
1388  {
1390  }
1391  }
1392  else
1393  {
1394  diag.UseDevice(true); // typically this is a large vector, so store on device
1395  diag = 0.0;
1396  for (int i = 0; i < iSz; ++i)
1397  {
1398  if (elem_restrict_trial)
1399  {
1400  integrators[i]->AssembleDiagonalPA_ADAt(localTrial, diag);
1401  }
1402  else
1403  {
1404  integrators[i]->AssembleDiagonalPA_ADAt(D, diag);
1405  }
1406  }
1407  }
1408 }
1409 
1411  DiscreteLinearOperator *linop) :
1413 {
1414 }
1415 
1416 const
1418 const
1419 {
1420  return a->GetOutputRestrictionTranspose();
1421 }
1422 
1424 {
1425  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1426  const int integratorCount = integrators.Size();
1427  for (int i = 0; i < integratorCount; ++i)
1428  {
1429  integrators[i]->AssemblePA(*trial_fes, *test_fes);
1430  }
1431 
1432  test_multiplicity.UseDevice(true);
1433  test_multiplicity.SetSize(elem_restrict_test->Width()); // l-vector
1434  Vector ones(elem_restrict_test->Height()); // e-vector
1435  ones = 1.0;
1436 
1437  const ElementRestriction* elem_restrict =
1438  dynamic_cast<const ElementRestriction*>(elem_restrict_test);
1439  if (elem_restrict)
1440  {
1441  elem_restrict->MultTransposeUnsigned(ones, test_multiplicity);
1442  }
1443  else
1444  {
1445  mfem_error("A real ElementRestriction is required in this setting!");
1446  }
1447 
1448  auto tm = test_multiplicity.ReadWrite();
1449  MFEM_FORALL(i, test_multiplicity.Size(),
1450  {
1451  tm[i] = 1.0 / tm[i];
1452  });
1453 }
1454 
1456  const Vector &x, Vector &y, const double c) const
1457 {
1458  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1459  const int iSz = integrators.Size();
1460 
1461  // * G operation
1463  elem_restrict_test, y, localTest, c);
1464 
1465  // * B^TDB operation
1466  for (int i = 0; i < iSz; ++i)
1467  {
1468  integrators[i]->AddMultPA(localTrial, localTest);
1469  }
1470 
1471  // do a kind of "set" rather than "add" in the below
1472  // operation as compared to the BilinearForm case
1473  // * G^T operation (kind of...)
1474  const ElementRestriction* elem_restrict =
1475  dynamic_cast<const ElementRestriction*>(elem_restrict_test);
1476  if (elem_restrict)
1477  {
1478  tempY.SetSize(y.Size());
1479  elem_restrict->MultLeftInverse(localTest, tempY);
1480  y += tempY;
1481  }
1482  else
1483  {
1484  mfem_error("In this setting you need a real ElementRestriction!");
1485  }
1486 }
1487 
1489  const Vector &x, Vector &y, const double c) const
1490 {
1491  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1492  const int iSz = integrators.Size();
1493 
1494  // do a kind of "set" rather than "add" in the below
1495  // operation as compared to the BilinearForm case
1496  // * G operation (kinda)
1497  Vector xscaled(x);
1498  MFEM_VERIFY(x.Size() == test_multiplicity.Size(), "Input vector of wrong size");
1499  auto xs = xscaled.ReadWrite();
1500  auto tm = test_multiplicity.Read();
1501  MFEM_FORALL(i, x.Size(),
1502  {
1503  xs[i] *= tm[i];
1504  });
1507 
1508  // * B^TD^TB operation
1509  for (int i = 0; i < iSz; ++i)
1510  {
1511  integrators[i]->AddMultTransposePA(localTest, localTrial);
1512  }
1513 
1514  // * G^T operation
1515  if (elem_restrict_trial)
1516  {
1517  tempY.SetSize(y.Size());
1519  y += tempY;
1520  }
1521  else
1522  {
1523  mfem_error("Trial ElementRestriction not defined");
1524  }
1525 }
1526 
1528  const Array<int>& ess1, const Array<int>& ess2, OperatorHandle &A)
1529 {
1530  const Operator *Pi = this->GetProlongation();
1531  const Operator *RoT = this->GetOutputRestrictionTranspose();
1532  Operator *rap = SetupRAP(Pi, RoT);
1533 
1535  = new RectangularConstrainedOperator(rap, ess1, ess2, rap != this);
1536 
1537  A.Reset(Arco);
1538 }
1539 
1540 } // namespace mfem
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...
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition: fespace.hpp:631
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
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:172
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:587
void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const
Assemble the diagonal of ADA^T for a diagonal vector D.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
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)
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 MultTransposeUnsigned(const Vector &x, Vector &y) const
Compute MultTranspose without applying signs based on DOF orientations.
void Mult(const Vector &x, Vector &y) const override=0
Extract the face degrees of freedom from x into y.
const FiniteElementSpace * test_fes
virtual const Operator * GetProlongation() const
Get the input finite element space prolongation matrix.
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
Array< BilinearFormIntegrator * > * GetDBFI()
Access all the integrators added with AddDomainIntegrator().
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
void DGMultTranspose(const Vector &x, Vector &y) const
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.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:73
void AddMultTranspose(const Vector &x, Vector &y, const double c=1.0) const
y += c*A^T*x
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
PAMixedBilinearFormExtension(MixedBilinearForm *form)
void MultTranspose(const Vector &x, Vector &y) const
y = A^T*x
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:923
void MultLeftInverse(const Vector &x, Vector &y) const
virtual const Operator * GetRestriction() const
Get the finite element space restriction matrix.
virtual const Operator * GetRestriction() const
Get the finite element space restriction matrix.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
Geometry::Type GetFaceBaseGeometry(int i) const
Definition: mesh.hpp:1066
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual const Operator * GetProlongation() const
Get the finite element space prolongation operator.
BilinearFormExtension(BilinearForm *form)
void SyncAliasMemory(const Vector &v) const
Update the alias memory location of the vector to match v.
Definition: vector.hpp:245
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Definition: fespace.hpp:957
Array< BilinearFormIntegrator * > * GetFBFI()
Access all integrators added with AddInteriorFaceIntegrator().
virtual const Operator * GetRestriction() const
Get the input finite element space restriction matrix.
void SortColumnIndices()
Sort the column indices corresponding to each row.
Definition: sparsemat.cpp:457
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:118
Memory< double > & GetMemoryData()
Definition: sparsemat.hpp:269
Memory< int > & GetMemoryI()
Definition: sparsemat.hpp:237
Data and methods for partially-assembled bilinear forms.
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...
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:94
Class extending the MixedBilinearForm class to support different AssemblyLevels.
Memory< int > & GetMemoryJ()
Definition: sparsemat.hpp:253
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:68
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:102
const Operator * GetOutputRestrictionTranspose() const
Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type...
virtual const Operator * GetOutputProlongation() const
Get the output 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 ...
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
PADiscreteLinearOperatorExtension(DiscreteLinearOperator *linop)
void Assemble()
Partial assembly of all internal integrators.
Data type sparse matrix.
Definition: sparsemat.hpp:50
Native ordering as defined by the FiniteElement.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:67
void FillI(SparseMatrix &mat) const
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:441
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:457
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...
void AddMult(const Vector &x, Vector &y, const double c=1.0) const
y += c*A*x
FiniteElementSpace * TestFESpace()
Return the test FE space associated with the BilinearForm.
virtual const Operator * GetOutputRestriction() const
Get the test finite element space restriction matrix.
void EliminateBC(const Array< int > &ess_dofs, OperatorHandle &A)
Array< BilinearFormIntegrator * > * GetBTFBFI()
Access all integrators added with AddBdrTraceFaceIntegrator().
Class extending the BilinearForm class to support different AssemblyLevels.
bool IsDGSpace() const
Return whether or not the space is discontinuous (L2)
Definition: fespace.hpp:925
void FormRectangularSystemOperator(const Array< int > &, const Array< int > &, OperatorHandle &A)
Setup OperatorHandle A to contain constrained linear operator.
void MultUnsigned(const Vector &x, Vector &y) const
Compute Mult without applying signs based on DOF orientations.
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
bool Conforming() const
Definition: fespace.hpp:447
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)
void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
virtual const FaceRestriction * GetFaceRestriction(ElementDofOrdering e_ordering, FaceType, L2FaceValues mul=L2FaceValues::DoubleValued) const
Return an Operator that converts L-vectors to E-vectors on each face.
Definition: fespace.cpp:1292
virtual const Operator * GetOutputRestrictionTranspose() const
Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type...
Definition: operator.hpp:130
Array< BilinearFormIntegrator * > * GetTFBFI()
Access all integrators added with AddTraceFaceIntegrator().
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 ...
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 ...
Array< BilinearFormIntegrator * > * GetBBFI()
Access all integrators added with AddBoundaryIntegrator().
void Assemble()
Assemble at the level given for the BilinearFormExtension subclass.
Vector & FaceNbrData()
Definition: pgridfunc.hpp:208
void SetupRestrictionOperators(const L2FaceValues m)
void Mult(const Vector &x, Vector &y) const
y = A*x
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Data and methods for element-assembled bilinear forms.
virtual const Operator * GetRestriction() const
Get the finite element space restriction operator.
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Operator that converts FiniteElementSpace L-vectors to E-vectors.
Definition: restriction.hpp:36
const FiniteElementSpace * test_fes
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
void SerialRAP(OperatorHandle &A)
Compute serial RAP operator and store it in A as a SparseMatrix.
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:729
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 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:51
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
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:267
EABilinearFormExtension(BilinearForm *form)
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
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.
double a
Definition: lissajous.cpp:41
MixedBilinearForm * a
Not owned.
A &quot;square matrix&quot; 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:465
void DGMult(const Vector &x, Vector &y) const
Rectangular Operator for imposing essential boundary conditions on the input space using only the act...
Definition: operator.hpp:911
void AddMultTranspose(const Vector &x, Vector &y, const double c=1.0) const
y += c*A^T*x
virtual const Operator * GetOutputProlongation() const
Get the test finite element space prolongation matrix.
virtual const Operator * GetOutputRestriction() const
Get the output finite element space restriction matrix.
void RAP(OperatorHandle &A)
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:164
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 MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:904
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:2781
void OverrideSize(int height_, int width_)
Sets the height and width of the matrix.
Definition: sparsemat.cpp:297
Lexicographic ordering for tensor-product FiniteElements.
Class for parallel bilinear form.
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:1259
Operator * SetupRAP(const Operator *Pi, const Operator *Po)
Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to &quot;P&quot;...
Definition: operator.cpp:105
int GetFaceNbrVSize() const
Definition: pfespace.hpp:394
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 mixed bilinear forms.
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Vector data type.
Definition: vector.hpp:60
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...
virtual void AddFaceMatricesToElementMatrices(const Vector &fea_data, Vector &ea_data) const
This methods adds the DG face matrices to the element matrices.
Array< BilinearFormIntegrator * > * GetDBFI()
Access all integrators added with AddDomainIntegrator().
int * HostReadWriteI()
Definition: sparsemat.hpp:249
void AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
Class for parallel grid function.
Definition: pgridfunc.hpp:32
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
FiniteElementSpace * TrialFESpace()
Return the trial FE space associated with the BilinearForm.
const FiniteElementSpace * trial_fes
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:449
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
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...
void FillSparseMatrix(const Vector &mat_ea, SparseMatrix &mat) const
Fill a Sparse Matrix with Element Matrices.
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
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:145
void AddMult(const Vector &x, Vector &y, const double c) const
y += c*A*x