Lösungsvorschlag

Zur ersten Frage: typename war notwendig, damit der Übersetzer bei der syntaktischen Analyse von GeMatrix::Index erkennt, dass es sich dabei um einen Typnamen handelt. Zwar weiß der Übersetzer bei dem Template, dass GeMatrix ein Typ ist, jedoch bleibt vor einer Instantiierung verborgen, um was es sich bei GeMatrix::Index handeln könnte.

Zur zweiten Frage: So könnte init_complex_value aussehen:

template<typename T, typename Index>
std::complex<T> init_complex_value(Index i, Index j, Index m, Index n) {
   return std::complex<T>(T(i), T(j));
}

Und der zugehörige Aufruf:

GeMatrix<std::complex<double>> B(2, 3, StorageOrder::ColMajor);
initGeMatrix(B, init_complex_value<double, std::size_t>);

Zum dritten Teil: Hier ist die angepasste Funktion in bench.h:

template <typename T, typename Index, typename InitValue>
void
initGeMatrix(Index m, Index n, T *A, Index incRowA, Index incColA,
             InitValue initValue)
{
    for (Index j=0; j<n; ++j) {
        for (Index i=0; i<m; ++i) {
            A[i*incRowA+j*incColA] = initValue(i, j, m, n);
        }
    }
}

Und in gematrix.h:

template <typename GeMatrix, typename InitValue>
void
initGeMatrix(GeMatrix &A, InitValue initValue)
{
    bench::initGeMatrix(A.m, A.n, A.data, A.incRow, A.incCol, initValue);
}

Hier sind alle Dateien zum Lösungsvorschlag:

#ifndef HPC_BENCH_H
#define HPC_BENCH_H 1

#include <chrono>
#include <cmath>
#include <complex>
#include <random>

namespace bench {

template <typename T, typename Index>
T
asumDiffGeMatrix(Index m, Index n,
                 const T *A, Index incRowA, Index incColA,
                 T *B, Index incRowB, Index incColB)
{

    T asum = 0;

    for (Index j=0; j<n; ++j) {
        for (Index i=0; i<m; ++i) {
            asum += std::abs(B[i*incRowB+j*incColB] - A[i*incRowA+j*incColA]);
        }
    }
    return asum;
}

template <typename T, typename Index, typename InitValue>
void
initGeMatrix(Index m, Index n, T *A, Index incRowA, Index incColA,
             InitValue initValue)
{
    for (Index j=0; j<n; ++j) {
        for (Index i=0; i<m; ++i) {
            A[i*incRowA+j*incColA] = initValue(i, j, m, n);
        }
    }
}

template <typename T>
struct WallTime
{
    void
    tic()
    {
        t0 = std::chrono::high_resolution_clock::now();
    }

    T
    toc()
    {
        using namespace std::chrono;

        elapsed = high_resolution_clock::now() - t0;
        return duration<T,seconds::period>(elapsed).count();
    }

    std::chrono::high_resolution_clock::time_point t0;
    std::chrono::high_resolution_clock::duration   elapsed;
};


} // namespace bench


#endif // HPC_BENCH_H
#ifndef HPC_GEMATRIX_H
#define HPC_GEMATRIX_H 1

#include <cstddef>
#include <cassert>
#include "bench.h"

namespace matvec {

enum class StorageOrder {
    ColMajor,
    RowMajor
};

template <typename T, typename I>
    struct GeMatrixView;

template <typename T, typename I=std::size_t>
struct GeMatrix
{
    typedef T                       ElementType;
    typedef I                       Index;
    typedef GeMatrix<T,Index>       NoView;
    typedef GeMatrixView<T,Index>   View;

    GeMatrix(Index m, Index n, StorageOrder order=StorageOrder::ColMajor)
        : m(m), n(n),
          incRow(order==StorageOrder::ColMajor ? 1: n),
          incCol(order==StorageOrder::RowMajor ? 1: m),
          data(new T[m*n])
    {
    }

