MFEM  v3.3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
bilinearform.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_BILINEARFORM
13 #define MFEM_BILINEARFORM
14 
15 #include "../config/config.hpp"
16 #include "../linalg/linalg.hpp"
17 #include "fespace.hpp"
18 #include "gridfunc.hpp"
19 #include "linearform.hpp"
20 #include "bilininteg.hpp"
21 #include "staticcond.hpp"
22 #include "hybridization.hpp"
23 
24 namespace mfem
25 {
26 
27 /** Class for bilinear form - "Matrix" with associated FE space and
28  BLFIntegrators. */
29 class BilinearForm : public Matrix
30 {
31 protected:
32  /// Sparse matrix to be associated with the form.
34 
35  /// Matrix used to eliminate b.c.
37 
38  /// FE space on which the form lives.
40 
41  /// Indicates the Mesh::sequence corresponding to the current state of the
42  /// BilinearForm.
43  long sequence;
44 
46 
47  /// Set of Domain Integrators to be applied.
49 
50  /// Set of Boundary Integrators to be applied.
52 
53  /// Set of interior face Integrators to be applied.
55 
56  /// Set of boundary face Integrators to be applied.
59 
62 
64 
67 
69  // Allocate appropriate SparseMatrix and assign it to mat
70  void AllocMat();
71 
72  void ConformingAssemble();
73 
74  // may be used in the construction of derived classes
76  {
77  fes = NULL; sequence = -1;
78  mat = mat_e = NULL; extern_bfs = 0; element_matrices = NULL;
79  static_cond = NULL; hybridization = NULL;
81  }
82 
83 public:
84  /// Creates bilinear form associated with FE space *f.
86 
87  BilinearForm(FiniteElementSpace *f, BilinearForm *bf, int ps = 0);
88 
89  /// Get the size of the BilinearForm as a square matrix.
90  int Size() const { return height; }
91 
92  /** Enable the use of static condensation. For details see the description
93  for class StaticCondensation in fem/staticcond.hpp This method should be
94  called before assembly. If the number of unknowns after static
95  condensation is not reduced, it is not enabled. */
97 
98  /** Check if static condensation was actually enabled by a previous call to
99  EnableStaticCondensation. */
100  bool StaticCondensationIsEnabled() const { return static_cond; }
101 
102  /// Return the trace FE space associated with static condensation.
104  { return static_cond ? static_cond->GetTraceFESpace() : NULL; }
105 
106  /** Enable hybridization; for details see the description for class
107  Hybridization in fem/hybridization.hpp. This method should be called
108  before assembly. */
109  void EnableHybridization(FiniteElementSpace *constr_space,
110  BilinearFormIntegrator *constr_integ,
111  const Array<int> &ess_tdof_list);
112 
113  /** For scalar FE spaces, precompute the sparsity pattern of the matrix
114  (assuming dense element matrices) based on the types of integrators
115  present in the bilinear form. */
116  void UsePrecomputedSparsity(int ps = 1) { precompute_sparsity = ps; }
117 
118  /** @brief Use the given CSR sparsity pattern to allocate the internal
119  SparseMatrix.
120 
121  - The @a I and @a J arrays must define a square graph with size equal to
122  GetVSize() of the associated FiniteElementSpace.
123  - This method should be called after enabling static condensation or
124  hybridization, if used.
125  - In the case of static condensation, @a I and @a J are not used.
126  - The ownership of the arrays @a I and @a J remains with the caller. */
127  void UseSparsity(int *I, int *J, bool isSorted);
128 
129  /// Use the sparsity of @a A to allocate the internal SparseMatrix.
130  void UseSparsity(SparseMatrix &A);
131 
132  /** Pre-allocate the internal SparseMatrix before assembly. If the flag
133  'precompute sparsity' is set, the matrix is allocated in CSR format (i.e.
134  finalized) and the entries are initialized with zeros. */
135  void AllocateMatrix() { if (mat == NULL) { AllocMat(); } }
136 
138 
140 
142 
144 
145  const double &operator()(int i, int j) { return (*mat)(i,j); }
146 
147  /// Returns reference to a_{ij}.
148  virtual double &Elem(int i, int j);
149 
150  /// Returns constant reference to a_{ij}.
151  virtual const double &Elem(int i, int j) const;
152 
153  /// Matrix vector multiplication.
154  virtual void Mult(const Vector &x, Vector &y) const { mat->Mult(x, y); }
155 
156  void FullMult(const Vector &x, Vector &y) const
157  { mat->Mult(x, y); mat_e->AddMult(x, y); }
158 
159  virtual void AddMult(const Vector &x, Vector &y, const double a = 1.0) const
160  { mat -> AddMult (x, y, a); }
161 
162  void FullAddMult(const Vector &x, Vector &y) const
163  { mat->AddMult(x, y); mat_e->AddMult(x, y); }
164 
165  virtual void AddMultTranspose(const Vector & x, Vector & y,
166  const double a = 1.0) const
167  { mat->AddMultTranspose(x, y, a); }
168 
169  void FullAddMultTranspose (const Vector & x, Vector & y) const
170  { mat->AddMultTranspose(x, y); mat_e->AddMultTranspose(x, y); }
171 
172  virtual void MultTranspose (const Vector & x, Vector & y) const
173  { y = 0.0; AddMultTranspose (x, y); }
174 
175  double InnerProduct(const Vector &x, const Vector &y) const
176  { return mat->InnerProduct (x, y); }
177 
178  /// Returns a pointer to (approximation) of the matrix inverse.
179  virtual MatrixInverse *Inverse() const;
180 
181  /// Finalizes the matrix initialization.
182  virtual void Finalize(int skip_zeros = 1);
183 
184  /// Returns a reference to the sparse matrix
185  const SparseMatrix &SpMat() const
186  {
187  MFEM_VERIFY(mat, "mat is NULL and can't be dereferenced");
188  return *mat;
189  }
191  {
192  MFEM_VERIFY(mat, "mat is NULL and can't be dereferenced");
193  return *mat;
194  }
195  SparseMatrix *LoseMat() { SparseMatrix *tmp = mat; mat = NULL; return tmp; }
196 
197  /// Returns a reference to the sparse matrix of eliminated b.c.
198  const SparseMatrix &SpMatElim() const
199  {
200  MFEM_VERIFY(mat_e, "mat_e is NULL and can't be dereferenced");
201  return *mat_e;
202  }
204  {
205  MFEM_VERIFY(mat_e, "mat_e is NULL and can't be dereferenced");
206  return *mat_e;
207  }
208 
209  /// Adds new Domain Integrator.
211 
212  /// Adds new Boundary Integrator.
214 
215  /// Adds new interior Face Integrator.
217 
218  /// Adds new boundary Face Integrator.
220 
221  /** @brief Adds new boundary Face Integrator, restricted to specific boundary
222  attributes. */
224  Array<int> &bdr_marker);
225 
226  void operator=(const double a)
227  {
228  if (mat != NULL) { *mat = a; }
229  if (mat_e != NULL) { *mat_e = a; }
230  }
231 
232  /// Assembles the form i.e. sums over all domain/bdr integrators.
233  void Assemble(int skip_zeros = 1);
234 
235  /// Get the finite element space prolongation matrix
236  virtual const Operator *GetProlongation() const
237  { return fes->GetConformingProlongation(); }
238  /// Get the finite element space restriction matrix
239  virtual const Operator *GetRestriction() const
240  { return fes->GetConformingRestriction(); }
241 
242  /** Form the linear system A X = B, corresponding to the current bilinear
243  form and b(.), by applying any necessary transformations such as:
244  eliminating boundary conditions; applying conforming constraints for
245  non-conforming AMR; static condensation; hybridization.
246 
247  The GridFunction-size vector x must contain the essential b.c. The
248  BilinearForm and the LinearForm-size vector b must be assembled.
249 
250  The vector X is initialized with a suitable initial guess: when using
251  hybridization, the vector X is set to zero; otherwise, the essential
252  entries of X are set to the corresponding b.c. and all other entries are
253  set to zero (copy_interior == 0) or copied from x (copy_interior != 0).
254 
255  This method can be called multiple times (with the same ess_tdof_list
256  array) to initialize different right-hand sides and boundary condition
257  values.
258 
259  After solving the linear system, the finite element solution x can be
260  recovered by calling RecoverFEMSolution (with the same vectors X, b, and
261  x).
262 
263  NOTE: If there are no transformations, X simply reuses the data of x. */
264  void FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x, Vector &b,
265  SparseMatrix &A, Vector &X, Vector &B,
266  int copy_interior = 0);
267 
268  /// Form the linear system matrix A, see FormLinearSystem for details.
269  void FormSystemMatrix(const Array<int> &ess_tdof_list, SparseMatrix &A);
270 
271  /** Call this method after solving a linear system constructed using the
272  FormLinearSystem method to recover the solution as a GridFunction-size
273  vector in x. Use the same arguments as in the FormLinearSystem call. */
274  virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x);
275 
276  /// Compute and store internally all element matrices.
277  void ComputeElementMatrices();
278 
279  /// Free the memory used by the element matrices.
281  { delete element_matrices; element_matrices = NULL; }
282 
283  void ComputeElementMatrix(int i, DenseMatrix &elmat);
284  void AssembleElementMatrix(int i, const DenseMatrix &elmat,
285  Array<int> &vdofs, int skip_zeros = 1);
286  void AssembleBdrElementMatrix(int i, const DenseMatrix &elmat,
287  Array<int> &vdofs, int skip_zeros = 1);
288 
289  /** Eliminate essential boundary DOFs from the system. The array
290  'bdr_attr_is_ess' marks boundary attributes that constitute the essential
291  part of the boundary. If d == 0, the diagonal at the essential DOFs is
292  set to 1.0, otherwise it is left the same. */
293  void EliminateEssentialBC(const Array<int> &bdr_attr_is_ess,
294  Vector &sol, Vector &rhs, int d = 0);
295 
296  void EliminateEssentialBC(const Array<int> &bdr_attr_is_ess, int d = 0);
297  /// Perform elimination and set the diagonal entry to the given value
298  void EliminateEssentialBCDiag(const Array<int> &bdr_attr_is_ess,
299  double value);
300 
301  /// Eliminate the given vdofs. NOTE: here, vdofs is a list of DOFs.
302  void EliminateVDofs(const Array<int> &vdofs, Vector &sol, Vector &rhs,
303  int d = 0);
304 
305  /** Eliminate the given vdofs storing the eliminated part internally; this
306  method works in conjunction with EliminateVDofsInRHS and allows
307  elimination of boundary conditions in multiple right-hand sides. In this
308  method, vdofs is a list of DOFs. */
309  void EliminateVDofs(const Array<int> &vdofs, int d = 0);
310 
311  /** Similar to EliminateVDofs but here ess_dofs is a marker
312  (boolean) array on all vdofs (ess_dofs[i] < 0 is true). */
313  void EliminateEssentialBCFromDofs(const Array<int> &ess_dofs, Vector &sol,
314  Vector &rhs, int d = 0);
315 
316  /** Similar to EliminateVDofs but here ess_dofs is a marker
317  (boolean) array on all vdofs (ess_dofs[i] < 0 is true). */
318  void EliminateEssentialBCFromDofs(const Array<int> &ess_dofs, int d = 0);
319  /// Perform elimination and set the diagonal entry to the given value
320  void EliminateEssentialBCFromDofsDiag(const Array<int> &ess_dofs,
321  double value);
322 
323  /** Use the stored eliminated part of the matrix (see EliminateVDofs) to
324  modify r.h.s.; vdofs is a list of DOFs (non-directional, i.e. >= 0). */
325  void EliminateVDofsInRHS(const Array<int> &vdofs, const Vector &x,
326  Vector &b);
327 
328  double FullInnerProduct(const Vector &x, const Vector &y) const
329  { return mat->InnerProduct(x, y) + mat_e->InnerProduct(x, y); }
330 
331  virtual void Update(FiniteElementSpace *nfes = NULL);
332 
333  /// (DEPRECATED) Return the FE space associated with the BilinearForm.
334  /** @deprecated Use FESpace() instead. */
336 
337  /// Return the FE space associated with the BilinearForm.
339  /// Read-only access to the associated FiniteElementSpace.
340  const FiniteElementSpace *FESpace() const { return fes; }
341 
342  /// Destroys bilinear form.
343  virtual ~BilinearForm();
344 };
345 
346 /**
347  Class for assembling of bilinear forms `a(u,v)` defined on different
348  trial and test spaces. The assembled matrix `A` is such that
349 
350  a(u,v) = V^t A U
351 
352  where `U` and `V` are the vectors representing the functions `u` and `v`,
353  respectively. The first argument, `u`, of `a(,)` is in the trial space
354  and the second argument, `v`, is in the test space. Thus,
355 
356  # of rows of A = dimension of the test space and
357  # of cols of A = dimension of the trial space.
358 
359  Both trial and test spaces should be defined on the same mesh.
360 */
361 class MixedBilinearForm : public Matrix
362 {
363 protected:
365 
367 
370  Array<BilinearFormIntegrator*> skt; // trace face integrators
371 
372 public:
374  FiniteElementSpace *te_fes);
375 
376  virtual double& Elem (int i, int j);
377 
378  virtual const double& Elem (int i, int j) const;
379 
380  virtual void Mult (const Vector & x, Vector & y) const;
381 
382  virtual void AddMult (const Vector & x, Vector & y,
383  const double a = 1.0) const;
384 
385  virtual void AddMultTranspose (const Vector & x, Vector & y,
386  const double a = 1.0) const;
387 
388  virtual void MultTranspose (const Vector & x, Vector & y) const
389  { y = 0.0; AddMultTranspose (x, y); }
390 
391  virtual MatrixInverse * Inverse() const;
392 
393  virtual void Finalize (int skip_zeros = 1);
394 
395  /** Extract the associated matrix as SparseMatrix blocks. The number of
396  block rows and columns is given by the vector dimensions (vdim) of the
397  test and trial spaces, respectively. */
398  void GetBlocks(Array2D<SparseMatrix *> &blocks) const;
399 
400  const SparseMatrix &SpMat() const { return *mat; }
401  SparseMatrix &SpMat() { return *mat; }
402  SparseMatrix *LoseMat() { SparseMatrix *tmp = mat; mat = NULL; return tmp; }
403 
405 
407 
408  /** Add a trace face integrator. This type of integrator assembles terms
409  over all faces of the mesh using the face FE from the trial space and the
410  two adjacent volume FEs from the test space. */
412 
414 
416 
418 
419  void operator= (const double a) { *mat = a; }
420 
421  void Assemble (int skip_zeros = 1);
422 
423  /** For partially conforming trial and/or test FE spaces, complete the
424  assembly process by performing A := P2^t A P1 where A is the internal
425  sparse matrix; P1 and P2 are the conforming prolongation matrices of the
426  trial and test FE spaces, respectively. After this call the
427  MixedBilinearForm becomes an operator on the conforming FE spaces. */
428  void ConformingAssemble();
429 
430  void EliminateTrialDofs(Array<int> &bdr_attr_is_ess,
431  Vector &sol, Vector &rhs);
432 
434  Vector &sol, Vector &rhs);
435 
436  virtual void EliminateTestDofs(Array<int> &bdr_attr_is_ess);
437 
438  void Update();
439 
440  virtual ~MixedBilinearForm();
441 };
442 
443 
444 /**
445  Class for constructing the matrix representation of a linear operator,
446  `v = L u`, from one FiniteElementSpace (domain) to another FiniteElementSpace
447  (range). The constructed matrix `A` is such that
448 
449  V = A U
450 
451  where `U` and `V` are the vectors of degrees of freedom representing the
452  functions `u` and `v`, respectively. The dimensions of `A` are
453 
454  number of rows of A = dimension of the range space and
455  number of cols of A = dimension of the domain space.
456 
457  This class is very similar to MixedBilinearForm. One difference is that
458  the linear operator `L` is defined using a special kind of
459  BilinearFormIntegrator (we reuse its functionality instead of defining a
460  new class). The other difference with the MixedBilinearForm class is that
461  the "assembly" process overwrites the global matrix entries using the
462  local element matrices instead of adding them.
463 
464  Note that if we define the bilinear form `b(u,v) := (Lu,v)` using an inner
465  product in the range space, then its matrix representation, `B`, is
466 
467  B = M A, (since V^t B U = b(u,v) = (Lu,v) = V^t M A U)
468 
469  where `M` denotes the mass matrix for the inner product in the range space:
470  `V1^t M V2 = (v1,v2)`. Similarly, if `c(u,w) := (Lu,Lw)` then
471 
472  C = A^t M A.
473 */
475 {
476 public:
478  FiniteElementSpace *range_fes)
479  : MixedBilinearForm(domain_fes, range_fes) { }
480 
482  { AddDomainIntegrator(di); }
483 
485  { AddTraceFaceIntegrator(di); }
486 
488 
489  virtual void Assemble(int skip_zeros = 1);
490 };
491 
492 }
493 
494 #endif
void FullAddMult(const Vector &x, Vector &y) const
SparseMatrix & SpMat()
Array< BilinearFormIntegrator * > * GetBBFI()
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
Array< BilinearFormIntegrator * > skt
const SparseMatrix * GetConformingProlongation() const
Definition: fespace.cpp:725
virtual MatrixInverse * Inverse() const
Returns a pointer to (approximation) of the matrix inverse.
Array< BilinearFormIntegrator * > fbfi
Set of interior face Integrators to be applied.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
Array< BilinearFormIntegrator * > * GetDBFI()
Array< BilinearFormIntegrator * > dbfi
Set of Domain Integrators to be applied.
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
Array< BilinearFormIntegrator * > * GetBFBFI()
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
void AddMult(const Vector &x, Vector &y, const double a=1.0) const
y += A * x (default) or y += a * A * x
Definition: sparsemat.cpp:481
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
y += At * x (default) or y += a * At * x
Definition: sparsemat.cpp:558
Abstract data type for matrix inverse.
Definition: matrix.hpp:58
void EliminateVDofs(const Array< int > &vdofs, Vector &sol, Vector &rhs, int d=0)
Eliminate the given vdofs. NOTE: here, vdofs is a list of DOFs.
void FormSystemMatrix(const Array< int > &ess_tdof_list, SparseMatrix &A)
Form the linear system matrix A, see FormLinearSystem for details.
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
void EliminateEssentialBCFromDofs(const Array< int > &ess_dofs, Vector &sol, Vector &rhs, int d=0)
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
Array< BilinearFormIntegrator * > * GetFBFI()
MixedBilinearForm(FiniteElementSpace *tr_fes, FiniteElementSpace *te_fes)
double FullInnerProduct(const Vector &x, const Vector &y) const
void AddDomainInterpolator(DiscreteInterpolator *di)
void EnableHybridization(FiniteElementSpace *constr_space, BilinearFormIntegrator *constr_integ, const Array< int > &ess_tdof_list)
void EliminateEssentialBCFromDofsDiag(const Array< int > &ess_dofs, double value)
Perform elimination and set the diagonal entry to the given value.
void EliminateEssentialBCFromTrialDofs(Array< int > &marked_vdofs, Vector &sol, Vector &rhs)
Array< BilinearFormIntegrator * > bdr
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
SparseMatrix & SpMat()
const SparseMatrix * GetConformingRestriction() const
Definition: fespace.cpp:732
void FullAddMultTranspose(const Vector &x, Vector &y) const
DiscreteLinearOperator(FiniteElementSpace *domain_fes, FiniteElementSpace *range_fes)
virtual void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
Data type sparse matrix.
Definition: sparsemat.hpp:38
StaticCondensation * static_cond
const FiniteElementSpace * FESpace() const
Read-only access to the associated FiniteElementSpace.
SparseMatrix * mat
Sparse matrix to be associated with the form.
int Size() const
Get the size of the BilinearForm as a square matrix.
Abstract data type matrix.
Definition: matrix.hpp:27
void Assemble(int skip_zeros=1)
void FullMult(const Vector &x, Vector &y) const
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
void EliminateEssentialBCDiag(const Array< int > &bdr_attr_is_ess, double value)
Perform elimination and set the diagonal entry to the given value.
virtual void Update(FiniteElementSpace *nfes=NULL)
void EliminateTrialDofs(Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs)
FiniteElementSpace * GetFES()
(DEPRECATED) Return the FE space associated with the BilinearForm.
void AddTraceFaceIntegrator(BilinearFormIntegrator *bfi)
void GetBlocks(Array2D< SparseMatrix * > &blocks) const
Array< BilinearFormIntegrator * > * GetTFBFI()
Array< BilinearFormIntegrator * > * GetBBFI()
void ComputeElementMatrix(int i, DenseMatrix &elmat)
SparseMatrix & SpMatElim()
SparseMatrix * mat_e
Matrix used to eliminate b.c.
bool StaticCondensationIsEnabled() const
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:22
virtual const Operator * GetRestriction() const
Get the finite element space restriction matrix.
void AssembleElementMatrix(int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1)
void UsePrecomputedSparsity(int ps=1)
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 ...
SparseMatrix * LoseMat()
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:475
Array< BilinearFormIntegrator * > bfbfi
Set of boundary face Integrators to be applied.
Array< BilinearFormIntegrator * > bbfi
Set of Boundary Integrators to be applied.
DenseTensor * element_matrices
SparseMatrix * LoseMat()
void operator=(const double a)
virtual void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
void AddTraceFaceInterpolator(DiscreteInterpolator *di)
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
Array< BilinearFormIntegrator * > dom
double InnerProduct(const Vector &x, const Vector &y) const
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:24
const SparseMatrix & SpMat() const
FiniteElementSpace * test_fes
const double & operator()(int i, int j)
FiniteElementSpace * GetTraceFESpace()
Return a pointer to the reduced/trace FE space.
Definition: staticcond.hpp:109
virtual void EliminateTestDofs(Array< int > &bdr_attr_is_ess)
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void EliminateVDofsInRHS(const Array< int > &vdofs, const Vector &x, Vector &b)
const SparseMatrix & SpMatElim() const
Returns a reference to the sparse matrix of eliminated b.c.
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds new Boundary Integrator.
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator.
FiniteElementSpace * fes
FE space on which the form lives.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
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 ...
void ComputeElementMatrices()
Compute and store internally all element matrices.
FiniteElementSpace * SCFESpace() const
Return the trace FE space associated with static condensation.
void AssembleBdrElementMatrix(int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1)
DenseMatrix elemmat
FiniteElementSpace * trial_fes
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
virtual ~BilinearForm()
Destroys bilinear form.
Vector data type.
Definition: vector.hpp:41
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
virtual MatrixInverse * Inverse() const
Returns a pointer to (an approximation) of the matrix inverse.
Hybridization * hybridization
void operator=(const double a)
Array< BilinearFormIntegrator * > * GetDBFI()
const SparseMatrix & SpMat() const
Returns a reference to the sparse matrix.
Array< Array< int > * > bfbfi_marker
Array< int > vdofs
void FreeElementMatrices()
Free the memory used by the element matrices.
Abstract operator.
Definition: operator.hpp:21
void EliminateEssentialBC(const Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs, int d=0)
virtual void Assemble(int skip_zeros=1)
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0)
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:634
Array< BilinearFormIntegrator * > * GetDI()
void UseSparsity(int *I, int *J, bool isSorted)
Use the given CSR sparsity pattern to allocate the internal SparseMatrix.
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition: sparsemat.cpp:679
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void EnableStaticCondensation()
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator.