MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
bilinearform_ext.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_EXT
13 #define MFEM_BILINEARFORM_EXT
14 
15 #include "../config/config.hpp"
16 #include "fespace.hpp"
17 #include "../general/device.hpp"
18 
19 namespace mfem
20 {
21 
22 class BilinearForm;
23 
24 
25 /** @brief Class extending the BilinearForm class to support the different
26  AssemblyLevel%s. */
28 {
29 protected:
30  BilinearForm *a; ///< Not owned
31 
32 public:
34 
35  virtual MemoryClass GetMemoryClass() const
36  { return Device::GetMemoryClass(); }
37 
38  /// Get the finite element space prolongation matrix
39  virtual const Operator *GetProlongation() const;
40 
41  /// Get the finite element space restriction matrix
42  virtual const Operator *GetRestriction() const;
43 
44  virtual void Assemble() = 0;
45  virtual void FormSystemMatrix(const Array<int> &ess_tdof_list,
46  OperatorHandle &A) = 0;
47  virtual void FormLinearSystem(const Array<int> &ess_tdof_list,
48  Vector &x, Vector &b,
49  OperatorHandle &A, Vector &X, Vector &B,
50  int copy_interior = 0) = 0;
51  virtual void Update() = 0;
52 };
53 
54 /// Data and methods for fully-assembled bilinear forms
56 {
57 public:
59  : BilinearFormExtension(form) { }
60 
61  /// TODO
62  void Assemble() {}
63  void FormSystemMatrix(const Array<int> &ess_tdof_list, OperatorHandle &A) {}
64  void FormLinearSystem(const Array<int> &ess_tdof_list,
65  Vector &x, Vector &b,
66  OperatorHandle &A, Vector &X, Vector &B,
67  int copy_interior = 0) {}
68  void Mult(const Vector &x, Vector &y) const {}
69  void MultTranspose(const Vector &x, Vector &y) const {}
70  void Update() {}
72 };
73 
74 /// Data and methods for element-assembled bilinear forms
76 {
77 public:
79  : BilinearFormExtension(form) { }
80 
81  /// TODO
82  void Assemble() {}
83  void FormSystemMatrix(const Array<int> &ess_tdof_list, OperatorHandle &A) {}
84  void FormLinearSystem(const Array<int> &ess_tdof_list,
85  Vector &x, Vector &b,
86  OperatorHandle &A, Vector &X, Vector &B,
87  int copy_interior = 0) {}
88  void Mult(const Vector &x, Vector &y) const {}
89  void MultTranspose(const Vector &x, Vector &y) const {}
90  void Update() {}
92 };
93 
94 /// Data and methods for partially-assembled bilinear forms
96 {
97 protected:
98  const FiniteElementSpace *trialFes, *testFes; // Not owned
99  mutable Vector localX, localY;
100  const Operator *elem_restrict_lex; // Not owned
101 
102 public:
104 
105  void Assemble();
106  void FormSystemMatrix(const Array<int> &ess_tdof_list, OperatorHandle &A);
107  void FormLinearSystem(const Array<int> &ess_tdof_list,
108  Vector &x, Vector &b,
109  OperatorHandle &A, Vector &X, Vector &B,
110  int copy_interior = 0);
111 
112  void Mult(const Vector &x, Vector &y) const;
113  void MultTranspose(const Vector &x, Vector &y) const;
114  void Update();
115 };
116 
117 /// Data and methods for matrix-free bilinear forms
119 {
120 public:
122  : BilinearFormExtension(form) { }
123 
124  /// TODO
125  void Assemble() {}
126  void FormSystemMatrix(const Array<int> &ess_tdof_list, OperatorHandle &A) {}
127  void FormLinearSystem(const Array<int> &ess_tdof_list,
128  Vector &x, Vector &b,
129  OperatorHandle &A, Vector &X, Vector &B,
130  int copy_interior = 0) {}
131  void Mult(const Vector &x, Vector &y) const {}
132  void MultTranspose(const Vector &x, Vector &y) const {}
133  void Update() {}
135 };
136 
137 }
138 
139 #endif
static MemoryClass GetMemoryClass()
Get the current Device MemoryClass. This is the MemoryClass used by most MFEM device kernels to acces...
Definition: device.hpp:214
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Data and methods for matrix-free bilinear forms.
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
MFBilinearFormExtension(BilinearForm *form)
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
virtual const Operator * GetRestriction() const
Get the finite element space restriction matrix.
BilinearFormExtension(BilinearForm *form)
Data and methods for partially-assembled bilinear forms.
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)=0
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 ...
Data and methods for fully-assembled bilinear forms.
const FiniteElementSpace * trialFes
BilinearForm * a
Not owned.
virtual void Assemble()=0
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Class extending the BilinearForm class to support the different AssemblyLevels.
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
FABilinearFormExtension(BilinearForm *form)
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
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 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 Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Data and methods for element-assembled bilinear forms.
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:85
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
EABilinearFormExtension(BilinearForm *form)
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
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 Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Vector data type.
Definition: vector.hpp:48
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)=0
const FiniteElementSpace * testFes
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
Abstract operator.
Definition: operator.hpp:21
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
MemoryClass
Memory classes identify subsets of memory types.
Definition: mem_manager.hpp:40