Reference documentation for deal.II version 8.5.1
Public Member Functions | Private Attributes | List of all members
ProductMatrix< VectorType > Class Template Reference

#include <deal.II/lac/matrix_lib.h>

Inheritance diagram for ProductMatrix< VectorType >:
[legend]

Public Member Functions

 ProductMatrix ()
 
 ProductMatrix (VectorMemory< VectorType > &mem)
 
template<typename MatrixType1 , typename MatrixType2 >
 ProductMatrix (const MatrixType1 &m1, const MatrixType2 &m2, VectorMemory< VectorType > &mem)
 
 ~ProductMatrix ()
 
template<typename MatrixType1 , typename MatrixType2 >
void reinit (const MatrixType1 &m1, const MatrixType2 &m2)
 
template<typename MatrixType1 , typename MatrixType2 >
void initialize (const MatrixType1 &m1, const MatrixType2 &m2, VectorMemory< VectorType > &mem)
 
void clear ()
 
virtual void vmult (VectorType &w, const VectorType &v) const
 
virtual void Tvmult (VectorType &w, const VectorType &v) const
 
virtual void vmult_add (VectorType &w, const VectorType &v) const
 
virtual void Tvmult_add (VectorType &w, const VectorType &v) const
 
- Public Member Functions inherited from PointerMatrixBase< VectorType >
virtual ~PointerMatrixBase ()
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&)
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&)
 
void subscribe (const char *identifier=0) const
 
void unsubscribe (const char *identifier=0) const
 
unsigned int n_subscriptions () const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Private Attributes

PointerMatrixBase< VectorType > * m1
 
PointerMatrixBase< VectorType > * m2
 
SmartPointer< VectorMemory< VectorType >, ProductMatrix< VectorType > > mem
 

Additional Inherited Members

- Public Types inherited from PointerMatrixBase< VectorType >
typedef VectorType::value_type value_type
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, char *arg2, std::string &arg3)
 
static ::ExceptionBaseExcNoSubscriber (char *arg1, char *arg2)
 

Detailed Description

template<typename VectorType>
class ProductMatrix< VectorType >

Poor man's matrix product of two quadratic matrices. Stores two quadratic matrices m1 and m2 of arbitrary types and implements matrix-vector multiplications for the product M1M2 by performing multiplication with both factors consecutively. Because the types of the matrices are opaque (i.e., they can be arbitrary), you can stack products of three or more matrices by making one of the two matrices an object of the current type handles be a ProductMatrix itself.

Here is an example multiplying two different FullMatrix objects:

// ---------------------------------------------------------------------
//
// Copyright (C) 2005 - 2014 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the deal.II distribution.
//
// ---------------------------------------------------------------------
// See documentation of ProductMatrix for documentation of this example
#include <deal.II/base/logstream.h>
#include <deal.II/lac/matrix_lib.h>
#include <deal.II/lac/full_matrix.h>
#include <deal.II/lac/vector.h>
using namespace dealii;
double Adata[] =
{
.5, .1,
.4, .2
};
double Bdata[] =
{
.866, .5,
-.5, .866
};
int main()
{
A.fill(Adata);
B.fill(Bdata);
u(0) = 1.;
u(1) = 2.;
AB.vmult(v,u);
deallog << v(0) << '\t' << v(1) << std::endl;
AB.Tvmult(v,u);
deallog << v(0) << '\t' << v(1) << std::endl;
}
Deprecated:
If deal.II was configured with C++11 support, use the LinearOperator class instead, see the module on linear operators for further details.
Author
Guido Kanschat, 2000, 2001, 2002, 2005

Definition at line 55 of file matrix_lib.h.

Constructor & Destructor Documentation

◆ ProductMatrix() [1/3]

template<typename VectorType >
ProductMatrix< VectorType >::ProductMatrix ( )

Standard constructor. Matrices and the memory pool must be added later using initialize().

Definition at line 591 of file matrix_lib.h.

