MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
pbilinearform.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#include "../config/config.hpp"
13
14#ifdef MFEM_USE_MPI
15
16#include "fem.hpp"
18
19namespace mfem
20{
21
23{
24 int nbr_size = pfes->GetFaceNbrVSize();
25
26 if (precompute_sparsity == 0 || fes->GetVDim() > 1)
27 {
29 {
30 mat = new SparseMatrix(height + nbr_size, width + nbr_size);
31 }
32 else
33 {
34 mat = new SparseMatrix(height, width + nbr_size);
35 }
36 return;
37 }
38
39 // the sparsity pattern is defined from the map: face->element->dof
40 const Table &lelem_ldof = fes->GetElementToDofTable(); // <-- dofs
41 const Table &nelem_ndof = pfes->face_nbr_element_dof; // <-- vdofs
42 Table elem_dof; // element + nbr-element <---> dof
43 if (nbr_size > 0)
44 {
45 // merge lelem_ldof and nelem_ndof into elem_dof
46 int s1 = lelem_ldof.Size(), s2 = nelem_ndof.Size();
47 const int *I1 = lelem_ldof.GetI(), *J1 = lelem_ldof.GetJ();
48 const int *I2 = nelem_ndof.GetI(), *J2 = nelem_ndof.GetJ();
49 const int nnz1 = I1[s1], nnz2 = I2[s2];
50
51 elem_dof.SetDims(s1 + s2, nnz1 + nnz2);
52
53 int *I = elem_dof.GetI(), *J = elem_dof.GetJ();
54 for (int i = 0; i <= s1; i++)
55 {
56 I[i] = I1[i];
57 }
58 for (int j = 0; j < nnz1; j++)
59 {
60 J[j] = J1[j];
61 }
62 for (int i = 0; i <= s2; i++)
63 {
64 I[s1+i] = I2[i] + nnz1;
65 }
66 for (int j = 0; j < nnz2; j++)
67 {
68 J[nnz1+j] = J2[j] + height;
69 }
70 }
71 // dof_elem x elem_face x face_elem x elem_dof (keep_nbr_block = true)
72 // ldof_lelem x lelem_face x face_elem x elem_dof (keep_nbr_block = false)
73 Table dof_dof;
74 {
75 Table face_dof; // face_elem x elem_dof
76 {
78 if (nbr_size > 0)
79 {
80 mfem::Mult(*face_elem, elem_dof, face_dof);
81 }
82 else
83 {
84 mfem::Mult(*face_elem, lelem_ldof, face_dof);
85 }
86 delete face_elem;
87 if (nbr_size > 0)
88 {
89 elem_dof.Clear();
90 }
91 }
92
94 {
95 Table dof_face;
96 Transpose(face_dof, dof_face, height + nbr_size);
97 mfem::Mult(dof_face, face_dof, dof_dof);
98 }
99 else
100 {
101 Table ldof_face;
102 {
103 Table face_ldof;
104 Table *face_lelem = fes->GetMesh()->GetFaceToElementTable();
105 mfem::Mult(*face_lelem, lelem_ldof, face_ldof);
106 delete face_lelem;
107 Transpose(face_ldof, ldof_face, height);
108 }
109 mfem::Mult(ldof_face, face_dof, dof_dof);
110 }
111 }
112
113 int *I = dof_dof.GetI();
114 int *J = dof_dof.GetJ();
115 int nrows = dof_dof.Size();
116 real_t *data = Memory<real_t>(I[nrows]);
117
118 mat = new SparseMatrix(I, J, data, nrows, height + nbr_size);
119 *mat = 0.0;
120
121 dof_dof.LoseData();
122}
123
125 bool steal_loc_A)
126{
127 ParFiniteElementSpace &pfespace = *ParFESpace();
128
129 // Create a block diagonal parallel matrix
131 A_diag.MakeSquareBlockDiag(pfespace.GetComm(),
132 pfespace.GlobalVSize(),
133 pfespace.GetDofOffsets(),
134 &loc_A);
135
136 // Parallel matrix assembly using P^t A P (if needed)
138 {
139 A_diag.SetOperatorOwner(false);
140 A.Reset(A_diag.As<HypreParMatrix>());
141 if (steal_loc_A)
142 {
144 }
145 }
146 else
147 {
149 P.ConvertFrom(pfespace.Dof_TrueDof_Matrix());
150 A.MakePtAP(A_diag, P);
151 }
152}
153
155{
156 A.Clear();
157
158 if (A_local == NULL) { return; }
159 MFEM_VERIFY(A_local->Finalized(), "the local matrix must be finalized");
160
161 OperatorHandle dA(A.Type()), Ph(A.Type()), hdA;
162
163 if (interior_face_integs.Size() == 0)
164 {
165 // construct a parallel block-diagonal matrix 'A' based on 'a'
167 pfes->GetDofOffsets(), A_local);
168 }
169 else
170 {
171 // handle the case when 'a' contains off-diagonal
172 int lvsize = pfes->GetVSize();
173 const HYPRE_BigInt *face_nbr_glob_ldof = pfes->GetFaceNbrGlobalDofMap();
174 HYPRE_BigInt ldof_offset = pfes->GetMyDofOffset();
175
176 Array<HYPRE_BigInt> glob_J(A_local->NumNonZeroElems());
177 int *J = A_local->GetJ();
178 for (int i = 0; i < glob_J.Size(); i++)
179 {
180 if (J[i] < lvsize)
181 {
182 glob_J[i] = J[i] + ldof_offset;
183 }
184 else
185 {
186 glob_J[i] = face_nbr_glob_ldof[J[i] - lvsize];
187 }
188 }
189
190 // TODO - construct dA directly in the A format
191 hdA.Reset(
192 new HypreParMatrix(pfes->GetComm(), lvsize, pfes->GlobalVSize(),
193 pfes->GlobalVSize(), A_local->GetI(), glob_J,
194 A_local->GetData(), pfes->GetDofOffsets(),
195 pfes->GetDofOffsets()));
196 // - hdA owns the new HypreParMatrix
197 // - the above constructor copies all input arrays
198 glob_J.DeleteAll();
199 dA.ConvertFrom(hdA);
200 }
201
202 // TODO - assemble the Dof_TrueDof_Matrix directly in the required format?
203 Ph.ConvertFrom(pfes->Dof_TrueDof_Matrix());
204 // TODO: When Ph.Type() == Operator::ANY_TYPE we want to use the Operator
205 // returned by pfes->GetProlongationMatrix(), however that Operator is a
206 // const Operator, so we cannot store it in OperatorHandle. We need a const
207 // version of class OperatorHandle, e.g. ConstOperatorHandle.
208
209 A.MakePtAP(dA, Ph);
210}
211
219
221{
222 ParMesh *pmesh = pfes->GetParMesh();
224 Array<int> vdofs1, vdofs2, vdofs_all;
226
227 int nfaces = pmesh->GetNSharedFaces();
228 for (int i = 0; i < nfaces; i++)
229 {
230 T = pmesh->GetSharedFaceTransformations(i);
231 int Elem2NbrNo = T->Elem2No - pmesh->GetNE();
232 pfes->GetElementVDofs(T->Elem1No, vdofs1);
233 pfes->GetFaceNbrElementVDofs(Elem2NbrNo, vdofs2);
234 vdofs1.Copy(vdofs_all);
235 for (int j = 0; j < vdofs2.Size(); j++)
236 {
237 if (vdofs2[j] >= 0)
238 {
239 vdofs2[j] += height;
240 }
241 else
242 {
243 vdofs2[j] -= height;
244 }
245 }
246 vdofs_all.Append(vdofs2);
247 for (int k = 0; k < interior_face_integs.Size(); k++)
248 {
250 AssembleFaceMatrix(*pfes->GetFE(T->Elem1No),
251 *pfes->GetFaceNbrFE(Elem2NbrNo),
252 *T, elemmat);
253 if (keep_nbr_block)
254 {
255 mat->AddSubMatrix(vdofs_all, vdofs_all, elemmat, skip_zeros);
256 }
257 else
258 {
259 mat->AddSubMatrix(vdofs1, vdofs_all, elemmat, skip_zeros);
260 }
261 }
262 }
263}
264
265void ParBilinearForm::Assemble(int skip_zeros)
266{
267 if (interior_face_integs.Size())
268 {
270 if (!ext && mat == NULL)
271 {
272 pAllocMat();
273 }
274 }
275
276 BilinearForm::Assemble(skip_zeros);
277
278 if (!ext && interior_face_integs.Size() > 0)
279 {
280 AssembleSharedFaces(skip_zeros);
281 }
282}
283
285{
286 MFEM_ASSERT(diag.Size() == fes->GetTrueVSize(),
287 "Vector for holding diagonal has wrong size!");
288 const Operator *P = fes->GetProlongationMatrix();
289 if (!ext)
290 {
291 MFEM_ASSERT(p_mat.Ptr(), "the ParBilinearForm is not assembled!");
292 p_mat->AssembleDiagonal(diag); // TODO: add support for PETSc matrices
293 return;
294 }
295 // Here, we have extension, ext.
297 {
298 ext->AssembleDiagonal(diag);
299 return;
300 }
301 // Here, we have extension, ext, and parallel/conforming prolongation, P.
302 Vector local_diag(P->Height());
303 ext->AssembleDiagonal(local_diag);
304 if (fes->Conforming())
305 {
306 P->MultTranspose(local_diag, diag);
307 return;
308 }
309 // For an AMR mesh, a convergent diagonal is assembled with |P^T| d_l,
310 // where |P^T| has the entry-wise absolute values of the conforming
311 // prolongation transpose operator.
312 const HypreParMatrix *HP = dynamic_cast<const HypreParMatrix*>(P);
313 if (HP)
314 {
315 HP->AbsMultTranspose(1.0, local_diag, 0.0, diag);
316 }
317 else
318 {
319 MFEM_ABORT("unsupported prolongation matrix type.");
320 }
321}
322
323void ParBilinearForm
324::ParallelEliminateEssentialBC(const Array<int> &bdr_attr_is_ess,
325 HypreParMatrix &A, const HypreParVector &X,
326 HypreParVector &B) const
327{
328 Array<int> dof_list;
329
330 pfes->GetEssentialTrueDofs(bdr_attr_is_ess, dof_list);
331
332 // do the parallel elimination
333 A.EliminateRowsCols(dof_list, X, B);
334}
335
337ParallelEliminateEssentialBC(const Array<int> &bdr_attr_is_ess,
338 HypreParMatrix &A) const
339{
340 Array<int> dof_list;
341
342 pfes->GetEssentialTrueDofs(bdr_attr_is_ess, dof_list);
343
344 return A.EliminateRowsCols(dof_list);
345}
346
348const
349{
350 const Operator *P = pfes->GetProlongationMatrix();
351 Xaux.SetSize(P->Height());
352 Yaux.SetSize(P->Height());
353 Ytmp.SetSize(P->Width());
354
355 P->Mult(x, Xaux);
356 if (ext)
357 {
358 ext->Mult(Xaux, Yaux);
359 }
360 else
361 {
362 MFEM_VERIFY(interior_face_integs.Size() == 0,
363 "the case of interior face integrators is not"
364 " implemented");
365 mat->Mult(Xaux, Yaux);
366 }
368 y.Add(a, Ytmp);
369}
370
372 const ParGridFunction &y) const
373{
374 MFEM_ASSERT(mat != NULL, "local matrix must be assembled");
375
376 real_t loc = InnerProduct(x, y);
377 real_t glob = 0.;
378
379 MPI_Allreduce(&loc, &glob, 1, MPITypeMap<real_t>::mpi_type, MPI_SUM,
380 pfes->GetComm());
381
382 return glob;
383}
384
386 const ParGridFunction &y) const
387{
388 MFEM_ASSERT(x.ParFESpace() == pfes, "the parallel spaces must match");
389 MFEM_ASSERT(y.ParFESpace() == pfes, "the parallel spaces must match");
390
393
394 real_t res = TrueInnerProduct(*x_p, *y_p);
395
396 delete x_p;
397 delete y_p;
398
399 return res;
400}
401
403 HypreParVector &y) const
404{
405 MFEM_VERIFY(p_mat.Ptr() != NULL, "parallel matrix must be assembled");
406
408 {
409 return TrueInnerProduct((const Vector&)x, (const Vector&)y);
410 }
411
414
415 A->Mult(x, *Ax);
416
417 real_t res = mfem::InnerProduct(y, *Ax);
418
419 delete Ax;
420
421 return res;
422}
423
425 const Vector &y) const
426{
427 MFEM_VERIFY(p_mat.Ptr() != NULL, "parallel matrix must be assembled");
428
429 Vector Ax(pfes->GetTrueVSize());
430 p_mat->Mult(x, Ax);
431
432 real_t res = mfem::InnerProduct(pfes->GetComm(), y, Ax);
433
434 return res;
435}
436
438 const Array<int> &ess_tdof_list, Vector &x, Vector &b,
439 OperatorHandle &A, Vector &X, Vector &B, int copy_interior)
440{
441 if (ext)
442 {
443 ext->FormLinearSystem(ess_tdof_list, x, b, A, X, B, copy_interior);
444 return;
445 }
446
447 // Finish the matrix assembly and perform BC elimination, storing the
448 // eliminated part of the matrix.
449 FormSystemMatrix(ess_tdof_list, A);
450
451 const Operator &P = *pfes->GetProlongationMatrix();
453
454 // Transform the system and perform the elimination in B, based on the
455 // essential BC values from x. Restrict the BC part of x in X, and set the
456 // non-BC part to zero. Since there is no good initial guess for the Lagrange
457 // multipliers, set X = 0.0 for hybridization.
458 if (static_cond)
459 {
460 // Schur complement reduction to the exposed dofs
461 static_cond->ReduceSystem(x, b, X, B, copy_interior);
462 }
463 else if (hybridization)
464 {
465 // Reduction to the Lagrange multipliers system
466 HypreParVector true_X(pfes), true_B(pfes);
467 P.MultTranspose(b, true_B);
468 R.Mult(x, true_X);
469 p_mat.EliminateBC(p_mat_e, ess_tdof_list, true_X, true_B);
470 R.MultTranspose(true_B, b);
471 hybridization->ReduceRHS(true_B, B);
472 X.SetSize(B.Size());
473 X = 0.0;
474 }
475 else
476 {
477 // Variational restriction with P
478 X.SetSize(P.Width());
479 B.SetSize(X.Size());
480 P.MultTranspose(b, B);
481 R.Mult(x, X);
482 p_mat.EliminateBC(p_mat_e, ess_tdof_list, X, B);
483 if (!copy_interior) { X.SetSubVectorComplement(ess_tdof_list, 0.0); }
484 }
485}
486
488 const Array<int> &vdofs, const Vector &x, Vector &b)
489{
491}
492
495{
496 if (ext)
497 {
498 ext->FormSystemMatrix(ess_tdof_list, A);
499 return;
500 }
501
502 // Finish the matrix assembly and perform BC elimination, storing the
503 // eliminated part of the matrix.
504 if (static_cond)
505 {
507 {
508 static_cond->SetEssentialTrueDofs(ess_tdof_list);
511 }
513 }
514 else
515 {
516 if (mat)
517 {
518 const int remove_zeros = 0;
519 Finalize(remove_zeros);
520 MFEM_VERIFY(p_mat.Ptr() == NULL && p_mat_e.Ptr() == NULL,
521 "The ParBilinearForm must be updated with Update() before "
522 "re-assembling the ParBilinearForm.");
524 delete mat;
525 mat = NULL;
526 delete mat_e;
527 mat_e = NULL;
528 p_mat_e.EliminateRowsCols(p_mat, ess_tdof_list);
529 }
530 if (hybridization)
531 {
533 }
534 else
535 {
536 A = p_mat;
537 }
538 }
539}
540
542 const Vector &X, const Vector &b, Vector &x)
543{
544 if (ext)
545 {
546 ext->RecoverFEMSolution(X, b, x);
547 return;
548 }
549
550 const Operator &P = *pfes->GetProlongationMatrix();
551
552 if (static_cond)
553 {
554 // Private dofs back solve
556 }
557 else if (hybridization)
558 {
559 // Primal unknowns recovery
560 HypreParVector true_X(pfes), true_B(pfes);
561 P.MultTranspose(b, true_B);
563 R.Mult(x, true_X); // get essential b.c. from x
564 hybridization->ComputeSolution(true_B, X, true_X);
565 x.SetSize(P.Height());
566 P.Mult(true_X, x);
567 }
568 else
569 {
570 // Apply conforming prolongation
572 P.Mult(X, x);
573 }
574}
575
577{
579
580 if (nfes)
581 {
582 pfes = dynamic_cast<ParFiniteElementSpace *>(nfes);
583 MFEM_VERIFY(pfes != NULL, "nfes must be a ParFiniteElementSpace!");
584 }
585
586 p_mat.Clear();
587 p_mat_e.Clear();
588}
589
590
592{
593 // construct the block-diagonal matrix A
594 HypreParMatrix *A =
600 mat);
601
604
605 delete A;
606
607 return rap;
608}
609
611{
612 // construct the rectangular block-diagonal matrix dA
613 OperatorHandle dA(A.Type());
619 mat);
620
621 OperatorHandle P_test(A.Type()), P_trial(A.Type());
622
623 // TODO - construct the Dof_TrueDof_Matrix directly in the required format.
624 P_test.ConvertFrom(test_pfes->Dof_TrueDof_Matrix());
626
627 A.MakeRAP(P_test, dA, P_trial);
628}
629
630/// Compute y += a (P^t A P) x, where x and y are vectors on the true dofs
632 const real_t a) const
633{
634 if (Xaux.ParFESpace() != trial_pfes)
635 {
638 }
639
640 Xaux.Distribute(&x);
641 mat->Mult(Xaux, Yaux);
643}
644
646 const Array<int>
647 &trial_tdof_list,
648 const Array<int> &test_tdof_list,
650{
651 if (ext)
652 {
653 ext->FormRectangularSystemOperator(trial_tdof_list, test_tdof_list, A);
654 return;
655 }
656
657 if (mat)
658 {
659 Finalize();
661 delete mat;
662 mat = NULL;
663 delete mat_e;
664 mat_e = NULL;
665 HypreParMatrix *temp =
666 p_mat.As<HypreParMatrix>()->EliminateCols(trial_tdof_list);
667 p_mat.As<HypreParMatrix>()->EliminateRows(test_tdof_list);
668 p_mat_e.Reset(temp, true);
669 }
670
671 A = p_mat;
672}
673
675 const Array<int>
676 &trial_tdof_list,
677 const Array<int> &test_tdof_list, Vector &x,
678 Vector &b, OperatorHandle &A, Vector &X,
679 Vector &B)
680{
681 if (ext)
682 {
683 ext->FormRectangularLinearSystem(trial_tdof_list, test_tdof_list,
684 x, b, A, X, B);
685 return;
686 }
687
688 FormRectangularSystemMatrix(trial_tdof_list, test_tdof_list, A);
689
690 const Operator *test_P = test_pfes->GetProlongationMatrix();
691 const SparseMatrix *trial_R = trial_pfes->GetRestrictionMatrix();
692
695 test_P->MultTranspose(b, B);
696 trial_R->Mult(x, X);
697
698 p_mat_e.As<HypreParMatrix>()->Mult(-1.0, X, 1.0, B);
699 B.SetSubVector(test_tdof_list, 0.0);
700}
701
703{
704 MFEM_ASSERT(mat, "Matrix is not assembled");
705 MFEM_ASSERT(mat->Finalized(), "Matrix is not finalized");
709 delete RA;
710 return RAP;
711}
712
714{
715 // construct the rectangular block-diagonal matrix dA
716 OperatorHandle dA(A.Type());
722 mat);
723
725 OperatorHandle R_test_transpose(A.Type());
726 R_test_transpose.MakeRectangularBlockDiag(range_fes->GetComm(),
731 Rt);
732
733 // TODO - construct the Dof_TrueDof_Matrix directly in the required format.
734 OperatorHandle P_trial(A.Type());
736
737 A.MakeRAP(R_test_transpose, dA, P_trial);
738 delete Rt;
739}
740
742{
743 if (ext)
744 {
745 Array<int> empty;
746 ext->FormRectangularSystemOperator(empty, empty, A);
747 return;
748 }
749
750 mfem_error("not implemented!");
751}
752
754const
755{
756 MFEM_VERIFY(mat->Finalized(), "Local matrix needs to be finalized for "
757 "GetParBlocks");
758
760
762
763 RLP->GetBlocks(blocks,
766
767 delete RLP;
768}
769
770}
771
772#endif
Dynamic 2D array using row-major layout.
Definition array.hpp:372
void SetSize(int m, int n)
Definition array.hpp:385
int Size() const
Return the logical size of the array.
Definition array.hpp:144
void DeleteAll()
Delete the whole array.
Definition array.hpp:864
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:769
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition array.hpp:874
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)=0
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)=0
virtual void AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY....
BilinearFormExtension * ext
Extension for supporting Full Assembly (FA), Element Assembly (EA),Partial Assembly (PA),...
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
Hybridization * hybridization
Owned.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
real_t InnerProduct(const Vector &x, const Vector &y) const
Compute .
Array< BilinearFormIntegrator * > interior_face_integs
Set of interior face Integrators to be applied.
SparseMatrix * mat
Sparse matrix to be associated with the form. Owned.
FiniteElementSpace * fes
FE space on which the form lives. Not owned.
StaticCondensation * static_cond
Owned.
SparseMatrix * mat_e
Sparse Matrix used to store the eliminations from the b.c. Owned. .
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
A specialized ElementTransformation class representing a face and its two neighboring elements.
Definition eltrans.hpp:484
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
const Table & GetElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element,...
Definition fespace.hpp:1136
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:716
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition fespace.hpp:597
bool Conforming() const
Definition fespace.hpp:565
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:287
Ordering::Type GetOrdering() const
Return the ordering method.
Definition fespace.hpp:725
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
int GetVDim() const
Returns vector dimension.
Definition fespace.hpp:706
void ComputeSolution(const Vector &b, const Vector &sol_r, Vector &sol) const
HypreParMatrix & GetParallelMatrix()
Return the parallel hybridized matrix.
void ReduceRHS(const Vector &b, Vector &b_r) const
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
void AbsMultTranspose(real_t a, const Vector &x, real_t b, Vector &y) const
Computes y = a * |At| * x + b * y, using entry-wise absolute values of the transpose of the matrix A.
Definition hypre.cpp:1977
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, real_t alpha=1.0, real_t beta=0.0) const
Computes y = alpha * A^t * x + beta * y.
Definition hypre.cpp:1951
HypreParMatrix * LeftDiagMult(const SparseMatrix &D, HYPRE_BigInt *row_starts=NULL) const
Multiply the HypreParMatrix on the left by a block-diagonal parallel matrix D and return the result a...
Definition hypre.cpp:1994
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition hypre.cpp:2356
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, real_t alpha=1.0, real_t beta=0.0) const
Computes y = alpha * A * x + beta * y.
Definition hypre.cpp:1815
void GetBlocks(Array2D< HypreParMatrix * > &blocks, bool interleaved_rows=false, bool interleaved_cols=false) const
Definition hypre.cpp:1660
Wrapper for hypre's parallel vector class.
Definition hypre.hpp:206
Class used by MFEM to store pointers to host and/or device memory.
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
Table * GetFaceToElementTable() const
Definition mesh.cpp:7167
virtual void FormRectangularSystemOperator(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)=0
virtual void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)=0
MixedBilinearFormExtension * ext
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY.
virtual void Mult(const Vector &x, Vector &y) const
Matrix multiplication: .
SparseMatrix * mat
Owned.
SparseMatrix * mat_e
Owned.
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 EliminateRowsCols(OperatorHandle &A, const Array< int > &ess_dof_list)
Reset the OperatorHandle to be the eliminated part of A after elimination of the essential dofs ess_d...
Definition handle.cpp:252
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition handle.hpp:120
void MakeRAP(OperatorHandle &Rt, OperatorHandle &A, OperatorHandle &P)
Reset the OperatorHandle to hold the product R A P, where R = Rt^t.
Definition handle.cpp:161
void ConvertFrom(OperatorHandle &A)
Convert the given OperatorHandle A to the currently set type id.
Definition handle.cpp:200
void MakePtAP(OperatorHandle &A, OperatorHandle &P)
Reset the OperatorHandle to hold the product P^t A P.
Definition handle.cpp:123
void MakeSquareBlockDiag(MPI_Comm comm, HYPRE_BigInt glob_size, HYPRE_BigInt *row_starts, SparseMatrix *diag)
Reset the OperatorHandle to hold a parallel square block-diagonal matrix using the currently set type...
Definition handle.cpp:60
Operator * Ptr() const
Access the underlying Operator pointer.
Definition handle.hpp:87
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged.
Definition handle.hpp:124
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition handle.hpp:145
void EliminateBC(const OperatorHandle &A_e, const Array< int > &ess_dof_list, const Vector &X, Vector &B) const
Eliminate essential dofs from the solution X into the r.h.s. B.
Definition handle.cpp:340
Operator::Type Type() const
Get the currently set operator type id.
Definition handle.hpp:99
void MakeRectangularBlockDiag(MPI_Comm comm, HYPRE_BigInt glob_num_rows, HYPRE_BigInt glob_num_cols, HYPRE_BigInt *row_starts, HYPRE_BigInt *col_starts, SparseMatrix *diag)
Reset the OperatorHandle to hold a parallel rectangular block-diagonal matrix using the currently set...
Definition handle.cpp:91
Abstract operator.
Definition operator.hpp:25
int width
Dimension of the input / number of columns in the matrix.
Definition operator.hpp:28
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
@ Hypre_ParCSR
ID for class HypreParMatrix.
Definition operator.hpp:287
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear sys...
Definition operator.cpp:148
Type GetType() const
Return the type ID of the Operator class.
Definition operator.hpp:307
virtual void AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
Definition operator.hpp:131
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
Vector Xaux
Auxiliary vectors used in TrueAddMult(): L-, L-, and T-vector, resp.
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(....
ParFiniteElementSpace * pfes
Points to the same object as fes.
void ParallelRAP(SparseMatrix &loc_A, OperatorHandle &A, bool steal_loc_A=false)
Compute parallel RAP operator and store it in A as a HypreParMatrix.
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void ParallelEliminateEssentialBC(const Array< int > &bdr_attr_is_ess, HypreParMatrix &A, const HypreParVector &X, HypreParVector &B) const
Eliminate essential boundary DOFs from a parallel assembled system.
real_t ParInnerProduct(const ParGridFunction &x, const ParGridFunction &y) const
Compute .
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
void AssembleSharedFaces(int skip_zeros=1)
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
virtual void AssembleDiagonal(Vector &diag) const
Assemble the diagonal of the bilinear form into diag. Note that diag is a true-dof Vector.
real_t TrueInnerProduct(const ParGridFunction &x, const ParGridFunction &y) const
Compute on true dofs (grid function version)
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
void TrueAddMult(const Vector &x, Vector &y, const real_t a=1.0) const
Compute y += a (P^t A P) x, where x and y are vectors on the true dofs.
void EliminateVDofsInRHS(const Array< int > &vdofs, const Vector &x, Vector &b)
void GetParBlocks(Array2D< HypreParMatrix * > &blocks) const
ParFiniteElementSpace * range_fes
Points to the same object as test_fes.
HypreParMatrix * ParallelAssemble() const
Returns the matrix "assembled" on the true dofs.
ParFiniteElementSpace * domain_fes
Points to the same object as trial_fes.
virtual void FormRectangularSystemMatrix(OperatorHandle &A)
Return in A a parallel (on truedofs) version of this operator.
Abstract parallel finite element space.
Definition pfespace.hpp:29
MPI_Comm GetComm() const
Definition pfespace.hpp:273
HYPRE_BigInt * GetTrueDofOffsets() const
Definition pfespace.hpp:282
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
HYPRE_BigInt GlobalVSize() const
Definition pfespace.hpp:283
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:285
int GetTrueVSize() const override
Return the number of local vector true dofs.
Definition pfespace.hpp:289
HYPRE_BigInt GetMyDofOffset() const
HYPRE_BigInt * GetDofOffsets() const
Definition pfespace.hpp:281
const Operator * GetProlongationMatrix() const override
The returned Operator is owned by the FiniteElementSpace.
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition pfespace.hpp:327
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs, DofTransformation &doftrans) const
const HYPRE_BigInt * GetFaceNbrGlobalDofMap()
Definition pfespace.hpp:401
const SparseMatrix * GetRestrictionMatrix() const override
Get the R matrix which restricts a local dof vector to true dof vector.
Definition pfespace.hpp:389
ParMesh * GetParMesh() const
Definition pfespace.hpp:277
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition pfespace.hpp:436
const FiniteElement * GetFaceNbrFE(int i) const
const FiniteElement * GetFE(int i) const override
Definition pfespace.cpp:534
Class for parallel grid function.
Definition pgridfunc.hpp:33
ParFiniteElementSpace * ParFESpace() const
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
void SetSpace(FiniteElementSpace *f) override
Associate a new FiniteElementSpace with the ParGridFunction.
Definition pgridfunc.cpp:96
void Distribute(const Vector *tv)
Class for parallel meshes.
Definition pmesh.hpp:34
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition pmesh.cpp:3153
Table * GetFaceToAllElementTable() const
Definition pmesh.cpp:2839
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Get the FaceElementTransformations for the given shared face (edge 2D) using the shared face index sf...
Definition pmesh.cpp:2923
ParGridFunction Xaux
Auxiliary objects used in TrueAddMult().
OperatorHandle p_mat
Matrix and eliminated matrix.
void TrueAddMult(const Vector &x, Vector &y, const real_t a=1.0) const
Compute y += a (P^t A P) x, where x and y are vectors on the true dofs.
virtual void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
Form the parallel linear system A X = B, corresponding to this mixed bilinear form and the linear for...
virtual void FormRectangularSystemMatrix(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)
Return in A a parallel (on truedofs) version of this operator.
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P_test^t A P_trial.
ParFiniteElementSpace * test_pfes
Points to the same object as test_fes.
ParFiniteElementSpace * trial_pfes
Points to the same object as trial_fes.
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.
bool Finalized() const
Returns whether or not CSR format has been finalized.
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
real_t * GetData()
Return the element data, i.e. the array A.
int * GetJ()
Return the array J.
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
int * GetI()
Return the array I.
void Finalize()
Finalize the construction of the Schur complement matrix.
void EliminateReducedTrueDofs(const Array< int > &ess_rtdof_list, Matrix::DiagonalPolicy dpolicy)
Eliminate the given reduced true dofs from the Schur complement matrix S.
void SetEssentialTrueDofs(const Array< int > &ess_tdof_list)
Determine and save internally essential reduced true dofs.
void ReduceSystem(Vector &x, Vector &b, Vector &X, Vector &B, int copy_interior=0) const
Set the reduced solution X and r.h.s B vectors from the full linear system solution x and r....
void ComputeSolution(const Vector &b, const Vector &sc_sol, Vector &sol) const
HypreParMatrix & GetParallelMatrix()
Return the parallel Schur complement matrix.
bool HasEliminatedBC() const
Return true if essential boundary conditions have been eliminated from the Schur complement matrix.
void LoseData()
Call this if data has been stolen.
Definition table.hpp:180
int * GetJ()
Definition table.hpp:114
void Clear()
Definition table.cpp:381
int Size() const
Returns the number of TYPE I elements.
Definition table.hpp:92
int * GetI()
Definition table.hpp:113
void SetDims(int rows, int nnz)
Definition table.cpp:140
Vector data type.
Definition vector.hpp:80
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
Definition vector.cpp:604
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
void SetSubVectorComplement(const Array< int > &dofs, const real_t val)
Set all vector entries NOT in the dofs Array to the given val.
Definition vector.cpp:739
Vector & Add(const real_t a, const Vector &Va)
(*this) += a * Va
Definition vector.cpp:247
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
MemoryType GetHypreMemoryType()
The MemoryType used by MFEM when allocating arrays for Hypre objects.
Definition hypre.hpp:179
void mfem_error(const char *msg)
Definition error.cpp:154
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition table.cpp:476
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
Definition table.cpp:414
real_t InnerProduct(HypreParVector *x, HypreParVector *y)
Definition hypre.cpp:439
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
bool IsIdentityProlongation(const Operator *P)
Definition operator.hpp:720
void HypreStealOwnership(HypreParMatrix &A_hyp, SparseMatrix &A_diag)
Make A_hyp steal ownership of its diagonal part A_diag.
Definition hypre.cpp:2840
float real_t
Definition config.hpp:43
Helper struct to convert a C++ type to an MPI type.