69 subroutine linear_fit_yerr(nDat, xDat,yDat, yErr,  nCoef,fCoef,freeCoef,  nCov,covMat, chiSq, myFunc)
 
   73    integer, 
parameter :: mmax = 50
 
   74    integer, 
intent(in) :: ndat, ncoef,freecoef(ncoef), ncov
 
   75    real(
double), 
intent(in) :: xdat(ndat),ydat(ndat), yerr(ndat)
 
   76    real(
double), 
intent(inout) :: fcoef(ncoef)
 
   77    real(
double), 
intent(out) :: chisq, covmat(ncov,ncov)
 
   79    integer :: ii,ij,ik,il,im, nfit
 
   80    real(
double) :: err2i,bf_e2i, tot,ym, basefunc(mmax),beta(mmax)
 
   83    if(minval(abs(yerr)).lt.10*tiny(yerr)) 
call quit_program_error(
'libSUFR linear_fit_yerr(): errors cannot be zero', 0)
 
   88       if(freecoef(ij).ne.0) nfit = nfit+1
 
   90    if(nfit.eq.0) 
call warn(
'libSUFR linear_fit_yerr():  No parameters to be fit for', 0)
 
   92    covmat(1:nfit,1:nfit) = 0.d0
 
   96       call myfunc(xdat(ii),ncoef, basefunc)
 
   98       if(nfit.lt.ncoef) 
then 
  100             if(freecoef(ij).eq.0) ym = ym - fcoef(ij)*basefunc(ij)
 
  103       err2i = 1.d0/yerr(ii)**2  
 
  107          if(freecoef(il).ne.0) 
then 
  109             bf_e2i = basefunc(il)*err2i
 
  113                if(freecoef(im).ne.0) 
then 
  115                   covmat(ij,ik) = covmat(ij,ik) + bf_e2i*basefunc(im)
 
  118             beta(ij) = beta(ij) + ym*bf_e2i
 
  125          covmat(ik,ij) = covmat(ij,ik)
 
  130    call solve_linear_equations_gauss_jordan(nfit, ncov,covmat,  1, 1,beta)
 
  134       if(freecoef(il).ne.0) 
then 
  142       call myfunc(xdat(ii),ncoef, basefunc)
 
  144       tot = sum(fcoef(1:ncoef)*basefunc(1:ncoef))
 
  145       chisq = chisq  +  ( (ydat(ii)-tot) / yerr(ii) )**2
 
  149    call sort_var_covar_matrix(ncov,covmat, ncoef,freecoef, nfit)
 
 
  191  subroutine nonlin_fit_yerr(nY,nDat, xDat,yDat, yErr,  nCoef,fCoef,freeCoef,  nMat,covMat,curvMat,  chiSq, iterStat, myFunc)
 
  196    integer, 
intent(in) :: ny,ndat, ncoef,freecoef(ncoef), nmat
 
  197    real(
double), 
intent(in) :: xdat(ndat), ydat(ny,ndat), yerr(ndat)
 
  198    real(
double), 
intent(inout) :: fcoef(ncoef), iterstat
 
  199    real(
double), 
intent(out) :: curvmat(nmat,nmat), covmat(nmat,nmat), chisq
 
  201    integer, 
parameter :: mmax = 20
 
  202    integer :: ifit,mfit, icoef,cntcoef
 
  203    real(
double) :: ochisq,trycoef(mmax),dchisq(mmax),dfcoef(mmax)
 
  205    save ochisq,trycoef,dchisq,dfcoef,mfit
 
  210    if(.false.) 
call myfunc()  
 
  214    if(iterstat.lt.0.d0) 
then 
  217          if(freecoef(icoef).ne.0) mfit = mfit + 1
 
  223       call nonlin_fit_eval(ny,ndat, xdat,ydat, yerr,  ncoef,fcoef,freecoef, nmat,curvmat,  dchisq,chisq, myfunc)
 
  226       trycoef(1:ncoef) = fcoef(1:ncoef)
 
  231    covmat(1:mfit,1:mfit) = curvmat(1:mfit,1:mfit)
 
  232    dfcoef(1:mfit)        = dchisq(1:mfit)
 
  235       covmat(ifit,ifit) = curvmat(ifit,ifit) * (1.d0+iterstat)
 
  240    call solve_linear_equations_gauss_jordan(mfit, nmat,covmat,  1, 1,dfcoef)
 
  244    if(
deq0(iterstat)) 
then 
  245       call sort_var_covar_matrix(nmat, covmat,   ncoef, freecoef,  mfit)
 
  246       call sort_var_covar_matrix(nmat, curvmat,  ncoef, freecoef,  mfit)
 
  254       if(freecoef(icoef).ne.0) 
