IFPACK  Development
 All Classes Namespaces Files Functions Variables Enumerations Friends Pages
Ifpack_DiagPreconditioner.h
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 #if defined(Ifpack_SHOW_DEPRECATED_WARNINGS)
47 #ifdef __GNUC__
48 #warning "The Ifpack package is deprecated"
49 #endif
50 #endif
51 
52 #include "Ifpack_ConfigDefs.h"
53 #include "Epetra_Operator.h"
54 #include "Epetra_Vector.h"
55 class Epetra_BlockMap;
56 class Epetra_Map;
57 class Epetra_MultiVector;
58 class Epetra_Comm;
59 
61 /*
62 Ifpack_DiagPreconditioner: a class to wrap a vector as diagonal preconditioner. The preconditioner is simply defined by
63 \f[
64 z_i = D_i r_i,
65 \f]
66 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.
67 
68 \author Marzio Sala, ETHZ/D-INFK
69 
70 \date Last updated on 17-Apr-06
71 
72  */
74 {
75  public:
76 
78  Ifpack_DiagPreconditioner(const Epetra_Map& DomainMap,
79  const Epetra_Map& RangeMap,
80  const Epetra_Vector& diag);
81 
84 
85  int SetUseTranspose(bool UseTranspose_in)
86  {
87  UseTranspose_ = UseTranspose_in;
88  return(0);
89  }
90 
91  int Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
92 
93  int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
94 
95  double NormInf() const
96  {
97  return(-1.0);
98  }
99 
100  const char* Label() const
101  {
102  return("Ifpack_DiagPreconditioner");
103  }
104 
105  bool UseTranspose() const
106  {
107  return(UseTranspose_);
108  }
109 
110  bool HasNormInf() const
111  {
112  return(false);
113  }
114 
115  const Epetra_Comm& Comm() const
116  {
117  return(diag_.Comm());
118  }
119 
120  const Epetra_Map& OperatorDomainMap() const
121  {
122  return(RangeMap_);
123  }
124 
125  const Epetra_Map& OperatorRangeMap() const
126  {
127  return(DomainMap_);
128  }
129 
130  const Epetra_BlockMap& Map() const
131  {
132  return(diag_.Map());
133  }
134 
135  private:
136  bool UseTranspose_;
137  const Epetra_Map& DomainMap_;
138  const Epetra_Map& RangeMap_;
139  const Epetra_Vector& diag_;
140 };
141 
142 #endif
Ifpack_DiagPreconditioner: a class for diagonal preconditioning.
int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Ifpack_DiagPreconditioner(const Epetra_Map &DomainMap, const Epetra_Map &RangeMap, const Epetra_Vector &diag)
ctor