MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AbstractLinAlgPack_MatrixOp.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 ABSTRACT_LINALG_PACK_MATRIX_WITH_OP_H
43 #define ABSTRACT_LINALG_PACK_MATRIX_WITH_OP_H
44 
45 #include <iosfwd>
46 
48 #include "Teuchos_RCP.hpp"
49 
50 namespace AbstractLinAlgPack {
51 
129 class MatrixOp : public virtual MatrixBase {
130 public:
131 
134 
136  friend
137  void Mt_S( MatrixOp* mwo_lhs, value_type alpha );
139  friend
140  void Mp_StM(
141  MatrixOp* mwo_lhs, value_type alpha, const MatrixOp& M_rhs
142  , BLAS_Cpp::Transp trans_rhs);
144  friend
145  void Mp_StMtP(
146  MatrixOp* mwo_lhs, value_type alpha
147  ,const MatrixOp& M_rhs, BLAS_Cpp::Transp M_trans
148  ,const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
149  );
151  friend
152  void Mp_StPtM(
153  MatrixOp* mwo_lhs, value_type alpha
154  ,const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
155  ,const MatrixOp& M_rhs, BLAS_Cpp::Transp M_trans
156  );
158  friend
159  void Mp_StPtMtP(
160  MatrixOp* mwo_lhs, value_type alpha
161  ,const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
162  ,const MatrixOp& M_rhs, BLAS_Cpp::Transp trans_rhs
163  ,const GenPermMatrixSlice& P_rhs2, BLAS_Cpp::Transp P_rhs2_trans
164  );
166  friend
167  void Vp_StMtV(
168  VectorMutable* v_lhs, value_type alpha, const MatrixOp& M_rhs1
169  ,BLAS_Cpp::Transp trans_rhs1, const Vector& v_rhs2, value_type beta
170  );
172  friend
173  void Vp_StMtV(
174  VectorMutable* v_lhs, value_type alpha, const MatrixOp& M_rhs1
175  ,BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice& sv_rhs2, value_type beta
176  );
178  friend
179  void Vp_StPtMtV(
180  VectorMutable* v_lhs, value_type alpha
181  ,const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
182  ,const MatrixOp& M_rhs2, BLAS_Cpp::Transp M_rhs2_trans
183  ,const Vector& v_rhs3, value_type beta
184  );
186  friend
187  void Vp_StPtMtV(
188  VectorMutable* v_lhs, value_type alpha
189  ,const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
190  ,const MatrixOp& M_rhs2, BLAS_Cpp::Transp M_rhs2_trans
191  ,const SpVectorSlice& sv_rhs3, value_type beta
192  );
194  friend
196  const Vector& v_rhs1, const MatrixOp& M_rhs2
197  ,BLAS_Cpp::Transp trans_rhs2, const Vector& v_rhs3
198  );
200  friend
202  const SpVectorSlice& sv_rhs1, const MatrixOp& M_rhs2
203  ,BLAS_Cpp::Transp trans_rhs2, const SpVectorSlice& sv_rhs3
204  );
206  friend
207  void syr2k(
208  const MatrixOp& M, BLAS_Cpp::Transp M_trans, value_type alpha
209  ,const GenPermMatrixSlice& P1, BLAS_Cpp::Transp P1_trans
210  ,const GenPermMatrixSlice& P2, BLAS_Cpp::Transp P2_trans
211  ,value_type beta, MatrixSymOp* symwo_lhs
212  );
214  friend
215  void Mp_StMtM(
216  MatrixOp* mwo_lhs, value_type alpha
217  ,const MatrixOp& mwo_rhs1, BLAS_Cpp::Transp trans_rhs1
218  ,const MatrixOp& mwo_rhs2, BLAS_Cpp::Transp trans_rhs2
219  ,value_type beta
220  );
222  friend
223  void syrk(
224  const MatrixOp &mwo_rhs
225  ,BLAS_Cpp::Transp M_trans
226  ,value_type alpha
227  ,value_type beta
228  ,MatrixSymOp *sym_lhs
229  );
230 
232 
235 
236 #ifndef DOXYGEN_COMPILE
237 
238  typedef Teuchos::RCP<const MatrixOp> mat_ptr_t;
240  typedef Teuchos::RCP<MatrixOp> mat_mut_ptr_t;
241 #endif
242 
249  };
250 
252  struct MatNorm {
253  MatNorm(value_type _value, EMatNormType _type) : value(_value), type(_type) {}
256  };
257 
259  class MethodNotImplemented : public std::runtime_error
260  {public: MethodNotImplemented(const std::string& what_arg) : std::runtime_error(what_arg) {}};
261 
263  class IncompatibleMatrices : public std::logic_error
264  {public: IncompatibleMatrices(const std::string& what_arg) : std::logic_error(what_arg) {}};
265 
267 
270 
280  virtual void zero_out();
281 
290  virtual void Mt_S( value_type alpha );
291 
299  virtual MatrixOp& operator=(const MatrixOp& mwo_rhs);
300 
302 
305 
318  virtual mat_mut_ptr_t clone();
319 
325  virtual mat_ptr_t clone() const;
326 
328 
331 
338  virtual std::ostream& output(std::ostream& out) const;
339 
341 
344 
375  const MatNorm calc_norm(
376  EMatNormType requested_norm_type = MAT_NORM_1
377  ,bool allow_replacement = false
378  ) const;
379 
381 
384 
400  virtual mat_ptr_t sub_view(const Range1D& row_rng, const Range1D& col_rng) const;
401 
404  mat_ptr_t sub_view(
405  const index_type& rl, const index_type& ru
406  ,const index_type& cl, const index_type& cu
407  ) const;
408 
410 
413 
454  virtual mat_ptr_t perm_view(
455  const Permutation *P_row
456  ,const index_type row_part[]
457  ,int num_row_part
458  ,const Permutation *P_col
459  ,const index_type col_part[]
460  ,int num_col_part
461  ) const;
462 
492  virtual mat_ptr_t perm_view_update(
493  const Permutation *P_row
494  ,const index_type row_part[]
495  ,int num_row_part
496  ,const Permutation *P_col
497  ,const index_type col_part[]
498  ,int num_col_part
499  ,const mat_ptr_t &perm_view
500  ) const;
501 
503 
504 #ifdef TEMPLATE_FRIENDS_NOT_SUPPORTED
505 public:
506 #else
507 protected:
508 #endif
509 
512 
520  virtual bool Mp_StM(
521  MatrixOp* mwo_lhs, value_type alpha
522  ,BLAS_Cpp::Transp trans_rhs
523  ) const;
524 
532  virtual bool Mp_StM(
533  value_type alpha,const MatrixOp& M_rhs, BLAS_Cpp::Transp trans_rhs
534  );
535 
543  virtual bool Mp_StMtP(
544  MatrixOp* mwo_lhs, value_type alpha
545  ,BLAS_Cpp::Transp M_trans
546  ,const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
547  ) const;
548 
556  virtual bool Mp_StMtP(
557  value_type alpha
558  ,const MatrixOp& mwo_rhs, BLAS_Cpp::Transp M_trans
559  ,const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
560  );
561 
569  virtual bool Mp_StPtM(
570  MatrixOp* mwo_lhs, value_type alpha
571  ,const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
572  ,BLAS_Cpp::Transp M_trans
573  ) const;
574 
582  virtual bool Mp_StPtM(
583  value_type alpha
584  ,const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
585  ,const MatrixOp& mwo_rhs, BLAS_Cpp::Transp M_trans
586  );
587 
595  virtual bool Mp_StPtMtP(
596  MatrixOp* mwo_lhs, value_type alpha
597  ,const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
598  ,BLAS_Cpp::Transp M_trans
599  ,const GenPermMatrixSlice& P_rhs2, BLAS_Cpp::Transp P_rhs2_trans
600  ) const;
601 
609  virtual bool Mp_StPtMtP(
610  value_type alpha
611  ,const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
612  ,const MatrixOp& mwo_rhs, BLAS_Cpp::Transp M_trans
613  ,const GenPermMatrixSlice& P_rhs2, BLAS_Cpp::Transp P_rhs2_trans
614  );
615 
616  // end Level-1 BLAS
618 
621 
623  virtual void Vp_StMtV(
624  VectorMutable* v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1
625  ,const Vector& v_rhs2, value_type beta
626  ) const = 0;
627 
629  virtual void Vp_StMtV(
630  VectorMutable* v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1
631  ,const SpVectorSlice& sv_rhs2, value_type beta
632  ) const;
633 
635  virtual void Vp_StPtMtV(
636  VectorMutable* v_lhs, value_type alpha
637  ,const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
638  ,BLAS_Cpp::Transp M_rhs2_trans
639  ,const Vector& v_rhs3, value_type beta
640  ) const;
641 
643  virtual void Vp_StPtMtV(
644  VectorMutable* v_lhs, value_type alpha
645  ,const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
646  ,BLAS_Cpp::Transp M_rhs2_trans
647  ,const SpVectorSlice& sv_rhs3, value_type beta
648  ) const;
649 
651  virtual value_type transVtMtV(
652  const Vector& v_rhs1, BLAS_Cpp::Transp trans_rhs2
653  ,const Vector& v_rhs3
654  ) const;
655 
657  virtual value_type transVtMtV(
658  const SpVectorSlice& sv_rhs1, BLAS_Cpp::Transp trans_rhs2
659  ,const SpVectorSlice& sv_rhs3
660  ) const;
661 
673  virtual void syr2k(
674  BLAS_Cpp::Transp M_trans, value_type alpha
675  ,const GenPermMatrixSlice& P1, BLAS_Cpp::Transp P1_trans
676  ,const GenPermMatrixSlice& P2, BLAS_Cpp::Transp P2_trans
677  ,value_type beta, MatrixSymOp* symwo_lhs
678  ) const;
679 
681 
684 
693  virtual bool Mp_StMtM(
694  MatrixOp* mwo_lhs, value_type alpha
695  ,BLAS_Cpp::Transp trans_rhs1
696  ,const MatrixOp& mwo_rhs2, BLAS_Cpp::Transp trans_rhs2
697  ,value_type beta
698  ) const;
699 
708  virtual bool Mp_StMtM(
709  MatrixOp* mwo_lhs, value_type alpha
710  ,const MatrixOp& mwo_rhs1, BLAS_Cpp::Transp trans_rhs1
711  ,BLAS_Cpp::Transp trans_rhs2
712  ,value_type beta
713  ) const;
714 
723  virtual bool Mp_StMtM(
724  value_type alpha
725  ,const MatrixOp& mwo_rhs1, BLAS_Cpp::Transp trans_rhs1
726  ,const MatrixOp& mwo_rhs2,BLAS_Cpp::Transp trans_rhs2
727  ,value_type beta
728  );
729 
741  virtual bool syrk(
742  BLAS_Cpp::Transp M_trans
743  ,value_type alpha
744  ,value_type beta
745  ,MatrixSymOp *sym_lhs
746  ) const;
747 
759  virtual bool syrk(
760  const MatrixOp &mwo_rhs
761  ,BLAS_Cpp::Transp M_trans
762  ,value_type alpha
763  ,value_type beta
764  );
765 
766  // end Level-3 BLAS
768 
769 }; // end class MatrixOp
770 
771 // ////////////////////////////////////////////////////////////////////////////////////////////////
781 
784 
785 //
792 void Mt_S( MatrixOp* mwo_lhs, value_type alpha );
793 
794 //
809 void Mp_StM(
810  MatrixOp* mwo_lhs, value_type alpha, const MatrixOp& M_rhs
811  , BLAS_Cpp::Transp trans_rhs);
812 
819 void Mp_StMtP(
820  MatrixOp* mwo_lhs, value_type alpha
821  ,const MatrixOp& M_rhs, BLAS_Cpp::Transp M_trans
822  ,const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
823  );
824 
831 void Mp_StPtM(
832  MatrixOp* mwo_lhs, value_type alpha
833  ,const GenPermMatrixSlice& P_rhs, BLAS_Cpp::Transp P_rhs_trans
834  ,const MatrixOp& M_rhs, BLAS_Cpp::Transp M_trans
835  );
836 
843 void Mp_StPtMtP(
844  MatrixOp* mwo_lhs, value_type alpha
845  ,const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
846  ,const MatrixOp& M_rhs, BLAS_Cpp::Transp trans_rhs
847  ,const GenPermMatrixSlice& P_rhs2, BLAS_Cpp::Transp P_rhs2_trans
848  );
849 
850 // end Level-1 BLAS
852 
855 
857 inline void Vp_StMtV(
858  VectorMutable* v_lhs, value_type alpha, const MatrixOp& M_rhs1
859  ,BLAS_Cpp::Transp trans_rhs1, const Vector& v_rhs2, value_type beta = 1.0
860  )
861 {
862  M_rhs1.Vp_StMtV(v_lhs,alpha,trans_rhs1,v_rhs2,beta);
863 }
864 
866 inline void Vp_StMtV(
867  VectorMutable* v_lhs, value_type alpha, const MatrixOp& M_rhs1
868  ,BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice& sv_rhs2, value_type beta = 1.0
869  )
870 {
871  M_rhs1.Vp_StMtV(v_lhs,alpha,trans_rhs1,sv_rhs2,beta);
872 }
873 
875 inline void Vp_StPtMtV(
876  VectorMutable* v_lhs, value_type alpha
877  ,const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
878  ,const MatrixOp& M_rhs2, BLAS_Cpp::Transp M_rhs2_trans
879  ,const Vector& v_rhs3, value_type beta = 1.0
880  )
881 {
882  M_rhs2.Vp_StPtMtV(v_lhs,alpha,P_rhs1,P_rhs1_trans,M_rhs2_trans,v_rhs3,beta);
883 }
884 
886 inline void Vp_StPtMtV(
887  VectorMutable* v_lhs, value_type alpha
888  ,const GenPermMatrixSlice& P_rhs1, BLAS_Cpp::Transp P_rhs1_trans
889  ,const MatrixOp& M_rhs2, BLAS_Cpp::Transp M_rhs2_trans
890  ,const SpVectorSlice& sv_rhs3, value_type beta = 1.0
891  )
892 {
893  M_rhs2.Vp_StPtMtV(v_lhs,alpha,P_rhs1,P_rhs1_trans,M_rhs2_trans,sv_rhs3,beta);
894 }
895 
898  const Vector& v_rhs1, const MatrixOp& M_rhs2
899  ,BLAS_Cpp::Transp trans_rhs2, const Vector& v_rhs3
900  )
901 {
902  return M_rhs2.transVtMtV(v_rhs1,trans_rhs2,v_rhs3);
903 }
904 
907  const SpVectorSlice& sv_rhs1, const MatrixOp& M_rhs2
908  ,BLAS_Cpp::Transp trans_rhs2, const SpVectorSlice& sv_rhs3
909  )
910 {
911  return M_rhs2.transVtMtV(sv_rhs1,trans_rhs2,sv_rhs3);
912 }
913 
915 inline void syr2k(
916  const MatrixOp& M, BLAS_Cpp::Transp M_trans, value_type alpha
917  ,const GenPermMatrixSlice& P1, BLAS_Cpp::Transp P1_trans
918  ,const GenPermMatrixSlice& P2, BLAS_Cpp::Transp P2_trans
919  ,value_type beta, MatrixSymOp* symwo_lhs
920  )
921 {
922  M.syr2k(M_trans,alpha,P1,P1_trans,P2,P2_trans,beta,symwo_lhs);
923 }
924 
925 // end Level-2 BLAS
927 
930 
943 void Mp_StMtM(
944  MatrixOp* mwo_lhs, value_type alpha
945  ,const MatrixOp& mwo_rhs1, BLAS_Cpp::Transp trans_rhs1
946  ,const MatrixOp& mwo_rhs2, BLAS_Cpp::Transp trans_rhs2
947  ,value_type beta = 1.0
948  );
949 
956 void syrk(
957  const MatrixOp &mwo_rhs
958  ,BLAS_Cpp::Transp M_trans
959  ,value_type alpha
960  ,value_type beta
961  ,MatrixSymOp *sym_lhs
962  );
963 
964 // end Level-3 BLAS
966 
967 // end non-member functions
969 
970 // //////////////////////////////////////////////////
971 // Inlined methods for MatrixOp
972 
973 inline
974 MatrixOp::mat_ptr_t
976  const index_type& rl, const index_type& ru
977  ,const index_type& cl, const index_type& cu
978  ) const
979 {
980  return this->sub_view(Range1D(rl,ru),Range1D(cl,cu));
981 }
982 
983 } // end namespace AbstractLinAlgPack
984 
985 #endif // ABSTRACT_LINALG_PACK_MATRIX_WITH_OP_H
virtual mat_ptr_t sub_view(const Range1D &row_rng, const Range1D &col_rng) const
Create a transient constant sub-matrix view of this matrix (if supported).
Base class for all polymorphic matrices.
friend void syr2k(const MatrixOp &M, 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, MatrixSymOp *symwo_lhs)
Abstract interface for immutable, finite dimensional, coordinate vectors {abstract}.
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).
virtual void zero_out()
M_lhs = 0 : Zero out the matrix.
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
friend void Mt_S(MatrixOp *mwo_lhs, value_type alpha)
const MatNorm calc_norm(EMatNormType requested_norm_type=MAT_NORM_1, bool allow_replacement=false) const
Compute a norm of this matrix.
The induced infinity norm ||M||inf, i.e. max abs row sum.
virtual mat_ptr_t perm_view_update(const Permutation *P_row, const index_type row_part[], int num_row_part, const Permutation *P_col, const index_type col_part[], int num_col_part, const mat_ptr_t &perm_view) const
Reinitialize a permuted view: M_perm = P_row' * M * P_col.
virtual mat_mut_ptr_t clone()
Clone the non-const matrix object (if supported).
friend value_type transVtMtV(const Vector &v_rhs1, const MatrixOp &M_rhs2, BLAS_Cpp::Transp trans_rhs2, const Vector &v_rhs3)
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)
Interface adding operations specific for a symmetric matrix {abstract}.
virtual std::ostream & output(std::ostream &out) const
Virtual output function.
void syr2k(const MatrixOp &M, 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, MatrixSymOp *symwo_lhs)
symwo_lhs += alpha*op(P1')*op(M)*op(P2) + alpha*op(P2')*op(M')*op(P1) + beta*symwo_lhs ...
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).
. One-based subregion index range class.
The induced two (i.e. Euclidean norm) norm ||M||2.
The Forbenious norm, i.e. max abs element.
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 syrk(const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta, MatrixSymOp *sym_lhs)
Perform a rank-k update of a symmetric matrix of the form:
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).
const LAPACK_C_Decl::f_int & M
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.
void Mt_S(MatrixOp *mwo_lhs, value_type alpha)
The induced one norm ||M||1, i.e. max abs col sum.
Abstract interface to permutation matrices.
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)
Abstract interface for mutable coordinate vectors {abstract}.
virtual mat_ptr_t perm_view(const Permutation *P_row, const index_type row_part[], int num_row_part, const Permutation *P_col, const index_type col_part[], int num_col_part) const
Create a permuted view: M_perm = P_row' * M * P_col.
friend void syrk(const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta, MatrixSymOp *sym_lhs)
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)
Transp
TRANS.
MatNorm(value_type _value, EMatNormType _type)
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)
virtual MatrixOp & operator=(const MatrixOp &mwo_rhs)
M_lhs = mwo_rhs : Virtual assignment operator.
friend 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)