int IsAdrInList(Adr adr, CellAdr *L, int *it_flag)
{
    /* teste si l'adresse adr est presente ds la liste L, si oui
       la fonction renvoie 1, sinon 0. On renvoit aussi it */

    if ( L == NULL )
    {
        return 0;
    }
    else if ( L->adr == adr )
    {
        *it_flag = L->it;
        return 1;
    }
    else
    {
        return ( IsAdrInList(adr, L->next, it_flag) );
    }
}
int sci_umf_luget(char* fname, void* pvApiCtx)
{
    /*
    *  LU_ptr is (a pointer to) a factorization of A, we have:
    *             -1
    *          P R  A Q = L U
    *
    *      A is n_row x n_col
    *      L is n_row x n
    *      U is n     x n_col     n = min(n_row, n_col)
    */

    SciErr sciErr;
    void* Numeric = NULL;
    int lnz = 0, unz = 0, n_row = 0, n_col = 0, n = 0, nz_udiag = 0, i = 0, stat = 0, do_recip = 0, it_flag = 0;
    int *L_mnel = NULL, *L_icol = NULL, *L_ptrow = NULL, *U_mnel = NULL, *U_icol = NULL, *U_ptrow = NULL, *V_irow = NULL, *V_ptcol = NULL;
    double *L_R = NULL, *L_I = NULL, *U_R = NULL, *U_I = NULL, *V_R = NULL, *V_I = NULL, *Rs = NULL;
    int *p = NULL, *q = NULL, pl_miss = 0, error_flag = 0 ;

    int* piAddr1 = NULL;
    int iType1   = 0;

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

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

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

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

    /* Check if the pointer is a valid ref to ... */
    if ( IsAdrInList(Numeric, ListNumeric, &it_flag) )
    {
        if (it_flag == 0 )
        {
            umfpack_di_get_lunz(&lnz, &unz, &n_row, &n_col, &nz_udiag, Numeric);
        }
        else
        {
            umfpack_zi_get_lunz(&lnz, &unz, &n_row, &n_col, &nz_udiag, Numeric);
        }
    }
    else
    {
        Scierror(999, _("%s: Wrong value for input argument #%d: Must be a valid reference to (umf) LU factors.\n"), fname, 1);
        return 1;
    }

    if (n_row <= n_col)
    {
        n = n_row;
    }
    else
    {
        n = n_col;
    }
    L_mnel  = (int*)MALLOC(n_row * sizeof(int));
    L_icol  = (int*)MALLOC(lnz * sizeof(int));
    L_ptrow = (int*)MALLOC((n_row + 1) * sizeof(int));
    L_R     = (double*)MALLOC( lnz * sizeof(double));
    U_mnel  = (int*)MALLOC(n * sizeof(int));
    U_icol  = (int*)MALLOC(unz * sizeof(int));
    U_ptrow = (int*)MALLOC((n + 1) * sizeof(int));
    U_R     = (double*)MALLOC( unz * sizeof(double));
    V_irow  = (int*)MALLOC(unz * sizeof(int));
    V_ptcol = (int*)MALLOC((n_col + 1) * sizeof(int));
    V_R     = (double*)MALLOC( unz * sizeof(double));
    p       = (int*)MALLOC(n_row * sizeof(int));
    q       = (int*)MALLOC(n_col * sizeof(int));
    Rs      = (double*)MALLOC(n_row * sizeof(double));

    if ( it_flag == 1 )
    {
        L_I = (double*)MALLOC(lnz * sizeof(double));
        U_I = (double*)MALLOC(unz * sizeof(double));
        V_I = (double*)MALLOC(unz * sizeof(double));
    }
    else
    {
        L_I = U_I = V_I = NULL;
    }

    if (    !(L_mnel && L_icol && L_R && L_ptrow  && p &&
              U_mnel && U_icol && U_R && U_ptrow  && q &&
              V_irow && V_R && V_ptcol  && Rs)
            || (it_flag && !(L_I && U_I && V_I))   )
    {
        error_flag = 1;
        goto the_end;
    }

    if ( it_flag == 0 )
    {
        stat = umfpack_di_get_numeric(L_ptrow, L_icol, L_R, V_ptcol, V_irow, V_R,
                                      p, q, (double *)NULL, &do_recip, Rs, Numeric);
    }
    else
    {
        stat = umfpack_zi_get_numeric(L_ptrow, L_icol, L_R, L_I, V_ptcol, V_irow, V_R, V_I,
                                      p, q, (double *)NULL, (double *)NULL, &do_recip, Rs, Numeric);
    }

    if ( stat != UMFPACK_OK )
    {
        error_flag = 2;
        goto the_end;
    };

    if ( do_recip )
    {
        for ( i = 0 ; i < n_row ; i++ )
        {
            Rs[i] = 1.0 / Rs[i];
        }
    }

    if ( it_flag == 0 )
    {
        stat = umfpack_di_transpose(n, n_col, V_ptcol, V_irow, V_R, (int *) NULL,
                                    (int*) NULL, U_ptrow, U_icol, U_R);
    }
    else
    {
        stat = umfpack_zi_transpose(n, n_col, V_ptcol, V_irow, V_R, V_I, (int *) NULL,
                                    (int*) NULL, U_ptrow, U_icol, U_R, U_I, 0);
    }

    if ( stat != UMFPACK_OK )
    {
        error_flag = 2;
        goto the_end;
    };

    for ( i = 0 ; i < n_row ; i++ )
    {
        L_mnel[i] = L_ptrow[i + 1] - L_ptrow[i];
    }
    for ( i = 0 ; i < n ; i++ )
    {
        U_mnel[i] = U_ptrow[i + 1] - U_ptrow[i];
    }

    for ( i = 0 ; i < lnz ; i++ )
    {
        L_icol[i]++;
    }
    for ( i = 0 ; i < unz ; i++ )
    {
        U_icol[i]++;
    }

    for ( i = 0 ; i < n_row ; i++ )
    {
        p[i]++;
    }
    for ( i = 0 ; i < n_col ; i++ )
    {
        q[i]++;
    }

    /* output L */
    if (it_flag) // complex
    {
        sciErr = createComplexSparseMatrix(pvApiCtx, 2, n_row, n, lnz, L_mnel, L_icol, L_R, L_I);
    }
    else
    {
        sciErr = createSparseMatrix(pvApiCtx, 2, n_row, n, lnz, L_mnel, L_icol, L_R);
    }

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

    /* output U */
    if (it_flag) // complex
    {
        sciErr = createComplexSparseMatrix(pvApiCtx, 3, n, n_col, unz, U_mnel, U_icol, U_R, U_I);
    }
    else
    {
        sciErr = createSparseMatrix(pvApiCtx, 3, n, n_col, unz, U_mnel, U_icol, U_R);
    }

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

    /* output p */
    sciErr = createMatrixOfDoubleAsInteger(pvApiCtx, 4, n_row, 1, p);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        FREE(L_mnel);
        FREE(U_mnel);
        return 1;
    }

    /* output q */
    sciErr = createMatrixOfDoubleAsInteger(pvApiCtx, 5, n_col, 1, q);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        FREE(L_mnel);
        FREE(U_mnel);
        return 1;
    }

    /* output res */
    sciErr = createMatrixOfDouble(pvApiCtx, 6, n_row, 1, Rs);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        FREE(L_mnel);
        FREE(U_mnel);
        return 1;
    }

