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