MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AbstractLinAlgPack_MatrixOpSerial.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 SPARSE_LINALG_PACK_MATRIX_WITH_OP_SERIAL_H
43 #define SPARSE_LINALG_PACK_MATRIX_WITH_OP_SERIAL_H
44 
45 #include <iosfwd>
46 
50 
51 namespace AbstractLinAlgPack {
52 
78  : virtual public AbstractLinAlgPack::MatrixOp // doxygen needs full name
79 {
80 public:
81 
83  using MatrixOp::Mp_StM;
85  using MatrixOp::Mp_StMtP;
87  using MatrixOp::Mp_StPtM;
91  using MatrixOp::Vp_StMtV;
93  using MatrixOp::Mp_StMtM;
94 
97 
99  virtual void Mp_StM(DMatrixSlice* gms_lhs, value_type alpha
100  , BLAS_Cpp::Transp trans_rhs) const;
101 
103  virtual void Mp_StMtP(DMatrixSlice* gms_lhs, value_type alpha
104  , BLAS_Cpp::Transp M_trans
105  , const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
106  ) const;
107 
109  virtual void Mp_StPtM(DMatrixSlice* gms_lhs, value_type alpha
110  , const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
111  , BLAS_Cpp::Transp M_trans
112  ) const;
113 
115  virtual void Mp_StPtMtP(DMatrixSlice* gms_lhs, value_type alpha
116  , const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
117  , BLAS_Cpp::Transp M_trans
118  , const GenPermMatrixSlice& P_rhs2, BLAS_Cpp::Transp P_rhs2_trans
119  ) const;
120 
122 
125 
127  virtual void Vp_StMtV(DVectorSlice* vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1
128  , const DVectorSlice& vs_rhs2, value_type beta) const = 0;
129 
131  virtual void Vp_StMtV(DVectorSlice* vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1
132  , const SpVectorSlice& sv_rhs2, value_type beta) const;
133 
135  virtual void Vp_StPtMtV(DVectorSlice* vs_lhs, value_type alpha
136  , const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
137  , BLAS_Cpp::Transp M_rhs2_trans
138  , const DVectorSlice& vs_rhs3, value_type beta) const;
139 
141  virtual void Vp_StPtMtV(DVectorSlice* vs_lhs, value_type alpha
142  , const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
143  , BLAS_Cpp::Transp M_rhs2_trans
144  , const SpVectorSlice& sv_rhs3, value_type beta) const;
145 
147  virtual value_type transVtMtV(const DVectorSlice& vs_rhs1, BLAS_Cpp::Transp trans_rhs2
148  , const DVectorSlice& vs_rhs3) const;
149 
161  virtual void syr2k(
162  BLAS_Cpp::Transp M_trans, value_type alpha
163  , const GenPermMatrixSlice& P1, BLAS_Cpp::Transp P1_trans
164  , const GenPermMatrixSlice& P2, BLAS_Cpp::Transp P2_trans
165  , value_type beta, DMatrixSliceSym* sym_lhs ) const;
166 
168 
171 
173  virtual void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha
174  , BLAS_Cpp::Transp trans_rhs1, const DMatrixSlice& gms_rhs2
175  , BLAS_Cpp::Transp trans_rhs2, value_type beta) const;
176 
178  virtual void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha, const DMatrixSlice& gms_rhs1
179  , BLAS_Cpp::Transp trans_rhs1, BLAS_Cpp::Transp trans_rhs2, value_type beta) const;
180 
182  virtual void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha
183  , BLAS_Cpp::Transp trans_rhs1, const MatrixOpSerial& mwo_rhs2
184  , BLAS_Cpp::Transp trans_rhs2, value_type beta) const;
185 
187  virtual void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha
188  , BLAS_Cpp::Transp trans_rhs1, const DMatrixSliceSym& sym_rhs2
189  , BLAS_Cpp::Transp trans_rhs2, value_type beta) const;
190 
192  virtual void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha, const DMatrixSliceSym& sym_rhs1
193  , BLAS_Cpp::Transp trans_rhs1, BLAS_Cpp::Transp trans_rhs2, value_type beta) const;
194 
196  virtual void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha
197  , BLAS_Cpp::Transp trans_rhs1, const DMatrixSliceTri& tri_rhs2
198  , BLAS_Cpp::Transp trans_rhs2, value_type beta) const;
199 
201  virtual void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha, const DMatrixSliceTri& tri_rhs1
202  , BLAS_Cpp::Transp trans_rhs1, BLAS_Cpp::Transp trans_rhs2, value_type beta) const;
203 
204 
213  virtual void syrk(
214  BLAS_Cpp::Transp M_trans, value_type alpha
215  , value_type beta, DMatrixSliceSym* sym_lhs ) const;
216 
218 
221 
223  const VectorSpace& space_cols() const;
225  const VectorSpace& space_rows() const;
227  std::ostream& output(std::ostream& out) const;
229  bool Mp_StM(
230  MatrixOp* mwo_lhs, value_type alpha
231  ,BLAS_Cpp::Transp trans_rhs
232  ) const;
234  bool Mp_StMtP(
235  MatrixOp* mwo_lhs, value_type alpha
236  , BLAS_Cpp::Transp M_trans
237  , const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
238  ) const;
240  bool Mp_StPtM(
241  MatrixOp* mwo_lhs, value_type alpha
242  , const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
243  , BLAS_Cpp::Transp M_trans
244  ) const;
246  bool Mp_StPtMtP(
247  MatrixOp* mwo_lhs, value_type alpha
248  ,const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
249  ,BLAS_Cpp::Transp M_trans
250  ,const GenPermMatrixSlice& P_rhs2, BLAS_Cpp::Transp P_rhs2_trans
251  ) const;
253  void Vp_StMtV(
254  VectorMutable* v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1
255  , const Vector& v_rhs2, value_type beta) const;
257  void Vp_StMtV(
258  VectorMutable* v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1
259  , const SpVectorSlice& sv_rhs2, value_type beta) const;
261  void Vp_StPtMtV(
262  VectorMutable* v_lhs, value_type alpha
263  , const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
264  , BLAS_Cpp::Transp M_rhs2_trans
265  , const Vector& v_rhs3, value_type beta) const;
267  void Vp_StPtMtV(
268  VectorMutable* v_lhs, value_type alpha
269  , const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
270  , BLAS_Cpp::Transp M_rhs2_trans
271  , const SpVectorSlice& sv_rhs3, value_type beta) const;
274  const Vector& v_rhs1, BLAS_Cpp::Transp trans_rhs2
275  , const Vector& v_rhs3) const;
278  const SpVectorSlice& sv_rhs1, BLAS_Cpp::Transp trans_rhs2
279  , const SpVectorSlice& sv_rhs3) const;
281  void syr2k(
282  BLAS_Cpp::Transp M_trans, value_type alpha
283  , const GenPermMatrixSlice& P1, BLAS_Cpp::Transp P1_trans
284  , const GenPermMatrixSlice& P2, BLAS_Cpp::Transp P2_trans
285  , value_type beta, MatrixSymOp* symwo_lhs ) const;
287  bool Mp_StMtM(
288  MatrixOp* mwo_lhs, value_type alpha
289  ,BLAS_Cpp::Transp trans_rhs1
290  ,const MatrixOp& mwo_rhs2, BLAS_Cpp::Transp trans_rhs2
291  ,value_type beta ) const;
293  bool syrk(
294  BLAS_Cpp::Transp M_trans, value_type alpha
295  ,value_type beta, MatrixSymOp* sym_lhs ) const;
296 
298 
299 private:
300 
301  // ////////////////////////////////////
302  // Private data members
303 
306 
307  // ////////////////////////////////////
308  // Private member functions
309 
310 }; // end class MatrixOpSerial
311 
318 
321 
323 inline void Mp_StM(DMatrixSlice* gms_lhs, value_type alpha, const MatrixOpSerial& M_rhs
324  , BLAS_Cpp::Transp trans_rhs)
325 {
326  M_rhs.Mp_StM(gms_lhs,alpha,trans_rhs);
327 }
328 
330 inline void Mp_StMtP(DMatrixSlice* gms_lhs, value_type alpha
331  , const MatrixOpSerial& M_rhs, BLAS_Cpp::Transp M_trans
332  , const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
333  )
334 {
335  M_rhs.Mp_StMtP(gms_lhs,alpha,M_trans,P_rhs,P_rhs_trans);
336 }
337 
339 inline void Mp_StPtM(DMatrixSlice* gms_lhs, value_type alpha
340  , const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
341  , const MatrixOpSerial& M_rhs, BLAS_Cpp::Transp M_trans
342  )
343 {
344  M_rhs.Mp_StPtM(gms_lhs,alpha,P_rhs,P_rhs_trans,M_trans);
345 }
346 
348 inline void Mp_StPtMtP(DMatrixSlice* gms_lhs, value_type alpha
349  , const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
350  , const MatrixOpSerial& M_rhs, BLAS_Cpp::Transp trans_rhs
351  , const GenPermMatrixSlice& P_rhs2, BLAS_Cpp::Transp P_rhs2_trans
352  )
353 {
354  M_rhs.Mp_StPtMtP(gms_lhs,alpha,P_rhs1,P_rhs1_trans,trans_rhs,P_rhs2,P_rhs2_trans);
355 }
356 
358 
361 
363 inline void Vp_StMtV(DVectorSlice* vs_lhs, value_type alpha, const MatrixOpSerial& M_rhs1
364  , BLAS_Cpp::Transp trans_rhs1, const DVectorSlice& vs_rhs2, value_type beta = 1.0)
365 {
366  M_rhs1.Vp_StMtV(vs_lhs,alpha,trans_rhs1,vs_rhs2,beta);
367 }
368 
370 inline void Vp_StMtV(DVectorSlice* vs_lhs, value_type alpha, const MatrixOpSerial& M_rhs1
371  , BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice& sv_rhs2, value_type beta = 1.0)
372 {
373  M_rhs1.Vp_StMtV(vs_lhs,alpha,trans_rhs1,sv_rhs2,beta);
374 }
375 
377 inline void Vp_StPtMtV(DVectorSlice* vs_lhs, value_type alpha
378  , const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
379  , const MatrixOpSerial& M_rhs2, BLAS_Cpp::Transp M_rhs2_trans
380  , const DVectorSlice& vs_rhs3, value_type beta = 1.0)
381 {
382  M_rhs2.Vp_StPtMtV(vs_lhs,alpha,P_rhs1,P_rhs1_trans,M_rhs2_trans,vs_rhs3,beta);
383 }
384 
386 inline void Vp_StPtMtV(DVectorSlice* vs_lhs, value_type alpha
387  , const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
388  , const MatrixOpSerial& M_rhs2, BLAS_Cpp::Transp M_rhs2_trans
389  , const SpVectorSlice& sv_rhs3, value_type beta = 1.0)
390 {
391  M_rhs2.Vp_StPtMtV(vs_lhs,alpha,P_rhs1,P_rhs1_trans,M_rhs2_trans,sv_rhs3,beta);
392 }
393 
395 inline value_type transVtMtV(const DVectorSlice& vs_rhs1, const MatrixOpSerial& M_rhs2
396  , BLAS_Cpp::Transp trans_rhs2, const DVectorSlice& vs_rhs3)
397 {
398  return M_rhs2.transVtMtV(vs_rhs1,trans_rhs2,vs_rhs3);
399 }
400 
402 inline value_type transVtMtV(const SpVectorSlice& sv_rhs1, const MatrixOpSerial& M_rhs2
403  , BLAS_Cpp::Transp trans_rhs2, const SpVectorSlice& sv_rhs3)
404 {
405  return M_rhs2.transVtMtV(sv_rhs1,trans_rhs2,sv_rhs3);
406 }
407 
409 
412 
414 inline void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha, const MatrixOpSerial& M_rhs1
415  , BLAS_Cpp::Transp trans_rhs1, const DMatrixSlice& gms_rhs2
416  , BLAS_Cpp::Transp trans_rhs2, value_type beta = 1.0)
417 {
418  M_rhs1.Mp_StMtM(gms_lhs,alpha,trans_rhs1,gms_rhs2,trans_rhs2,beta);
419 }
420 
422 inline void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha, const DMatrixSlice& gms_rhs1
423  , BLAS_Cpp::Transp trans_rhs1, const MatrixOpSerial& M_rhs2, BLAS_Cpp::Transp trans_rhs2
424  , value_type beta = 1.0)
425 {
426  M_rhs2.Mp_StMtM(gms_lhs,alpha,gms_rhs1,trans_rhs1,trans_rhs2,beta);
427 }
428 
430 inline void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha, const MatrixOpSerial& mwo_rhs1
431  , BLAS_Cpp::Transp trans_rhs1, const MatrixOpSerial& gms_rhs2
432  , BLAS_Cpp::Transp trans_rhs2, value_type beta = 1.0)
433 {
434  mwo_rhs1.Mp_StMtM(gms_lhs,alpha,trans_rhs1,gms_rhs2,trans_rhs2,beta);
435 }
436 
438 inline void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha, const MatrixOpSerial& M_rhs1
439  , BLAS_Cpp::Transp trans_rhs1, const DMatrixSliceSym& sym_rhs2
440  , BLAS_Cpp::Transp trans_rhs2, value_type beta = 1.0)
441 {
442  M_rhs1.Mp_StMtM(gms_lhs,alpha,trans_rhs1,sym_rhs2,trans_rhs2,beta);
443 }
444 
446 inline void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha, const DMatrixSliceSym& sym_rhs1
447  , BLAS_Cpp::Transp trans_rhs1, const MatrixOpSerial& M_rhs2, BLAS_Cpp::Transp trans_rhs2
448  , value_type beta = 1.0)
449 {
450  M_rhs2.Mp_StMtM(gms_lhs,alpha,sym_rhs1,trans_rhs1,trans_rhs2,beta);
451 }
452 
454 inline void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha, const MatrixOpSerial& M_rhs1
455  , BLAS_Cpp::Transp trans_rhs1, const DMatrixSliceTri& tri_rhs2, BLAS_Cpp::Transp trans_rhs2
456  , value_type beta = 1.0)
457 {
458  M_rhs1.Mp_StMtM(gms_lhs,alpha,trans_rhs1,tri_rhs2,trans_rhs2,beta);
459 }
460 
462 inline void Mp_StMtM(DMatrixSlice* gms_lhs, value_type alpha, const DMatrixSliceTri& tri_rhs1
463  , BLAS_Cpp::Transp trans_rhs1, const MatrixOpSerial& M_rhs2, BLAS_Cpp::Transp trans_rhs2
464  , value_type beta = 1.0)
465 {
466  M_rhs2.Mp_StMtM(gms_lhs,alpha,tri_rhs1,trans_rhs1,trans_rhs2,beta);
467 }
468 
470 
472 
473 } // end namespace AbstractLinAlgPack
474 
475 #endif // SPARSE_LINALG_PACK_MATRIX_WITH_OP_SERIAL_H
Abstract interface for immutable, finite dimensional, coordinate vectors {abstract}.
std::ostream & output(std::ostream &out) const
virtual value_type transVtMtV(const DVectorSlice &vs_rhs1, BLAS_Cpp::Transp trans_rhs2, const DVectorSlice &vs_rhs3) const
result = vs_rhs1' * op(M_rhs2) * vs_rhs3
virtual void Vp_StPtMtV(DVectorSlice *vs_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_rhs2_trans, const DVectorSlice &vs_rhs3, value_type beta) const
vs_lhs = alpha * op(P_rhs1) * op(M_rhs2) * vs_rhs3 + beta * vs_rhs
void Mp_StMtM(MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, const MatrixOp &mwo_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta=1.0)
mwo_lhs = alpha * op(mwo_rhs1) * op(mwo_rhs2) + beta * mwo_lhs (right) (xGEMM).
void Vp_StPtMtV(VectorMutable *v_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, const MatrixOp &M_rhs2, BLAS_Cpp::Transp M_rhs2_trans, const Vector &v_rhs3, value_type beta=1.0)
v_lhs = alpha * op(P_rhs1) * op(M_rhs2) * v_rhs3 + beta * v_rhs
value_type transVtMtV(const Vector &v_rhs1, const MatrixOp &M_rhs2, BLAS_Cpp::Transp trans_rhs2, const Vector &v_rhs3)
result = v_rhs1' * op(M_rhs2) * v_rhs3
friend void Mp_StPtMtP(MatrixOp *mwo_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs, const GenPermMatrixSlice &P_rhs2, BLAS_Cpp::Transp P_rhs2_trans)
Base class for all matrices implemented in a shared memory address space.
Interface adding operations specific for a symmetric matrix {abstract}.
void Mp_StPtMtP(MatrixOp *mwo_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs, const GenPermMatrixSlice &P_rhs2, BLAS_Cpp::Transp P_rhs2_trans)
mwo_lhs += alpha * op(P_rhs1) * op(M_rhs) * op(P_rhs2).
Abstract interface for objects that represent a space for mutable coordinate vectors.
virtual void syrk(BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta, DMatrixSliceSym *sym_lhs) const
Perform a rank-k update of a dense symmetric matrix of the form:
friend void Mp_StMtM(MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, const MatrixOp &mwo_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta)
void Mp_StMtP(MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans)
mwo_lhs += alpha * op(M_rhs) * op(P_rhs).
std::ostream * out
void Mp_StPtM(MatrixOp *mwo_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans, const MatrixOp &M_rhs, BLAS_Cpp::Transp M_trans)
mwo_lhs += alpha * op(P) * op(M_rhs).
virtual void Mp_StMtM(DMatrixSlice *gms_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const DMatrixSlice &gms_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta) const
gms_lhs = alpha * op(M_rhs1) * op(gms_rhs2) + beta * gms_lhs (right) (xGEMM)
virtual void Mp_StPtMtP(DMatrixSlice *gms_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs2, BLAS_Cpp::Transp P_rhs2_trans) const
gms_lhs += alpha * op(P_rhs1) * op(M_rhs) * op(P_rhs2)
void Mp_StM(MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs)
void Vp_StMtV(VectorMutable *v_lhs, value_type alpha, const MatrixOp &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta=1.0)
v_lhs = alpha * op(M_rhs1) * v_rhs2 + beta * v_lhs (BLAS xGEMV)
friend void Mp_StM(MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs)
Base class for all matrices that support basic matrix operations.
virtual void Mp_StMtP(DMatrixSlice *gms_lhs, value_type alpha, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans) const
gms_lhs += alpha * op(M_rhs) * op(P_rhs)
virtual void Vp_StMtV(DVectorSlice *vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const DVectorSlice &vs_rhs2, value_type beta) const =0
vs_lhs = alpha * op(M_rhs1) * vs_rhs2 + beta * vs_lhs (BLAS xGEMV)
friend void Mp_StMtP(MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans)
virtual void Mp_StM(DMatrixSlice *gms_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs) const
gms_lhs += alpha * op(M_rhs) (BLAS xAXPY)
Abstract interface for mutable coordinate vectors {abstract}.
virtual void syr2k(BLAS_Cpp::Transp M_trans, value_type alpha, const GenPermMatrixSlice &P1, BLAS_Cpp::Transp P1_trans, const GenPermMatrixSlice &P2, BLAS_Cpp::Transp P2_trans, value_type beta, DMatrixSliceSym *sym_lhs) const
Perform a specialized rank-2k update of a dense symmetric matrix of the form:
Subclass for serial vector space objects that create VectorMutableDense vector and MultiVectorMutable...
friend void Vp_StMtV(VectorMutable *v_lhs, value_type alpha, const MatrixOp &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta)
virtual void Mp_StPtM(DMatrixSlice *gms_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans, BLAS_Cpp::Transp M_trans) const
gms_lhs += alpha * op(P) * op(M_rhs)
Transp
TRANS.
Concrete matrix type to represent general permutation (mapping) matrices.
friend void Mp_StPtM(MatrixOp *mwo_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans, const MatrixOp &M_rhs, BLAS_Cpp::Transp M_trans)