Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
NOX_Epetra_LinearSystem_SGGS.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef NOX_EPETRA_LINEARSYSTEMSGGS_H
11 #define NOX_EPETRA_LINEARSYSTEMSGGS_H
12 
13 #include "Stokhos_ConfigDefs.h"
14 
15 #ifdef HAVE_STOKHOS_NOX
16 
17 #include "NOX_Common.H"
18 
19 #include "NOX_Epetra_LinearSystem.H" // base class
20 #include "NOX_Utils.H" // class data element
21 
22 #include "Stokhos_ParallelData.hpp"
27 #include "Stokhos_SGOperator.hpp"
29 #include "Epetra_Export.h"
30 
31 namespace NOX {
32  namespace Epetra {
33  namespace Interface {
34  class Required;
35  class Jacobian;
36  class Preconditioner;
37  }
38  }
39 }
40 
41 namespace NOX {
42 
43  namespace Epetra {
44 
49  class LinearSystemSGGS : public virtual NOX::Epetra::LinearSystem {
50  public:
51 
53  LinearSystemSGGS(
54  Teuchos::ParameterList& printingParams,
55  Teuchos::ParameterList& linearSolverParams,
59  const Teuchos::RCP<const Stokhos::OrthogPolyBasis<int,double> >& sg_basis,
60  const Teuchos::RCP<const Stokhos::ParallelData>& sg_parallel_data,
62  const Teuchos::RCP<const Epetra_Map>& base_map,
63  const Teuchos::RCP<const Epetra_Map>& sg_map,
64  const Teuchos::RCP<NOX::Epetra::Scaling> scalingObject =
66 
68  virtual ~LinearSystemSGGS();
69 
74  virtual bool applyJacobian(const NOX::Epetra::Vector& input,
75  NOX::Epetra::Vector& result) const;
76 
81  virtual bool applyJacobianTranspose(const NOX::Epetra::Vector& input,
82  NOX::Epetra::Vector& result) const;
83 
88  virtual bool applyJacobianInverse(Teuchos::ParameterList &params,
89  const NOX::Epetra::Vector &input,
90  NOX::Epetra::Vector &result);
91 
93  virtual bool applyRightPreconditioning(bool useTranspose,
94  Teuchos::ParameterList& params,
95  const NOX::Epetra::Vector& input,
96  NOX::Epetra::Vector& result) const;
97 
99  virtual Teuchos::RCP<NOX::Epetra::Scaling> getScaling();
100 
102  virtual void resetScaling(const Teuchos::RCP<NOX::Epetra::Scaling>& s);
103 
105  virtual bool computeJacobian(const NOX::Epetra::Vector& x);
106 
108  virtual bool createPreconditioner(const NOX::Epetra::Vector& x,
110  bool recomputeGraph) const;
111 
113  virtual bool destroyPreconditioner() const;
114 
116  virtual bool recomputePreconditioner(const NOX::Epetra::Vector& x,
117  Teuchos::ParameterList& linearSolverParams) const;
118 
120  virtual PreconditionerReusePolicyType
121  getPreconditionerPolicy(bool advanceReuseCounter=true);
122 
124  virtual bool isPreconditionerConstructed() const;
125 
127  virtual bool hasPreconditioner() const;
128 
131  getJacobianOperator() const;
132 
134  virtual Teuchos::RCP<Epetra_Operator> getJacobianOperator();
135 
138  getGeneratedPrecOperator() const;
139 
141  virtual Teuchos::RCP<Epetra_Operator> getGeneratedPrecOperator();
142 
144  virtual void setJacobianOperatorForSolve(const Teuchos::RCP<const Epetra_Operator>& solveJacOp);
145 
147  virtual void setPrecOperatorForSolve(const Teuchos::RCP<const Epetra_Operator>& solvePrecOp);
148 
149  protected:
150 
153 
156 
159 
161  bool is_stoch_parallel;
162 
165 
168 
171 
174 
176  mutable Teuchos::RCP<Stokhos::SGOperator> sg_op;
177 
180 
183 
186 
189 
191  NOX::Utils utils;
192 
194  mutable Teuchos::RCP<EpetraExt::BlockVector> sg_df_block;
195 
197  mutable Teuchos::RCP<EpetraExt::BlockVector> sg_y_block;
198 
200  mutable Teuchos::RCP<Epetra_Vector> kx;
201 
203  bool is_parallel;
204 
207 
209  Teuchos::RCP<Epetra_Export> col_exporter;
210 
213 
216 
218  bool only_use_linear;
219 
221  int k_limit;
222  };
223 
224  } // namespace Epetra
225 } // namespace NOX
226 
227 #endif
228 
229 #endif
Stokhos::Sparse3Tensor< int, double > Cijk_type