Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Ifpack_DiagPreconditioner.h
Go to the documentation of this file.
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack: Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2002) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 */
42 
43 #ifndef IFPACK_DIAG_PRECONDITIONER_H
44 #define IFPACK_DIAG_PRECONDITIONER_H
45 
46 #include "Ifpack_ConfigDefs.h"
47 #include "Epetra_Operator.h"
48 #include "Epetra_Vector.h"
49 class Epetra_BlockMap;
50 class Epetra_Map;
51 class Epetra_MultiVector;
52 class Epetra_Comm;
53 
55 /*
56 Ifpack_DiagPreconditioner: a class to wrap a vector as diagonal preconditioner. The preconditioner is simply defined by
57 \f[
58 z_i = D_i r_i,
59 \f]
60 where \f$r,z\f$ are the vector to be preconditioned and the preconditioned vector, and \f$D_i\f$ is the i-th element of the scaling vector.
61 
62 \author Marzio Sala, ETHZ/D-INFK
63 
64 \date Last updated on 17-Apr-06
65 
66  */
68 {
69  public:
70 
72  Ifpack_DiagPreconditioner(const Epetra_Map& DomainMap,
73  const Epetra_Map& RangeMap,
74  const Epetra_Vector& diag);
75 
78 
79  int SetUseTranspose(bool UseTranspose_in)
80  {
81  UseTranspose_ = UseTranspose_in;
82  return(0);
83  }
84 
85  int Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
86 
87  int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
88 
89  double NormInf() const
90  {
91  return(-1.0);
92  }
93 
94  const char* Label() const
95  {
96  return("Ifpack_DiagPreconditioner");
97  }
98 
99  bool UseTranspose() const
100  {
101  return(UseTranspose_);
102  }
103 
104  bool HasNormInf() const
105  {
106  return(false);
107  }
108 
109  const Epetra_Comm& Comm() const
110  {
111  return(diag_.Comm());
112  }
113 
115  {
116  return(RangeMap_);
117  }
118 
120  {
121  return(DomainMap_);
122  }
123 
124  const Epetra_BlockMap& Map() const
125  {
126  return(diag_.Map());
127  }
128 
129  private:
134 };
135 
136 #endif
137 
138 #if defined(Ifpack_SHOW_DEPRECATED_WARNINGS)
139 #ifdef __GNUC__
140 #warning "The Ifpack package is deprecated"
141 #endif
142 #endif
143 
const Epetra_Comm & Comm() const
const Epetra_Map & OperatorRangeMap() const
int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Ifpack_DiagPreconditioner: a class for diagonal preconditioning.
int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
const Epetra_Map & OperatorDomainMap() const
Ifpack_DiagPreconditioner(const Epetra_Map &DomainMap, const Epetra_Map &RangeMap, const Epetra_Vector &diag)
ctor
virtual const Epetra_BlockMap & Map() const =0
const Epetra_BlockMap & Map() const
int SetUseTranspose(bool UseTranspose_in)