/************************************************************************* Testing *************************************************************************/ public static bool testfbls(bool silent) { bool result = new bool(); int n = 0; int m = 0; int mx = 0; int i = 0; int j = 0; bool waserrors = new bool(); bool cgerrors = new bool(); double v = 0; double v1 = 0; double v2 = 0; double[] tmp1 = new double[0]; double[] tmp2 = new double[0]; double[,] a = new double[0,0]; double[] b = new double[0]; double[] x = new double[0]; double[] xe = new double[0]; double[] buf = new double[0]; double alpha = 0; double e1 = 0; double e2 = 0; fbls.fblslincgstate cgstate = new fbls.fblslincgstate(); int i_ = 0; mx = 10; waserrors = false; cgerrors = false; // // Test CG solver: // * generate problem (A, B, Alpha, XE - exact solution) and initial approximation X // * E1 = ||A'A*x-b|| // * solve // * E2 = ||A'A*x-b|| // * test that E2<0.001*E1 // for(n=1; n<=mx; n++) { for(m=1; m<=mx; m++) { a = new double[m, n]; b = new double[n]; x = new double[n]; xe = new double[n]; tmp1 = new double[m]; tmp2 = new double[n]; // // init A, alpha, B, X (initial approximation), XE (exact solution) // X is initialized in such way that is has no chances to be equal to XE. // for(i=0; i<=m-1; i++) { for(j=0; j<=n-1; j++) { a[i,j] = 2*math.randomreal()-1; } } alpha = math.randomreal()+0.1; for(i=0; i<=n-1; i++) { b[i] = 2*math.randomreal()-1; xe[i] = 2*math.randomreal()-1; x[i] = (2*math.randominteger(2)-1)*(2+math.randomreal()); } // // Test dense CG (which solves A'A*x=b and accepts dense A) // for(i=0; i<=n-1; i++) { x[i] = (2*math.randominteger(2)-1)*(2+math.randomreal()); } ablas.rmatrixmv(m, n, a, 0, 0, 0, x, 0, ref tmp1, 0); ablas.rmatrixmv(n, m, a, 0, 0, 1, tmp1, 0, ref tmp2, 0); for(i_=0; i_<=n-1;i_++) { tmp2[i_] = tmp2[i_] + alpha*x[i_]; } for(i_=0; i_<=n-1;i_++) { tmp2[i_] = tmp2[i_] - b[i_]; } v = 0.0; for(i_=0; i_<=n-1;i_++) { v += tmp2[i_]*tmp2[i_]; } e1 = Math.Sqrt(v); fbls.fblssolvecgx(a, m, n, alpha, b, ref x, ref buf); ablas.rmatrixmv(m, n, a, 0, 0, 0, x, 0, ref tmp1, 0); ablas.rmatrixmv(n, m, a, 0, 0, 1, tmp1, 0, ref tmp2, 0); for(i_=0; i_<=n-1;i_++) { tmp2[i_] = tmp2[i_] + alpha*x[i_]; } for(i_=0; i_<=n-1;i_++) { tmp2[i_] = tmp2[i_] - b[i_]; } v = 0.0; for(i_=0; i_<=n-1;i_++) { v += tmp2[i_]*tmp2[i_]; } e2 = Math.Sqrt(v); cgerrors = cgerrors | (double)(e2)>(double)(0.001*e1); // // Test sparse CG (which relies on reverse communication) // for(i=0; i<=n-1; i++) { x[i] = (2*math.randominteger(2)-1)*(2+math.randomreal()); } ablas.rmatrixmv(m, n, a, 0, 0, 0, x, 0, ref tmp1, 0); ablas.rmatrixmv(n, m, a, 0, 0, 1, tmp1, 0, ref tmp2, 0); for(i_=0; i_<=n-1;i_++) { tmp2[i_] = tmp2[i_] + alpha*x[i_]; } for(i_=0; i_<=n-1;i_++) { tmp2[i_] = tmp2[i_] - b[i_]; } v = 0.0; for(i_=0; i_<=n-1;i_++) { v += tmp2[i_]*tmp2[i_]; } e1 = Math.Sqrt(v); fbls.fblscgcreate(x, b, n, cgstate); while( fbls.fblscgiteration(cgstate) ) { ablas.rmatrixmv(m, n, a, 0, 0, 0, cgstate.x, 0, ref tmp1, 0); ablas.rmatrixmv(n, m, a, 0, 0, 1, tmp1, 0, ref cgstate.ax, 0); for(i_=0; i_<=n-1;i_++) { cgstate.ax[i_] = cgstate.ax[i_] + alpha*cgstate.x[i_]; } v1 = 0.0; for(i_=0; i_<=m-1;i_++) { v1 += tmp1[i_]*tmp1[i_]; } v2 = 0.0; for(i_=0; i_<=n-1;i_++) { v2 += cgstate.x[i_]*cgstate.x[i_]; } cgstate.xax = v1+alpha*v2; } ablas.rmatrixmv(m, n, a, 0, 0, 0, cgstate.xk, 0, ref tmp1, 0); ablas.rmatrixmv(n, m, a, 0, 0, 1, tmp1, 0, ref tmp2, 0); for(i_=0; i_<=n-1;i_++) { tmp2[i_] = tmp2[i_] + alpha*cgstate.xk[i_]; } for(i_=0; i_<=n-1;i_++) { tmp2[i_] = tmp2[i_] - b[i_]; } v = 0.0; for(i_=0; i_<=n-1;i_++) { v += tmp2[i_]*tmp2[i_]; } e2 = Math.Sqrt(v); cgerrors = cgerrors | (double)(Math.Abs(e1-cgstate.e1))>(double)(100*math.machineepsilon*e1); cgerrors = cgerrors | (double)(Math.Abs(e2-cgstate.e2))>(double)(100*math.machineepsilon*e1); cgerrors = cgerrors | (double)(e2)>(double)(0.001*e1); } } // // report // waserrors = cgerrors; if( !silent ) { System.Console.Write("TESTING FBLS"); System.Console.WriteLine(); System.Console.Write("CG ERRORS: "); if( cgerrors ) { System.Console.Write("FAILED"); System.Console.WriteLine(); } else { System.Console.Write("OK"); System.Console.WriteLine(); } if( waserrors ) { System.Console.Write("TEST FAILED"); System.Console.WriteLine(); } else { System.Console.Write("TEST PASSED"); System.Console.WriteLine(); } System.Console.WriteLine(); System.Console.WriteLine(); } result = !waserrors; return result; }
/************************************************************************* Weighted fitting by penalized cubic spline. Equidistant grid with M nodes on [min(x,xc),max(x,xc)] is used to build basis functions. Basis functions are cubic splines with natural boundary conditions. Problem is regularized by adding non-linearity penalty to the usual least squares penalty function: S(x) = arg min { LS + P }, where LS = SUM { w[i]^2*(y[i] - S(x[i]))^2 } - least squares penalty P = C*10^rho*integral{ S''(x)^2*dx } - non-linearity penalty rho - tunable constant given by user C - automatically determined scale parameter, makes penalty invariant with respect to scaling of X, Y, W. INPUT PARAMETERS: X - points, array[0..N-1]. Y - function values, array[0..N-1]. W - weights, array[0..N-1] Each summand in square sum of approximation deviations from given values is multiplied by the square of corresponding weight. Fill it by 1's if you don't want to solve weighted problem. N - number of points (optional): * N>0 * if given, only first N elements of X/Y/W are processed * if not given, automatically determined from X/Y/W sizes M - number of basis functions ( = number_of_nodes), M>=4. Rho - regularization constant passed by user. It penalizes nonlinearity in the regression spline. It is logarithmically scaled, i.e. actual value of regularization constant is calculated as 10^Rho. It is automatically scaled so that: * Rho=2.0 corresponds to moderate amount of nonlinearity * generally, it should be somewhere in the [-8.0,+8.0] If you do not want to penalize nonlineary, pass small Rho. Values as low as -15 should work. OUTPUT PARAMETERS: Info- same format as in LSFitLinearWC() subroutine. * Info>0 task is solved * Info<=0 an error occured: -4 means inconvergence of internal SVD or Cholesky decomposition; problem may be too ill-conditioned (very rare) S - spline interpolant. Rep - Following fields are set: * RMSError rms error on the (X,Y). * AvgError average error on the (X,Y). * AvgRelError average relative error on the non-zero Y * MaxError maximum error NON-WEIGHTED ERRORS ARE CALCULATED IMPORTANT: this subroitine doesn't calculate task's condition number for K<>0. NOTE 1: additional nodes are added to the spline outside of the fitting interval to force linearity when x<min(x,xc) or x>max(x,xc). It is done for consistency - we penalize non-linearity at [min(x,xc),max(x,xc)], so it is natural to force linearity outside of this interval. NOTE 2: function automatically sorts points, so caller may pass unsorted array. -- ALGLIB PROJECT -- Copyright 19.10.2010 by Bochkanov Sergey *************************************************************************/ public static void spline1dfitpenalizedw(double[] x, double[] y, double[] w, int n, int m, double rho, ref int info, spline1d.spline1dinterpolant s, spline1dfitreport rep) { int i = 0; int j = 0; int b = 0; double v = 0; double relcnt = 0; double xa = 0; double xb = 0; double sa = 0; double sb = 0; double[] xoriginal = new double[0]; double[] yoriginal = new double[0]; double pdecay = 0; double tdecay = 0; double[,] fmatrix = new double[0,0]; double[] fcolumn = new double[0]; double[] y2 = new double[0]; double[] w2 = new double[0]; double[] xc = new double[0]; double[] yc = new double[0]; int[] dc = new int[0]; double fdmax = 0; double admax = 0; double[,] amatrix = new double[0,0]; double[,] d2matrix = new double[0,0]; double fa = 0; double ga = 0; double fb = 0; double gb = 0; double lambdav = 0; double[] bx = new double[0]; double[] by = new double[0]; double[] bd1 = new double[0]; double[] bd2 = new double[0]; double[] tx = new double[0]; double[] ty = new double[0]; double[] td = new double[0]; spline1d.spline1dinterpolant bs = new spline1d.spline1dinterpolant(); double[,] nmatrix = new double[0,0]; double[] rightpart = new double[0]; fbls.fblslincgstate cgstate = new fbls.fblslincgstate(); double[] c = new double[0]; double[] tmp0 = new double[0]; int i_ = 0; int i1_ = 0; x = (double[])x.Clone(); y = (double[])y.Clone(); w = (double[])w.Clone(); info = 0; alglib.ap.assert(n>=1, "Spline1DFitPenalizedW: N<1!"); alglib.ap.assert(m>=4, "Spline1DFitPenalizedW: M<4!"); alglib.ap.assert(alglib.ap.len(x)>=n, "Spline1DFitPenalizedW: Length(X)<N!"); alglib.ap.assert(alglib.ap.len(y)>=n, "Spline1DFitPenalizedW: Length(Y)<N!"); alglib.ap.assert(alglib.ap.len(w)>=n, "Spline1DFitPenalizedW: Length(W)<N!"); alglib.ap.assert(apserv.isfinitevector(x, n), "Spline1DFitPenalizedW: X contains infinite or NAN values!"); alglib.ap.assert(apserv.isfinitevector(y, n), "Spline1DFitPenalizedW: Y contains infinite or NAN values!"); alglib.ap.assert(apserv.isfinitevector(w, n), "Spline1DFitPenalizedW: Y contains infinite or NAN values!"); alglib.ap.assert(math.isfinite(rho), "Spline1DFitPenalizedW: Rho is infinite!"); // // Prepare LambdaV // v = -(Math.Log(math.machineepsilon)/Math.Log(10)); if( (double)(rho)<(double)(-v) ) { rho = -v; } if( (double)(rho)>(double)(v) ) { rho = v; } lambdav = Math.Pow(10, rho); // // Sort X, Y, W // spline1d.heapsortdpoints(ref x, ref y, ref w, n); // // Scale X, Y, XC, YC // lsfitscalexy(ref x, ref y, ref w, n, ref xc, ref yc, dc, 0, ref xa, ref xb, ref sa, ref sb, ref xoriginal, ref yoriginal); // // Allocate space // fmatrix = new double[n, m]; amatrix = new double[m, m]; d2matrix = new double[m, m]; bx = new double[m]; by = new double[m]; fcolumn = new double[n]; nmatrix = new double[m, m]; rightpart = new double[m]; tmp0 = new double[Math.Max(m, n)]; c = new double[m]; // // Fill: // * FMatrix by values of basis functions // * TmpAMatrix by second derivatives of I-th function at J-th point // * CMatrix by constraints // fdmax = 0; for(b=0; b<=m-1; b++) { // // Prepare I-th basis function // for(j=0; j<=m-1; j++) { bx[j] = (double)(2*j)/(double)(m-1)-1; by[j] = 0; } by[b] = 1; spline1d.spline1dgriddiff2cubic(bx, by, m, 2, 0.0, 2, 0.0, ref bd1, ref bd2); spline1d.spline1dbuildcubic(bx, by, m, 2, 0.0, 2, 0.0, bs); // // Calculate B-th column of FMatrix // Update FDMax (maximum column norm) // spline1d.spline1dconvcubic(bx, by, m, 2, 0.0, 2, 0.0, x, n, ref fcolumn); for(i_=0; i_<=n-1;i_++) { fmatrix[i_,b] = fcolumn[i_]; } v = 0; for(i=0; i<=n-1; i++) { v = v+math.sqr(w[i]*fcolumn[i]); } fdmax = Math.Max(fdmax, v); // // Fill temporary with second derivatives of basis function // for(i_=0; i_<=m-1;i_++) { d2matrix[b,i_] = bd2[i_]; } } // // * calculate penalty matrix A // * calculate max of diagonal elements of A // * calculate PDecay - coefficient before penalty matrix // for(i=0; i<=m-1; i++) { for(j=i; j<=m-1; j++) { // // calculate integral(B_i''*B_j'') where B_i and B_j are // i-th and j-th basis splines. // B_i and B_j are piecewise linear functions. // v = 0; for(b=0; b<=m-2; b++) { fa = d2matrix[i,b]; fb = d2matrix[i,b+1]; ga = d2matrix[j,b]; gb = d2matrix[j,b+1]; v = v+(bx[b+1]-bx[b])*(fa*ga+(fa*(gb-ga)+ga*(fb-fa))/2+(fb-fa)*(gb-ga)/3); } amatrix[i,j] = v; amatrix[j,i] = v; } } admax = 0; for(i=0; i<=m-1; i++) { admax = Math.Max(admax, Math.Abs(amatrix[i,i])); } pdecay = lambdav*fdmax/admax; // // Calculate TDecay for Tikhonov regularization // tdecay = fdmax*(1+pdecay)*10*math.machineepsilon; // // Prepare system // // NOTE: FMatrix is spoiled during this process // for(i=0; i<=n-1; i++) { v = w[i]; for(i_=0; i_<=m-1;i_++) { fmatrix[i,i_] = v*fmatrix[i,i_]; } } ablas.rmatrixgemm(m, m, n, 1.0, fmatrix, 0, 0, 1, fmatrix, 0, 0, 0, 0.0, nmatrix, 0, 0); for(i=0; i<=m-1; i++) { for(j=0; j<=m-1; j++) { nmatrix[i,j] = nmatrix[i,j]+pdecay*amatrix[i,j]; } } for(i=0; i<=m-1; i++) { nmatrix[i,i] = nmatrix[i,i]+tdecay; } for(i=0; i<=m-1; i++) { rightpart[i] = 0; } for(i=0; i<=n-1; i++) { v = y[i]*w[i]; for(i_=0; i_<=m-1;i_++) { rightpart[i_] = rightpart[i_] + v*fmatrix[i,i_]; } } // // Solve system // if( !trfac.spdmatrixcholesky(ref nmatrix, m, true) ) { info = -4; return; } fbls.fblscholeskysolve(nmatrix, 1.0, m, true, ref rightpart, ref tmp0); for(i_=0; i_<=m-1;i_++) { c[i_] = rightpart[i_]; } // // add nodes to force linearity outside of the fitting interval // spline1d.spline1dgriddiffcubic(bx, c, m, 2, 0.0, 2, 0.0, ref bd1); tx = new double[m+2]; ty = new double[m+2]; td = new double[m+2]; i1_ = (0) - (1); for(i_=1; i_<=m;i_++) { tx[i_] = bx[i_+i1_]; } i1_ = (0) - (1); for(i_=1; i_<=m;i_++) { ty[i_] = rightpart[i_+i1_]; } i1_ = (0) - (1); for(i_=1; i_<=m;i_++) { td[i_] = bd1[i_+i1_]; } tx[0] = tx[1]-(tx[2]-tx[1]); ty[0] = ty[1]-td[1]*(tx[2]-tx[1]); td[0] = td[1]; tx[m+1] = tx[m]+(tx[m]-tx[m-1]); ty[m+1] = ty[m]+td[m]*(tx[m]-tx[m-1]); td[m+1] = td[m]; spline1d.spline1dbuildhermite(tx, ty, td, m+2, s); spline1d.spline1dlintransx(s, 2/(xb-xa), -((xa+xb)/(xb-xa))); spline1d.spline1dlintransy(s, sb-sa, sa); info = 1; // // Fill report // rep.rmserror = 0; rep.avgerror = 0; rep.avgrelerror = 0; rep.maxerror = 0; relcnt = 0; spline1d.spline1dconvcubic(bx, rightpart, m, 2, 0.0, 2, 0.0, x, n, ref fcolumn); for(i=0; i<=n-1; i++) { v = (sb-sa)*fcolumn[i]+sa; rep.rmserror = rep.rmserror+math.sqr(v-yoriginal[i]); rep.avgerror = rep.avgerror+Math.Abs(v-yoriginal[i]); if( (double)(yoriginal[i])!=(double)(0) ) { rep.avgrelerror = rep.avgrelerror+Math.Abs(v-yoriginal[i])/Math.Abs(yoriginal[i]); relcnt = relcnt+1; } rep.maxerror = Math.Max(rep.maxerror, Math.Abs(v-yoriginal[i])); } rep.rmserror = Math.Sqrt(rep.rmserror/n); rep.avgerror = rep.avgerror/n; if( (double)(relcnt)!=(double)(0) ) { rep.avgrelerror = rep.avgrelerror/relcnt; } }