then 
  256          trycoef(icoef) = fcoef(icoef) + dfcoef(cntcoef)
 
  262    call nonlin_fit_eval(ny,ndat, xdat,ydat, yerr,  ncoef,trycoef,freecoef,  nmat,covmat,  dfcoef, chisq,  myfunc)
 
  266    if(chisq.lt.ochisq) 
then 
  267       iterstat = 0.1d0*iterstat
 
  270       curvmat(1:mfit,1:mfit) = covmat(1:mfit,1:mfit)
 
  271       dchisq(1:mfit) = dfcoef(1:mfit)
 
  273       fcoef(1:ncoef) = trycoef(1:ncoef)
 
  277       if(chisq.gt.ochisq) iterstat = 10*iterstat
 
 
  319  subroutine nonlin_fit_eval(nY,nDat, xDat,yDat, yErr,  nCoef,fCoef,freeCoef, nMat,curvMat,  dChiSq,chiSq, myFunc)
 
  322    integer, 
intent(in) :: ny,ndat, ncoef,freecoef(ncoef), nmat
 
  323    real(
double), 
intent(in) :: xdat(ndat),ydat(ny,ndat),yerr(ndat), fcoef(ncoef)
 
  324    real(
double), 
intent(out) :: chisq,curvmat(nmat,nmat),dchisq(ncoef)
 
  326    integer, 
parameter :: mmax = 20
 
  327    integer :: ifit,jfit,mfit, iy,idat, icoef,jcoef, cntcoefi,cntcoefj
 
  328    real(
double) :: err2i,dyda_e2i, dy,ymod(ny),dyda(ny,mmax)
 
  335       if(freecoef(icoef).ne.0) mfit = mfit + 1
 
  347       call myfunc(ny, xdat(idat),  ncoef,fcoef,  ymod,dyda)
 
  349       err2i = 1.d0/(yerr(idat)**2)  
 
  352          dy = ydat(iy,idat) - ymod(iy)
 
  356             if(freecoef(icoef).ne.0) 
then 
  357                cntcoefi = cntcoefi + 1
 
  358                dyda_e2i = dyda(iy,icoef)*err2i  
 
  362                   if(freecoef(jcoef).ne.0) 
then 
  363                      cntcoefj = cntcoefj + 1
 
  364                      curvmat(cntcoefi,cntcoefj) = curvmat(cntcoefi,cntcoefj) + &
 
  365                           dyda_e2i*dyda(iy,jcoef) 
 
  369                dchisq(cntcoefi) = dchisq(cntcoefi) + dy*dyda_e2i
 
  373          chisq = chisq + dy**2 * err2i  
 
  383          curvmat(jfit,ifit) = curvmat(ifit,jfit)
 
  387  end subroutine nonlin_fit_eval
 
  411    integer, 
intent(in) :: ny, ncoef
 
  412    real(
double), 
intent(in) :: xdat, fcoef(ncoef)
 
  413    real(
double), 
intent(out) :: ydat(ny), dyda(ny,ncoef)
 
  415    if(ny.ne.1)    
call quit_program_error(
'nonlin_fit_eval_example_myFunc():  nY must be equal to 1 ',    1)
 
  416    if(ncoef.ne.3) 
call quit_program_error(
'nonlin_fit_eval_example_myFunc():  nCoef must be equal to 3 ', 1)
 
  418    ydat(1)   = fcoef(1)*xdat**2 + fcoef(2)   
 
  419    dyda(1,1) = fcoef(2) * xdat               
 
  421    dyda(1,ncoef) = fcoef(1)                  
 
 
  441  pure subroutine sort_var_covar_matrix(nCov,covMat, nCoef,freeCoef, nFit)
 
  446    integer, 
intent(in) :: ncov,ncoef, freecoef(ncoef),nfit
 
  447    real(
double), 
intent(inout) :: covmat(ncov,ncov)
 
  461       if(freecoef(ij).ne.0) 
