MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AbstractLinAlgPack_DirectSparseSolverMA28.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Moocho: Multi-functional Object-Oriented arCHitecture for Optimization
5 // Copyright (2003) 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 Roscoe A. Bartlett (rabartl@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef ALAP_DIRECT_SPARSE_SOLVER_MA28_H
43 #define ALAP_DIRECT_SPARSE_SOLVER_MA28_H
44 
45 #include <valarray>
46 #include <vector>
47 #include <string>
48 
54 
55 namespace AbstractLinAlgPack {
56 
62 public:
63 
66 
69 
71 
74 
77 
80 
83 
86 
89 
92 
94  STANDARD_MEMBER_COMPOSITION_MEMBERS( std::string, output_file_name );
95 
97 
100 
104  ,value_type u = 0.1
105  ,bool grow = false
106  ,value_type tol = 0.0
107  ,index_type nsrch = 4
108  ,bool lbig = false
109  ,bool print_ma28_outputs = false
110  ,const std::string& output_file_name = ""
111  );
112 
114 
117 
122 
124 
125 protected:
126 
129 
132  class BasisMatrixMA28 : public BasisMatrixImp {
133  public:
134 
137 
141  void V_InvMtV(
142  VectorMutable* v_lhs, BLAS_Cpp::Transp trans_rhs1
143  ,const Vector& v_rhs2) const ;
144 
146 
147  }; // end class BasisMatrixMA28
148 
151  class FactorizationStructureMA28 : public FactorizationStructure {
152  public:
156  friend class BasisMatrixMA28;
157  private:
158  // /////////////////////////////////////////
159  // Private types
161  // /////////////////////////////////////////
162  // Private data members
163  mutable MA28_Cpp::MA28Solver ma28_; // Management of common block data
164  // Keep a memory of the size of the system to check for consistent usage.
165  index_type m_; // number of rows (keep for checks on consistancy)
166  index_type n_; // number of columns ("")
167  index_type max_n_; // max(m_,n_)
168  index_type nz_; // numner of non-zero elements in unfactorized matrix ("")
169  index_type licn_; // size of icn_ and a_ (default = 3 * nz_)
170  index_type lirn_; // size of irn_ (default = 3 * nz_)
171  // Control parameters
172  value_type u_; // fill-in vs. stability ratio (default = 0.1)
173  EScalingMethod scaling_; // Scaling method
174  // Matrix scaling
176  // Memory for factorization structure
177  std::valarray<index_type> ivect_;
178  std::valarray<index_type> jvect_;
179  std::valarray<index_type> icn_;
180  std::valarray<index_type> ikeep_;
181  // Basis matrix selection
182  IVector row_perm_;
183  IVector col_perm_;
185  // /////////////////////////////////////////
186  // Private member functions
189  }; // end class FactorizationStructureMA28
190 
194  public:
198  friend class BasisMatrixMA28;
199  private:
200  std::valarray<value_type> a_; // holds the non-zeros of the factorized matrix 'a'
201  }; // end class FactorizationNonzerosMA28
202 
204 
207 
216  ,FactorizationNonzeros *fact_nonzeros
217  ,DenseLinAlgPack::IVector *row_perm
218  ,DenseLinAlgPack::IVector *col_perm
219  ,size_type *rank
220  ,std::ostream *out
221  );
223  void imp_factor(
226  ,FactorizationNonzeros *fact_nonzeros
227  ,std::ostream *out
228  );
229 
231 
232 private:
233 
234  // /////////////////////////////////
235  // Private types
236 
238  enum E_IFlag {
245  NZ_LE_ZERO = -10,
257  };
258 
259  // /////////////////////////////////
260  // Private data members
261 
265 
266  // ////////////////////////////////
267  // Private member functions
268 
269  // Set MA28 control parameters
271 
272  // Print MA28 return parameters
273  void print_ma28_outputs(
274  bool ma28ad_bd
276  ,const FactorizationStructureMA28 &fs
277  ,const value_type w[]
278  ,std::ostream *out
279  );
280 
281  // Throw an exception for an iflag error
283 
284 }; // end class DirectSparseSolverMA28
285 
286 // ////////////////////////////////////////
287 // Inline members
288 
289 } // end namespace AbstractLinAlgPack
290 
291 #endif // ALAP_DIRECT_SPARSE_SOLVER_MA28_H
DirectSparseSolverMA28(value_type estimated_fillin_ratio=10.0, value_type u=0.1, bool grow=false, value_type tol=0.0, index_type nsrch=4, bool lbig=false, bool print_ma28_outputs=false, const std::string &output_file_name="")
Default constructor.
Teuchos::RCP< const Teuchos::AbstractFactory< BasisMatrix > > basis_matrix_factory_ptr_t
Abstract interface for immutable, finite dimensional, coordinate vectors {abstract}.
const Teuchos::RCP< FactorizationNonzeros > create_fact_nonzeros() const
Concreate sparse solver subclass that uses MA28.
void print_ma28_outputs(bool ma28ad_bd, index_type iflag, const FactorizationStructureMA28 &fs, const value_type w[], std::ostream *out)
Teuchos::RCP< BasisMatrixImp > create_matrix() const
const Teuchos::RCP< FactorizationStructure > create_fact_struc() const
void V_InvMtV(VectorMutable *v_lhs, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2) const
const basis_matrix_factory_ptr_t basis_matrix_factory() const
const f_int f_dbl_prec const f_int f_int const f_int f_int const f_dbl_prec f_int f_int f_dbl_prec w[]
void imp_factor(const AbstractLinAlgPack::MatrixConvertToSparse &A, const FactorizationStructure &fact_struc, FactorizationNonzeros *fact_nonzeros, std::ostream *out)
Mix-in interface for extracing explicit elements from a sparse matrix in one of several Fortran compa...
const f_int f_dbl_prec const f_int f_int const f_int f_int const f_dbl_prec f_int f_int f_dbl_prec f_int * iflag
std::ostream * out
Implementation node class for DirectSparseSolver that takes care of the memory management details...
Abstract class for objects that represent the factorization structure of a particular matrix...
void set_ma28_parameters(FactorizationStructureMA28 *fs)
Abstract interface for mutable coordinate vectors {abstract}.
void estimated_fillin_ratio(value_type estimated_fillin_ratio)
Abstract class for objects that represent the factorization nonzeros of a particular matrix...
void ThrowIFlagException(index_type iflag)
STANDARD_MEMBER_COMPOSITION_MEMBERS(value_type, u)
Pivot tolerance versus sparsity.
Transp
TRANS.
const f_int f_dbl_prec const f_int f_int const f_int f_int const f_dbl_prec & u
MA28 Basic Encapsulation Class.
void imp_analyze_and_factor(const AbstractLinAlgPack::MatrixConvertToSparse &A, FactorizationStructure *fact_struc, FactorizationNonzeros *fact_nonzeros, DenseLinAlgPack::IVector *row_perm, DenseLinAlgPack::IVector *col_perm, size_type *rank, std::ostream *out)