MFEM  v3.1
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 
29 class BilinearForm : public Matrix
30 {
31 protected:
34 
35  // Matrix used to eliminate b.c.
37 
40 
42 
45 
48 
51 
54 
57 
59 
62 
64  // Allocate appropriate SparseMatrix and assign it to mat
65  void AllocMat();
66 
67  // may be used in the construction of derived classes
69  {
70  fes = NULL; mat = mat_e = NULL; extern_bfs = 0; element_matrices = NULL;
71  static_cond = NULL; hybridization = NULL;
73  }
74 
75 public:
78 
79  BilinearForm(FiniteElementSpace *f, BilinearForm *bf, int ps = 0);
80 
82  int Size() const { return height; }
83 
89 
92  bool StaticCondensationIsEnabled() const { return static_cond; }
93 
96  { return static_cond->GetTraceFESpace(); }
97 
101  void EnableHybridization(FiniteElementSpace *constr_space,
102  BilinearFormIntegrator *constr_integ,
103  const Array<int> &ess_tdof_list);
104 
108  void UsePrecomputedSparsity(int ps = 1) { precompute_sparsity = ps; }
109 
113  void AllocateMatrix() { if (mat == NULL) { AllocMat(); } }
114 
116 
118 
120 
122 
123  const double &operator()(int i, int j) { return (*mat)(i,j); }
124 
126  virtual double &Elem(int i, int j);
127 
129  virtual const double &Elem(int i, int j) const;
130 
132  virtual void Mult(const Vector &x, Vector &y) const;
133 
134  void FullMult(const Vector &x, Vector &y) const
135  { mat->Mult(x, y); mat_e->AddMult(x, y); }
136 
137  virtual void AddMult(const Vector &x, Vector &y, const double a = 1.0) const
138  { mat -> AddMult (x, y, a); }
139 
140  void FullAddMult(const Vector &x, Vector &y) const
141  { mat->AddMult(x, y); mat_e->AddMult(x, y); }
142 
143  double InnerProduct(const Vector &x, const Vector &y) const
144  { return mat->InnerProduct (x, y); }
145 
147  virtual MatrixInverse *Inverse() const;
148 
150  virtual void Finalize(int skip_zeros = 1);
151 
153  const SparseMatrix &SpMat() const { return *mat; }
154  SparseMatrix &SpMat() { return *mat; }
155  SparseMatrix *LoseMat() { SparseMatrix *tmp = mat; mat = NULL; return tmp; }
156 
159 
162 
165 
168 
169  void operator=(const double a)
170  {
171  if (mat != NULL) { *mat = a; }
172  if (mat_e != NULL) { *mat_e = a; }
173  }
174 
176  void Assemble(int skip_zeros = 1);
177 
182  void ConformingAssemble();
183 
186  {
188  rhs.ConformingAssemble();
189  sol.ConformingProject();
190  }
191 
212  void FormLinearSystem(Array<int> &ess_tdof_list, Vector &x, Vector &b,
213  SparseMatrix &A, Vector &X, Vector &B,
214  int copy_interior = 0);
215 
219  void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x);
220 
222  void ComputeElementMatrices();
223 
226  { delete element_matrices; element_matrices = NULL; }
227 
228  void ComputeElementMatrix(int i, DenseMatrix &elmat);
229  void AssembleElementMatrix(int i, const DenseMatrix &elmat,
230  Array<int> &vdofs, int skip_zeros = 1);
231 
236  void EliminateEssentialBC(Array<int> &bdr_attr_is_ess,
237  Vector &sol, Vector &rhs, int d = 0);
238 
239  void EliminateEssentialBC(Array<int> &bdr_attr_is_ess, int d = 0);
241  void EliminateEssentialBCDiag(Array<int> &bdr_attr_is_ess, double value);
242 
244  void EliminateVDofs(Array<int> &vdofs, Vector &sol, Vector &rhs, int d = 0);
245 
250  void EliminateVDofs(Array<int> &vdofs, int d = 0);
251 
254  void EliminateEssentialBCFromDofs(Array<int> &ess_dofs, Vector &sol,
255  Vector &rhs, int d = 0);
256 
259  void EliminateEssentialBCFromDofs(Array<int> &ess_dofs, int d = 0);
261  void EliminateEssentialBCFromDofsDiag(Array<int> &ess_dofs, double value);
262 
265  void EliminateVDofsInRHS(Array<int> &vdofs, const Vector &x, Vector &b);
266 
267  double FullInnerProduct(const Vector &x, const Vector &y) const
268  { return mat->InnerProduct(x, y) + mat_e->InnerProduct(x, y); }
269 
270  virtual void Update(FiniteElementSpace *nfes = NULL);
271 
274 
276  virtual ~BilinearForm();
277 };
278 
294 class MixedBilinearForm : public Matrix
295 {
296 protected:
298 
300 
303  Array<BilinearFormIntegrator*> skt; // trace face integrators
304 
305 public:
307  FiniteElementSpace *te_fes);
308 
309  virtual double& Elem (int i, int j);
310 
311  virtual const double& Elem (int i, int j) const;
312 
313  virtual void Mult (const Vector & x, Vector & y) const;
314 
315  virtual void AddMult (const Vector & x, Vector & y,
316  const double a = 1.0) const;
317 
318  virtual void AddMultTranspose (const Vector & x, Vector & y,
319  const double a = 1.0) const;
320 
321  virtual void MultTranspose (const Vector & x, Vector & y) const
322  { y = 0.0; AddMultTranspose (x, y); }
323 
324  virtual MatrixInverse * Inverse() const;
325 
326  virtual void Finalize (int skip_zeros = 1);
327 
331  void GetBlocks(Array2D<SparseMatrix *> &blocks) const;
332 
333  const SparseMatrix &SpMat() const { return *mat; }
334  SparseMatrix &SpMat() { return *mat; }
335  SparseMatrix *LoseMat() { SparseMatrix *tmp = mat; mat = NULL; return tmp; }
336 
338 
340 
345 
347 
349 
351 
352  void operator= (const double a) { *mat = a; }
353 
354  void Assemble (int skip_zeros = 1);
355 
361  void ConformingAssemble();
362 
363  void EliminateTrialDofs(Array<int> &bdr_attr_is_ess,
364  Vector &sol, Vector &rhs);
365 
367  Vector &sol, Vector &rhs);
368 
369  virtual void EliminateTestDofs(Array<int> &bdr_attr_is_ess);
370 
371  void Update();
372 
373  virtual ~MixedBilinearForm();
374 };
375 
376 
408 {
409 public:
411  FiniteElementSpace *range_fes)
412  : MixedBilinearForm(domain_fes, range_fes) { }
413 
415  { AddDomainIntegrator(di); }
416 
418  { AddTraceFaceIntegrator(di); }
419 
421 
422  virtual void Assemble(int skip_zeros = 1);
423 };
424 
425 }
426 
427 #endif
void EliminateEssentialBCFromDofs(Array< int > &ess_dofs, Vector &sol, Vector &rhs, int d=0)
void FullAddMult(const Vector &x, Vector &y) const
SparseMatrix & SpMat()
Array< BilinearFormIntegrator * > * GetBBFI()
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
Array< BilinearFormIntegrator * > skt
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.
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()
void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
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:423
Abstract data type for matrix inverse.
Definition: matrix.hpp:58
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
Array< BilinearFormIntegrator * > * GetFBFI()
MixedBilinearForm(FiniteElementSpace *tr_fes, FiniteElementSpace *te_fes)
void FormLinearSystem(Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0)
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 EliminateEssentialBCFromTrialDofs(Array< int > &marked_vdofs, Vector &sol, Vector &rhs)
Array< BilinearFormIntegrator * > bdr
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
SparseMatrix & SpMat()
void EliminateVDofs(Array< int > &vdofs, Vector &sol, Vector &rhs, int d=0)
Eliminate the given vdofs. NOTE: here, vdofs is a list of DOFs.
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
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 Update(FiniteElementSpace *nfes=NULL)
void EliminateTrialDofs(Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs)
FiniteElementSpace * GetFES()
Return the FE space associated with the BilinearForm.
void ConformingAssemble(GridFunction &sol, LinearForm &rhs)
void AddTraceFaceIntegrator(BilinearFormIntegrator *bfi)
void GetBlocks(Array2D< SparseMatrix * > &blocks) const
Array< BilinearFormIntegrator * > * GetTFBFI()
Array< BilinearFormIntegrator * > * GetBBFI()
void ComputeElementMatrix(int i, DenseMatrix &elmat)
void EliminateEssentialBCDiag(Array< int > &bdr_attr_is_ess, double value)
Perform elimination and set the diagonal entry to the given value.
SparseMatrix * mat_e
bool StaticCondensationIsEnabled() const
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:22
void AssembleElementMatrix(int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1)
void UsePrecomputedSparsity(int ps=1)
Abstract finite element space.
Definition: fespace.hpp:62
SparseMatrix * LoseMat()
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:417
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)
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
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:107
virtual void EliminateTestDofs(Array< int > &bdr_attr_is_ess)
void EliminateVDofsInRHS(Array< int > &vdofs, const Vector &x, Vector &b)
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void ConformingProject(Vector &x) const
Definition: gridfunc.cpp:2114
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.
void ComputeElementMatrices()
Compute and store internally all element matrices.
FiniteElementSpace * SCFESpace() const
Return the trace FE space associated with static condensation.
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:33
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)
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
Array< BilinearFormIntegrator * > * GetDBFI()
const SparseMatrix & SpMat() const
Returns a reference to the sparse matrix.
Array< int > vdofs
void FreeElementMatrices()
Free the memory used by the element matrices.
void ConformingAssemble(Vector &b) const
Apply the conforming interpolation matrix and return &#39;b&#39;: b = P&#39;*this.
Definition: linearform.cpp:89
void EliminateEssentialBCFromDofsDiag(Array< int > &ess_dofs, double value)
Perform elimination and set the diagonal entry to the given value.
virtual void Assemble(int skip_zeros=1)
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:565
Array< BilinearFormIntegrator * > * GetDI()
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition: sparsemat.cpp:619
void EliminateEssentialBC(Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs, int d=0)
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
void EnableStaticCondensation()
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator.