void mexFunction
(
    int nargout,
    mxArray *pargout [ ],
    int nargin,
    const mxArray *pargin [ ]
)
{
    if (nargout > 1 || nargin < 2 || nargin > 4)
    {
        mexErrMsgTxt ("Usage: C = cs_add(A,B,alpha,beta)") ;
    }
    if (mxIsComplex (pargin [0]) || mxIsComplex (pargin [1])
        || (nargin > 2 && mxIsComplex (pargin [2]))
        || (nargin > 3 && mxIsComplex (pargin [3])))
    {
#ifndef NCOMPLEX
        cs_complex_t alpha, beta ;
        cs_cl Amatrix, Bmatrix, *A, *B, *C, *D ;
        A = cs_cl_mex_get_sparse (&Amatrix, 0, pargin [0]) ;    /* get A */
        B = cs_cl_mex_get_sparse (&Bmatrix, 0, pargin [1]) ;    /* get B */
        alpha = (nargin < 3) ? 1 : get_complex (pargin [2]) ;   /* get alpha */
        beta  = (nargin < 4) ? 1 : get_complex (pargin [3]) ;   /* get beta */
        C = cs_cl_add (A,B,alpha,beta) ;    /* C = alpha*A + beta *B */
        cs_cl_dropzeros (C) ;           /* drop zeros */
        D = cs_cl_transpose (C, 1) ;    /* sort result via double transpose */
        cs_cl_spfree (C) ;
        C = cs_cl_transpose (D, 1) ;
        cs_cl_spfree (D) ;
        pargout [0] = cs_cl_mex_put_sparse (&C) ;       /* return C */
#else
        mexErrMsgTxt ("complex matrices not supported") ;
#endif
    }
    else
    {
        double alpha, beta ;
        cs_dl Amatrix, Bmatrix, *A, *B, *C, *D ;
        A = cs_dl_mex_get_sparse (&Amatrix, 0, 1, pargin [0]) ;    /* get A */
        B = cs_dl_mex_get_sparse (&Bmatrix, 0, 1, pargin [1]) ;    /* get B */
        alpha = (nargin < 3) ? 1 : mxGetScalar (pargin [2]) ;   /* get alpha */
        beta  = (nargin < 4) ? 1 : mxGetScalar (pargin [3]) ;   /* get beta */
        C = cs_dl_add (A,B,alpha,beta) ;        /* C = alpha*A + beta *B */
        cs_dl_dropzeros (C) ;           /* drop zeros */
        D = cs_dl_transpose (C, 1) ;    /* sort result via double transpose */
        cs_dl_spfree (C) ;
        C = cs_dl_transpose (D, 1) ;
        cs_dl_spfree (D) ;
        pargout [0] = cs_dl_mex_put_sparse (&C) ;       /* return C */
    }
}
Exemple #2
0
int main (void)
{
    cs_dl *T, *A, *Eye, *AT, *C, *D ;
    UF_long i, m ;
    T = cs_dl_load (stdin) ;               /* load triplet matrix T from stdin */
    printf ("T:\n") ; cs_dl_print (T, 0) ; /* print T */
    A = cs_dl_compress (T) ;               /* A = compressed-column form of T */
    printf ("A:\n") ; cs_dl_print (A, 0) ; /* print A */
    cs_dl_spfree (T) ;                     /* clear T */
    AT = cs_dl_transpose (A, 1) ;          /* AT = A' */
    printf ("AT:\n") ; cs_dl_print (AT, 0) ; /* print AT */
    m = A ? A->m : 0 ;                  /* m = # of rows of A */
    T = cs_dl_spalloc (m, m, m, 1, 1) ;    /* create triplet identity matrix */
    for (i = 0 ; i < m ; i++) cs_dl_entry (T, i, i, 1) ;
    Eye = cs_dl_compress (T) ;             /* Eye = speye (m) */
    cs_dl_spfree (T) ;
    C = cs_dl_multiply (A, AT) ;           /* C = A*A' */
    D = cs_dl_add (C, Eye, 1, cs_dl_norm (C)) ;   /* D = C + Eye*norm (C,1) */
    printf ("D:\n") ; cs_dl_print (D, 0) ; /* print D */
    cs_dl_spfree (A) ;                     /* clear A AT C D Eye */
    cs_dl_spfree (AT) ;
    cs_dl_spfree (C) ;
    cs_dl_spfree (D) ;
    cs_dl_spfree (Eye) ;
    return (0) ;
}
Exemple #3
0
void mexFunction
(
    int nargout,
    mxArray *pargout [ ],
    int nargin,
    const mxArray *pargin [ ]
)
{
    CS_INT n, nel, s ;
    cs *A, *AT ;
    if (nargout > 1 || nargin != 3)
    {
        mexErrMsgTxt ("Usage: C = cs_frand(n,nel,s)") ;
    }
    n = mxGetScalar (pargin [0]) ;
    nel = mxGetScalar (pargin [1]) ;
    s = mxGetScalar (pargin [2]) ;

    n = CS_MAX (1,n) ;
    nel = CS_MAX (1,nel) ;
    s = CS_MAX (1,s) ;

    AT = cs_dl_frand (n, nel, s) ;
    A = cs_dl_transpose (AT, 1) ;
    cs_dl_spfree (AT) ;
    cs_dl_dropzeros (A) ;

    pargout [0] = cs_dl_mex_put_sparse (&A) ;
}
Exemple #4
0
/* C = cs_transpose (A), computes C=A', where A must be sparse.
   C = cs_transpose (A,kind) computes C=A.' if kind <= 0, C=A' if kind > 0 */
