int read_sparse(char *fname, unsigned long fname_len)
{
    SciErr sciErr;
    int i, j, k;
    int* piAddr			= NULL;
    int iRows			= 0;
    int iCols			= 0;
    int iNbItem			= 0;
    int* piNbItemRow	= NULL;
    int* piColPos		= NULL;
    double* pdblReal	= NULL;
    double* pdblImg		= NULL;

    CheckInputArgument(pvApiCtx, 1, 1);

    sciErr = getVarAddressFromPosition(pvApiCtx, 1, &piAddr);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 0;
    }

    if (isVarComplex(pvApiCtx, piAddr))
    {
        sciErr = getComplexSparseMatrix(pvApiCtx, piAddr, &iRows, &iCols, &iNbItem, &piNbItemRow, &piColPos, &pdblReal, &pdblImg);
    }
    else
    {
        sciErr = getSparseMatrix(pvApiCtx, piAddr, &iRows, &iCols, &iNbItem, &piNbItemRow, &piColPos, &pdblReal);
    }

    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 0;
    }

    sciprint("Sparse %d item(s)\n", iNbItem);
    k = 0;

    for (i = 0 ; i < iRows ; i++)
    {
        for (j = 0 ; j < piNbItemRow[i] ; j++)
        {
            sciprint("(%d,%d) = %f", i + 1, piColPos[k], pdblReal[k]);
            if (isVarComplex(pvApiCtx, piAddr))
            {
                sciprint(" %+fi", pdblImg[k]);
            }

            sciprint("\n");
            k++;
        }
    }

    //assign allocated variables to Lhs position
    AssignOutputVariable(pvApiCtx, 1) = 0;
    return 0;
}
static bool export_sparse(int _iH5File, int *_piVar, char* _pstName)
{
    int iRet = 0;
    int iNbItem = 0;
    int* piNbItemRow = NULL;
    int* piColPos = NULL;
    double* pdblReal = NULL;
    double* pdblImg = NULL;
    int piDims[2];
    SciErr sciErr;

    if (isVarComplex(pvApiCtx, _piVar))
    {
        sciErr = getComplexSparseMatrix(pvApiCtx, _piVar, &piDims[0], &piDims[1], &iNbItem, &piNbItemRow, &piColPos, &pdblReal, &pdblImg);
        if (sciErr.iErr)
        {
            printError(&sciErr, 0);
            return false;
        }

        iRet = writeSparseComplexMatrix(_iH5File, _pstName, piDims[0], piDims[1], iNbItem, piNbItemRow, piColPos, pdblReal, pdblImg);
    }
    else
    {
        sciErr = getSparseMatrix(pvApiCtx, _piVar, &piDims[0], &piDims[1], &iNbItem, &piNbItemRow, &piColPos, &pdblReal);
        if (sciErr.iErr)
        {
            printError(&sciErr, 0);
            return false;
        }

        iRet = writeSparseMatrix(_iH5File, _pstName, piDims[0], piDims[1], iNbItem, piNbItemRow, piColPos, pdblReal);
    }

    if (iRet)
    {
        return false;
    }

    char pstMsg[512];
    sprintf(pstMsg, "sparse (%d x %d)", piDims[0], piDims[1]);
    print_type(pstMsg);
    return true;
}
int get_sparse_info(void* _pvCtx, int _iRhs, int* _piParent, int *_piAddr, int _iItemPos)
{
    SciErr sciErr;
    int iRows           = 0;
    int iCols           = 0;
    int iItem           = 0;
    int* piNbRow        = NULL;
    int* piColPos       = NULL;
    double* pdblReal    = NULL;
    double* pdblImg     = NULL;

    if (_iItemPos == 0)
    {
        //Not in list
        if (isVarComplex(_pvCtx, _piAddr))
        {
            sciErr = getComplexSparseMatrix(_pvCtx, _piAddr, &iRows, &iCols, &iItem, &piNbRow, &piColPos, &pdblReal, &pdblImg);
        }
        else
        {
            sciErr = getSparseMatrix(_pvCtx, _piAddr, &iRows, &iCols, &iItem, &piNbRow, &piColPos, &pdblReal);
        }
    }
    else
    {
        if (isVarComplex(_pvCtx, _piAddr))
        {
            sciErr = getComplexSparseMatrixInList(_pvCtx, _piParent, _iItemPos, &iRows, &iCols, &iItem, &piNbRow, &piColPos, &pdblReal, &pdblImg);
        }
        else
        {
            sciErr = getSparseMatrixInList(_pvCtx, _piParent, _iItemPos, &iRows, &iCols, &iItem, &piNbRow, &piColPos, &pdblReal);
        }
    }

    FREE(piNbRow);
    FREE(piColPos);
    FREE(pdblReal);
    FREE(pdblImg);

    insert_indent();
    sciprint("Sparse (%d x %d), Item(s) : %d \n", iRows, iCols, iItem);
    return 0;;
}
int sci_umf_lufact(char* fname, void* pvApiCtx)
{
    SciErr sciErr;
    int stat = 0;
    SciSparse AA;
    CcsSparse A;

    int mA              = 0; // rows
    int nA              = 0; // cols
    int iNbItem         = 0;
    int* piNbItemRow    = NULL;
    int* piColPos       = NULL;
    double* pdblSpReal  = NULL;
    double* pdblSpImg   = NULL;

    /* umfpack stuff */
    double* Control = NULL;
    double* Info    = NULL;
    void* Symbolic  = NULL;
    void* Numeric   = NULL;

    int* piAddr1 = NULL;
    int iComplex = 0;
    int iType1   = 0;

    /* Check numbers of input/output arguments */
    CheckInputArgument(pvApiCtx, 1, 1);
    CheckOutputArgument(pvApiCtx, 1, 1);

    /* get A the sparse matrix to factorize */
    sciErr = getVarAddressFromPosition(pvApiCtx, 1, &piAddr1);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    /* check if the first argument is a sparse matrix */
    sciErr = getVarType(pvApiCtx, piAddr1, &iType1);
    if (sciErr.iErr || iType1 != sci_sparse)
    {
        printError(&sciErr, 0);
        Scierror(999, _("%s: Wrong type for input argument #%d: A sparse matrix expected.\n"), fname, 1);
        return 1;
    }

    if (isVarComplex(pvApiCtx, piAddr1))
    {
        iComplex = 1;
        sciErr = getComplexSparseMatrix(pvApiCtx, piAddr1, &mA, &nA, &iNbItem, &piNbItemRow, &piColPos, &pdblSpReal, &pdblSpImg);
    }
    else
    {
        sciErr = getSparseMatrix(pvApiCtx, piAddr1, &mA, &nA, &iNbItem, &piNbItemRow, &piColPos, &pdblSpReal);
    }

    if (sciErr.iErr)
    {
        FREE(piNbItemRow);
        FREE(piColPos);
        FREE(pdblSpReal);
        if (pdblSpImg)
        {
            FREE(pdblSpImg);
        }
        printError(&sciErr, 0);
        return 1;
    }

    // fill struct sparse
    AA.m     = mA;
    AA.n     = nA;
    AA.it    = iComplex;
    AA.nel   = iNbItem;
    AA.mnel  = piNbItemRow;
    AA.icol  = piColPos;
    AA.R     = pdblSpReal;
    AA.I     = pdblSpImg;

    if (nA <= 0 || mA <= 0)
    {
        FREE(piNbItemRow);
        FREE(piColPos);
        FREE(pdblSpReal);
        if (pdblSpImg)
        {
            FREE(pdblSpImg);
        }
        Scierror(999, _("%s: Wrong size for input argument #%d.\n"), fname, 1);
        return 1;
    }

    SciSparseToCcsSparse(&AA, &A);

    FREE(piNbItemRow);
    FREE(piColPos);
    FREE(pdblSpReal);
    if (pdblSpImg)
    {
        FREE(pdblSpImg);
    }

    /* symbolic factorization */
    if (A.it == 1)
    {
        stat = umfpack_zi_symbolic(nA, mA, A.p, A.irow, A.R, A.I, &Symbolic, Control, Info);
    }
    else
    {
        stat = umfpack_di_symbolic(nA, mA, A.p, A.irow, A.R, &Symbolic, Control, Info);
    }

    if (stat != UMFPACK_OK)
    {
        freeCcsSparse(A);
        Scierror(999, _("%s: An error occurred: %s: %s\n"), fname, _("symbolic factorization"), UmfErrorMes(stat));
        return 1;
    }

    /* numeric factorization */
    if (A.it == 1)
    {
        stat = umfpack_zi_numeric(A.p, A.irow, A.R, A.I, Symbolic, &Numeric, Control, Info);
    }
    else
    {
        stat = umfpack_di_numeric(A.p, A.irow, A.R, Symbolic, &Numeric, Control, Info);
    }

    if (A.it == 1)
    {
        umfpack_zi_free_symbolic(&Symbolic);
    }
    else
    {
        umfpack_di_free_symbolic(&Symbolic);
    }

    if ( stat != UMFPACK_OK  &&  stat != UMFPACK_WARNING_singular_matrix )
    {
        freeCcsSparse(A);
        Scierror(999, _("%s: An error occurred: %s: %s\n"), fname, _("symbolic factorization"), UmfErrorMes(stat));
        return 1;
    }

    if ( stat == UMFPACK_WARNING_singular_matrix  &&  mA == nA )
    {
        if (getWarningMode())
        {
            Sciwarning("\n%s:%s\n", _("Warning"), _("The (square) matrix appears to be singular."));
        }
    }

    /*  add the pointer in the list ListNumeric  */
    if (! AddAdrToList(Numeric, A.it, &ListNumeric))
    {
        /* AddAdrToList return 0 if malloc have failed : as it is just
        for storing 2 pointers this is unlikely to occurs but ... */
        if (A.it == 1)
        {
            umfpack_zi_free_numeric(&Numeric);
        }
        else
        {
            umfpack_di_free_numeric(&Numeric);
        }

        freeCcsSparse(A);
        Scierror(999, _("%s: An error occurred: %s\n"), fname, _("no place to store the LU pointer in ListNumeric."));
        return 1;
    }

    freeCcsSparse(A);

    /* create the scilab object to store the pointer onto the LU factors */
    sciErr = createPointer(pvApiCtx, 2, Numeric);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    /* return the pointer */
    AssignOutputVariable(pvApiCtx, 1) = 2;
    ReturnArguments(pvApiCtx);
    return 0;
}
int sci_taucs_chfact(char* fname, void* pvApiCtx)
{
    SciErr sciErr;
    int stat = 0;
    int* perm = NULL;
    int* invperm = NULL;
    taucs_ccs_matrix *PAPT;
    taucs_ccs_matrix B;
    void *C = NULL;
    taucs_handle_factors *pC;

    SciSparse A;
    int mA              = 0; // rows
    int nA              = 0; // cols
    int iNbItem         = 0;
    int* piNbItemRow    = NULL;
    int* piColPos       = NULL;
    double* pdblSpReal  = NULL;
    double* pdblSpImg   = NULL;
    int iComplex        = 0;
    int* piAddr1        = NULL;

    /* Check numbers of input/output arguments */
    CheckInputArgument(pvApiCtx, 1, 1);
    CheckOutputArgument(pvApiCtx, 1, 1);

    /* get A the sparse matrix to factorize */
    sciErr = getVarAddressFromPosition(pvApiCtx, 1, &piAddr1);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    if (isVarComplex(pvApiCtx, piAddr1))
    {
        iComplex = 1;
        sciErr = getComplexSparseMatrix(pvApiCtx, piAddr1, &mA, &nA, &iNbItem, &piNbItemRow, &piColPos, &pdblSpReal, &pdblSpImg);
    }
    else
    {
        sciErr = getSparseMatrix(pvApiCtx, piAddr1, &mA, &nA, &iNbItem, &piNbItemRow, &piColPos, &pdblSpReal);
    }

    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    // fill struct sparse
    A.m     = mA;
    A.n     = nA;
    A.it    = iComplex;
    A.nel   = iNbItem;
    A.mnel  = piNbItemRow;
    A.icol  = piColPos;
    A.R     = pdblSpReal;
    A.I     = pdblSpImg;

    stat = spd_sci_sparse_to_taucs_sparse(&A, &B);
    if ( stat != A_PRIORI_OK )
    {
        if ( stat == MAT_IS_NOT_SPD )
        {
            freeTaucsSparse(B);
            Scierror(999, _("%s: Wrong value for input argument #%d: Must be symmetric positive definite matrix."), fname, 1);
        }
        /* the message for the other problem (not enough memory in stk) is treated automaticaly */
        return 1;
    }

    /* find the permutation */
    taucs_ccs_genmmd(&B, &perm, &invperm);
    if ( !perm )
    {
        freeTaucsSparse(B);
        Scierror(999, _("%s: No more memory.\n") , fname);
        return 1;
    }

    /* apply permutation */
    PAPT = taucs_ccs_permute_symmetrically(&B, perm, invperm);
    FREE(invperm);
    freeTaucsSparse(B);

    /* factor */
    C = taucs_ccs_factor_llt_mf(PAPT);
    taucs_ccs_free(PAPT);

    if (C == NULL)
    {
        /*   Note : an error indicator is given in the main scilab window
         *          (out of memory, no positive definite matrix , etc ...)
         */
        Scierror(999, _("%s: An error occurred: %s\n"), fname, _("factorization"));
        return 1;
    }

    /* put in an handle (Chol fact + perm + size) */
    pC = (taucs_handle_factors*)MALLOC( sizeof(taucs_handle_factors) );
    pC->p = perm;
    pC->C = C;
    pC->n = A.n;

    /* add in the list of Chol Factors  */
    AddAdrToList((Adr) pC, 0, &ListCholFactors);  /* FIXME add a test here .. */

    /* create the scilab object to store the pointer onto the Chol handle */
    sciErr = createPointer(pvApiCtx, 2, (void *)pC);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    /* return the pointer */
    AssignOutputVariable(pvApiCtx, 1) = 2;
    ReturnArguments(pvApiCtx);
    return 0;
}
int sci_umfpack(char* fname, void* pvApiCtx)
{
    SciErr sciErr;

    int mb      = 0;
    int nb      = 0;
    int i       = 0;
    int num_A   = 0;
    int num_b   = 0;
    int mW      = 0;
    int Case    = 0;
    int stat    = 0;

    SciSparse AA;
    CcsSparse A;

    int* piAddrA = NULL;
    int* piAddr2 = NULL;
    int* piAddrB = NULL;

    double* pdblBR = NULL;
    double* pdblBI = NULL;
    double* pdblXR = NULL;
    double* pdblXI = NULL;

    int iComplex = 0;
    int freepdblBI = 0;

    int mA              = 0; // rows
    int nA              = 0; // cols
    int iNbItem         = 0;
    int* piNbItemRow    = NULL;
    int* piColPos       = NULL;
    double* pdblSpReal  = NULL;
    double* pdblSpImg   = NULL;

    /* umfpack stuff */
    double Info[UMFPACK_INFO];
    double* Control = NULL;
    void* Symbolic  = NULL;
    void* Numeric   = NULL;
    int* Wi         = NULL;
    double* W       = NULL;
    char* pStr      = NULL;
    int iType2      = 0;
    int iTypeA      = 0;
    int iTypeB      = 0;

    /* Check numbers of input/output arguments */
    CheckInputArgument(pvApiCtx, 3, 3);
    CheckOutputArgument(pvApiCtx, 1, 1);

    /* First get arg #2 : a string of length 1 */
    sciErr = getVarAddressFromPosition(pvApiCtx, 2, &piAddr2);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    sciErr = getVarType(pvApiCtx, piAddr2, &iType2);
    if (sciErr.iErr || iType2 != sci_strings)
    {
        printError(&sciErr, 0);
        Scierror(999, _("%s: Wrong type for input argument #%d: string expected.\n"), fname, 2);
        return 1;
    }

    if (getAllocatedSingleString(pvApiCtx, piAddr2, &pStr))
    {
        return 1;
    }

    /* select Case 1 or 2 depending (of the first char of) the string ... */
    if (pStr[0] == '\\') // compare pStr[0] with '\'
    {
        Case  = 1;
        num_A = 1;
        num_b = 3;
    }
    else if (pStr[0] == '/')
    {
        Case  = 2;
        num_A = 3;
        num_b = 1;
    }
    else
    {
        Scierror(999, _("%s: Wrong input argument #%d: '%s' or '%s' expected.\n"), fname, 2, "\\", "/");
        FREE(pStr);
        return 1;
    }
    FREE(pStr);

    /* get A */
    sciErr = getVarAddressFromPosition(pvApiCtx, num_A, &piAddrA);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    sciErr = getVarType(pvApiCtx, piAddrA, &iTypeA);
    if (sciErr.iErr || iTypeA != sci_sparse)
    {
        printError(&sciErr, 0);
        Scierror(999, _("%s: Wrong type for input argument #%d: A sparse matrix expected.\n"), fname, 1);
        return 1;
    }

    if (isVarComplex(pvApiCtx, piAddrA))
    {
        AA.it = 1;
        iComplex = 1;
        sciErr = getComplexSparseMatrix(pvApiCtx, piAddrA, &mA, &nA, &iNbItem, &piNbItemRow, &piColPos, &pdblSpReal, &pdblSpImg);
    }
    else
    {
        AA.it = 0;
        sciErr = getSparseMatrix(pvApiCtx, piAddrA, &mA, &nA, &iNbItem, &piNbItemRow, &piColPos, &pdblSpReal);
    }

    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    // fill struct sparse
    AA.m     = mA;
    AA.n     = nA;
    AA.nel   = iNbItem;
    AA.mnel  = piNbItemRow;
    AA.icol  = piColPos;
    AA.R     = pdblSpReal;
    AA.I     = pdblSpImg;

    if ( mA != nA || mA < 1 )
    {
        Scierror(999, _("%s: Wrong size for input argument #%d.\n"), fname, num_A);
        return 1;
    }

    /* get B*/
    sciErr = getVarAddressFromPosition(pvApiCtx, num_b, &piAddrB);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    sciErr = getVarType(pvApiCtx, piAddrB, &iTypeB);
    if (sciErr.iErr || iTypeB != sci_matrix)
    {
        printError(&sciErr, 0);
        Scierror(999, _("%s: Wrong type for input argument #%d: A matrix expected.\n"), fname, 3);
        return 1;
    }

    if (isVarComplex(pvApiCtx, piAddrB))
    {
        iComplex = 1;
        sciErr = getComplexMatrixOfDouble(pvApiCtx, piAddrB, &mb, &nb, &pdblBR, &pdblBI);
    }
    else
    {
        sciErr = getMatrixOfDouble(pvApiCtx, piAddrB, &mb, &nb, &pdblBR);
    }

    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    if ( (Case == 1 && ( mb != mA || nb < 1 )) || (Case == 2 && ( nb != mA || mb < 1 )) )
    {
        Scierror(999, _("%s: Wrong size for input argument #%d.\n"), fname, num_b);
        return 1;
    }

    SciSparseToCcsSparse(&AA, &A);

    /* allocate memory for the solution x */
    if (iComplex)
    {
        sciErr = allocComplexMatrixOfDouble(pvApiCtx, 4, mb, nb, &pdblXR, &pdblXI);
    }
    else
    {
        sciErr = allocMatrixOfDouble(pvApiCtx, 4, mb, nb, &pdblXR);
    }

    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        freeCcsSparse(A);
        return 1;
    }

    if (A.it == 1)
    {
        mW = 10 * mA;
    }
    else
    {
        mW = 5 * mA;
    }

    if (A.it == 1  &&  pdblBI == NULL)
    {
        int iSize = mb * nb * sizeof(double);
        pdblBI = (double*)MALLOC(iSize);
        memset(pdblBI, 0x00, iSize);
        freepdblBI = 1;
    }

    /* Now calling umfpack routines */
    if (A.it == 1)
    {
        stat = umfpack_zi_symbolic(mA, nA, A.p, A.irow, A.R, A.I, &Symbolic, Control, Info);
    }
    else
    {
        stat = umfpack_di_symbolic(mA, nA, A.p, A.irow, A.R, &Symbolic, Control, Info);
    }

    if ( stat  != UMFPACK_OK )
    {
        Scierror(999, _("%s: An error occurred: %s: %s\n"), fname, _("symbolic factorization"), UmfErrorMes(stat));
        freeCcsSparse(A);
        if (freepdblBI)
        {
            FREE(pdblBI);
        }
        return 1;
    }

    if (A.it == 1)
    {
        stat = umfpack_zi_numeric(A.p, A.irow, A.R, A.I, Symbolic, &Numeric, Control, Info);
    }
    else
    {
        stat = umfpack_di_numeric(A.p, A.irow, A.R, Symbolic, &Numeric, Control, Info);
    }

    if (A.it == 1)
    {
        umfpack_zi_free_symbolic(&Symbolic);
    }
    else
    {
        umfpack_di_free_symbolic(&Symbolic);
    }

    if ( stat  != UMFPACK_OK )
    {
        Scierror(999, _("%s: An error occurred: %s: %s\n"), fname, _("numeric factorization"), UmfErrorMes(stat));
        if (A.it == 1)
        {
            umfpack_zi_free_numeric(&Numeric);
        }
        else
        {
            umfpack_di_free_numeric(&Numeric);
        }
        freeCcsSparse(A);
        if (freepdblBI)
        {
            FREE(pdblBI);
        }
        return 1;
    }
 
    /* allocate memory for umfpack_di_wsolve usage or umfpack_zi_wsolve usage*/
    Wi = (int*)MALLOC(mA * sizeof(int));
    W = (double*)MALLOC(mW * sizeof(double));

    if ( Case == 1 )   /*  x = A\b  <=> Ax = b */
    {
        if (A.it == 0)
        {
            for ( i = 0 ; i < nb ; i++ )
            {
                umfpack_di_wsolve(UMFPACK_A, A.p, A.irow, A.R, &pdblXR[i * mb], &pdblBR[i * mb],
                                  Numeric, Control, Info, Wi, W);
            }

            if (isVarComplex(pvApiCtx, piAddrB))
            {
                for ( i = 0 ; i < nb ; i++ )
                {
                    umfpack_di_wsolve(UMFPACK_A, A.p, A.irow, A.R, &pdblXI[i * mb], &pdblBI[i * mb],
                                      Numeric, Control, Info, Wi, W);
                }
            }
        }
        else /*  A.it == 1  */
        {
            for ( i = 0 ; i < nb ; i++ )
            {
                umfpack_zi_wsolve(UMFPACK_A, A.p, A.irow, A.R, A.I, &pdblXR[i * mb], &pdblXI[i * mb],
                                  &pdblBR[i * mb], &pdblBI[i * mb], Numeric, Control, Info, Wi, W);
            }
        }
    }
    else  /* Case == 2,   x = b/A  <=> x A = b <=> A.'x.' = b.' */
    {
        if (A.it == 0)
        {
            TransposeMatrix(pdblBR, mb, nb, pdblXR);    /* put b in x (with transposition) */
            for ( i = 0 ; i < mb ; i++ )
            {
                umfpack_di_wsolve(UMFPACK_At, A.p, A.irow, A.R, &pdblBR[i * nb], &pdblXR[i * nb],
                                  Numeric, Control, Info, Wi, W);      /* the solutions are in br */
            }

            TransposeMatrix(pdblBR, nb, mb, pdblXR);         /* put now br in xr with transposition */

            if (isVarComplex(pvApiCtx, piAddrB))
            {
                TransposeMatrix(pdblBI, mb, nb, pdblXI);    /* put b in x (with transposition) */
                for ( i = 0 ; i < mb ; i++ )
                {
                    umfpack_di_wsolve(UMFPACK_At, A.p, A.irow, A.R, &pdblBI[i * nb], &pdblXI[i * nb],
                                      Numeric, Control, Info, Wi, W);      /* the solutions are in bi */
                }
                TransposeMatrix(pdblBI, nb, mb, pdblXI);         /* put now bi in xi with transposition */
            }
        }
        else /*  A.it==1  */
        {
            TransposeMatrix(pdblBR, mb, nb, pdblXR);
            TransposeMatrix(pdblBI, mb, nb, pdblXI);
            for ( i = 0 ; i < mb ; i++ )
            {
                umfpack_zi_wsolve(UMFPACK_Aat, A.p, A.irow, A.R, A.I, &pdblBR[i * nb], &pdblBI[i * nb],
                                  &pdblXR[i * nb], &pdblXI[i * nb], Numeric, Control, Info, Wi, W);
            }
            TransposeMatrix(pdblBR, nb, mb, pdblXR);
            TransposeMatrix(pdblBI, nb, mb, pdblXI);
        }
    }

    if (A.it == 1)
    {
        umfpack_zi_free_numeric(&Numeric);
    }
    else
    {
        umfpack_di_free_numeric(&Numeric);
    }

    if (piNbItemRow != NULL)
    {
        FREE(piNbItemRow);
    }
    if (piColPos != NULL)
    {
        FREE(piColPos);
    }
    if (pdblSpReal != NULL)
    {
        FREE(pdblSpReal);
    }
    if (pdblSpImg != NULL)
    {
        FREE(pdblSpImg);
    }
    FREE(W);
    FREE(Wi);
    if (freepdblBI)
    {
        FREE(pdblBI);
    }
    freeCcsSparse(A);

    AssignOutputVariable(pvApiCtx, 1) = 4;
    ReturnArguments(pvApiCtx);
    return 0;
}
SciErr readCommonNamedSparseMatrix(void* _pvCtx, const char* _pstName, int _iComplex, int* _piRows, int* _piCols, int* _piNbItem, int* _piNbItemRow, int* _piColPos, double* _pdblReal, double* _pdblImg)
{
	SciErr sciErr; sciErr.iErr = 0; sciErr.iMsgCount = 0;
	int* piAddr		= NULL;
	int* piNbItemRow	= 0;
	int* piColPos		= 0;
	int iOne		= 1;

	double* pdblReal	= NULL;
	double* pdblImg		= NULL;

	sciErr = getVarAddressFromName(_pvCtx, _pstName, &piAddr);
	if(sciErr.iErr)
	{
		addErrorMessage(&sciErr, API_ERROR_READ_NAMED_SPARSE, _("%s: Unable to get variable \"%s\""), _iComplex ? "readNamedComplexSparseMatrix" : "readNamedSparseMatrix", _pstName);
		return sciErr;
	}
	
	if(_iComplex == 1)
	{
		sciErr = getComplexSparseMatrix(_pvCtx, piAddr, _piRows, _piCols, _piNbItem, &piNbItemRow, &piColPos, &pdblReal, &pdblImg);
	}
	else
	{
		sciErr = getSparseMatrix(_pvCtx, piAddr, _piRows, _piCols, _piNbItem, &piNbItemRow, &piColPos, &pdblReal);
	}

	if(sciErr.iErr)
	{
		addErrorMessage(&sciErr, API_ERROR_READ_NAMED_SPARSE, _("%s: Unable to get variable \"%s\""), _iComplex ? "readNamedComplexSparseMatrix" : "readNamedSparseMatrix", _pstName);
		return sciErr;
	}

	if(_piNbItemRow == NULL)
	{
		return sciErr;
	}

	memcpy(_piNbItemRow, piNbItemRow, *_piRows * sizeof(int));

	if(_piColPos == NULL)
	{
		return sciErr;
	}

	memcpy(_piColPos, piColPos, *_piNbItem * sizeof(int));


	if(_pdblReal == NULL)
	{
		return sciErr;
	}

	C2F(dcopy)(_piNbItem, pdblReal, &iOne, _pdblReal, &iOne);

	if(_iComplex && _pdblImg)
	{
		C2F(dcopy)(_piNbItem, pdblImg, &iOne, _pdblImg, &iOne);
	}

	return sciErr;
}
Exemple #8
0
/*--------------------------------------------------------------------------*/
int sci_qp_solve(char *fname, unsigned long fname_len)
{
    SciErr sciErr;
    static int un = 1, deux = 2;
    // n : first dimension of Q
    // nbis : second dimension of Q (nbis is expected to be equal to n)
    static int n = 0, nbis = 0;
    static int unbis = 0;
    static int m = 0, next = 0;
    static int mbis = 0;
    static int pipo = 0;
    static int nact = 0;
    int r = 0;
    static int lw = 0,  k = 0;
    static SciSparse Sp;
    static int issparse = 0;
    double *work = NULL;

    int* piAddr1 = NULL;
    int* piAddr2 = NULL;
    int* piAddr3 = NULL;
    int* piAddr4 = NULL;
    int* piAddr5 = NULL;

    double* Q = NULL;
    double* C = NULL;
    double* p = NULL;
    double* b = NULL;
    int* me   = NULL;

    double* x = NULL;
    int* iact = NULL;
    int* iter = NULL;
    double* crval = NULL;
    int *ierr = NULL;


    /*   Check rhs and lhs   */
    CheckInputArgument(pvApiCtx, 5, 5) ;
    CheckOutputArgument(pvApiCtx, 1, 5) ;

    /*Warning this interface does not support arguments passed by reference */

    /* RhsVar: qp_solve(Q,p,C,b,me) */
    /*                1,2,3,4,5   */

    next = nbInputArgument(pvApiCtx) + 1;
    /*   Variable 1 (Q)   */
    //get variable address
    sciErr = getVarAddressFromPosition(pvApiCtx, 1, &piAddr1);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    // Retrieve a matrix of double at position 1.
    sciErr = getMatrixOfDouble(pvApiCtx, piAddr1, &n, &nbis, &Q);
    if (sciErr.iErr)
    {
        Scierror(202, _("%s: Wrong type for argument #%d: A real expected.\n"), fname, 1);
        printError(&sciErr, 0);
        return 1;
    }

    //CheckSquare
    if (n != nbis)
    {
        Scierror(999, _("%s: Wrong size for input argument #%d: A square matrix expected.\n"), fname, 1);
        return 1;
    }

    /*   Variable 2 (p)   */
    //get variable address
    sciErr = getVarAddressFromPosition(pvApiCtx, 2, &piAddr2);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    // Retrieve a matrix of double at position 2.
    sciErr = getMatrixOfDouble(pvApiCtx, piAddr2, &nbis, &unbis, &p);
    if (sciErr.iErr)
    {
        Scierror(202, _("%s: Wrong type for argument #%d: A real expected.\n"), fname, 2);
        printError(&sciErr, 0);
        return 1;
    }

    //CheckLength
    if (nbis * unbis != n)
    {
        Scierror(999, _("%s: Wrong size for input argument #%d: %d expected.\n"), fname, 2, nbis * unbis);
        return 1;
    }


    /*   Variable 3 (C)   */
    issparse =  (checkInputArgumentType(pvApiCtx, 3, 5));
    //get variable address
    sciErr = getVarAddressFromPosition(pvApiCtx, 3, &piAddr3);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    if (!issparse)
    {
        // Retrieve a matrix of double at position 3.
        sciErr = getMatrixOfDouble(pvApiCtx, piAddr3, &nbis, &m, &C);
        if (sciErr.iErr)
        {
            Scierror(202, _("%s: Wrong type for argument #%d: A real expected.\n"), fname, 3);
            printError(&sciErr, 0);
            return 1;
        }
    }
    else
    {
        if (isVarComplex(pvApiCtx, piAddr3))
        {
            Sp.it = 1;
            sciErr = getComplexSparseMatrix(pvApiCtx, piAddr3, &(Sp.m), &(Sp.n), &(Sp.nel), &(Sp.mnel), &(Sp.icol), &(Sp.R), &(Sp.I));
        }
        else
        {
            Sp.it = 0;
            sciErr = getSparseMatrix(pvApiCtx, piAddr3, &(Sp.m), &(Sp.n), &(Sp.nel), &(Sp.mnel), &(Sp.icol), &(Sp.R));
        }

        if (sciErr.iErr)
        {
            printError(&sciErr, 0);
            return 1;
        }

        nbis = Sp.m;
        m = Sp.n;
    }

    if ( nbis != n ) // car C est passee en transposee dans la macro qpsolve
    {
        Scierror(999, _("%s: Wrong size for input argument #%d: %d column(s) expected for matrix %s.\n"), fname, 3, n, "C");
        return 0;
    }

    /*   Variable 4 (b)   */
    //get variable address
    sciErr = getVarAddressFromPosition(pvApiCtx, 4, &piAddr4);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    // Retrieve a matrix of double at position 4.
    sciErr = getMatrixOfDouble(pvApiCtx, piAddr4, &mbis, &unbis, &b);
    if (sciErr.iErr)
    {
        Scierror(202, _("%s: Wrong type for argument #%d: A real expected.\n"), fname, 4);
        printError(&sciErr, 0);
        return 1;
    }

    //CheckLength
    if (mbis * unbis != m)
    {
        Scierror(999, _("%s: Wrong size for input argument #%d: %d expected.\n"), fname, 4, mbis * unbis);
        return 1;
    }


    /*   Variable 5 (me)   */
    //get variable address
    sciErr = getVarAddressFromPosition(pvApiCtx, 5, &piAddr5);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    // Retrieve a matrix of double at position 5.
    sciErr = getMatrixOfDoubleAsInteger(pvApiCtx, piAddr5, &pipo, &unbis, &me);
    if (sciErr.iErr)
    {
        Scierror(202, _("%s: Wrong type for argument #%d: A real expected.\n"), fname, 5);
        printError(&sciErr, 0);
        return 1;
    }

    //CheckScalar
    if (pipo != 1 || unbis != 1)
    {
        Scierror(999, _("%s: Wrong size for input argument #%d: A real scalar expected.\n"), fname, 5);
        return 1;
    }

    if ((*(me) < 0) || (*(me) > n))
    {
        Scierror(999, _("%s: Wrong value for input argument #%d: %s must be an integer in the range 0 to %d.\n"), fname, 5, "me", n);
        return 0;
    }

    /* nbOutputArgument(pvApiCtx) variables: x, iact, iter, crval */
    next = Rhs;
    sciErr = allocMatrixOfDouble(pvApiCtx, next + 1, n, un, &x);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        Scierror(999, _("%s: Memory allocation error.\n"), fname);
        return 1;
    }

    sciErr = allocMatrixOfDoubleAsInteger(pvApiCtx, next + 2, m, un, &iact);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        Scierror(999, _("%s: Memory allocation error.\n"), fname);
        return 1;
    }

    sciErr = allocMatrixOfDoubleAsInteger(pvApiCtx, next + 3, deux, un, &iter);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        Scierror(999, _("%s: Memory allocation error.\n"), fname);
        return 1;
    }

    sciErr = allocMatrixOfDouble(pvApiCtx, next + 4, un, un, &crval);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        Scierror(999, _("%s: Memory allocation error.\n"), fname);
        return 1;
    }

    sciErr = allocMatrixOfDoubleAsInteger(pvApiCtx, next + 5, un, un, &ierr);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        Scierror(999, _("%s: Memory allocation error.\n"), fname);
        return 1;
    }


    r = Min(n, m);
    lw =  2 * n + r * (r + 5) / 2 + 2 * m + 1;
    if ((work = (double *)MALLOC(lw * sizeof(double))) == NULL)
    {
        Scierror(999, _("%s: Cannot allocate more memory.\n"), fname);
        return 1;
    }
    /* change the sign of  C and b.*/
    *ierr = 0;
    if (!issparse)
    {
        /* linear constraints matrix is stored full */
        C2F(qpgen2)((Q), (p), &n, &n,  (x), (crval), (C),
                    (b), &n, &m, (me), (iact), &nact, (iter), work,
                    (ierr));
    }
    else
    {
        /* linear constraints matrix is a sparse matrix */
        /* Change the linear constraints matrix representation:
        qpgen1sci requires column-compressed sparse matrix internal
        representation while Scilab sparse matrices are row-compressed */
        double *R = NULL, *I = NULL;
        int *ind = NULL;

        if ((R = (double *)MALLOC(Sp.nel * sizeof(double))) == NULL)
        {
            FREE(work);
            work = NULL;
            Scierror(999, _("%s: Cannot allocate more memory.\n"), fname);
            return 1;
        }
        if ((ind = (int *)MALLOC((m + Sp.nel) * sizeof(int))) == NULL)
        {
            FREE(work);
            work = NULL;
            FREE(R);
            R = NULL;
            Scierror(999, _("%s: Cannot allocate more memory.\n"), fname);
            return 1;
        }

        // Transpose the sparse matrix A
        C2F(spt)(&n, &m, &(Sp.nel) ,  &(Sp.it), (int *)work,
                 Sp.R,  Sp.I,  Sp.mnel,  Sp.icol, R, I, ind, ind + m);

        C2F(qpgen1sci)((Q), (p), &n, &n,  (x), (crval),
                       ind, ind + m,  R,
                       (b), &m, (me), (iact), &nact, (iter),
                       work, (ierr));
        FREE(work);
        work = NULL;
        FREE(R);
        R = NULL;
        FREE(ind);
        ind = NULL;
    }

    for (k = nact; k < m; k++)
    {
        (iact)[k] = 0;
    }
    /* LhsVar: [x, iact, iter, f] = qp_solve(...) */
    if (Lhs != 5)
    {
        if (*ierr == 0)
        {
            for (k = 0; k < Lhs; k++)
            {
                AssignOutputVariable(pvApiCtx, 1 + k) = next + 1 + k;
            }
            ReturnArguments(pvApiCtx);
        }
        else if (*ierr == 1)
        {
            Scierror(999, _("%s: The minimization problem has no solution.\n"), fname);
        }
        else if (*ierr == 2)
        {
            Scierror(999, _("%s: Q is not symmetric positive definite.\n"), fname);
        }
    }
    else
    {
        for (k = 0; k < Lhs; k++)
        {
            AssignOutputVariable(pvApiCtx, 1 + k) = next + 1 + k;
        }
        if (*ierr == 1)
        {
            if (getWarningMode())
            {
                sciprint(_("\n%s: Warning: The minimization problem has no solution. The results may be inaccurate.\n\n"), fname);
            }
        }
        else if (*ierr == 2)
        {
            if (getWarningMode())
            {
                sciprint(_("\n%s: Warning: Q is not symmetric positive definite. The results may be inaccurate.\n\n"), fname);
            }
        }
        ReturnArguments(pvApiCtx);
    }
    return 0;
}
Exemple #9
0
int sci_res_with_prec(char* fname, void* pvApiCtx)
{
    SciErr sciErr;

    int mx = 0, nx = 0, mb = 0, nb = 0, i = 0;

    SciSparse A;
    int mA              = 0; // rows
    int nA              = 0; // cols
    int iNbItem         = 0;
    int* piNbItemRow    = NULL;
    int* piColPos       = NULL;
    double* pdblSpReal  = NULL;
    double* pdblSpImg   = NULL;

    int iComplex = 0;

    int* piAddr1 = NULL;
    int* piAddr2 = NULL;
    int* piAddr3 = NULL;

    double* pdblXR = NULL;
    double* pdblXI = NULL;
    double* pdblBR = NULL;
    double* pdblBI = NULL;
    double* pdblNR = NULL;
    double* pdblNI = NULL;
    double* pdblRR = NULL;
    double* pdblRI = NULL;

    int nbInputArg = nbInputArgument(pvApiCtx);

    /* Check numbers of input/output arguments */
    CheckInputArgument(pvApiCtx, 3, 3);
    CheckOutputArgument(pvApiCtx, 1, 2);

    /* get A the sparse matrix */
    sciErr = getVarAddressFromPosition(pvApiCtx, 1, &piAddr1);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    if (isVarComplex(pvApiCtx, piAddr1))
    {
        iComplex = 1;
        sciErr = getComplexSparseMatrix(pvApiCtx, piAddr1, &mA, &nA, &iNbItem, &piNbItemRow, &piColPos, &pdblSpReal, &pdblSpImg);
    }
    else
    {
        sciErr = getSparseMatrix(pvApiCtx, piAddr1, &mA, &nA, &iNbItem, &piNbItemRow, &piColPos, &pdblSpReal);
    }

    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    // fill struct sparse
    A.m     = mA;
    A.n     = nA;
    A.it    = iComplex;
    A.nel   = iNbItem;
    A.mnel  = piNbItemRow;
    A.icol  = piColPos;
    A.R     = pdblSpReal;
    A.I     = pdblSpImg;

    /* get x */
    sciErr = getVarAddressFromPosition(pvApiCtx, 2, &piAddr2);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    if (isVarComplex(pvApiCtx, piAddr2))
    {
        iComplex = 1;
        sciErr = getComplexMatrixOfDouble(pvApiCtx, piAddr2, &mx, &nx, &pdblXR, &pdblXI);
    }
    else
    {
        sciErr = getMatrixOfDouble(pvApiCtx, piAddr2, &mx, &nx, &pdblXR);
    }

    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    /* get b */
    sciErr = getVarAddressFromPosition(pvApiCtx, 3, &piAddr3);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    if (isVarComplex(pvApiCtx, piAddr3))
    {
        iComplex = 1;
        sciErr = getComplexMatrixOfDouble(pvApiCtx, piAddr3, &mb, &nb, &pdblBR, &pdblBI);
    }
    else
    {
        sciErr = getMatrixOfDouble(pvApiCtx, piAddr3, &mb, &nb, &pdblBR);
    }

    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    /* check size of inputs */
    if ( nx < 1 || nb != nx || mx != nA || mb != mA )
    {
        Scierror(999, _("%s: Wrong size for input arguments: Same sizes expected.\n"), fname);
        return 1;
    }

    /* Create the matrix as return of the function */
    if (iComplex)
    {
        sciErr = allocComplexMatrixOfDouble(pvApiCtx, 4, mb, nb, &pdblRR, &pdblRI);
    }
    else
    {
        sciErr = allocMatrixOfDouble(pvApiCtx, 4, mb, nb, &pdblRR);
    }

    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    /* Create the matrix as return of the function */
    sciErr = allocMatrixOfDouble(pvApiCtx, 5, 1, nb, &pdblNR);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    /* perform operations */
    if (iComplex == 0)
    {
        for ( i = 0 ; i < nb ; i++ )
        {
            residu_with_prec(&A, pdblXR + i * mx, pdblBR + i * mb, pdblRR + i * mb, pdblNR + i);
        }
    }
    else
    {
        if (pdblXI == NULL)
        {
            int iSize = mx * nx * sizeof(double);
            pdblXI = (double*)MALLOC(iSize);
            memset(pdblXI, 0x00, iSize);
        }

        if (pdblBI == NULL)
        {
            int iSize = mb * nb * sizeof(double);
            pdblBI = (double*)MALLOC(iSize);
            memset(pdblBI, 0x00, iSize);
        }

        if (pdblSpImg == NULL)
        {
            /* Create the matrix as return of the function */
            int iSize = nb * sizeof(double);
            pdblNI = (double*)MALLOC(iSize);
            memset(pdblNI, 0x00, iSize);

            for ( i = 0 ; i < nb ; i++ )
            {
                residu_with_prec(&A, pdblXR + i * mx, pdblBR + i * mb, pdblRR + i * mb, pdblNR + i);
            }

            for ( i = 0 ; i < nb ; i++ )
            {
                residu_with_prec(&A, pdblXI + i * mx, pdblBI + i * mb, pdblRI + i * mb, pdblNI + i);
            }

            for ( i = 0 ; i < nb ; i++ )
            {
                pdblNR[i] = sqrt(pdblNR[i] * pdblNR[i] + pdblNI[i] * pdblNI[i]);
            }
        }
        else
        {
            for ( i = 0 ; i < nb ; i++ )
            {
                cmplx_residu_with_prec(&A,  pdblXR + i * mx, pdblXI + i * mx,
                                       pdblBR + i * mb, pdblBI + i * mb,
                                       pdblRR + i * mb, pdblRI + i * mb,
                                       pdblNR + i);
            }
        }
    }

    if (isVarComplex(pvApiCtx, piAddr1) == 0)
    {
        FREE(pdblNI);
    }

    if (isVarComplex(pvApiCtx, piAddr2) == 0)
    {
        FREE(pdblXI);
    }

    if (isVarComplex(pvApiCtx, piAddr3) == 0)
    {
        FREE(pdblBI);
    }

    AssignOutputVariable(pvApiCtx, 1) = 4;
    AssignOutputVariable(pvApiCtx, 2) = 5;
    ReturnArguments(pvApiCtx);
    return 0;
}
Exemple #10
0
int sci_umf_lusolve(char* fname, unsigned long l)
{
    SciErr sciErr;

    int mb      = 0;
    int nb      = 0;
    int it_flag = 0;
    int i       = 0;
    int j       = 0;

    int NoTranspose = 0;
    int NoRaffinement = 0;
    SciSparse AA;
    CcsSparse A;

    /* umfpack stuff */
    double Info[UMFPACK_INFO]; // double *Info = (double *) NULL;
    double Control[UMFPACK_CONTROL];
    void* Numeric = NULL;
    int lnz = 0, unz = 0, n = 0, n_col = 0, nz_udiag = 0, umf_flag = 0;
    int* Wi = NULL;
    int mW = 0;
    double *W = NULL;

    int iComplex = 0;

    int* piAddr1 = NULL;
    int* piAddr2 = NULL;
    int* piAddr3 = NULL;
    int* piAddr4 = NULL;

    double* pdblBR = NULL;
    double* pdblBI = NULL;
    double* pdblXR = NULL;
    double* pdblXI = NULL;

    int mA              = 0; // rows
    int nA              = 0; // cols
    int iNbItem         = 0;
    int* piNbItemRow    = NULL;
    int* piColPos       = NULL;
    double* pdblSpReal  = NULL;
    double* pdblSpImg   = NULL;

    /* Check numbers of input/output arguments */
    CheckInputArgument(pvApiCtx, 2, 4);
    CheckOutputArgument(pvApiCtx, 1, 1);

    /* First get arg #1 : the pointer to the LU factors */
    sciErr = getVarAddressFromPosition(pvApiCtx, 1, &piAddr1);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    sciErr = getPointer(pvApiCtx, piAddr1, &Numeric);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    /* Check if this pointer is a valid ref to a umfpack LU numeric object */
    if ( ! IsAdrInList(Numeric, ListNumeric, &it_flag) )
    {
        Scierror(999, _("%s: Wrong value for input argument #%d: Must be a valid reference to (umf) LU factors.\n"), fname, 1);
        return 1;
    }

    /*  get some parameters of the factorization (for some checking) */
    if ( it_flag == 0 )
    {
        umfpack_di_get_lunz(&lnz, &unz, &n, &n_col, &nz_udiag, Numeric);
    }
    else
    {
        iComplex = 1;
        umfpack_zi_get_lunz(&lnz, &unz, &n, &n_col, &nz_udiag, Numeric);
    }

    if ( n != n_col )
    {
        Scierror(999, _("%s: An error occurred: %s.\n"), fname, _("This is not a factorization of a square matrix"));
        return 1;
    }

    if ( nz_udiag < n )
    {
        Scierror(999, _("%s: An error occurred: %s.\n"), fname, _("This is a factorization of a singular matrix"));
        return 1;
    }

    /* Get now arg #2 : the vector b */
    sciErr = getVarAddressFromPosition(pvApiCtx, 2, &piAddr2);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    if (isVarComplex(pvApiCtx, piAddr2))
    {
        iComplex = 1;
        sciErr = getComplexMatrixOfDouble(pvApiCtx, piAddr2, &mb, &nb, &pdblBR, &pdblBI);
    }
    else
    {
        sciErr = getMatrixOfDouble(pvApiCtx, piAddr2, &mb, &nb, &pdblBR);
    }

    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    if (mb != n || nb < 1)    /* test if the right hand side is compatible */
    {
        Scierror(999, _("%s: Wrong size for input argument #%d.\n"), fname, 2);
        return 1;
    }

    /* allocate memory for the solution x */
    if (iComplex)
    {
        sciErr = allocComplexMatrixOfDouble(pvApiCtx, nbInputArgument(pvApiCtx) + 1, mb, nb, &pdblXR, &pdblXI);
    }
    else
    {
        sciErr = allocMatrixOfDouble(pvApiCtx, nbInputArgument(pvApiCtx) + 1, mb, nb, &pdblXR);
    }

    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    /*  selection between the different options :
     *   -- solving Ax=b or A'x=b (Note: we could add  A.'x=b)
     *   -- with or without raffinement
     */

    if (nbInputArgument(pvApiCtx) == 2)
    {
        NoTranspose = 1;
        NoRaffinement = 1;
    }
    else  /* 3 or 4 input arguments but the third must be a string */
    {
        char* pStr = NULL;
        sciErr = getVarAddressFromPosition(pvApiCtx, 3, &piAddr3);
        if (sciErr.iErr)
        {
            printError(&sciErr, 0);
            return 1;
        }

        getAllocatedSingleString(pvApiCtx, piAddr3, &pStr);
        if (strcmp(pStr, "Ax=b") == 0)
        {
            NoTranspose = 1;
        }
        else if ( strcmp(pStr, "A'x=b") == 0 )
        {
            NoTranspose = 0;
        }
        else
        {
            Scierror(999, _("%s: Wrong input argument #%d: '%s' or '%s' expected.\n"), fname, 3, "Ax=b", "A'x=b");
            return 1;
        }

        if (nbInputArgument(pvApiCtx) == 4)
        {
            sciErr = getVarAddressFromPosition(pvApiCtx, 4, &piAddr4);
            if (sciErr.iErr)
            {
                printError(&sciErr, 0);
                return 1;
            }

            if (isVarComplex(pvApiCtx, piAddr4))
            {
                AA.it = 1;
                sciErr = getComplexSparseMatrix(pvApiCtx, piAddr4, &mA, &nA, &iNbItem, &piNbItemRow, &piColPos, &pdblSpReal, &pdblSpImg);
            }
            else
            {
                AA.it = 0;
                sciErr = getSparseMatrix(pvApiCtx, piAddr4, &mA, &nA, &iNbItem, &piNbItemRow, &piColPos, &pdblSpReal);
            }

            if (sciErr.iErr)
            {
                printError(&sciErr, 0);
                return 1;
            }

            // fill struct sparse
            AA.m     = mA;
            AA.n     = nA;
            AA.nel   = iNbItem;
            AA.mnel  = piNbItemRow;
            AA.icol  = piColPos;
            AA.R     = pdblSpReal;
            AA.I     = pdblSpImg;

            /*  some check... but we can't be sure that the matrix corresponds to the LU factors */
            if ( mA != nA || mA != n || AA.it != it_flag )
            {
                Scierror(999, _("%s: Wrong size for input argument #%d: %s.\n"), fname, 4, _("Matrix is not compatible with the given LU factors"));
                return 1;
            }

            NoRaffinement = 0;
        }
        else
        {
            NoRaffinement = 1;   /* only 3 input var => no raffinement */
        }
    }

    /* allocate memory for umfpack_di_wsolve usage or umfpack_zi_wsolve usage*/
    Wi = (int*)MALLOC(n * sizeof(int));

    if (it_flag == 1)
    {
        if (NoRaffinement)
        {
            mW = 4 * n;
        }
        else
        {
            mW = 10 * n;
        }
    }
    else
    {
        if (NoRaffinement)
        {
            mW = n;
        }
        else
        {
            mW = 5 * n;
        }
    }

    W = (double*)MALLOC(mW * sizeof(double));

    if (NoRaffinement == 0)
    {
        SciSparseToCcsSparse(&AA, &A);
    }
    else
    {
        A.p = NULL;
        A.irow = NULL;
        A.R = NULL;
        A.I = NULL;
    }

    /* get the pointer for b */
    if (it_flag == 1  &&  pdblBI == NULL)
    {
        int iSize = mb * nb * sizeof(double);
        pdblBI = (double*)MALLOC(iSize);
        memset(pdblBI, 0x00, iSize);
    }

    /* init Control */
    if (it_flag == 0)
    {
        umfpack_di_defaults(Control);
    }
    else
    {
        umfpack_zi_defaults(Control);
    }

    if (NoRaffinement)
    {
        Control[UMFPACK_IRSTEP] = 0;
    }

    if (NoTranspose)
    {
        umf_flag = UMFPACK_A;
    }
    else
    {
        umf_flag = UMFPACK_At;
    }

    if (it_flag == 0)
    {
        for (j = 0; j < nb ; j++)
        {
            umfpack_di_wsolve(umf_flag, A.p, A.irow, A.R, &pdblXR[j * mb], &pdblBR[j * mb], Numeric, Control, Info, Wi, W);
        }

        if (iComplex == 1)
        {
            for (j = 0; j < nb ; j++)
            {
                umfpack_di_wsolve(umf_flag, A.p, A.irow, A.R, &pdblXI[j * mb], &pdblBI[j * mb], Numeric, Control, Info, Wi, W);
            }
        }
    }
    else
    {
        for (j = 0; j < nb ; j++)
        {
            umfpack_zi_wsolve(umf_flag, A.p, A.irow, A.R, A.I, &pdblXR[j * mb], &pdblXI[j * mb], &pdblBR[j * mb], &pdblBI[j * mb], Numeric, Control, Info, Wi, W);
        }
    }

    if (isVarComplex(pvApiCtx, piAddr2) == 0)
    {
        FREE(pdblBI);
    }

    freeCcsSparse(A);

    FREE(W);
    FREE(Wi);

    AssignOutputVariable(pvApiCtx, 1) = nbInputArgument(pvApiCtx) + 1;
    ReturnArguments(pvApiCtx);
    return 0;
}