◆ ProductMatrix() [2/3]

template<typename VectorType >
ProductMatrix< VectorType >::ProductMatrix ( VectorMemory< VectorType > &  mem)

Constructor only assigning the memory pool. Matrices must be added by reinit() later.

Definition at line 597 of file matrix_lib.h.

◆ ProductMatrix() [3/3]

template<typename VectorType >
template<typename MatrixType1 , typename MatrixType2 >
ProductMatrix< VectorType >::ProductMatrix ( const MatrixType1 &  m1,
const MatrixType2 &  m2,
VectorMemory< VectorType > &  mem 
)

Constructor. Additionally to the two constituting matrices, a memory pool for the auxiliary vector must be provided.

Definition at line 604 of file matrix_lib.h.

◆ ~ProductMatrix()

template<typename VectorType >
ProductMatrix< VectorType >::~ProductMatrix ( )

Destructor.

Definition at line 642 of file matrix_lib.h.

Member Function Documentation

◆ reinit()

template<typename VectorType >
template<typename MatrixType1 , typename MatrixType2 >
void ProductMatrix< VectorType >::reinit ( const MatrixType1 &  m1,
const MatrixType2 &  m2 
)

Change the matrices.

Definition at line 617 of file matrix_lib.h.

◆ initialize()

template<typename VectorType >
template<typename MatrixType1 , typename MatrixType2 >
void ProductMatrix< VectorType >::initialize ( const MatrixType1 &  m1,
const MatrixType2 &  m2,
VectorMemory< VectorType > &  mem 
)

Change the matrices and memory pool.

Definition at line 629 of file matrix_lib.h.

◆ clear()

template<typename VectorType >
void ProductMatrix< VectorType >::clear ( )
virtual

Reset the object to its original state.

Implements PointerMatrixBase< VectorType >.

Definition at line 651 of file matrix_lib.h.

◆ vmult()

template<typename VectorType >
void ProductMatrix< VectorType >::vmult ( VectorType &  w,
const VectorType &  v 
) const
virtual

Matrix-vector product w = m1 * m2 * v.

Implements PointerMatrixBase< VectorType >.

Definition at line 662 of file matrix_lib.h.

◆ Tvmult()

template<typename VectorType >
void ProductMatrix< VectorType >::Tvmult ( VectorType &  w,
const VectorType &  v 
) const
virtual

Transposed matrix-vector product w = m2T * m1T * v.

Implements PointerMatrixBase< VectorType >.

Definition at line 694 of file matrix_lib.h.

◆ vmult_add()

template<typename VectorType >
void ProductMatrix< VectorType >::vmult_add ( VectorType &  w,
const VectorType &  v 
) const
virtual

Adding matrix-vector product w += m1 * m2 * v

Implements PointerMatrixBase< VectorType >.

Definition at line 678 of file matrix_lib.h.

◆ Tvmult_add()

template<typename VectorType >
void ProductMatrix< VectorType >::Tvmult_add ( VectorType &  w,
const VectorType &  v 
) const
virtual

Adding, transposed matrix-vector product w += m2T * m1T * v.

Implements PointerMatrixBase< VectorType >.

Definition at line 710 of file matrix_lib.h.

Member Data Documentation

◆ m1

template<typename VectorType>
PointerMatrixBase<VectorType>* ProductMatrix< VectorType >::m1
private

The left matrix of the product.

Definition at line 131 of file matrix_lib.h.

◆ m2

template<typename VectorType>
PointerMatrixBase<VectorType>* ProductMatrix< VectorType >::m2
private

The right matrix of the product.

Definition at line 136 of file matrix_lib.h.

◆ mem

template<typename VectorType>
SmartPointer<VectorMemory<VectorType>,ProductMatrix<VectorType> > ProductMatrix< VectorType >::mem
private

Memory for auxiliary vector.

Definition at line 141 of file matrix_lib.h.


The documentation for this class was generated from the following file: