This is the mail archive of the gsl-discuss@sources.redhat.com mailing list for the GSL project.


Index Nav: [Date Index] [Subject Index] [Author Index] [Thread Index]
Message Nav: [Date Prev] [Date Next] [Thread Prev] [Thread Next]

GSL matrix transpose view


GSL function

    int gsl_matrix_transpose(gsl_matrix*)

transposes a square matrix in place.

But, more often, it is useful to create
a transpose view of a matrix object.
All that is required, is a second stride
to create a transpose view

    gsl_matrix T = gsl_matrix_t(&M);

of matrix M.

GSL functions must be modified slightly
so that transpose views can be used
just like any other submatrix views.


        typedef struct {
          gsl_handle H;
          gsl_offset O;
          gsl_extent N1;
          gsl_stride S1;
          gsl_extent N2;
          gsl_stride S2;
          } gsl_matrix;

        inline void
        gsl_matrix_init(gsl_matrix* pM,
          gsl_handle h, gsl_offset o,
          gsl_extent n2, gsl_stride s2,
          gsl_extent n1, gsl_stride s1) {
          pM->H = h; pM->O = o;
          pM->N1 = n1; pM->S1 = s1;
          pM->N2 = n2; pM->s2 = s2; }

        inline gsl_handle
        gsl_matrix_handle(const gsl_matrix* pM) {
          return pM->H; }

        inline gsl_offset
        gsl_matrix_offset(const gsl_matrix* pM) {
          return pM->O; }

        inline gsl_extent
        gsl_matrix_extent1(const gsl_matrix* pM) {
          return pM->N1; }

        inline gsl_stride
        gsl_matrix_stride1(const gsl_matrix* pM) {
          return pM->S1; }

        inline gsl_extent
        gsl_matrix_extent2(const gsl_matrix* pM) {
          return pM->N2; }

        inline gsl_stride
        gsl_matrix_extent1(const gsl_matrix* pM) {
          return pM->S2; }

        inline gsl_matrix
        gsl_matrix_t(const gsl_matrix* pM) {
          gsl_matrix T;
          gsl_matrix_init(&T,
          gsl_matrix_handle(pM), gsl_matrix_offset(pM),
          gsl_matrix_extent1(pM), gsl_matrix_stride1(pM),
          gsl_matrix_extent2(pM), gsl_matrix_stride2(pM));
          return T; }


Index Nav: [Date Index] [Subject Index] [Author Index] [Thread Index]
Message Nav: [Date Prev] [Date Next] [Thread Prev] [Thread Next]