MFEM  v3.0
 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.googlecode.com.
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 
22 namespace mfem
23 {
24 
27 class BilinearForm : public Matrix
28 {
29 protected:
32 
33  // Matrix used to eliminate b.c.
35 
38 
40 
43 
46 
49 
52 
55 
57 
59  // Allocate appropriate SparseMatrix and assign it to mat
60  void AllocMat();
61 
62  // may be used in the construction of derived classes
64  { fes = NULL; mat = mat_e = NULL; extern_bfs = 0; element_matrices = NULL;
65  precompute_sparsity = 0; }
66 
67 public:
70 
71  BilinearForm(FiniteElementSpace *f, BilinearForm *bf, int ps = 0);
72 
74  int Size() const { return height; }
75 
79  void UsePrecomputedSparsity(int ps = 1) { precompute_sparsity = ps; }
80 
84  void AllocateMatrix() { if (mat == NULL) AllocMat(); }
85 
87 
89 
91 
93 
94  const double &operator()(int i, int j) { return (*mat)(i,j); }
95 
97  virtual double &Elem(int i, int j);
98 
100  virtual const double &Elem(int i, int j) const;
101 
103  virtual void Mult(const Vector &x, Vector &y) const;
104 
105  void FullMult(const Vector &x, Vector &y) const
106  { mat->Mult(x, y); mat_e->AddMult(x, y); }
107 
108  virtual void AddMult(const Vector &x, Vector &y, const double a = 1.0) const
109  { mat -> AddMult (x, y, a); }
110 
111  void FullAddMult(const Vector &x, Vector &y) const
112  { mat->AddMult(x, y); mat_e->AddMult(x, y); }
113 
114  double InnerProduct(const Vector &x, const Vector &y) const
115  { return mat->InnerProduct (x, y); }
116 
118  virtual MatrixInverse *Inverse() const;
119 
121  virtual void Finalize(int skip_zeros = 1);
122 
124  const SparseMatrix &SpMat() const { return *mat; }
125  SparseMatrix &SpMat() { return *mat; }
126  SparseMatrix *LoseMat() { SparseMatrix *tmp = mat; mat = NULL; return tmp; }
127 
130 
133 
136 
139 
140  void operator=(const double a)
141  { if (mat != NULL) *mat = a; if (mat_e != NULL) *mat_e = a; }
142 
144  void Assemble(int skip_zeros = 1);
145 
150  void ConformingAssemble();
151 
154  {
156  rhs.ConformingAssemble();
157  sol.ConformingProject();
158  }
159 
161  void ComputeElementMatrices();
162 
165  { delete element_matrices; element_matrices = NULL; }
166 
167  void ComputeElementMatrix(int i, DenseMatrix &elmat);
168  void AssembleElementMatrix(int i, const DenseMatrix &elmat,
169  Array<int> &vdofs, int skip_zeros = 1);
170 
173  void EliminateEssentialBC(Array<int> &bdr_attr_is_ess,
174  Vector &sol, Vector &rhs, int d = 0);
175 
177  void EliminateVDofs(Array<int> &vdofs, Vector &sol, Vector &rhs, int d = 0);
178 
181  void EliminateVDofs(Array<int> &vdofs, int d = 0);
182 
185  void EliminateVDofsInRHS(Array<int> &vdofs, const Vector &x, Vector &b);
186 
187  double FullInnerProduct(const Vector &x, const Vector &y) const
188  { return mat->InnerProduct(x, y) + mat_e->InnerProduct(x, y); }
189 
190  void EliminateEssentialBC(Array<int> &bdr_attr_is_ess, int d = 0);
191 
194  void EliminateEssentialBCFromDofs(Array<int> &ess_dofs, Vector &sol,
195  Vector &rhs, int d = 0);
196 
199  void EliminateEssentialBCFromDofs(Array<int> &ess_dofs, int d = 0);
200 
201  void Update(FiniteElementSpace *nfes = NULL);
202 
204 
206  virtual ~BilinearForm();
207 };
208 
224 class MixedBilinearForm : public Matrix
225 {
226 protected:
228 
230 
233  Array<BilinearFormIntegrator*> skt; // trace face integrators
234 
235 public:
237  FiniteElementSpace *te_fes);
238 
239  virtual double& Elem (int i, int j);
240 
241  virtual const double& Elem (int i, int j) const;
242 
243  virtual void Mult (const Vector & x, Vector & y) const;
244 
245  virtual void AddMult (const Vector & x, Vector & y,
246  const double a = 1.0) const;
247 
248  virtual void AddMultTranspose (const Vector & x, Vector & y,
249  const double a = 1.0) const;
250 
251  virtual void MultTranspose (const Vector & x, Vector & y) const
252  { y = 0.0; AddMultTranspose (x, y); }
253 
254  virtual MatrixInverse * Inverse() const;
255 
256  virtual void Finalize (int skip_zeros = 1);
257 
261  void GetBlocks(Array2D<SparseMatrix *> &blocks) const;
262 
263  const SparseMatrix &SpMat() const { return *mat; }
264  SparseMatrix &SpMat() { return *mat; }
265  SparseMatrix *LoseMat() { SparseMatrix *tmp = mat; mat = NULL; return tmp; }
266 
268 
270 
275 
277 
279 
281 
282  void operator= (const double a) { *mat = a; }
283 
284  void Assemble (int skip_zeros = 1);
285 
291  void ConformingAssemble();
292 
293  void EliminateTrialDofs(Array<int> &bdr_attr_is_ess,
294  Vector &sol, Vector &rhs);
295 
297  Vector &sol, Vector &rhs);
298 
299  virtual void EliminateTestDofs(Array<int> &bdr_attr_is_ess);
300 
301  void Update();
302 
303  virtual ~MixedBilinearForm();
304 };
305 
306 
338 {
339 public:
341  FiniteElementSpace *range_fes)
342  : MixedBilinearForm(domain_fes, range_fes) { }
343 
345  { AddDomainIntegrator(di); }
346 
348 
349  virtual void Assemble(int skip_zeros = 1);
350 };
351 
352 }
353 
354 #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:26
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()
Data type dense matrix.
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:318
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)
double FullInnerProduct(const Vector &x, const Vector &y) const
void AddDomainInterpolator(DiscreteInterpolator *di)
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)
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
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
void Update(FiniteElementSpace *nfes=NULL)
void EliminateTrialDofs(Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs)
FiniteElementSpace * GetFES()
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)
SparseMatrix * mat_e
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:61
SparseMatrix * LoseMat()
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:312
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 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)
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:1773
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.
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:29
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
virtual MatrixInverse * Inverse() const
Returns a pointer to (an approximation) of the matrix inverse.
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 martix.
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
virtual void Assemble(int skip_zeros=1)
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:428
Array< BilinearFormIntegrator * > * GetDI()
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition: sparsemat.cpp:468
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 AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator.