DOLFIN
DOLFIN C++ interface
VectorSpaceBasis.h
1 // Copyright (C) 2013-2017 Patrick E. Farrell and Garth N. Wells
2 //
3 // This file is part of DOLFIN.
4 //
5 // DOLFIN is free software: you can redistribute it and/or modify
6 // it under the terms of the GNU Lesser General Public License as published by
7 // the Free Software Foundation, either version 3 of the License, or
8 // (at your option) any later version.
9 //
10 // DOLFIN is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 // GNU Lesser General Public License for more details.
14 //
15 // You should have received a copy of the GNU Lesser General Public License
16 // along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
17 
18 #ifndef __VECTOR_SPACE_BASIS_H
19 #define __VECTOR_SPACE_BASIS_H
20 
21 #include <memory>
22 #include <vector>
23 
24 namespace dolfin
25 {
26 
27  class GenericVector;
28 
32 
34  {
35  public:
36 
38  VectorSpaceBasis(const std::vector<std::shared_ptr<GenericVector>> basis);
39 
42 
46  void orthonormalize(double tol=1.0e-10);
47 
49  bool is_orthonormal(double tol=1.0e-10) const;
50 
52  bool is_orthogonal(double tol=1.0e-10) const;
53 
55  void orthogonalize(GenericVector& x) const;
56 
58  std::size_t dim() const;
59 
61  std::shared_ptr<const GenericVector> operator[] (std::size_t i) const;
62 
63  private:
64 
65  // Basis vectors
66  const std::vector<std::shared_ptr<GenericVector>> _basis;
67 
68  };
69 }
70 
71 #endif
~VectorSpaceBasis()
Destructor.
Definition: VectorSpaceBasis.h:41
VectorSpaceBasis(const std::vector< std::shared_ptr< GenericVector >> basis)
Constructor.
Definition: VectorSpaceBasis.cpp:29
Definition: adapt.h:29
bool is_orthonormal(double tol=1.0e-10) const
Test if basis is orthonormal.
Definition: VectorSpaceBasis.cpp:60
void orthonormalize(double tol=1.0e-10)
Definition: VectorSpaceBasis.cpp:35
std::size_t dim() const
Number of vectors in the basis.
Definition: VectorSpaceBasis.cpp:108
Definition: VectorSpaceBasis.h:33
bool is_orthogonal(double tol=1.0e-10) const
Test if basis is orthogonal.
Definition: VectorSpaceBasis.cpp:78
std::shared_ptr< const GenericVector > operator[](std::size_t i) const
Get a particular basis vector.
Definition: VectorSpaceBasis.cpp:114
void orthogonalize(GenericVector &x) const
Orthogonalize x with respect to basis.
Definition: VectorSpaceBasis.cpp:98
This class defines a common interface for vectors.
Definition: GenericVector.h:47