void mexFunction
(
    int nargout,
    mxArray *pargout [ ],
    int nargin,
    const mxArray *pargin [ ]
)
{
    CS_INT values ;
    if (nargout > 1 || nargin < 1 || nargin > 2)
    {
        mexErrMsgTxt ("Usage: C = cs_transpose(A,kind)") ;
    }
    values = (nargin > 1) ? mxGetScalar (pargin [1]) : 1 ;
    values = (values <= 0) ? -1 : 1 ;
    if (mxIsComplex (pargin [0]))
    {
#ifndef NCOMPLEX
        cs_cl Amatrix, *A, *C ;
        A = cs_cl_mex_get_sparse (&Amatrix, 0, pargin [0]) ;    /* get A */
        C = cs_cl_transpose (A, values) ;                       /* C = A' */
        pargout [0] = cs_cl_mex_put_sparse (&C) ;               /* return C */
#else
        mexErrMsgTxt ("complex matrices not supported") ;
#endif
    }
    else
    {
        cs_dl Amatrix, *A, *C ;
        A = cs_dl_mex_get_sparse (&Amatrix, 0, 1, pargin [0]) ; /* get A */
        C = cs_dl_transpose (A, values) ;                       /* C = A' */
        pargout [0] = cs_dl_mex_put_sparse (&C) ;               /* return C */
    }
}
Exemple #5
0
/* cs_sparse: convert triplet form into compress-column form sparse matrix */
void mexFunction
(
    int nargout,
    mxArray *pargout [ ],
    int nargin,
    const mxArray *pargin [ ]
)
{
    if (nargout > 1 || nargin != 3)
    {
        mexErrMsgTxt ("Usage: A = cs_sparse(i,j,x)") ;
    }
    if (mxIsComplex (pargin [2]))
    {
#ifndef NCOMPLEX
        cs_cl *A, *C, *T, Tmatrix ;
        T = &Tmatrix ;                  /* get i,j,x and copy to triplet form */
        T->nz = mxGetM (pargin [0]) ;
        T->p = cs_dl_mex_get_int (T->nz, pargin [0], &(T->n), 1) ;
        T->i = cs_dl_mex_get_int (T->nz, pargin [1], &(T->m), 1) ;
        cs_mex_check (1, T->nz, 1, 0, 0, 1, pargin [2]) ;
        T->x = cs_cl_mex_get_double (T->nz, pargin [2]) ;
        T->nzmax = T->nz ;
        C = cs_cl_compress (T) ;                /* create sparse matrix C */
        cs_cl_dupl (C) ;                        /* remove duplicates from C */
        cs_cl_dropzeros (C) ;                   /* remove zeros from C */
        A = cs_cl_transpose (C, -1) ;           /* A=C.' */
        cs_cl_spfree (C) ;
        pargout [0] = cs_cl_mex_put_sparse (&A) ;       /* return A */
        cs_free (T->p) ;
        cs_free (T->i) ;
        cs_free (T->x) ;                        /* free copy of complex values*/
#else
        mexErrMsgTxt ("complex matrices not supported") ;
#endif
    }
    else
    {
        cs_dl *A, *C, *T, Tmatrix ;
        T = &Tmatrix ;                  /* get i,j,x and copy to triplet form */
        T->nz = mxGetM (pargin [0]) ;
        T->p = cs_dl_mex_get_int (T->nz, pargin [0], &(T->n), 1) ;
        T->i = cs_dl_mex_get_int (T->nz, pargin [1], &(T->m), 1) ;
        cs_mex_check (1, T->nz, 1, 0, 0, 1, pargin [2]) ;
        T->x = mxGetPr (pargin [2]) ;
        T->nzmax = T->nz ;
        C = cs_dl_compress (T) ;                /* create sparse matrix C */
        cs_dl_dupl (C) ;                        /* remove duplicates from C */
        cs_dl_dropzeros (C) ;                   /* remove zeros from C */
        A = cs_dl_transpose (C, 1) ;            /* A=C' */
        cs_dl_spfree (C) ;
        pargout [0] = cs_dl_mex_put_sparse (&A) ;       /* return A */
        cs_free (T->p) ;
        cs_free (T->i) ;
    }
}
Exemple #6
0
/* cs_lu: sparse LU factorization, with optional fill-reducing ordering */
void mexFunction
(
    int nargout,
    mxArray *pargout [ ],
    int nargin,
    const mxArray *pargin [ ]
)
{
    CS_INT n, order, *p ;
    double tol ;
    if (nargout > 4 || nargin > 3 || nargin < 1)
    {
        mexErrMsgTxt ("Usage: [L,U,p,q] = cs_lu (A,tol)") ;
    }
    if (nargin == 2)                        /* determine tol and ordering */
    {
        tol = mxGetScalar (pargin [1]) ;
        order = (nargout == 4) ? 1 : 0 ;    /* amd (A+A'), or natural */
    }
    else
    {
        tol = 1 ;
        order = (nargout == 4) ? 2 : 0 ;    /* amd(S'*S) w/dense rows or I */
    }
    if (mxIsComplex (pargin [0]))
    {
#ifndef NCOMPLEX
        cs_cls *S ;
        cs_cln *N ;
        cs_cl Amatrix, *A, *D ;
        A = cs_cl_mex_get_sparse (&Amatrix, 1, pargin [0]) ;    /* get A */
        n = A->n ;
        S = cs_cl_sqr (order, A, 0) ;       /* symbolic ordering, no QR bound */
        N = cs_cl_lu (A, S, tol) ;          /* numeric factorization */
        if (!N) mexErrMsgTxt ("cs_lu failed (singular, or out of memory)") ;
        cs_cl_free (A->x) ;                 /* complex copy no longer needed */
        cs_cl_dropzeros (N->L) ;            /* drop zeros from L and sort it */
        D = cs_cl_transpose (N->L, 1) ;
        cs_cl_spfree (N->L) ;
        N->L = cs_cl_transpose (D, 1) ;
        cs_cl_spfree (D) ;
        cs_cl_dropzeros (N->U) ;            /* drop zeros from U and sort it */
        D = cs_cl_transpose (N->U, 1) ;
        cs_cl_spfree (N->U) ;
        N->U = cs_cl_transpose (D, 1) ;
        cs_cl_spfree (D) ;
        p = cs_cl_pinv (N->pinv, n) ;                       /* p=pinv' */
        pargout [0] = cs_cl_mex_put_sparse (&(N->L)) ;      /* return L */
        pargout [1] = cs_cl_mex_put_sparse (&(N->U)) ;      /* return U */
        pargout [2] = cs_dl_mex_put_int (p, n, 1, 1) ;      /* return p */
        /* return Q */
        if (nargout == 4) pargout [3] = cs_dl_mex_put_int (S->q, n, 1, 0) ;
        cs_cl_nfree (N) ;
        cs_cl_sfree (S) ;
#else
        mexErrMsgTxt ("complex matrices not supported") ;
#endif
    }
    else
    {
        cs_dls *S ;
        cs_dln *N ;
        cs_dl Amatrix, *A, *D ;
        A = cs_dl_mex_get_sparse (&Amatrix, 1, 1, pargin [0]) ; /* get A */
        n = A->n ;
        S = cs_dl_sqr (order, A, 0) ;       /* symbolic ordering, no QR bound */
        N = cs_dl_lu (A, S, tol) ;          /* numeric factorization */
        if (!N) mexErrMsgTxt ("cs_lu failed (singular, or out of memory)") ;
        cs_dl_dropzeros (N->L) ;            /* drop zeros from L and sort it */
        D = cs_dl_transpose (N->L, 1) ;
        cs_dl_spfree (N->L) ;
        N->L = cs_dl_transpose (D, 1) ;
        cs_dl_spfree (D) ;
        cs_dl_dropzeros (N->U) ;            /* drop zeros from U and sort it */
        D = cs_dl_transpose (N->U, 1) ;
        cs_dl_spfree (N->U) ;
        N->U = cs_dl_transpose (D, 1) ;
        cs_dl_spfree (D) ;
        p = cs_dl_pinv (N->pinv, n) ;                       /* p=pinv' */
        pargout [0] = cs_dl_mex_put_sparse (&(N->L)) ;      /* return L */
        pargout [1] = cs_dl_mex_put_sparse (&(N->U)) ;      /* return U */
        pargout [2] = cs_dl_mex_put_int (p, n, 1, 1) ;      /* return p */
        /* return Q */
        if (nargout == 4) pargout [3] = cs_dl_mex_put_int (S->q, n, 1, 0) ;
        cs_dl_nfree (N) ;
        cs_dl_sfree (S) ;
    }
}
/* cs_qr: sparse QR factorization */
void mexFunction
(
    int nargout,
    mxArray *pargout [ ],
    int nargin,
    const mxArray *pargin [ ]
)
{
    CS_INT m, n, order, *p ;
    if (nargout > 5 || nargin != 1)
    {
        mexErrMsgTxt ("Usage: [V,beta,p,R,q] = cs_qr(A)") ;
    }
    order = (nargout == 5) ? 3 : 0 ;        /* determine ordering */
    m = mxGetM (pargin [0]) ;
    n = mxGetN (pargin [0]) ;
    if (m < n) mexErrMsgTxt ("A must have # rows >= # columns") ;
    if (mxIsComplex (pargin [0]))
    {
#ifndef NCOMPLEX
        cs_cls *S ;
        cs_cln *N ;
        cs_cl Amatrix, *A, *D ;
        A = cs_cl_mex_get_sparse (&Amatrix, 0, pargin [0]) ;    /* get A */
        S = cs_cl_sqr (order, A, 1) ;       /* symbolic QR ordering & analysis*/
        N = cs_cl_qr (A, S) ;               /* numeric QR factorization */
        cs_free (A->x) ;
        if (!N) mexErrMsgTxt ("qr failed") ;
        cs_cl_dropzeros (N->L) ;            /* drop zeros from V and sort */
        D = cs_cl_transpose (N->L, 1) ;
        cs_cl_spfree (N->L) ;
        N->L = cs_cl_transpose (D, 1) ;
        cs_cl_spfree (D) ;
        cs_cl_dropzeros (N->U) ;            /* drop zeros from R and sort */
        D = cs_cl_transpose (N->U, 1) ;
        cs_cl_spfree (N->U) ;
        N->U = cs_cl_transpose (D, 1) ;
        cs_cl_spfree (D) ;
        m = N->L->m ;                               /* m may be larger now */
        p = cs_cl_pinv (S->pinv, m) ;                   /* p = pinv' */
        pargout [0] = cs_cl_mex_put_sparse (&(N->L)) ;  /* return V */
        cs_dl_mex_put_double (n, N->B, &(pargout [1])) ;   /* return beta */
        pargout [2] = cs_dl_mex_put_int (p, m, 1, 1) ;  /* return p */
        pargout [3] = cs_cl_mex_put_sparse (&(N->U)) ;  /* return R */
        pargout [4] = cs_dl_mex_put_int (S->q, n, 1, 0) ;  /* return q */
        cs_cl_nfree (N) ;
        cs_cl_sfree (S) ;
#else
        mexErrMsgTxt ("complex matrices not supported") ;
#endif
    }
    else
    {
        cs_dls *S ;
        cs_dln *N ;
        cs_dl Amatrix, *A, *D ;
        A = cs_dl_mex_get_sparse (&Amatrix, 0, 1, pargin [0]) ; /* get A */
        S = cs_dl_sqr (order, A, 1) ;       /* symbolic QR ordering & analysis*/
        N = cs_dl_qr (A, S) ;               /* numeric QR factorization */
        if (!N) mexErrMsgTxt ("qr failed") ;
        cs_dl_dropzeros (N->L) ;            /* drop zeros from V and sort */
        D = cs_dl_transpose (N->L, 1) ;
        cs_dl_spfree (N->L) ;
        N->L = cs_dl_transpose (D, 1) ;
        cs_dl_spfree (D) ;
        cs_dl_dropzeros (N->U) ;            /* drop zeros from R and sort */
        D = cs_dl_transpose (N->U, 1) ;
        cs_dl_spfree (N->U) ;
        N->U = cs_dl_transpose (D, 1) ;
        cs_dl_spfree (D) ;
        m = N->L->m ;                               /* m may be larger now */
        p = cs_dl_pinv (S->pinv, m) ;                   /* p = pinv' */
        pargout [0] = cs_dl_mex_put_sparse (&(N->L)) ;  /* return V */
        cs_dl_mex_put_double (n, N->B, &(pargout [1])) ;   /* return beta */
        pargout [2] = cs_dl_mex_put_int (p, m, 1, 1) ;  /* return p */
        pargout [3] = cs_dl_mex_put_sparse (&(N->U)) ;  /* return R */
        pargout [4] = cs_dl_mex_put_int (S->q, n, 1, 0) ;  /* return q */
        cs_dl_nfree (N) ;
        cs_dl_sfree (S) ;
    }
}
void mexFunction
(
    int nargout,
    mxArray *pargout [ ],
    int nargin,
    const mxArray *pargin [ ]
)
{
    int iterations;
    double mu, lrate,lrateB, lambda,lambdaB;
    double *buin, *biin, *qin, *xin, *yin, *pin, *zin;
    double *bu,   *bi,   *q,   *x,   *y,   *p,   *z;
    int ls;
    int itemsNum, usersNum;
    cs_dl Amatrix, *A;
    bool useruser = false;
    FILE *verbose;
    
    if (nargout < 5 || nargout > 7 || nargin < 10 || nargin > 12)
    {
        mexErrMsgTxt ("Usage: [bu,bi,q,x,y,[p,z]]=learnFactorModel(urm,mu,bu,bi,iterations,lrate,lrateB,lambda,lambdaB,q,x,y,[p,z]) \n urm must be sparse, bu and bi must be dense") ;
    }
    
    if (nargout>5 && nargin>11)
    {
        useruser=true; /* enable the integration of the user-user model into the item-item model */
        if (log) fprintf(stdout,"user-user model enabled\n");
    }    
    
    int indexArgin = 0, indexArgout = 0;;
    
    if (log>1)
    {
        verbose = fopen("/home/roby/verbose.txt", "w+"); 
        fprintf(stdout, "verbose mode");   
    }
    
    if (log) fprintf(stdout,"--input started\n");
    A = cs_dl_mex_get_sparse (&Amatrix, 0, 1, pargin [indexArgin]) ;  usersNum = mxGetM(pargin[indexArgin++]);  /* get A=urm */
    mu = mxGetScalar (pargin [indexArgin++]);
    buin = mxGetPr(pargin[indexArgin++]); 
    biin = mxGetPr(pargin[indexArgin++]);
    iterations = mxGetScalar (pargin [indexArgin++]);
    lrate = mxGetScalar (pargin [indexArgin++]);
    lrateB = mxGetScalar (pargin [indexArgin++]);
    lambda = mxGetScalar (pargin [indexArgin++]);
    lambdaB = mxGetScalar (pargin [indexArgin++]);
    qin = mxGetPr(pargin[indexArgin]); ls = mxGetM(pargin[indexArgin]); itemsNum = mxGetN(pargin[indexArgin++]); 
    xin = mxGetPr(pargin[indexArgin++]);
    yin = mxGetPr(pargin[indexArgin++]);
    if (useruser)
    {
        pin = mxGetPr(pargin[indexArgin++]);
        zin = mxGetPr(pargin[indexArgin++]);   
        if (log) fprintf(stdout," read: p and z \n");            
    }  
        
    if (log) fprintf(stdout,"users=%d, items=%d\n",usersNum,itemsNum);
    if (log) fprintf(stdout,"number of factors=%d, number of iterations=%d\n",ls,iterations);
    if (log) fprintf(stdout," bu =[%d x %d] \n", mxGetM(pargin[2]), mxGetN(pargin[2]));
    if (log) fprintf(stdout," bi =[%d x %d] \n", mxGetM(pargin[3]), mxGetN(pargin[3]));       
    
    if (log) fprintf(stdout,"--input completed\n");
    
    UF_long Am, An, Anzmax, Anz, *Ap, *Ai ;
    double *Ax ;
    if (!A) { printf ("(null)\n") ; return  ; }
    
    Am = A->m ; An = A->n ; 	    Ap = A->p ;   Ai = A->i ; 
    Ax = A->x ;	Anzmax = A->nzmax ; Anz = A->nz ;
    
    /* 
        B = transpose of URM
    */
    cs_dl *B;
    UF_long Bm, Bn, Bnzmax, Bnz, *Bp, *Bi ;
    double *Bx ;
    B = cs_dl_transpose (A, 1) ;                       /* B = A' = urm' */
    Bm = B->m ; Bn = B->n ; 	    Bp = B->p ;   Bi = B->i ; 
    Bx = B->x ;	Bnzmax = B->nzmax ; Bnz = B->nz ;
    
    int ii;
    pargout[indexArgout]=mxCreateDoubleMatrix(usersNum,1,mxREAL);
    bu=mxGetPr(pargout[indexArgout++]);
    for (ii=0;ii<usersNum;ii++) bu[ii] = buin[ii];

    pargout[indexArgout]=mxCreateDoubleMatrix(itemsNum,1,mxREAL);
    bi=mxGetPr(pargout[indexArgout++]);
    for (ii=0;ii<itemsNum;ii++) bi[ii] = biin[ii];

    pargout[indexArgout]=mxCreateDoubleMatrix(ls,itemsNum,mxREAL);
    q=mxGetPr(pargout[indexArgout++]);
    for (ii=0;ii<ls*itemsNum;ii++) q[ii] = qin[ii];

    pargout[indexArgout]=mxCreateDoubleMatrix(ls,itemsNum,mxREAL);
    x=mxGetPr(pargout[indexArgout++]);
    for (ii=0;ii<ls*itemsNum;ii++) x[ii] = xin[ii];

    pargout[indexArgout]=mxCreateDoubleMatrix(ls,itemsNum,mxREAL);
    y=mxGetPr(pargout[indexArgout++]);
    for (ii=0;ii<ls*itemsNum;ii++) y[ii] = yin[ii];
    
    if (useruser)    
    {
        pargout[indexArgout]=mxCreateDoubleMatrix(ls,usersNum,mxREAL);
        p=mxGetPr(pargout[indexArgout++]);
        for (ii=0;ii<ls*usersNum;ii++) p[ii] = pin[ii];
        if (log) fprintf(stdout," p initialized \n"); 

        pargout[indexArgout]=mxCreateDoubleMatrix(ls,usersNum,mxREAL);
        z=mxGetPr(pargout[indexArgout++]);
        for (ii=0;ii<ls*usersNum;ii++) z[ii] = zin[ii];
        if (log) fprintf(stdout," z initialized \n");         
    }            
    if (log) fprintf(stdout," \n allocated %d outputs \n",indexArgout);    
    
    int count, u, i, j, k;
    int iii;
    int itemj;
    int numRatedItems, numRatingUsers;
    double puCoeff, zvCoeff;
    double *pu=calloc(ls, sizeof(double));
    double *zv=calloc(ls, sizeof(double)); 
    double sum[ls];
    double *sumUU=calloc(ls, sizeof(double));
    double cumerror=0.0, pseudoRMSE=0.0; 
    long numTests=0;
    int anomaliesCounter=0;
    const int maxAnomalies=100000;
    time_t t1,t2, t0;
    (void) time(&t0);
    for (count=0;count<iterations;count++) /* FOR count=1,#iterations DO */
    {   
        (void) time(&t1);
        cumerror=0.0;
        numTests = 0; 
        anomaliesCounter = 0;
        for (u=0;u<usersNum;u++) /* FOR u=1,..,m DO */
        {
            double errorUU=0;
            for (iii=0;iii<ls;iii++) pu[iii]=0.0; /* pu[ls] <- 0 */
            numRatedItems = (int) (Bp [u+1] - Bp [u]);
            if (numRatedItems==0) continue;
           
            /* compute the component independent of i*/
            for (j=Bp[u]; j<Bp[u+1]; j++)
            {   
                double ruj, biasi, xjk, yjk;
                itemj = Bi[j];
                ruj=Bx[j]; /* r_uj */
                
                biasi=bi[itemj];
                puCoeff = (ruj - (mu + bu[u] + biasi)) / (sqrt((double) numRatedItems)); /* |R(u)|^-1/2 * (r_uj - b_uj) */
                for (k=0;k<ls;k++) /* compute pu for each feature k*/
                {   
                    xjk = x[ls*itemj+k]; /* x[ls*itemj+k] = x[ls,itemj] */ 
                    yjk = y[ls*itemj+k]; /* y[ls*itemj+k] = y[ls,itemj] */
                    pu[k] = pu[k] + (puCoeff * xjk);
                    pu[k] = pu[k] + (yjk / (sqrt((double) numRatedItems)));
                }
            }
            for (iii=0;iii<ls;iii++) sum[iii]=0.0; /* sum[ls] <- 0 */
            if (useruser)
            {
                for (iii=0;iii<ls;iii++) sumUU[iii]=0.0; /* sumUU[ls] <- 0 */           
            }
            /* FOR ALL i (itemj) IN R(u) DO */
            for (j=Bp[u]; j<Bp[u+1]; j++)
            {   
                double r_hat_ui=0.0, e_ui;
                double biasi;
                itemj = Bi[j];
                if (useruser) /* if user-user is enabled */
                {
                    int useri;
                    numRatingUsers = Ap [itemj+1] - Ap [itemj]; /* |R(i)| */ 
                    for (iii=0;iii<ls;iii++) zv[iii]=0.0; /* zv[ls] <- 0 */
                    biasi = bi[itemj]; /* bias item j */
                    for (i=Ap[itemj]; i<Ap[itemj+1]; i++)
                    {   
                        double rvj, biasvj, zjk;
                        useri = Ai[i];
                        rvj=Ax[i];
                        
                        biasvj=bu[useri];
                        zvCoeff = (rvj - (mu + biasvj + biasi)) / (sqrt((double) numRatingUsers)); /* |R(i)|^-1/2 * (r_vj - b_vj) */                    
                        for (k=0;k<ls;k++) /* compute pu for each feature k*/
                        {   
                            zjk = z[ls*useri+k]; /* z[ls*useri+k] = z[ls,itemj] */                                                          
                            zv[k] += (zvCoeff * zjk);                         
                        }
                    }
                    for (k=0; k<ls; k++) /* p_u' * z_v */
                    {                    
                        r_hat_ui += ( p[ls*u+k] * zv[k] ); /* p_u(k) * zv(k) */                          
                    }  
                }

                for (k=0; k<ls; k++) /* q_i' * pu */
                {                
                    r_hat_ui += ( q[ls*itemj+k] * pu[k] ); /* q_i(k) * pu(k) */                  
                }
                r_hat_ui += (mu + bu[u] + bi[itemj]); /* r_hat_ui = mu + bu + bi + q_i'*pu */
                e_ui = Bx[j] - r_hat_ui; /* e_ui = r_ui - r_hat_ui */
                if (log>1) fprintf(verbose,"(%d,%d)=%f\n",u,itemj,e_ui);
                if (mxIsNaN(e_ui) || mxIsInf(e_ui) || absolute(e_ui)>100000) 
                {
                    fprintf(stdout, " -!- [%d] -!- user=%d, item=%d, r=%f, r_hat=%f, e_ui=%f  gradient error too large \n",anomaliesCounter,u,itemj,Bx[j],r_hat_ui,e_ui);
                    fprintf(stdout, "  biasu=%f, biasi=%f \n", bu[u], bi[itemj]);
                    if (log>1) fclose(verbose);
                    return;    
                }                
                if (absolute(e_ui)>7) 
                {
                    fprintf(stdout, " -!- [%d] -!- user=%d, item=%d, r=%f, r_hat=%f, e_ui=%f  gradient error too large \n",anomaliesCounter,u,itemj,Bx[j],r_hat_ui,e_ui);
                    fprintf(stdout, "  biasu=%f, biasi=%f \n", bu[u], bi[itemj]);
                    if (log>1) fclose(verbose);
                    anomaliesCounter++;
                    if (anomaliesCounter>maxAnomalies) return;    
                }
                cumerror += (e_ui*e_ui);
                numTests++;
                for (k=0; k<ls; k++) /* sum <- sum + e_ui * q_i */
                {
                    sum[k] += (e_ui*q[ls*itemj+k]); /* sum(k) = e_ui * q_i(k) */
                }
                if (useruser) /* sumUU <- sumUU + e_ui * p_u */
                {
                    for (k=0; k<ls; k++)
                    {
                        sumUU[k] += (e_ui*p[ls*u+k]); /* sumUU(k) = e_ui * p_u(k) */     
                    }
                }
                /* perform gradient step on qi, bu, bi  AND on "pu" if user-user model is enabled */
                bu[u] += (lrateB * (e_ui - lambdaB*bu[u])); /* bu <- bu + gamma*(e_ui-lambdaB*bu) */        
                bi[itemj] += (lrateB * (e_ui - lambdaB*bi[itemj])); /* bi <- bi + gamma*(e_ui-lambdaB*bi) */                                   
                for (k=0; k<ls; k++) 
                {
                    q[ls*itemj+k] += (lrate * (e_ui*pu[k] - lambda* q[ls*itemj+k])); /* q_i <- q_i + gamma*(e_ui*pu-lambda*q_i) */
                    if (useruser)
                    {
                        p[ls*u+k] += (lrate * (e_ui*zv[k] - lambda* p[ls*u+k])); /* p_u <- p_u + gamma*(e_ui*zv-lambda*p_u) */
                    }
                }                
            }
            /* FOR ALL i IN R(u) DO */
            for (j=Bp[u]; j<Bp[u+1]; j++) /* perform gradient step on xi*/
            {   
                itemj = Bi[j];
                double ruj = Bx[j];
                double biasi = bi[itemj];
                puCoeff = (1/sqrt(numRatedItems)) * (ruj - (mu + bu[u] + biasi)); /* |R(u)|^-1/2 * (r_ui - b_ui) */
                for (k=0; k<ls; k++)
                {
                    x[ls*itemj+k] += (lrate * ( puCoeff*sum[k] - lambda*x[ls*itemj+k] ));   /* update of every feature of xi */
                }   
                if (useruser)
                {
                    int useri;
                    numRatingUsers = Ap [itemj+1] - Ap [itemj]; /* |R(i)| */ 
                    /* FOR ALL useri IN R(i) DO */
                    for (i=Ap[itemj]; i<Ap[itemj+1]; i++) /* perform gradient step on zv*/
                    {   
                        double rvj = Ax[i];
                        double biasvj;
                        
                        useri = Ai[i]; /* we are looping on all useri who rated itemj*/
                        biasvj = bu[useri]; /* bias user v*/
    
                        zvCoeff = (1/sqrt(numRatingUsers)) * (rvj - (mu + biasvj + biasi)); /* |R(i)|^-1/2 * (r_vj - b_vj) */
                        for (k=0; k<ls; k++)
                        {
                            z[ls*useri+k] += (lrate * ( zvCoeff*sumUU[k] - lambda*z[ls*useri+k] ) );   /* update of every feature of zv */
                        }          
                    }            
                }                       
            }
            /* FOR ALL i IN N(u) DO */
            for (j=Bp[u]; j<Bp[u+1]; j++) /* perform gradient step on yi*/
            {   
                itemj = Bi[j];
                double ruj = Bx[j];
                puCoeff = (1/sqrt(numRatedItems)); /* |N(u)|^-1/2 */
                for (k=0; k<ls; k++)
                {
                    y[ls*itemj+k] += (lrate * ( puCoeff*sum[k] - lambda*y[ls*itemj+k] ) );   /* update of every feature of yi */
                }          
            }                   
            if (log && ( ((u-1) % 10000 == 0)) )
            {
                (void) time(&t2);
                fprintf(stdout, "[%d] time for last group (up to user %d) is %d secs - remaining Time =%d secs\n", (int) t2-t0,u, (int) t2-t1, (int) ( ( (t2-t0)/((double) u))*((double)((usersNum-u)*(count+1)))));
                if ((((double) (t2-t1))/10000.0)>1.0) fprintf (stdout, "warning... high computing time");
                if (file_exists("~/stopnow")) return;
                (void) time(&t1);
            }


        }
        pseudoRMSE = cumerror / ((double) numTests);
        if (log) 
        {   
            fprintf(stdout, " cumulative error iteration %d = %g (%d tests) \n", count, pseudoRMSE,numTests);
            if (file_exists("~/stopiter")) return;
        }
        if (mxIsNaN(pseudoRMSE))
        {
            fprintf(stdout, " -!- STOPPED -!- ");
            return;        
        }
    }
    if (log>1) fclose(verbose);
    
}