the_end:
    FREE(L_mnel);
    FREE(L_icol);
    FREE(L_R);
    FREE(L_ptrow);
    FREE(p);
    FREE(U_mnel);
    FREE(U_icol);
    FREE(U_R);
    FREE(U_ptrow);
    FREE(q);
    FREE(V_irow);
    FREE(V_R);
    FREE(V_ptcol);
    FREE(Rs);

    if ( it_flag == 1 )
    {
        FREE(L_I);
        FREE(V_I);
        FREE(U_I);
    }

    switch (error_flag)
    {
        case 0:   /* no error */
            AssignOutputVariable(pvApiCtx, 1) = 2;
            AssignOutputVariable(pvApiCtx, 2) = 3;
            AssignOutputVariable(pvApiCtx, 3) = 4;
            AssignOutputVariable(pvApiCtx, 4) = 5;
            AssignOutputVariable(pvApiCtx, 5) = 6;
            ReturnArguments(pvApiCtx);
            return 0;

        case 1:   /* enough memory (with malloc) */
            Scierror(999, _("%s: No more memory.\n"), fname);
            break;

        case 2:   /* a problem with one umfpack routine */
            Scierror(999, "%s: %s\n", fname, UmfErrorMes(stat));
            break;
    }

    return 1;
}
Exemple #3
0
void DialogDasm::RefreshDasm()
{

    if ( pPC->pCPU->pDEBUG->debugged)
    {
        QString	text;

        if (! IsAdrInList(pPC->pCPU->pDEBUG->DasmAdr))
        {
            if (MaxAdr > pPC->pCPU->pDEBUG->DasmAdr)
            {
                // effacer tout et recommencer au debut
                ui->codelistWidget->clear();

                Index		= 0;
                MaxAdr		= pPC->pCPU->pDEBUG->DasmAdr;
                NextMaxAdr	= pPC->pCPU->pDEBUG->NextDasmAdr;
            }
            else
            {
                if (pPC->pCPU->pDEBUG->DasmAdr > NextMaxAdr)
                {
                    // Insert a separator
                    text = QString("%1:---------------------").arg(pPC->pCPU->pDEBUG->DasmAdr-1,5,16,QChar('0'));
                    QListWidgetItem *item = new QListWidgetItem(text);
                    int adr = pPC->pCPU->pDEBUG->DasmAdr;
                    item->setData(Qt::UserRole,QVariant(adr));
                    ui->codelistWidget->addItem(item);
                    selectRow(ui->codelistWidget->count()-1);

                    NextMaxAdr = pPC->pCPU->pDEBUG->NextDasmAdr;
                    Index++;
                }
                MaxAdr		= pPC->pCPU->pDEBUG->DasmAdr;
                NextMaxAdr	= pPC->pCPU->pDEBUG->NextDasmAdr;
            }

            QListWidgetItem *item = new QListWidgetItem(pPC->pCPU->pDEBUG->Buffer);
            int adr = pPC->pCPU->pDEBUG->DasmAdr;
            item->setData(Qt::UserRole,adr);
            ui->codelistWidget->addItem(item);
            selectRow(ui->codelistWidget->count()-1);
            Index++;
// full until 15 lines
            for (int j=Index;j<15;j++)
            {
                pPC->pCPU->pDEBUG->DisAsm_1(NextMaxAdr);
                MaxAdr		= pPC->pCPU->pDEBUG->DasmAdr;
                NextMaxAdr	= pPC->pCPU->pDEBUG->NextDasmAdr;
                QListWidgetItem *item = new QListWidgetItem(pPC->pCPU->pDEBUG->Buffer);
                int adr = pPC->pCPU->pDEBUG->DasmAdr;
                item->setData(Qt::UserRole,adr);

                ui->codelistWidget->addItem(item);
                Index++;
            }
        }

    }
    if (regwidget) regwidget->refresh();
    loadImem();
    loadMem();

}
int sci_taucs_chsolve(char* fname, void* pvApiCtx)
{
    SciErr sciErr;

    int mb = 0, nb = 0;
    int i = 0, j = 0, n = 0, it_flag = 0, Refinement = 0;
    double norm_res = 0., norm_res_bis = 0.;
    long double *wk = NULL;
    int A_is_upper_triangular = 0;
    taucs_handle_factors * pC = NULL;

    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;

    void* pvPtr     = NULL;
    double* pdblB   = NULL;
    double* pdblX   = NULL;
    double* pdblV   = NULL;
    double* pdblRes = NULL;

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

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

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

    pC = (taucs_handle_factors *)pvPtr;

    /* Check if this pointer is a valid ref to a Cholesky factor object */
    if ( ! IsAdrInList( (Adr)pC, ListCholFactors, &it_flag) )
    {
        Scierror(999, _("%s: Wrong value for input argument #%d: not a valid reference to Cholesky factors"), fname, 1);
        return 1;
    }

    /*  the number of rows/lines of the matrix  */
    n = pC->n;
    /* Get now arg #2 : the vector b */
    sciErr = getVarAddressFromPosition(pvApiCtx, 2, &piAddr2);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    sciErr = getMatrixOfDouble(pvApiCtx, piAddr2, &mb, &nb, &pdblB);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

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

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

        if (isVarComplex(pvApiCtx, piAddr3))
        {
            Scierror(999, _("%s: Wrong type for input argument #%d: not compatible with the Cholesky factorization.\n"), fname, 3);
            return 1;
        }

        sciErr = getSparseMatrix(pvApiCtx, piAddr3, &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;

        if (mA != nA || mA != n)
        {
            Scierror(999, _("%s: Wrong size for input argument #%d: not compatible with the Cholesky factorization.\n"), fname, 3);
            return 1;
        }

        Refinement = 1;
        A_is_upper_triangular = is_sparse_upper_triangular(&A);
    }
    else
    {
        Refinement = 0;
    }

    /* allocate memory for the solution x */
    sciErr = allocMatrixOfDouble(pvApiCtx, nbInputArgument(pvApiCtx) + 1, mb, nb, &pdblX);
    if (sciErr.iErr)
    {
        printError(&sciErr, 0);
        return 1;
    }

    if (Refinement)
    {
        pdblRes = (double*)MALLOC(mb * sizeof(double));
        if ( A_is_upper_triangular )
        {
            if ( (wk = (long double*)MALLOC( n * sizeof(long double))) == NULL )
            {
                if (pdblRes)
                {
                    FREE(pdblRes);
                }
                Scierror(999, _("%s: not enough memory.\n"), fname);
                return 1;
            }
        }
    }
    
    /* allocate memory for a temporary vector v */
    pdblV = (double*)MALLOC(mb * sizeof(double));

    for (j = 0; j < nb ; j++)
    {
        taucs_vec_permute(n, &pdblB[j * mb], &pdblX[j * mb], pC->p);
        taucs_supernodal_solve_llt(pC->C, pdblV, &pdblX[j * mb]); /* FIXME : add a test here */
        taucs_vec_ipermute(n, pdblV, &pdblX[j * mb], pC->p);
        if (Refinement)
        {
            /* do one iterative refinement */
            residu_with_prec_for_chol(&A, &pdblX[j * mb], &pdblV[j * mb], pdblRes, &norm_res, A_is_upper_triangular, wk);
            /*  FIXME: do a test if the norm_res has an anormal value and send a warning
             *         (the user has certainly not give the good matrix A
             */
            taucs_vec_permute(n, pdblRes, pdblV, pC->p);
            taucs_supernodal_solve_llt(pC->C, pdblRes, pdblV);  /* FIXME : add a test here */
            taucs_vec_ipermute(n, pdblRes, pdblV, pC->p);
            for ( i = 0 ; i < n ; i++ )
            {
                pdblV[i] = pdblX[j * mb + i] - pdblV[i]; /* v is the refined solution */
            }

            residu_with_prec_for_chol(&A, pdblV, &pdblB[j * mb], pdblRes, &norm_res_bis, A_is_upper_triangular, wk);
            /* accept it if the 2 norm of the residual is improved */
            if ( norm_res_bis < norm_res )
            {
                for ( i = 0 ; i < n ; i++ )
                {
                    pdblX[j * mb + i] = pdblV[i];
                }
            }
        }
    }

    FREE(wk);
    FREE(pdblV);
    FREE(pdblRes);

    AssignOutputVariable(pvApiCtx, 1) = nbInputArgument(pvApiCtx) + 1;
    ReturnArguments(pvApiCtx);
    return 0;
}
Exemple #5
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;
}