    ~GeMatrix()
    {
        delete[] data;
    }

    const ElementType &
    operator()(Index i, Index j) const
    {
        assert(i<m && j<n);
        return data[i*incRow + j*incCol];
    }

    ElementType &
    operator()(Index i, Index j)
    {
        assert(i<m && j<n);
        return data[i*incRow + j*incCol];
    }

    View
    operator()(Index i, Index j, Index numRows, Index numCols)
    {
        assert(i+numRows<=m);
        assert(j+numCols<=n);
        return View(numRows, numCols, &(operator()(i,j)), incRow, incCol);
    }

    const Index     m, n, incRow, incCol;
    ElementType*    data;
};

template <typename T, typename I=std::size_t>
struct GeMatrixView
{
    typedef T                       ElementType;
    typedef I                       Index;
    typedef GeMatrix<T,Index>       NoView;
    typedef GeMatrixView<T,Index>   View;

    GeMatrixView(Index m, Index n, T *data, Index incRow, Index incCol)
        : m(m), n(n), incRow(incRow), incCol(incCol), data(data)
    {
    }

    GeMatrixView(const GeMatrixView &rhs)
        : m(rhs.m), n(rhs.n), incRow(rhs.incRow), incCol(rhs.incCol),
          data(rhs.data)
    {
    }

    const ElementType &
    operator()(Index i, Index j) const
    {
        assert(i<m && j<n);
        return data[i*incRow + j*incCol];
    }

    ElementType &
    operator()(Index i, Index j)
    {
        assert(i<m && j<n);
        return data[i*incRow + j*incCol];
    }

    View
    operator()(Index i, Index j, Index numRows, Index numCols)
    {
        assert(i+numRows<=m);
        assert(j+numCols<=n);
        return View(numRows, numCols, &(operator()(i,j)), incRow, incCol);
    }

    const Index     m, n, incRow, incCol;
    ElementType*    data;
};

//
//  Interface for bench
//
template <typename GeMatrix, typename InitValue>
void
initGeMatrix(GeMatrix &A, InitValue initValue)
{
    bench::initGeMatrix(A.m, A.n, A.data, A.incRow, A.incCol, initValue);
}

template <typename GeMatrix>
typename GeMatrix::ElementType
asumDiffGeMatrix(const GeMatrix &A, const GeMatrix &B)
{
    return bench::asumDiffGeMatrix(A.m, A.n,
                                   A.data, A.incRow, A.incCol,
                                   B.data, B.incRow, B.incCol);
}

} // namespace matvec

#endif // HPC_GEMATRIX_H
#include <cstdlib>
#include <cstdio>
#include "gematrix.h"

template<typename T, typename Index>
T init_value(Index i, Index j, Index m, Index n) {
   return T(j)*T(n) + T(i) + T(1);
}

template<typename T, typename Index>
std::complex<T> init_complex_value(Index i, Index j, Index m, Index n) {
   return std::complex<T>(T(i), T(j));
}

void print_value(double value) {
   std::printf(" %4.1lf", value);
}

void print_value(std::complex<double> value) {
   std::printf(" (%4.1lf, %4.1lf)", value.real(), value.imag());
}

template<typename Matrix>
void print_matrix(const Matrix& A) {
   for (std::size_t i = 0; i < A.m; ++i) {
      std::printf("  ");
      for (std::size_t j = 0; j < A.n; ++j) {
         print_value(A(i, j));
      }
      std::printf("\n");
   }
}

int main() {
   using namespace matvec;
   GeMatrix<double> A(3, 7, StorageOrder::ColMajor);
   initGeMatrix(A, init_value<double, std::size_t>);
   printf("A:\n");
   print_matrix(A);
   GeMatrix<std::complex<double>> B(2, 3, StorageOrder::ColMajor);
   initGeMatrix(B, init_complex_value<double, std::size_t>);
   printf("B:\n");
   print_matrix(B);
}