MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pbilinearform.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "../config/config.hpp"
13 
14 #ifdef MFEM_USE_MPI
15 
16 #include "fem.hpp"
17 #include "../general/sort_pairs.hpp"
18 
19 namespace mfem
20 {
21 
23 {
24  int nbr_size = pfes->GetFaceNbrVSize();
25 
26  if (precompute_sparsity == 0 || fes->GetVDim() > 1)
27  {
28  if (keep_nbr_block)
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  {
77  Table *face_elem = pfes->GetParMesh()->GetFaceToAllElementTable();
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 
93  if (keep_nbr_block)
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  double *data = Memory<double>(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 {
126  A.Clear();
127 
128  if (A_local == NULL) { return; }
129  MFEM_VERIFY(A_local->Finalized(), "the local matrix must be finalized");
130 
131  OperatorHandle dA(A.Type()), Ph(A.Type()), hdA;
132 
133  if (interior_face_integs.Size() == 0)
134  {
135  // construct a parallel block-diagonal matrix 'A' based on 'a'
137  pfes->GetDofOffsets(), A_local);
138  }
139  else
140  {
141  // handle the case when 'a' contains off-diagonal
142  int lvsize = pfes->GetVSize();
143  const HYPRE_BigInt *face_nbr_glob_ldof = pfes->GetFaceNbrGlobalDofMap();
144  HYPRE_BigInt ldof_offset = pfes->GetMyDofOffset();
145 
146  Array<HYPRE_BigInt> glob_J(A_local->NumNonZeroElems());
147  int *J = A_local->GetJ();
148  for (int i = 0; i < glob_J.Size(); i++)
149  {
150  if (J[i] < lvsize)
151  {
152  glob_J[i] = J[i] + ldof_offset;
153  }
154  else
155  {
156  glob_J[i] = face_nbr_glob_ldof[J[i] - lvsize];
157  }
158  }
159 
160  // TODO - construct dA directly in the A format
161  hdA.Reset(
162  new HypreParMatrix(pfes->GetComm(), lvsize, pfes->GlobalVSize(),
163  pfes->GlobalVSize(), A_local->GetI(), glob_J,
164  A_local->GetData(), pfes->GetDofOffsets(),
165  pfes->GetDofOffsets()));
166  // - hdA owns the new HypreParMatrix
167  // - the above constructor copies all input arrays
168  glob_J.DeleteAll();
169  dA.ConvertFrom(hdA);
170  }
171 
172  // TODO - assemble the Dof_TrueDof_Matrix directly in the required format?
173  Ph.ConvertFrom(pfes->Dof_TrueDof_Matrix());
174  // TODO: When Ph.Type() == Operator::ANY_TYPE we want to use the Operator
175  // returned by pfes->GetProlongationMatrix(), however that Operator is a
176  // const Operator, so we cannot store it in OperatorHandle. We need a const
177  // version of class OperatorHandle, e.g. ConstOperatorHandle.
178 
179  A.MakePtAP(dA, Ph);
180 }
181 
183 {
185  ParallelAssemble(Mh, m);
186  Mh.SetOperatorOwner(false);
187  return Mh.As<HypreParMatrix>();
188 }
189 
191 {
192  ParMesh *pmesh = pfes->GetParMesh();
194  Array<int> vdofs1, vdofs2, vdofs_all;
196 
197  int nfaces = pmesh->GetNSharedFaces();
198  for (int i = 0; i < nfaces; i++)
199  {
200  T = pmesh->GetSharedFaceTransformations(i);
201  int Elem2NbrNo = T->Elem2No - pmesh->GetNE();
202  pfes->GetElementVDofs(T->Elem1No, vdofs1);
203  pfes->GetFaceNbrElementVDofs(Elem2NbrNo, vdofs2);
204  vdofs1.Copy(vdofs_all);
205  for (int j = 0; j < vdofs2.Size(); j++)
206  {
207  if (vdofs2[j] >= 0)
208  {
209  vdofs2[j] += height;
210  }
211  else
212  {
213  vdofs2[j] -= height;
214  }
215  }
216  vdofs_all.Append(vdofs2);
217  for (int k = 0; k < interior_face_integs.Size(); k++)
218  {
220  AssembleFaceMatrix(*pfes->GetFE(T->Elem1No),
221  *pfes->GetFaceNbrFE(Elem2NbrNo),
222  *T, elemmat);
223  if (keep_nbr_block)
224  {
225  mat->AddSubMatrix(vdofs_all, vdofs_all, elemmat, skip_zeros);
226  }
227  else
228  {
229  mat->AddSubMatrix(vdofs1, vdofs_all, elemmat, skip_zeros);
230  }
231  }
232  }
233 }
234 
235 void ParBilinearForm::Assemble(int skip_zeros)
236 {
237  if (interior_face_integs.Size())
238  {
240  if (!ext && mat == NULL)
241  {
242  pAllocMat();
243  }
244  }
245 
246  BilinearForm::Assemble(skip_zeros);
247 
248  if (!ext && interior_face_integs.Size() > 0)
249  {
250  AssembleSharedFaces(skip_zeros);
251  }
252 }
253 
255 {
256  MFEM_ASSERT(diag.Size() == fes->GetTrueVSize(),
257  "Vector for holding diagonal has wrong size!");
258  const Operator *P = fes->GetProlongationMatrix();
259  if (!ext)
260  {
261  MFEM_ASSERT(p_mat.Ptr(), "the ParBilinearForm is not assembled!");
262  p_mat->AssembleDiagonal(diag); // TODO: add support for PETSc matrices
263  return;
264  }
265  // Here, we have extension, ext.
266  if (IsIdentityProlongation(P))
267  {
268  ext->AssembleDiagonal(diag);
269  return;
270  }
271  // Here, we have extension, ext, and parallel/conforming prolongation, P.
272  Vector local_diag(P->Height());
273  ext->AssembleDiagonal(local_diag);
274  if (fes->Conforming())
275  {
276  P->MultTranspose(local_diag, diag);
277  return;
278  }
279  // For an AMR mesh, a convergent diagonal is assembled with |P^T| d_l,
280  // where |P^T| has the entry-wise absolute values of the conforming
281  // prolongation transpose operator.
282  const HypreParMatrix *HP = dynamic_cast<const HypreParMatrix*>(P);
283  if (HP)
284  {
285  HP->AbsMultTranspose(1.0, local_diag, 0.0, diag);
286  }
287  else
288  {
289  MFEM_ABORT("unsupported prolongation matrix type.");
290  }
291 }
292 
293 void ParBilinearForm
295  HypreParMatrix &A, const HypreParVector &X,
296  HypreParVector &B) const
297 {
298  Array<int> dof_list;
299 
300  pfes->GetEssentialTrueDofs(bdr_attr_is_ess, dof_list);
301 
302  // do the parallel elimination
303  A.EliminateRowsCols(dof_list, X, B);
304 }
305 
308  HypreParMatrix &A) const
309 {
310  Array<int> dof_list;
311 
312  pfes->GetEssentialTrueDofs(bdr_attr_is_ess, dof_list);
313 
314  return A.EliminateRowsCols(dof_list);
315 }
316 
317 void ParBilinearForm::TrueAddMult(const Vector &x, Vector &y, const double a)
318 const
319 {
320  if (Xaux.ParFESpace() != pfes)
321  {
322  Xaux.SetSpace(pfes);
323  Yaux.SetSpace(pfes);
325  }
326 
327  Xaux.Distribute(&x);
328  if (ext)
329  {
330  ext->Mult(Xaux, Yaux);
331  }
332  else
333  {
334  MFEM_VERIFY(interior_face_integs.Size() == 0,
335  "the case of interior face integrators is not"
336  " implemented");
337  mat->Mult(Xaux, Yaux);
338  }
340  y.Add(a,Ytmp);
341 }
342 
344  const Array<int> &ess_tdof_list, Vector &x, Vector &b,
345  OperatorHandle &A, Vector &X, Vector &B, int copy_interior)
346 {
347  if (ext)
348  {
349  ext->FormLinearSystem(ess_tdof_list, x, b, A, X, B, copy_interior);
350  return;
351  }
352 
353  // Finish the matrix assembly and perform BC elimination, storing the
354  // eliminated part of the matrix.
355  FormSystemMatrix(ess_tdof_list, A);
356 
357  const Operator &P = *pfes->GetProlongationMatrix();
358  const SparseMatrix &R = *pfes->GetRestrictionMatrix();
359 
360  // Transform the system and perform the elimination in B, based on the
361  // essential BC values from x. Restrict the BC part of x in X, and set the
362  // non-BC part to zero. Since there is no good initial guess for the Lagrange
363  // multipliers, set X = 0.0 for hybridization.
364  if (static_cond)
365  {
366  // Schur complement reduction to the exposed dofs
367  static_cond->ReduceSystem(x, b, X, B, copy_interior);
368  }
369  else if (hybridization)
370  {
371  // Reduction to the Lagrange multipliers system
372  HypreParVector true_X(pfes), true_B(pfes);
373  P.MultTranspose(b, true_B);
374  R.Mult(x, true_X);
375  p_mat.EliminateBC(p_mat_e, ess_tdof_list, true_X, true_B);
377  R.MultTranspose(true_B, b);
378  hybridization->ReduceRHS(true_B, B);
379  X.SetSize(B.Size());
380  X = 0.0;
381  }
382  else
383  {
384  // Variational restriction with P
385  X.SetSize(P.Width());
386  B.SetSize(X.Size());
387  P.MultTranspose(b, B);
388  R.Mult(x, X);
389  p_mat.EliminateBC(p_mat_e, ess_tdof_list, X, B);
390  if (!copy_interior) { X.SetSubVectorComplement(ess_tdof_list, 0.0); }
391  }
392 }
393 
395  const Array<int> &vdofs, const Vector &x, Vector &b)
396 {
397  p_mat.EliminateBC(p_mat_e, vdofs, x, b);
398 }
399 
401  OperatorHandle &A)
402 {
403  if (ext)
404  {
405  ext->FormSystemMatrix(ess_tdof_list, A);
406  return;
407  }
408 
409  // Finish the matrix assembly and perform BC elimination, storing the
410  // eliminated part of the matrix.
411  if (static_cond)
412  {
414  {
415  static_cond->SetEssentialTrueDofs(ess_tdof_list);
418  }
420  }
421  else
422  {
423  if (mat)
424  {
425  const int remove_zeros = 0;
426  Finalize(remove_zeros);
427  MFEM_VERIFY(p_mat.Ptr() == NULL && p_mat_e.Ptr() == NULL,
428  "The ParBilinearForm must be updated with Update() before "
429  "re-assembling the ParBilinearForm.");
431  delete mat;
432  mat = NULL;
433  delete mat_e;
434  mat_e = NULL;
435  p_mat_e.EliminateRowsCols(p_mat, ess_tdof_list);
436  }
437  if (hybridization)
438  {
440  }
441  else
442  {
443  A = p_mat;
444  }
445  }
446 }
447 
449  const Vector &X, const Vector &b, Vector &x)
450 {
451  if (ext)
452  {
453  ext->RecoverFEMSolution(X, b, x);
454  return;
455  }
456 
457  const Operator &P = *pfes->GetProlongationMatrix();
458 
459  if (static_cond)
460  {
461  // Private dofs back solve
462  static_cond->ComputeSolution(b, X, x);
463  }
464  else if (hybridization)
465  {
466  // Primal unknowns recovery
467  HypreParVector true_X(pfes), true_B(pfes);
468  P.MultTranspose(b, true_B);
469  const SparseMatrix &R = *pfes->GetRestrictionMatrix();
470  R.Mult(x, true_X); // get essential b.c. from x
471  hybridization->ComputeSolution(true_B, X, true_X);
472  x.SetSize(P.Height());
473  P.Mult(true_X, x);
474  }
475  else
476  {
477  // Apply conforming prolongation
479  P.Mult(X, x);
480  }
481 }
482 
484 {
485  BilinearForm::Update(nfes);
486 
487  if (nfes)
488  {
489  pfes = dynamic_cast<ParFiniteElementSpace *>(nfes);
490  MFEM_VERIFY(pfes != NULL, "nfes must be a ParFiniteElementSpace!");
491  }
492 
493  p_mat.Clear();
494  p_mat_e.Clear();
495 }
496 
497 
499 {
500  // construct the block-diagonal matrix A
501  HypreParMatrix *A =
507  mat);
508 
511 
512  delete A;
513 
514  return rap;
515 }
516 
518 {
519  // construct the rectangular block-diagonal matrix dA
520  OperatorHandle dA(A.Type());
526  mat);
527 
528  OperatorHandle P_test(A.Type()), P_trial(A.Type());
529 
530  // TODO - construct the Dof_TrueDof_Matrix directly in the required format.
532  P_trial.ConvertFrom(trial_pfes->Dof_TrueDof_Matrix());
533 
534  A.MakeRAP(P_test, dA, P_trial);
535 }
536 
537 /// Compute y += a (P^t A P) x, where x and y are vectors on the true dofs
539  const double a) const
540 {
541  if (Xaux.ParFESpace() != trial_pfes)
542  {
545  }
546 
547  Xaux.Distribute(&x);
548  mat->Mult(Xaux, Yaux);
550 }
551 
553  const Array<int>
554  &trial_tdof_list,
555  const Array<int> &test_tdof_list,
556  OperatorHandle &A)
557 {
558  if (ext)
559  {
560  ext->FormRectangularSystemOperator(trial_tdof_list, test_tdof_list, A);
561  return;
562  }
563 
564  if (mat)
565  {
566  Finalize();
568  delete mat;
569  mat = NULL;
570  delete mat_e;
571  mat_e = NULL;
572  HypreParMatrix *temp =
573  p_mat.As<HypreParMatrix>()->EliminateCols(trial_tdof_list);
574  p_mat.As<HypreParMatrix>()->EliminateRows(test_tdof_list);
575  p_mat_e.Reset(temp, true);
576  }
577 
578  A = p_mat;
579 }
580 
582  const Array<int>
583  &trial_tdof_list,
584  const Array<int> &test_tdof_list, Vector &x,
585  Vector &b, OperatorHandle &A, Vector &X,
586  Vector &B)
587 {
588  if (ext)
589  {
590  ext->FormRectangularLinearSystem(trial_tdof_list, test_tdof_list,
591  x, b, A, X, B);
592  return;
593  }
594 
595  FormRectangularSystemMatrix(trial_tdof_list, test_tdof_list, A);
596 
597  const Operator *test_P = test_pfes->GetProlongationMatrix();
598  const SparseMatrix *trial_R = trial_pfes->GetRestrictionMatrix();
599 
602  test_P->MultTranspose(b, B);
603  trial_R->Mult(x, X);
604 
605  p_mat_e.As<HypreParMatrix>()->Mult(-1.0, X, 1.0, B);
606  B.SetSubVector(test_tdof_list, 0.0);
607 }
608 
610 {
611  MFEM_ASSERT(mat, "Matrix is not assembled");
612  MFEM_ASSERT(mat->Finalized(), "Matrix is not finalized");
616  delete RA;
617  return RAP;
618 }
619 
621 {
622  // construct the rectangular block-diagonal matrix dA
623  OperatorHandle dA(A.Type());
629  mat);
630 
631  OperatorHandle R_test_transpose(A.Type()), P_trial(A.Type());
632 
633  // TODO - construct the Dof_TrueDof_Matrix directly in the required format.
634  R_test_transpose.ConvertFrom(range_fes->Dof_TrueDof_Matrix());
635  P_trial.ConvertFrom(domain_fes->Dof_TrueDof_Matrix());
636 
637  A.MakeRAP(R_test_transpose, dA, P_trial);
638 }
639 
641 {
642  if (ext)
643  {
644  Array<int> empty;
645  ext->FormRectangularSystemOperator(empty, empty, A);
646  return;
647  }
648 
649  mfem_error("not implemented!");
650 }
651 
653 const
654 {
655  MFEM_VERIFY(mat->Finalized(), "Local matrix needs to be finalized for "
656  "GetParBlocks");
657 
659 
660  blocks.SetSize(range_fes->GetVDim(), domain_fes->GetVDim());
661 
662  RLP->GetBlocks(blocks,
665 
666  delete RLP;
667 }
668 
669 }
670 
671 #endif
MemoryType GetHypreMemoryType()
The MemoryType used by MFEM when allocating arrays for Hypre objects.
Definition: hypre.hpp:136
virtual void FormRectangularSystemMatrix(OperatorHandle &A)
Return in A a parallel (on truedofs) version of this operator.
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:575
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition: hypre.cpp:2240
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:560
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:563
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
Definition: sparsemat.cpp:1497
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1032
void ConvertFrom(OperatorHandle &A)
Convert the given OperatorHandle A to the currently set type id.
Definition: handle.cpp:200
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
int * GetJ()
Definition: table.hpp:114
virtual void AssembleDiagonal(Vector &diag) const
Assemble the diagonal of the bilinear form into diag. Note that diag is a true-dof Vector...
ParMesh * GetParMesh() const
Definition: pfespace.hpp:277
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:2673
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(...
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:1153
BilinearFormExtension * ext
Extension for supporting Full Assembly (FA), Element Assembly (EA), Partial Assembly (PA)...
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition: pfespace.hpp:434
ParGridFunction Xaux
Auxiliary objects used in TrueAddMult().
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:472
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
HypreParMatrix & GetParallelMatrix()
Return the parallel hybridized matrix.
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition: pmesh.cpp:3025
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:214
void SetDims(int rows, int nnz)
Definition: table.cpp:140
Operator::Type Type() const
Get the currently set operator type id.
Definition: handle.hpp:99
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition: array.hpp:856
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:480
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
Array< BilinearFormIntegrator * > interior_face_integs
Set of interior face Integrators to be applied.
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
void TrueAddMult(const Vector &x, Vector &y, const double a=1.0) const
Compute y += a (P^t A P) x, where x and y are vectors on the 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...
Definition: staticcond.cpp:416
int * GetI()
Return the array I.
Definition: sparsemat.hpp:209
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
void EliminateReducedTrueDofs(const Array< int > &ess_rtdof_list, Matrix::DiagonalPolicy dpolicy)
Eliminate the given reduced true dofs from the Schur complement matrix S.
Definition: staticcond.cpp:286
OperatorHandle p_mat_e
void LoseData()
Call this if data has been stolen.
Definition: table.hpp:154
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the ParGridFunction.
Definition: pgridfunc.cpp:87
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AssembleSharedFaces(int skip_zeros=1)
void ComputeSolution(const Vector &b, const Vector &sc_sol, Vector &sol) const
Definition: staticcond.cpp:476
void DeleteAll()
Delete the whole array.
Definition: array.hpp:846
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: operator.hpp:93
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)=0
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
double * GetData()
Return the element data, i.e. the array A.
Definition: sparsemat.hpp:219
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:289
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P_test^t A P_trial.
void AbsMultTranspose(double a, const Vector &x, double 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:1862
void SetSize(int m, int n)
Definition: array.hpp:370
HYPRE_BigInt GlobalVSize() const
Definition: pfespace.hpp:283
ParFiniteElementSpace * pfes
Points to the same object as fes.
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1517
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:281
MPI_Comm GetComm() const
Definition: pfespace.hpp:273
virtual const FiniteElement * GetFE(int i) const
Definition: pfespace.cpp:535
virtual void FormRectangularSystemOperator(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)=0
Data type sparse matrix.
Definition: sparsemat.hpp:46
virtual void AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
StaticCondensation * static_cond
Owned.
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:751
void EnsureMultTranspose() const
Ensures that the matrix is capable of performing MultTranspose(), AddMultTranspose(), and AbsMultTranspose().
Definition: sparsemat.cpp:881
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:433
void Clear()
Definition: table.cpp:377
SparseMatrix * mat
Sparse matrix to be associated with the form. Owned.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:154
double b
Definition: lissajous.cpp:42
const HYPRE_BigInt * GetFaceNbrGlobalDofMap()
Definition: pfespace.hpp:399
void ComputeSolution(const Vector &b, const Vector &sol_r, Vector &sol) const
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
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Definition: pmesh.cpp:2909
void EliminateVDofsInRHS(const Array< int > &vdofs, const Vector &x, Vector &b)
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:321
HypreParMatrix & GetParallelMatrix()
Return the parallel Schur complement matrix.
Definition: staticcond.hpp:159
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:471
void Assemble(int skip_zeros=1)
Assemble the local matrix.
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.
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
HYPRE_BigInt GetMyDofOffset() const
Definition: pfespace.cpp:1143
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:566
bool Conforming() const
Definition: fespace.hpp:439
ParFiniteElementSpace * test_pfes
Points to the same object as test_fes.
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:92
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:557
SparseMatrix * mat
Owned.
Set the diagonal value to one.
Definition: operator.hpp:49
HYPRE_BigInt * GetDofOffsets() const
Definition: pfespace.hpp:281
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2584
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:87
Dynamic 2D array using row-major layout.
Definition: array.hpp:351
void SetSubVectorComplement(const Array< int > &dofs, const double val)
Set all vector entries NOT in the dofs Array to the given val.
Definition: vector.cpp:695
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
bool Finalized() const
Returns whether or not CSR format has been finalized.
Definition: sparsemat.hpp:521
HypreParMatrix * ParallelAssemble() const
Returns the matrix &quot;assembled&quot; on the true dofs.
SparseMatrix * mat_e
Sparse Matrix used to store the eliminations from the b.c. Owned. .
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:148
ParFiniteElementSpace * domain_fes
Points to the same object as trial_fes.
HYPRE_BigInt * GetTrueDofOffsets() const
Definition: pfespace.hpp:282
bool IsIdentityProlongation(const Operator *P)
Definition: operator.hpp:689
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
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
MixedBilinearFormExtension * ext
ParGridFunction Xaux
Auxiliary objects used in TrueAddMult().
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:638
void Finalize()
Finalize the construction of the Schur complement matrix.
Definition: staticcond.cpp:233
HYPRE_Int HYPRE_BigInt
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
Definition: table.cpp:410
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
Definition: handle.hpp:124
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:132
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A^t * x + beta * y.
Definition: hypre.cpp:1837
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:85
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
Table * GetFaceToElementTable() const
Definition: mesh.cpp:6083
double a
Definition: lissajous.cpp:41
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:252
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:120
virtual void AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
Definition: operator.hpp:107
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:389
SparseMatrix * mat_e
Owned.
ParFiniteElementSpace * range_fes
Points to the same object as test_fes.
bool HasEliminatedBC() const
Return true if essential boundary conditions have been eliminated from the Schur complement matrix...
Definition: staticcond.hpp:142
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:813
FiniteElementSpace * fes
FE space on which the form lives. Not owned.
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.
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...
ParFiniteElementSpace * trial_pfes
Points to the same object as trial_fes.
int GetFaceNbrVSize() const
Definition: pfespace.hpp:394
OperatorHandle p_mat
Vector data type.
Definition: vector.hpp:60
void GetParBlocks(Array2D< HypreParMatrix * > &blocks) const
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)=0
ID for class HypreParMatrix.
Definition: operator.hpp:259
Hybridization * hybridization
Owned.
void ReduceRHS(const Vector &b, Vector &b_r) const
int * GetI()
Definition: table.hpp:113
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:726
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:1879
void SetEssentialTrueDofs(const Array< int > &ess_tdof_list)
Determine and save internally essential reduced true dofs.
Definition: staticcond.hpp:128
Table * GetFaceToAllElementTable() const
Definition: pmesh.cpp:2803
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:337
void TrueAddMult(const Vector &x, Vector &y, const double a=1.0) const
Compute y += a (P^t A P) x, where x and y are vectors on the true dofs.
Class for parallel meshes.
Definition: pmesh.hpp:32
OperatorHandle p_mat
Matrix and eliminated matrix.
virtual void Mult(const Vector &x, Vector &y) const
Matrix multiplication: .
void MakePtAP(OperatorHandle &A, OperatorHandle &P)
Reset the OperatorHandle to hold the product P^t A P.
Definition: handle.cpp:123
ParGridFunction Yaux
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
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 Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:145
DofTransformation * GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1464
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
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:113