MFEM  v4.4.0
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/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();
66  if (elem_restrict && !DeviceCanUseCeed())
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();
131  if (DeviceCanUseCeed() || !elem_restrict)
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 {
258  if (elem_restrict)
259  {
260  localX.SetSize(elem_restrict->Height(), Device::GetDeviceMemoryType());
261  localY.SetSize(elem_restrict->Height(), Device::GetDeviceMemoryType());
262  localY.UseDevice(true); // ensure 'localY = 0.0' is done on device
263  }
264 
265  // Construct face restriction operators only if the bilinear form has
266  // interior or boundary face integrators
267  if (int_face_restrict_lex == NULL && a->GetFBFI()->Size() > 0)
268  {
274  int_face_Y.UseDevice(true); // ensure 'int_face_Y = 0.0' is done on device
275  }
276 
277  if (bdr_face_restrict_lex == NULL && a->GetBFBFI()->Size() > 0)
278  {
282  m);
285  bdr_face_Y.UseDevice(true); // ensure 'faceBoundY = 0.0' is done on device
286  }
287 }
288 
290 {
292 
293  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
294  const int integratorCount = integrators.Size();
295  for (int i = 0; i < integratorCount; ++i)
296  {
297  integrators[i]->AssemblePA(*a->FESpace());
298  }
299 
300  MFEM_VERIFY(a->GetBBFI()->Size() == 0,
301  "Partial assembly does not support AddBoundaryIntegrator yet.");
302 
303  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
304  const int intFaceIntegratorCount = intFaceIntegrators.Size();
305  for (int i = 0; i < intFaceIntegratorCount; ++i)
306  {
307  intFaceIntegrators[i]->AssemblePAInteriorFaces(*a->FESpace());
308  }
309 
310  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
311  const int boundFaceIntegratorCount = bdrFaceIntegrators.Size();
312  for (int i = 0; i < boundFaceIntegratorCount; ++i)
313  {
314  bdrFaceIntegrators[i]->AssemblePABoundaryFaces(*a->FESpace());
315  }
316 }
317 
319 {
320  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
321 
322  const int iSz = integrators.Size();
323  if (elem_restrict && !DeviceCanUseCeed())
324  {
325  localY = 0.0;
326  for (int i = 0; i < iSz; ++i)
327  {
328  integrators[i]->AssembleDiagonalPA(localY);
329  }
330  const ElementRestriction* H1elem_restrict =
331  dynamic_cast<const ElementRestriction*>(elem_restrict);
332  if (H1elem_restrict)
333  {
334  H1elem_restrict->MultTransposeUnsigned(localY, y);
335  }
336  else
337  {
339  }
340  }
341  else
342  {
343  y.UseDevice(true); // typically this is a large vector, so store on device
344  y = 0.0;
345  for (int i = 0; i < iSz; ++i)
346  {
347  integrators[i]->AssembleDiagonalPA(y);
348  }
349  }
350 }
351 
353 {
354  FiniteElementSpace *fes = a->FESpace();
355  height = width = fes->GetVSize();
356  trial_fes = fes;
357  test_fes = fes;
358 
359  elem_restrict = nullptr;
360  int_face_restrict_lex = nullptr;
361  bdr_face_restrict_lex = nullptr;
362 }
363 
365  OperatorHandle &A)
366 {
367  Operator *oper;
368  Operator::FormSystemOperator(ess_tdof_list, oper);
369  A.Reset(oper); // A will own oper
370 }
371 
373  Vector &x, Vector &b,
374  OperatorHandle &A,
375  Vector &X, Vector &B,
376  int copy_interior)
377 {
378  Operator *oper;
379  Operator::FormLinearSystem(ess_tdof_list, x, b, oper, X, B, copy_interior);
380  A.Reset(oper); // A will own oper
381 }
382 
384 {
385  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
386 
387  const int iSz = integrators.Size();
388  if (DeviceCanUseCeed() || !elem_restrict)
389  {
390  y.UseDevice(true); // typically this is a large vector, so store on device
391  y = 0.0;
392  for (int i = 0; i < iSz; ++i)
393  {
394  integrators[i]->AddMultPA(x, y);
395  }
396  }
397  else
398  {
400  localY = 0.0;
401  for (int i = 0; i < iSz; ++i)
402  {
403  integrators[i]->AddMultPA(localX, localY);
404  }
406  }
407 
408  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
409  const int iFISz = intFaceIntegrators.Size();
410  if (int_face_restrict_lex && iFISz>0)
411  {
413  if (int_face_X.Size()>0)
414  {
415  int_face_Y = 0.0;
416  for (int i = 0; i < iFISz; ++i)
417  {
418  intFaceIntegrators[i]->AddMultPA(int_face_X, int_face_Y);
419  }
421  }
422  }
423 
424  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
425  const int bFISz = bdrFaceIntegrators.Size();
426  if (bdr_face_restrict_lex && bFISz>0)
427  {
429  if (bdr_face_X.Size()>0)
430  {
431  bdr_face_Y = 0.0;
432  for (int i = 0; i < bFISz; ++i)
433  {
434  bdrFaceIntegrators[i]->AddMultPA(bdr_face_X, bdr_face_Y);
435  }
437  }
438  }
439 }
440 
442 {
443  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
444  const int iSz = integrators.Size();
445  if (elem_restrict)
446  {
448  localY = 0.0;
449  for (int i = 0; i < iSz; ++i)
450  {
451  integrators[i]->AddMultTransposePA(localX, localY);
452  }
454  }
455  else
456  {
457  y.UseDevice(true);
458  y = 0.0;
459  for (int i = 0; i < iSz; ++i)
460  {
461  integrators[i]->AddMultTransposePA(x, y);
462  }
463  }
464 
465  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
466  const int iFISz = intFaceIntegrators.Size();
467  if (int_face_restrict_lex && iFISz>0)
468  {
470  if (int_face_X.Size()>0)
471  {
472  int_face_Y = 0.0;
473  for (int i = 0; i < iFISz; ++i)
474  {
475  intFaceIntegrators[i]->AddMultTransposePA(int_face_X, int_face_Y);
476  }
478  }
479  }
480 
481  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
482  const int bFISz = bdrFaceIntegrators.Size();
483  if (bdr_face_restrict_lex && bFISz>0)
484  {
486  if (bdr_face_X.Size()>0)
487  {
488  bdr_face_Y = 0.0;
489  for (int i = 0; i < bFISz; ++i)
490  {
491  bdrFaceIntegrators[i]->AddMultTransposePA(bdr_face_X, bdr_face_Y);
492  }
494  }
495  }
496 }
497 
498 // Data and methods for element-assembled bilinear forms
500  : PABilinearFormExtension(form),
501  factorize_face_terms(false)
502 {
503  if (form->FESpace()->IsDGSpace() && form->FESpace()->Conforming())
504  {
505  factorize_face_terms = true;
506  }
507 }
508 
510 {
512 
513  ne = trial_fes->GetMesh()->GetNE();
514  elemDofs = trial_fes->GetFE(0)->GetDof();
515 
517  ea_data.UseDevice(true);
518 
519  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
520  const int integratorCount = integrators.Size();
521  if ( integratorCount == 0 )
522  {
523  ea_data = 0.0;
524  }
525  for (int i = 0; i < integratorCount; ++i)
526  {
527  integrators[i]->AssembleEA(*a->FESpace(), ea_data, i);
528  }
529 
530  faceDofs = trial_fes ->
531  GetTraceElement(0, trial_fes->GetMesh()->GetFaceBaseGeometry(0)) ->
532  GetDof();
533 
534  MFEM_VERIFY(a->GetBBFI()->Size() == 0,
535  "Element assembly does not support AddBoundaryIntegrator yet.");
536 
537  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
538  const int intFaceIntegratorCount = intFaceIntegrators.Size();
539  if (intFaceIntegratorCount>0)
540  {
543  ea_data_ext.SetSize(2*nf_int*faceDofs*faceDofs, Device::GetMemoryType());
544  }
545  for (int i = 0; i < intFaceIntegratorCount; ++i)
546  {
547  intFaceIntegrators[i]->AssembleEAInteriorFaces(*a->FESpace(),
548  ea_data_int,
549  ea_data_ext,
550  i);
551  }
552 
553  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
554  const int boundFaceIntegratorCount = bdrFaceIntegrators.Size();
555  if (boundFaceIntegratorCount>0)
556  {
559  ea_data_bdr = 0.0;
560  }
561  for (int i = 0; i < boundFaceIntegratorCount; ++i)
562  {
563  bdrFaceIntegrators[i]->AssembleEABoundaryFaces(*a->FESpace(),ea_data_bdr,i);
564  }
565 
567  {
568  auto restFint = dynamic_cast<const L2FaceRestriction*>(int_face_restrict_lex);
570  }
572  {
573  auto restFbdr = dynamic_cast<const L2FaceRestriction*>(bdr_face_restrict_lex);
575  }
576 }
577 
579 {
580  // Apply the Element Restriction
581  const bool useRestrict = !DeviceCanUseCeed() && elem_restrict;
582  if (!useRestrict)
583  {
584  y.UseDevice(true); // typically this is a large vector, so store on device
585  y = 0.0;
586  }
587  else
588  {
590  localY = 0.0;
591  }
592  // Apply the Element Matrices
593  {
594  const int NDOFS = elemDofs;
595  auto X = Reshape(useRestrict?localX.Read():x.Read(), NDOFS, ne);
596  auto Y = Reshape(useRestrict?localY.ReadWrite():y.ReadWrite(), NDOFS, ne);
597  auto A = Reshape(ea_data.Read(), NDOFS, NDOFS, ne);
598  MFEM_FORALL(glob_j, ne*NDOFS,
599  {
600  const int e = glob_j/NDOFS;
601  const int j = glob_j%NDOFS;
602  double res = 0.0;
603  for (int i = 0; i < NDOFS; i++)
604  {
605  res += A(i, j, e)*X(i, e);
606  }
607  Y(j, e) += res;
608  });
609  // Apply the Element Restriction transposed
610  if (useRestrict)
611  {
613  }
614  }
615 
616  // Treatment of interior faces
617  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
618  const int iFISz = intFaceIntegrators.Size();
619  if (int_face_restrict_lex && iFISz>0)
620  {
621  // Apply the Interior Face Restriction
623  if (int_face_X.Size()>0)
624  {
625  int_face_Y = 0.0;
626  // Apply the interior face matrices
627  const int NDOFS = faceDofs;
628  auto X = Reshape(int_face_X.Read(), NDOFS, 2, nf_int);
629  auto Y = Reshape(int_face_Y.ReadWrite(), NDOFS, 2, nf_int);
631  {
632  auto A_int = Reshape(ea_data_int.Read(), NDOFS, NDOFS, 2, nf_int);
633  MFEM_FORALL(glob_j, nf_int*NDOFS,
634  {
635  const int f = glob_j/NDOFS;
636  const int j = glob_j%NDOFS;
637  double res = 0.0;
638  for (int i = 0; i < NDOFS; i++)
639  {
640  res += A_int(i, j, 0, f)*X(i, 0, f);
641  }
642  Y(j, 0, f) += res;
643  res = 0.0;
644  for (int i = 0; i < NDOFS; i++)
645  {
646  res += A_int(i, j, 1, f)*X(i, 1, f);
647  }
648  Y(j, 1, f) += res;
649  });
650  }
651  auto A_ext = Reshape(ea_data_ext.Read(), NDOFS, NDOFS, 2, nf_int);
652  MFEM_FORALL(glob_j, nf_int*NDOFS,
653  {
654  const int f = glob_j/NDOFS;
655  const int j = glob_j%NDOFS;
656  double res = 0.0;
657  for (int i = 0; i < NDOFS; i++)
658  {
659  res += A_ext(i, j, 0, f)*X(i, 0, f);
660  }
661  Y(j, 1, f) += res;
662  res = 0.0;
663  for (int i = 0; i < NDOFS; i++)
664  {
665  res += A_ext(i, j, 1, f)*X(i, 1, f);
666  }
667  Y(j, 0, f) += res;
668  });
669  // Apply the Interior Face Restriction transposed
671  }
672  }
673 
674  // Treatment of boundary faces
675  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
676  const int bFISz = bdrFaceIntegrators.Size();
677  if (!factorize_face_terms && bdr_face_restrict_lex && bFISz>0)
678  {
679  // Apply the Boundary Face Restriction
681  if (bdr_face_X.Size()>0)
682  {
683  bdr_face_Y = 0.0;
684  // Apply the boundary face matrices
685  const int NDOFS = faceDofs;
686  auto X = Reshape(bdr_face_X.Read(), NDOFS, nf_bdr);
687  auto Y = Reshape(bdr_face_Y.ReadWrite(), NDOFS, nf_bdr);
688  auto A = Reshape(ea_data_bdr.Read(), NDOFS, NDOFS, nf_bdr);
689  MFEM_FORALL(glob_j, nf_bdr*NDOFS,
690  {
691  const int f = glob_j/NDOFS;
692  const int j = glob_j%NDOFS;
693  double res = 0.0;
694  for (int i = 0; i < NDOFS; i++)
695  {
696  res += A(i, j, f)*X(i, f);
697  }
698  Y(j, f) += res;
699  });
700  // Apply the Boundary Face Restriction transposed
702  }
703  }
704 }
705 
707 {
708  // Apply the Element Restriction
709  const bool useRestrict = !DeviceCanUseCeed() && elem_restrict;
710  if (!useRestrict)
711  {
712  y.UseDevice(true); // typically this is a large vector, so store on device
713  y = 0.0;
714  }
715  else
716  {
718  localY = 0.0;
719  }
720  // Apply the Element Matrices transposed
721  {
722  const int NDOFS = elemDofs;
723  auto X = Reshape(useRestrict?localX.Read():x.Read(), NDOFS, ne);
724  auto Y = Reshape(useRestrict?localY.ReadWrite():y.ReadWrite(), NDOFS, ne);
725  auto A = Reshape(ea_data.Read(), NDOFS, NDOFS, ne);
726  MFEM_FORALL(glob_j, ne*NDOFS,
727  {
728  const int e = glob_j/NDOFS;
729  const int j = glob_j%NDOFS;
730  double res = 0.0;
731  for (int i = 0; i < NDOFS; i++)
732  {
733  res += A(j, i, e)*X(i, e);
734  }
735  Y(j, e) += res;
736  });
737  // Apply the Element Restriction transposed
738  if (useRestrict)
739  {
741  }
742  }
743 
744  // Treatment of interior faces
745  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
746  const int iFISz = intFaceIntegrators.Size();
747  if (int_face_restrict_lex && iFISz>0)
748  {
749  // Apply the Interior Face Restriction
751  if (int_face_X.Size()>0)
752  {
753  int_face_Y = 0.0;
754  // Apply the interior face matrices transposed
755  const int NDOFS = faceDofs;
756  auto X = Reshape(int_face_X.Read(), NDOFS, 2, nf_int);
757  auto Y = Reshape(int_face_Y.ReadWrite(), NDOFS, 2, nf_int);
759  {
760  auto A_int = Reshape(ea_data_int.Read(), NDOFS, NDOFS, 2, nf_int);
761  MFEM_FORALL(glob_j, nf_int*NDOFS,
762  {
763  const int f = glob_j/NDOFS;
764  const int j = glob_j%NDOFS;
765  double res = 0.0;
766  for (int i = 0; i < NDOFS; i++)
767  {
768  res += A_int(j, i, 0, f)*X(i, 0, f);
769  }
770  Y(j, 0, f) += res;
771  res = 0.0;
772  for (int i = 0; i < NDOFS; i++)
773  {
774  res += A_int(j, i, 1, f)*X(i, 1, f);
775  }
776  Y(j, 1, f) += res;
777  });
778  }
779  auto A_ext = Reshape(ea_data_ext.Read(), NDOFS, NDOFS, 2, nf_int);
780  MFEM_FORALL(glob_j, nf_int*NDOFS,
781  {
782  const int f = glob_j/NDOFS;
783  const int j = glob_j%NDOFS;
784  double res = 0.0;
785  for (int i = 0; i < NDOFS; i++)
786  {
787  res += A_ext(j, i, 1, f)*X(i, 0, f);
788  }
789  Y(j, 1, f) += res;
790  res = 0.0;
791  for (int i = 0; i < NDOFS; i++)
792  {
793  res += A_ext(j, i, 0, f)*X(i, 1, f);
794  }
795  Y(j, 0, f) += res;
796  });
797  // Apply the Interior Face Restriction transposed
799  }
800  }
801 
802  // Treatment of boundary faces
803  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
804  const int bFISz = bdrFaceIntegrators.Size();
805  if (!factorize_face_terms && bdr_face_restrict_lex && bFISz>0)
806  {
807  // Apply the Boundary Face Restriction
809  if (bdr_face_X.Size()>0)
810  {
811  bdr_face_Y = 0.0;
812  // Apply the boundary face matrices transposed
813  const int NDOFS = faceDofs;
814  auto X = Reshape(bdr_face_X.Read(), NDOFS, nf_bdr);
815  auto Y = Reshape(bdr_face_Y.ReadWrite(), NDOFS, nf_bdr);
816  auto A = Reshape(ea_data_bdr.Read(), NDOFS, NDOFS, nf_bdr);
817  MFEM_FORALL(glob_j, nf_bdr*NDOFS,
818  {
819  const int f = glob_j/NDOFS;
820  const int j = glob_j%NDOFS;
821  double res = 0.0;
822  for (int i = 0; i < NDOFS; i++)
823  {
824  res += A(j, i, f)*X(i, f);
825  }
826  Y(j, f) += res;
827  });
828  // Apply the Boundary Face Restriction transposed
830  }
831  }
832 }
833 
834 // Data and methods for fully-assembled bilinear forms
836  : EABilinearFormExtension(form),
837  mat(a->mat)
838 {
839 #ifdef MFEM_USE_MPI
840  ParFiniteElementSpace *pfes = nullptr;
841  if ( a->GetFBFI()->Size()>0 &&
842  (pfes = dynamic_cast<ParFiniteElementSpace*>(form->FESpace())) )
843  {
844  pfes->ExchangeFaceNbrData();
845  }
846 #endif
847 }
848 
850 {
852  FiniteElementSpace &fes = *a->FESpace();
853  int width = fes.GetVSize();
854  int height = fes.GetVSize();
855  bool keep_nbr_block = false;
856 #ifdef MFEM_USE_MPI
857  ParFiniteElementSpace *pfes = nullptr;
858  if ( a->GetFBFI()->Size()>0 &&
859  (pfes = dynamic_cast<ParFiniteElementSpace*>(&fes)) )
860  {
861  pfes->ExchangeFaceNbrData();
862  width += pfes->GetFaceNbrVSize();
863  dg_x.SetSize(width);
864  ParBilinearForm *pb = nullptr;
865  if ((pb = dynamic_cast<ParBilinearForm*>(a)) && (pb->keep_nbr_block))
866  {
867  height += pfes->GetFaceNbrVSize();
868  dg_y.SetSize(height);
869  keep_nbr_block = true;
870  }
871  }
872 #endif
873  if (a->mat) // We reuse the sparse matrix memory
874  {
875  if (fes.IsDGSpace())
876  {
877  const L2ElementRestriction *restE =
878  static_cast<const L2ElementRestriction*>(elem_restrict);
879  const L2FaceRestriction *restF =
880  static_cast<const L2FaceRestriction*>(int_face_restrict_lex);
881  MFEM_VERIFY(
882  fes.Conforming(),
883  "Full Assembly not yet supported on NCMesh.");
884  // 1. Fill J and Data
885  // 1.1 Fill J and Data with Elem ea_data
886  restE->FillJAndData(ea_data, *mat);
887  // 1.2 Fill J and Data with Face ea_data_ext
888  if (restF) { restF->FillJAndData(ea_data_ext, *mat, keep_nbr_block); }
889  // 1.3 Shift indirections in I back to original
890  auto I = mat->HostReadWriteI();
891  for (int i = height; i > 0; i--)
892  {
893  I[i] = I[i-1];
894  }
895  I[0] = 0;
896  }
897  else
898  {
899  const ElementRestriction &rest =
900  static_cast<const ElementRestriction&>(*elem_restrict);
901  rest.FillJAndData(ea_data, *mat);
902  }
903  }
904  else // We create, compute the sparsity, and fill the sparse matrix
905  {
906  mat = new SparseMatrix(height, width, 0);
907  if (fes.IsDGSpace())
908  {
909  const L2ElementRestriction *restE =
910  static_cast<const L2ElementRestriction*>(elem_restrict);
911  const L2FaceRestriction *restF =
912  static_cast<const L2FaceRestriction*>(int_face_restrict_lex);
913  MFEM_VERIFY(
914  fes.Conforming(),
915  "Full Assembly not yet supported on NCMesh.");
916  // 1. Fill I
917  mat->GetMemoryI().New(height+1, mat->GetMemoryI().GetMemoryType());
918  // 1.1 Increment with restE
919  restE->FillI(*mat);
920  // 1.2 Increment with restF
921  if (restF) { restF->FillI(*mat, keep_nbr_block); }
922  // 1.3 Sum the non-zeros in I
923  auto h_I = mat->HostReadWriteI();
924  int cpt = 0;
925  for (int i = 0; i < height; i++)
926  {
927  const int nnz = h_I[i];
928  h_I[i] = cpt;
929  cpt += nnz;
930  }
931  const int nnz = cpt;
932  h_I[height] = nnz;
933  mat->GetMemoryJ().New(nnz, mat->GetMemoryJ().GetMemoryType());
934  mat->GetMemoryData().New(nnz, mat->GetMemoryData().GetMemoryType());
935  // 2. Fill J and Data
936  // 2.1 Fill J and Data with Elem ea_data
937  restE->FillJAndData(ea_data, *mat);
938  // 2.2 Fill J and Data with Face ea_data_ext
939  if (restF) { restF->FillJAndData(ea_data_ext, *mat, keep_nbr_block); }
940  // 2.3 Shift indirections in I back to original
941  auto I = mat->HostReadWriteI();
942  for (int i = height; i > 0; i--)
943  {
944  I[i] = I[i-1];
945  }
946  I[0] = 0;
947  }
948  else // continuous Galerkin case
949  {
950  const ElementRestriction &rest =
951  static_cast<const ElementRestriction&>(*elem_restrict);
952  rest.FillSparseMatrix(ea_data, *mat);
953  }
954  a->mat = mat;
955  }
956 }
957 
959 {
960 #ifdef MFEM_USE_MPI
961  const ParFiniteElementSpace *pfes;
962  if ( (pfes = dynamic_cast<const ParFiniteElementSpace*>(test_fes)) )
963  {
964  // DG Prolongation
965  ParGridFunction x_gf;
966  x_gf.MakeRef(const_cast<ParFiniteElementSpace*>(pfes),
967  const_cast<Vector&>(x),0);
968  x_gf.ExchangeFaceNbrData();
969  Vector &shared_x = x_gf.FaceNbrData();
970  const int local_size = a->FESpace()->GetVSize();
971  auto dg_x_ptr = dg_x.Write();
972  auto x_ptr = x.Read();
973  MFEM_FORALL(i,local_size,
974  {
975  dg_x_ptr[i] = x_ptr[i];
976  });
977  const int shared_size = shared_x.Size();
978  auto shared_x_ptr = shared_x.Read();
979  MFEM_FORALL(i,shared_size,
980  {
981  dg_x_ptr[local_size+i] = shared_x_ptr[i];
982  });
983  ParBilinearForm *pform = nullptr;
984  if ((pform = dynamic_cast<ParBilinearForm*>(a)) && (pform->keep_nbr_block))
985  {
986  mat->Mult(dg_x, dg_y);
987  // DG Restriction
988  auto dg_y_ptr = dg_y.Read();
989  auto y_ptr = y.ReadWrite();
990  MFEM_FORALL(i,local_size,
991  {
992  y_ptr[i] += dg_y_ptr[i];
993  });
994  }
995  else
996  {
997  mat->Mult(dg_x, y);
998  }
999  }
1000  else
1001 #endif
1002  {
1003  mat->Mult(x, y);
1004  }
1005 }
1006 
1008 {
1009  if ( a->GetFBFI()->Size()>0 )
1010  {
1011  DGMult(x, y);
1012  }
1013  else
1014  {
1015  mat->Mult(x, y);
1016  }
1017 }
1018 
1020 {
1021 #ifdef MFEM_USE_MPI
1022  const ParFiniteElementSpace *pfes;
1023  if ( (pfes = dynamic_cast<const ParFiniteElementSpace*>(test_fes)) )
1024  {
1025  // DG Prolongation
1026  ParGridFunction x_gf;
1027  x_gf.MakeRef(const_cast<ParFiniteElementSpace*>(pfes),
1028  const_cast<Vector&>(x),0);
1029  x_gf.ExchangeFaceNbrData();
1030  Vector &shared_x = x_gf.FaceNbrData();
1031  const int local_size = a->FESpace()->GetVSize();
1032  auto dg_x_ptr = dg_x.Write();
1033  auto x_ptr = x.Read();
1034  MFEM_FORALL(i,local_size,
1035  {
1036  dg_x_ptr[i] = x_ptr[i];
1037  });
1038  const int shared_size = shared_x.Size();
1039  auto shared_x_ptr = shared_x.Read();
1040  MFEM_FORALL(i,shared_size,
1041  {
1042  dg_x_ptr[local_size+i] = shared_x_ptr[i];
1043  });
1044  ParBilinearForm *pb = nullptr;
1045  if ((pb = dynamic_cast<ParBilinearForm*>(a)) && (pb->keep_nbr_block))
1046  {
1047  mat->MultTranspose(dg_x, dg_y);
1048  // DG Restriction
1049  auto dg_y_ptr = dg_y.Read();
1050  auto y_ptr = y.ReadWrite();
1051  MFEM_FORALL(i,local_size,
1052  {
1053  y_ptr[i] += dg_y_ptr[i];
1054  });
1055  }
1056  else
1057  {
1058  mat->MultTranspose(dg_x, y);
1059  }
1060  }
1061  else
1062 #endif
1063  {
1064  mat->MultTranspose(x, y);
1065  }
1066 }
1067 
1069 {
1070  if ( a->GetFBFI()->Size()>0 )
1071  {
1072  DGMultTranspose(x, y);
1073  }
1074  else
1075  {
1076  mat->MultTranspose(x, y);
1077  }
1078 }
1079 
1080 
1082  : Operator(form->Height(), form->Width()), a(form)
1083 {
1084  // empty
1085 }
1086 
1088 {
1089  return a->GetProlongation();
1090 }
1091 
1093 {
1094  return a->GetRestriction();
1095 }
1096 
1098 {
1099  return a->GetOutputProlongation();
1100 }
1101 
1103 {
1104  return a->GetOutputRestriction();
1105 }
1106 
1107 // Data and methods for partially-assembled bilinear forms
1108 
1110  MixedBilinearForm *form)
1112  trial_fes(form->TrialFESpace()),
1113  test_fes(form->TestFESpace()),
1114  elem_restrict_trial(NULL),
1115  elem_restrict_test(NULL)
1116 {
1117  Update();
1118 }
1119 
1121 {
1122  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1123  const int integratorCount = integrators.Size();
1124  for (int i = 0; i < integratorCount; ++i)
1125  {
1126  integrators[i]->AssemblePA(*trial_fes, *test_fes);
1127  }
1128  MFEM_VERIFY(a->GetBBFI()->Size() == 0,
1129  "Partial assembly does not support AddBoundaryIntegrator yet.");
1130  MFEM_VERIFY(a->GetTFBFI()->Size() == 0,
1131  "Partial assembly does not support AddTraceFaceIntegrator yet.");
1132  MFEM_VERIFY(a->GetBTFBFI()->Size() == 0,
1133  "Partial assembly does not support AddBdrTraceFaceIntegrator yet.");
1134 }
1135 
1137 {
1138  trial_fes = a->TrialFESpace();
1139  test_fes = a->TestFESpace();
1140  height = test_fes->GetVSize();
1141  width = trial_fes->GetVSize();
1146  if (elem_restrict_trial)
1147  {
1148  localTrial.UseDevice(true);
1151  }
1152  if (elem_restrict_test)
1153  {
1154  localTest.UseDevice(true); // ensure 'localY = 0.0' is done on device
1156  }
1157 }
1158 
1160  const Array<int> &trial_tdof_list,
1161  const Array<int> &test_tdof_list,
1162  OperatorHandle &A)
1163 {
1164  Operator * oper;
1165  Operator::FormRectangularSystemOperator(trial_tdof_list, test_tdof_list,
1166  oper);
1167  A.Reset(oper); // A will own oper
1168 }
1169 
1171  const Array<int> &trial_tdof_list,
1172  const Array<int> &test_tdof_list,
1173  Vector &x, Vector &b,
1174  OperatorHandle &A,
1175  Vector &X, Vector &B)
1176 {
1177  Operator *oper;
1178  Operator::FormRectangularLinearSystem(trial_tdof_list, test_tdof_list, x, b,
1179  oper, X, B);
1180  A.Reset(oper); // A will own oper
1181 }
1182 
1184  const Operator *elem_restrict_x,
1185  const Vector &x,
1186  Vector &localX,
1187  const Operator *elem_restrict_y,
1188  Vector &y,
1189  Vector &localY,
1190  const double c) const
1191 {
1192  // * G operation: localX = c*local(x)
1193  if (elem_restrict_x)
1194  {
1195  elem_restrict_x->Mult(x, localX);
1196  if (c != 1.0)
1197  {
1198  localX *= c;
1199  }
1200  }
1201  else
1202  {
1203  if (c == 1.0)
1204  {
1205  localX.SyncAliasMemory(x);
1206  }
1207  else
1208  {
1209  localX.Set(c, x);
1210  }
1211  }
1212  if (elem_restrict_y)
1213  {
1214  localY = 0.0;
1215  }
1216  else
1217  {
1218  y.UseDevice(true);
1219  localY.SyncAliasMemory(y);
1220  }
1221 }
1222 
1224 {
1225  y = 0.0;
1226  AddMult(x, y);
1227 }
1228 
1230  const double c) const
1231 {
1232  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1233  const int iSz = integrators.Size();
1234 
1235  // * G operation
1237  elem_restrict_test, y, localTest, c);
1238 
1239  // * B^TDB operation
1240  for (int i = 0; i < iSz; ++i)
1241  {
1242  integrators[i]->AddMultPA(localTrial, localTest);
1243  }
1244 
1245  // * G^T operation
1246  if (elem_restrict_test)
1247  {
1248  tempY.SetSize(y.Size());
1250  y += tempY;
1251  }
1252 }
1253 
1255  Vector &y) const
1256 {
1257  y = 0.0;
1258  AddMultTranspose(x, y);
1259 }
1260 
1262  const double c) const
1263 {
1264  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1265  const int iSz = integrators.Size();
1266 
1267  // * G operation
1270 
1271  // * B^TD^TB operation
1272  for (int i = 0; i < iSz; ++i)
1273  {
1274  integrators[i]->AddMultTransposePA(localTest, localTrial);
1275  }
1276 
1277  // * G^T operation
1278  if (elem_restrict_trial)
1279  {
1280  tempY.SetSize(y.Size());
1282  y += tempY;
1283  }
1284 }
1285 
1287  Vector &diag) const
1288 {
1289  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1290 
1291  const int iSz = integrators.Size();
1292 
1293  if (elem_restrict_trial)
1294  {
1295  const ElementRestriction* H1elem_restrict_trial =
1296  dynamic_cast<const ElementRestriction*>(elem_restrict_trial);
1297  if (H1elem_restrict_trial)
1298  {
1299  H1elem_restrict_trial->MultUnsigned(D, localTrial);
1300  }
1301  else
1302  {
1304  }
1305  }
1306 
1307  if (elem_restrict_test)
1308  {
1309  localTest = 0.0;
1310  for (int i = 0; i < iSz; ++i)
1311  {
1312  if (elem_restrict_trial)
1313  {
1314  integrators[i]->AssembleDiagonalPA_ADAt(localTrial, localTest);
1315  }
1316  else
1317  {
1318  integrators[i]->AssembleDiagonalPA_ADAt(D, localTest);
1319  }
1320  }
1321  const ElementRestriction* H1elem_restrict_test =
1322  dynamic_cast<const ElementRestriction*>(elem_restrict_test);
1323  if (H1elem_restrict_test)
1324  {
1325  H1elem_restrict_test->MultTransposeUnsigned(localTest, diag);
1326  }
1327  else
1328  {
1330  }
1331  }
1332  else
1333  {
1334  diag.UseDevice(true); // typically this is a large vector, so store on device
1335  diag = 0.0;
1336  for (int i = 0; i < iSz; ++i)
1337  {
1338  if (elem_restrict_trial)
1339  {
1340  integrators[i]->AssembleDiagonalPA_ADAt(localTrial, diag);
1341  }
1342  else
1343  {
1344  integrators[i]->AssembleDiagonalPA_ADAt(D, diag);
1345  }
1346  }
1347  }
1348 }
1349 
1351  DiscreteLinearOperator *linop) :
1353 {
1354 }
1355 
1356 const
1358 const
1359 {
1360  return a->GetOutputRestrictionTranspose();
1361 }
1362 
1364 {
1365  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1366  const int integratorCount = integrators.Size();
1367  for (int i = 0; i < integratorCount; ++i)
1368  {
1369  integrators[i]->AssemblePA(*trial_fes, *test_fes);
1370  }
1371 
1372  test_multiplicity.UseDevice(true);
1373  test_multiplicity.SetSize(elem_restrict_test->Width()); // l-vector
1374  Vector ones(elem_restrict_test->Height()); // e-vector
1375  ones = 1.0;
1376 
1377  const ElementRestriction* elem_restrict =
1378  dynamic_cast<const ElementRestriction*>(elem_restrict_test);
1379  if (elem_restrict)
1380  {
1381  elem_restrict->MultTransposeUnsigned(ones, test_multiplicity);
1382  }
1383  else
1384  {
1385  mfem_error("A real ElementRestriction is required in this setting!");
1386  }
1387 
1388  auto tm = test_multiplicity.ReadWrite();
1389  MFEM_FORALL(i, test_multiplicity.Size(),
1390  {
1391  tm[i] = 1.0 / tm[i];
1392  });
1393 }
1394 
1396  const Vector &x, Vector &y, const double c) const
1397 {
1398  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1399  const int iSz = integrators.Size();
1400 
1401  // * G operation
1403  elem_restrict_test, y, localTest, c);
1404 
1405  // * B^TDB operation
1406  for (int i = 0; i < iSz; ++i)
1407  {
1408  integrators[i]->AddMultPA(localTrial, localTest);
1409  }
1410 
1411  // do a kind of "set" rather than "add" in the below
1412  // operation as compared to the BilinearForm case
1413  // * G^T operation (kind of...)
1414  const ElementRestriction* elem_restrict =
1415  dynamic_cast<const ElementRestriction*>(elem_restrict_test);
1416  if (elem_restrict)
1417  {
1418  tempY.SetSize(y.Size());
1419  elem_restrict->MultLeftInverse(localTest, tempY);
1420  y += tempY;
1421  }
1422  else
1423  {
1424  mfem_error("In this setting you need a real ElementRestriction!");
1425  }
1426 }
1427 
1429  const Vector &x, Vector &y, const double c) const
1430 {
1431  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
1432  const int iSz = integrators.Size();
1433 
1434  // do a kind of "set" rather than "add" in the below
1435  // operation as compared to the BilinearForm case
1436  // * G operation (kinda)
1437  Vector xscaled(x);
1438  MFEM_VERIFY(x.Size() == test_multiplicity.Size(), "Input vector of wrong size");
1439  auto xs = xscaled.ReadWrite();
1440  auto tm = test_multiplicity.Read();
1441  MFEM_FORALL(i, x.Size(),
1442  {
1443  xs[i] *= tm[i];
1444  });
1447 
1448  // * B^TD^TB operation
1449  for (int i = 0; i < iSz; ++i)
1450  {
1451  integrators[i]->AddMultTransposePA(localTest, localTrial);
1452  }
1453 
1454  // * G^T operation
1455  if (elem_restrict_trial)
1456  {
1457  tempY.SetSize(y.Size());
1459  y += tempY;
1460  }
1461  else
1462  {
1463  mfem_error("Trial ElementRestriction not defined");
1464  }
1465 }
1466 
1468  const Array<int>& ess1, const Array<int>& ess2, OperatorHandle &A)
1469 {
1470  const Operator *Pi = this->GetProlongation();
1471  const Operator *RoT = this->GetOutputRestrictionTranspose();
1472  Operator *rap = SetupRAP(Pi, RoT);
1473 
1475  = new RectangularConstrainedOperator(rap, ess1, ess2, rap != this);
1476 
1477  A.Reset(Arco);
1478 }
1479 
1480 } // namespace mfem
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition: fespace.hpp:607
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:563
void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const
Assemble the diagonal of ADA^T for a diagonal vector D.
Array< BilinearFormIntegrator * > * GetBBFI()
Access all the integrators added with AddBoundaryIntegrator().
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:521
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:72
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:199
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
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:1059
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:244
bool UsesTensorBasis(const FiniteElementSpace &fes)
Definition: fespace.hpp:983
Array< BilinearFormIntegrator * > * GetFBFI()
Access all integrators added with AddInteriorFaceIntegrator().
virtual const Operator * GetRestriction() const
Get the input finite element space restriction matrix.
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:117
Memory< double > & GetMemoryData()
Definition: sparsemat.hpp:256
Memory< int > & GetMemoryI()
Definition: sparsemat.hpp:224
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:93
Class extending the MixedBilinearForm class to support different AssemblyLevels.
Memory< int > & GetMemoryJ()
Definition: sparsemat.hpp:240
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
const Operator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:1261
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition: dtensor.hpp:131
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:46
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:66
void FillI(SparseMatrix &mat) const
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:433
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:446
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.
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:901
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:439
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:1294
virtual const Operator * GetOutputRestrictionTranspose() const
Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type...
Definition: operator.hpp:129
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).
virtual void AddMultTranspose(const Vector &x, Vector &y) const =0
Add the face degrees of freedom x to the element degrees of freedom y.
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:27
const FiniteElementSpace * test_fes
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:638
void Assemble()
Assemble at the level given for the BilinearFormExtension subclass.
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
const FaceRestriction * int_face_restrict_lex
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:454
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:905
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 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:66
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:813
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2783
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 &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:236
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
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:438
MixedBilinearFormExtension(MixedBilinearForm *form)
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Operator that converts L2 FiniteElementSpace L-vectors to E-vectors.
Definition: restriction.hpp:84
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.
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