then 
  463             call swapdbl(covmat(ii,ij), covmat(ii,ik))  
 
  467             call swapdbl(covmat(ij,ii), covmat(ik,ii))  
 
  474  end subroutine sort_var_covar_matrix
 
  493  subroutine solve_linear_equations_gauss_jordan(nMat, nMatArr,matArr,  nVec, nVecArr,vecArr)
 
  498    integer, 
intent(in) :: nmat,nmatarr, nvec,nvecarr
 
  499    real(
double), 
intent(inout) :: matarr(nmatarr,nmatarr),vecarr(nmatarr,nvecarr)
 
  501    integer, 
parameter :: nmax = 50  
 
  502    integer :: ii,icol,irow,ij,ik,il,ll,indxc(nmax),indxr(nmax),ipiv(nmax)
 
  503    real(
double) :: arrmax,tmpdbl,pivinv
 
  515          if(ipiv(ij).ne.1) 
then 
  517                if(ipiv(ik).eq.0) 
then 
  518                   if(abs(matarr(ij,ik)).ge.arrmax) 
then 
  519                      arrmax = abs(matarr(ij,ik))
 
  523                else if(ipiv(ik).gt.1) 
then 
  524                   call warn(
'libSUFR solve_linear_equations_Gauss_Jordan():  singular matrix',0)
 
  530       ipiv(icol) = ipiv(icol) + 1
 
  532       if(irow.ne.icol) 
then 
  534             call swapdbl(matarr(irow,il), matarr(icol,il))
 
  538             call swapdbl(vecarr(irow,il), vecarr(icol,il))
 
  545       if(abs(matarr(icol,icol)).lt.10*tiny(matarr)) 
call warn(
'libSUFR solve_linear_equations_Gauss_Jordan():  singular matrix', 0)
 
  547       pivinv = 1.d0/matarr(icol,icol)
 
  548       matarr(icol,icol) = 1.d0
 
  550       matarr(icol,1:nmat) = matarr(icol,1:nmat)*pivinv
 
  551       vecarr(icol,1:nvec) = vecarr(icol,1:nvec)*pivinv
 
  555             tmpdbl = matarr(ll,icol)
 
  556             matarr(ll,icol) = 0.d0
 
  557             matarr(ll,1:nmat) = matarr(ll,1:nmat) - matarr(icol,1:nmat)*tmpdbl
 
  558             vecarr(ll,1:nvec) = vecarr(ll,1:nvec) - vecarr(icol,1:nvec)*tmpdbl
 
  564       if(indxr(il).ne.indxc(il)) 
then 
  566             call swapdbl(matarr(ik,indxr(il)), matarr(ik,indxc(il)))
 
  571  end subroutine solve_linear_equations_gauss_jordan
 
  588    integer, 
intent(in) :: ncoef
 
  589    real(
double), 
intent(in) :: xval
 
  590    real(
double), 
intent(out) :: basefunc(ncoef)
 
  596       basefunc(coefi) = basefunc(coefi-1) * xval
 
 
Procedures to fit functions to data.
 
subroutine, public linear_fit_yerr(ndat, xdat, ydat, yerr, ncoef, fcoef, freecoef, ncov, covmat, chisq, myfunc)
Fit the model coefficients using a chi^2 method.
 
subroutine, public nonlin_fit_example_myfunc(ny, xdat, ncoef, fcoef, ydat, dyda)
Dummy example function myFunc for nonlin_fit_yerr(): return the value and partial derivatives.
 
subroutine, public nonlin_fit_yerr(ny, ndat, xdat, ydat, yerr, ncoef, fcoef, freecoef, nmat, covmat, curvmat, chisq, iterstat, myfunc)
Levenberg-Marquardt method to reduce the value of chi-squared of a fit of a set of data points with e...
 
pure subroutine, public basefunc_polynomial(xval, ncoef, basefunc)
User-provided base function to fit a polynomial to.
 
Provides kinds and related constants/routines.
 
integer, parameter double
Double-precision float. Precision = 15, range = 307.
 
Procedures for numerical operations.
 
elemental logical function deq0(x0, eps)
Test whether a double-precision variable is equal to zero better than a given value (default: 2x mach...
 
System-related procedures.
 
subroutine quit_program_error(message, status)
Print an error message to StdErr and stop the execution of the current program.
 
subroutine warn(message, unit)
Print a warning to StdOut or StErr.
 
pure subroutine swapdbl(dbl1, dbl2)
Swap two double-precision real variables.