MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AbstractLinAlgPack_LinAlgOpPackHack.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 
43 #ifndef LIN_ALG_OP_PACK_HACK_H
44 #define LIN_ALG_OP_PACK_HACK_H
45 
47 
50 
51 namespace LinAlgOpPack {
52 
61 
63 void assign(DMatrixSlice* gms_lhs, const MatrixOp& M_rhs, BLAS_Cpp::Transp trans_rhs);
64 
67 void Mp_StM(
68  DMatrixSlice* vs_lhs, value_type alpha
69  ,const MatrixOp& mwo_rhs1, BLAS_Cpp::Transp trans_rhs1
70  );
71 
74 void Vp_StMtV(
75  DVectorSlice* vs_lhs, value_type alpha, const MatrixOp& mwo_rhs1
76  ,BLAS_Cpp::Transp trans_rhs1, const DVectorSlice& vs_rhs2
77  ,value_type beta = 1.0 );
78 
81 void Vp_MtV(
82  DVectorSlice* vs_lhs, const MatrixOp& mwo_rhs1
83  ,BLAS_Cpp::Transp trans_rhs1, const DVectorSlice& vs_rhs2
84  ,value_type beta = 1.0 );
85 
88 void Vp_StMtV(
89  DVectorSlice* vs_lhs, value_type alpha, const MatrixOp& mwo_rhs1
90  ,BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice& sv_rhs2
91  ,value_type beta = 1.0 );
92 
95 void V_MtV(
96  DVectorSlice* vs_lhs, const MatrixOp& mwo_rhs1
97  ,BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice& sv_rhs2 );
98 
101 void V_InvMtV(
102  DVectorSlice* vs_lhs, const MatrixOpNonsing& mwo_rhs1
103  ,BLAS_Cpp::Transp trans_rhs1, const DVectorSlice& vs_rhs2 );
104 
107 void V_InvMtV(
108  DVector* v_lhs, const MatrixOpNonsing& mwo_rhs1
109  ,BLAS_Cpp::Transp trans_rhs1, const DVectorSlice& vs_rhs2 );
110 
113 void V_InvMtV(
114  DVectorSlice* vs_lhs, const MatrixOpNonsing& mwo_rhs1
115  ,BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice& sv_rhs2 );
116 
119 void V_InvMtV(
120  DVector* v_lhs, const MatrixOpNonsing& mwo_rhs1
121  ,BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice& sv_rhs2 );
122 
125 void Vp_StPtMtV(
126  DVectorSlice* vs_lhs, value_type alpha
127  ,const GenPermMatrixSlice& gpms_rhs1, BLAS_Cpp::Transp trans_rhs1
128  ,const MatrixOp& mwo_rhs2, BLAS_Cpp::Transp trans_rhs2
129  ,const DVectorSlice& vs_rhs3, value_type beta = 1.0 );
130 
133 void Vp_StPtMtV(
134  DVectorSlice* vs_lhs, value_type alpha
135  ,const GenPermMatrixSlice& gpms_rhs1, BLAS_Cpp::Transp trans_rhs1
136  ,const MatrixOp& mwo_rhs2, BLAS_Cpp::Transp trans_rhs2
137  ,const SpVectorSlice& sv_rhs3, value_type beta = 1.0 );
138 
139 } // end namespace LinAlgOpPack
140 
141 // //////////////////////
142 // Inline functions
143 
144 inline
146  DVectorSlice* vs_lhs, const MatrixOp& mwo_rhs1
147  ,BLAS_Cpp::Transp trans_rhs1, const DVectorSlice& vs_rhs2
148  ,value_type beta )
149 {
150  Vp_StMtV(vs_lhs,1.0,mwo_rhs1,trans_rhs1,vs_rhs2,beta);
151 }
152 
153 inline
155  DVectorSlice* vs_lhs, const MatrixOp& mwo_rhs1
156  ,BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice& sv_rhs2 )
157 {
158  Vp_StMtV(vs_lhs,1.0,mwo_rhs1,trans_rhs1,sv_rhs2);
159 }
160 
161 
162 
163 
164 #endif // LIN_ALG_OP_PACK_HACK_H
void Vp_StMtV(DVectorSlice *vs_lhs, value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, const DVectorSlice &vs_rhs2, value_type beta=1.0)
vs_lhs = alpha * op(mwo_rhs1) * vs_rhs2 + beta * vs_lhs.
void assign(VectorMutable *v_lhs, const V &V_rhs)
v_lhs = V_rhs.
VectorSliceTmpl< value_type > DVectorSlice
void Mp_StM(DMatrixSlice *vs_lhs, value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1)
m_lhs += alpha * op(mwo_rhs1).
void V_InvMtV(DVectorSlice *vs_lhs, const MatrixOpNonsing &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, const DVectorSlice &vs_rhs2)
vs_lhs = inv(op(mwo_rhs1)) * vs_rhs2.
void Vp_MtV(VectorMutable *v_lhs, const MatrixOp &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const V &V_rhs2)
v_lhs += op(M_rhs1) * V_rhs2.
void Vp_StPtMtV(DVectorSlice *vs_lhs, value_type alpha, const GenPermMatrixSlice &gpms_rhs1, BLAS_Cpp::Transp trans_rhs1, const MatrixOp &mwo_rhs2, BLAS_Cpp::Transp trans_rhs2, const DVectorSlice &vs_rhs3, value_type beta=1.0)
vs_lhs = alpha * op(gpms_rhs1) * op(mwo_rhs2) * vs_rhs3 + beta * vs_lhs.
Base class for all matrices that support basic matrix operations.
SparseVectorSlice< SparseElement< index_type, value_type > > SpVectorSlice
DenseLinAlgPack::VectorSliceTmpl< value_type > DVectorSlice
void V_MtV(VectorMutable *v_lhs, const MatrixOp &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const V &V_rhs2)
v_lhs = op(M_rhs1) * V_rhs2.
AbstractLinAlgPack::value_type value_type
VectorTmpl< value_type > DVector
Abstract base class for all nonsingular polymorphic matrices that can be used to compute matrix-vecto...
DenseLinAlgPack::DMatrixSlice DMatrixSlice
Abstract base class for all nonsingular polymorphic matrices that can solve for linear system with bu...
Transp
TRANS.
void Vp_StMtV(DVectorSlice *vs_lhs, value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2, value_type beta=1.0)
vs_lhs = alpha * op(mwo_rhs1) * vs_rhs2 + beta * sv_lhs.
Concrete matrix type to represent general permutation (mapping) matrices.