SparseVector BSplineBasis1D::evaluateDerivative(double x, int r) const
{
    // Evaluate rth derivative of basis functions at x
    // Returns vector [D^(r)B_(u-p,p)(x) ... D^(r)B_(u,p)(x)]
    // where u is the knot index and p is the degree
    int p = degree;

    // Continuity requirement
    //assert(p >= r+1);
    if(!(p >= r+1))
    {
        // Return zero-gradient
        SparseVector DB(numBasisFunctions());
        return DB;
    }

    // Check for knot multiplicity here!

    supportHack(x);

    int knotIndex = indexHalfopenInterval(x);

    // Algorithm 3.18 from Lyche and Moerken 2011
    SparseMatrix B(1,1);
    B.insert(0,0) = 1;

    for(int i = 1; i <= p-r; i++)
    {
        SparseMatrix R = buildBasisMatrix(x, knotIndex, i);
        B = B*R;
    }

    for(int i = p-r+1; i <= p; i++)
    {
        SparseMatrix DR = buildBasisMatrix(x, knotIndex, i, true);
        B = B*DR;
    }
    double factorial = std::tgamma(p+1)/std::tgamma(p-r+1);
    B = B*factorial;

    assert(B.cols() == p+1);

    // From row vector to extended column vector
    SparseVector DB(numBasisFunctions());
    DB.reserve(p+1);
    int i = knotIndex-p; // First insertion index
    for(int k = 0; k < B.outerSize(); ++k)
    for(SparseMatrix::InnerIterator it(B,k); it; ++it)
    {
        DB.insert(i+it.col()) = it.value();
    }

    return DB;
}
// Old implementation of Jacobian
DenseMatrix Basis::evalBasisJacobianOld(DenseVector &x) const
{
    // Jacobian basis matrix
    DenseMatrix J; J.setZero(numBasisFunctions(), numVariables);

    // Calculate partial derivatives
    for (unsigned int i = 0; i < numVariables; i++)
    {
        // One column in basis jacobian
        DenseVector bi; bi.setOnes(1);
        for (unsigned int j = 0; j < numVariables; j++)
        {
            DenseVector temp = bi;
            DenseVector xi;
            if (j == i)
            {
                // Differentiated basis 
                xi = bases.at(j).evaluateFirstDerivative(x(j));
            }
            else
            {
                // Normal basis
                xi = bases.at(j).evaluate(x(j));
            }

            bi = kroneckerProduct(temp, xi);
            //bi = kroneckerProduct(bi, xi).eval();
        }

        // Fill out column
        J.block(0,i,bi.rows(),1) = bi.block(0,0,bi.rows(),1);
    }

    return J;
}
SparseVector Basis::eval(const DenseVector &x) const
{
    // Set correct starting size and values: one element of 1.
    SparseMatrix tv1(1,1);  //tv1.resize(1);
    tv1.reserve(1);
    tv1.insert(0,0) = 1;
    SparseMatrix tv2 = tv1;

    // Evaluate all basisfunctions in each dimension i and generate the tensor product of the function values
    for (int dim = 0; dim < x.size(); dim++)
    {
        SparseVector xi = bases.at(dim).evaluate(x(dim));

        // Avoid matrix copy
        if (dim % 2 == 0)
        {
            tv1 = kroneckerProduct(tv2, xi); // tv1 = tv1 x xi
        }
        else
        {
            tv2 = kroneckerProduct(tv1, xi); // tv2 = tv1 x xi
        }
    }

    // Return correct vector
    if (tv1.size() == numBasisFunctions())
    {
        return tv1;
    }
    return tv2;
}
SparseVector BSplineBasis1D::evaluate(double x) const
{
    SparseVector basisvalues(numBasisFunctions());

    supportHack(x);

    std::vector<int> indexSupported = indexSupportedBasisfunctions(x);

    basisvalues.reserve(indexSupported.size());

    // Iterate through the nonzero basisfunctions and store functionvalues
    for(auto it = indexSupported.begin(); it != indexSupported.end(); it++)
    {
        basisvalues.insert(*it) = deBoorCox(x, *it, degree);
    }

    // Alternative evaluation using basis matrix
//    int knotIndex = indexHalfopenInterval(x); // knot index

//    SparseMatrix basisvalues2 = buildBsplineMatrix(x, knotIndex, 1);
//    for (int i = 2; i <= basisDegree; i++)
//    {
//        SparseMatrix Ri = buildBsplineMatrix(x, knotIndex, i);
//        basisvalues2 = basisvalues2*Ri;
//    }
//    basisvalues2.makeCompressed();

//    assert(basisvalues2.rows() == 1);
//    assert(basisvalues2.cols() == basisDegree + 1);

    return basisvalues;
}
// Old implementation of first derivative of basis functions
DenseVector BSplineBasis1D::evaluateFirstDerivative(double x) const
{
    DenseVector values;
    values.setZero(numBasisFunctions());

    supportHack(x);

    std::vector<int> supportedBasisFunctions = indexSupportedBasisfunctions(x);

    for(int i : supportedBasisFunctions)
    {
        // Differentiate basis function
        // Equation 3.35 in Lyche & Moerken (2011)
        double b1 = deBoorCox(x, i, degree-1);
        double b2 = deBoorCox(x, i+1, degree-1);

        double t11 = knots.at(i);
        double t12 = knots.at(i+degree);
        double t21 = knots.at(i+1);
        double t22 = knots.at(i+degree+1);

        (t12 == t11) ? b1 = 0 : b1 = b1/(t12-t11);
        (t22 == t21) ? b2 = 0 : b2 = b2/(t22-t21);

        values(i) = degree*(b1 - b2);
    }

    return values;
}
SparseMatrix Basis::evalBasisJacobian(DenseVector &x) const
{
    // Jacobian basis matrix
    SparseMatrix J(numBasisFunctions(), numVariables);
    //J.setZero(numBasisFunctions(), numInputs);

    // Calculate partial derivatives
    for (unsigned int i = 0; i < numVariables; i++)
    {
        // One column in basis jacobian
        SparseMatrix Ji(1,1);
        Ji.insert(0,0) = 1;
        for (unsigned int j = 0; j < numVariables; j++)
        {
            SparseMatrix temp = Ji;
            SparseMatrix xi;
            if (j == i)
            {
                // Differentiated basis
                xi = bases.at(j).evaluateDerivative(x(j),1);
            }
            else
            {
                // Normal basis
                xi = bases.at(j).evaluate(x(j));
            }
            Ji = kroneckerProduct(temp, xi);
            //myKroneckerProduct(temp,xi,Ji);
            //Ji = kroneckerProduct(Ji, xi).eval();
        }

        // Fill out column
        for (int k = 0; k < Ji.outerSize(); ++k)
        for (SparseMatrix::InnerIterator it(Ji,k); it; ++it)
        {
            if (it.value() != 0)
                J.insert(it.row(),i) = it.value();
        }
        //J.block(0,i,Ji.rows(),1) = bi.block(0,0,Ji.rows(),1);
    }

    J.makeCompressed();

    return J;
}
SparseMatrix Basis::evalBasisHessian(DenseVector &x) const
{
    // Hessian basis matrix
    /* Hij = B1 x ... x DBi x ... x DBj x ... x Bn
     * (Hii = B1 x ... x DDBi x ... x Bn)
     * Where B are basis functions evaluated at x,
     * DB are the derivative of the basis functions,
     * and x is the kronecker product.
     * Hij is in R^(numBasisFunctions x 1)
     * so that basis hessian H is in R^(numBasisFunctions*numInputs x numInputs)
     * The real B-spline Hessian is calculated as (c^T x 1^(numInputs x 1))*H
     */
    SparseMatrix H(numBasisFunctions()*numVariables, numVariables);
    //H.setZero(numBasisFunctions()*numInputs, numInputs);

    // Calculate partial derivatives
    // Utilizing that Hessian is symmetric
    // Filling out lower left triangular
    for (unsigned int i = 0; i < numVariables; i++) // row
    {
        for (unsigned int j = 0; j <= i; j++) // col
        {
            // One column in basis jacobian
            SparseMatrix Hi(1,1);
            Hi.insert(0,0) = 1;

            for (unsigned int k = 0; k < numVariables; k++)
            {
                SparseMatrix temp = Hi;
                SparseMatrix Bk;
                if (i == j && k == i)
                {
                    // Diagonal element
                    Bk = bases.at(k).evaluateDerivative(x(k),2);
                }
                else if (k == i || k == j)
                {
                    Bk = bases.at(k).evaluateDerivative(x(k),1);
                }
                else
                {
                    Bk = bases.at(k).evaluate(x(k));
                }
                Hi = kroneckerProduct(temp, Bk);
            }

            // Fill out column
            for (int k = 0; k < Hi.outerSize(); ++k)
            for (SparseMatrix::InnerIterator it(Hi,k); it; ++it)
            {
                if (it.value() != 0)
                {
                    int row = i*numBasisFunctions()+it.row();
                    int col = j;
                    H.insert(row,col) = it.value();
                }
            }
        }
    }

    H.makeCompressed();

    return H;
}