Contents

SparseSolve(_:_:_:_:_:)

Solves the equation AX = B for matrices of double-precision values, treating A as an operator and using the specified iterative method.

Declaration

func SparseSolve(_ method: SparseIterativeMethod, _ ApplyOperator: @escaping (Bool, CBLAS_TRANSPOSE, DenseMatrix_Double, DenseMatrix_Double) -> Void, _ B: DenseMatrix_Double, _ X: DenseMatrix_Double, _ Preconditioner: SparseOpaquePreconditioner_Double) -> SparseIterativeStatus_t

Parameters

  • method:

    The iterative method.

  • ApplyOperator:

    The apply operator block to run. The block takes the following parameters:

    accumulate

    Indicates whether to perform Y += op(A)X (if true), or Y = op(A)X (if false).

    trans

    Indicates whether op(A) is the application of A if CblasNoTrans, or Aᵀ if CblasTrans

    X

    The matrix to multiply.

    Y

    The matrix for accumulating or storing the result.

  • B:

    The matrix B.

  • X:

    The matrix X.

  • Preconditioner:

    The preconditioner to apply.

Return Value

A SparseIterativeStatus_t enumeration that represents the status of the iterative solve.

Discussion

Use this function to solve a system of linear equations using a factored coefficient matrix. Preconditioning the coefficient matrix can reduce the number of iterations the function requires to converge the system. In cases where the matrix A isn’t explicitly available or you need control over the multiplication, this function allows you to provide an apply block.

The following figure shows two systems of equations where the coefficient matrix is sparse:

[Image]

The following code solves this system by applying a diagonal scaling preconditioner and using the least squares minimum residual method:

/// Create the coefficient matrix _A_
let rowIndices: [Int32] =    [ 0,  1, 1,  2]
let columnIndices: [Int32] = [ 2,  0, 2,  1]
let aValues: [Double] =       [10, 20, 5, 50]

let A = SparseConvertFromCoordinate(3, 3,
                                    4, 1,
                                    SparseAttributes_t(),
                                    rowIndices, columnIndices,
                                    aValues)
let preconditioner = SparseCreatePreconditioner(SparsePreconditionerDiagScaling,
                                                A)

defer {
    SparseCleanup(A)
    SparseCleanup(preconditioner)
}

/// Create the right-hand-side matrix, _B_.
var bValues: [Double] = [30, 35, 100,
                         300, 350, 1000]

var xValues: [Double] = [0, 0, 0,
                         0, 0, 0]

/// Create the apply operator block.
func applyOperator(accumulate: Bool,
                   trans: CBLAS_TRANSPOSE,
                   X: DenseMatrix_Double,
                   Y: DenseMatrix_Double) {
    switch(accumulate, trans == CblasTrans) {
        case (false, false):
            SparseMultiply(A, X, Y)
        case (false, true):
            SparseMultiply(SparseGetTranspose(A), X, Y)
        case (true, false):
            SparseMultiplyAdd(A, X, Y)
        case (true, true):
            SparseMultiplyAdd(SparseGetTranspose(A), X, Y)
    }
}

xValues.withUnsafeMutableBufferPointer { xPtr in
    bValues.withUnsafeMutableBufferPointer { bPtr in
        let B = DenseMatrix_Double(rowCount: 3,
                                   columnCount: 2,
                                   columnStride: 3,
                                   attributes: SparseAttributes_t(),
                                   data: bPtr.baseAddress!)
        
        let X = DenseMatrix_Double(rowCount: 3,
                                   columnCount: 2,
                                   columnStride: 3,
                                   attributes: SparseAttributes_t(),
                                   data: xPtr.baseAddress!)
        
        SparseSolve(SparseLSMR(),
                    applyOperator,
                    B, X,
                    preconditioner)
    }
}

On return, xValues contains the values [1.0, 2.0, 3.0, 10.0, 20.0, 30.0].

See Also

Iterative sparse solve functions with preconditioner