2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
38 * Defines a driver routine for lmfit, and a callback for it to use.
40 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
41 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \ingroup module_correlationfunctions
46 #include "gmx_lmcurve.h"
57 #include "gromacs/correlationfunctions/expfit.h"
58 #include "gromacs/math/functions.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/smalloc.h"
68 double (*f)(const double t, const double* par);
69 } lmcurve_data_struct;
71 //! Callback function used by lmmin
72 static void lmcurve_evaluate(
73 const double* par, const int m_dat, const void* data, double* fvec,
76 const lmcurve_data_struct* D = reinterpret_cast<const lmcurve_data_struct*>(data);
77 for (int i = 0; i < m_dat; i++)
84 fvec[i] = (D->y[i] - D->f(D->t[i], par))/dy;
89 //! Calls lmmin with the given data, with callback function \c f.
90 static void gmx_lmcurve(
91 const int n_par, double* par, const int m_dat,
92 const double* t, const double* y, const double *dy,
93 double (*f)(double t, const double* par),
94 const lm_control_struct* control, lm_status_struct* status)
96 lmcurve_data_struct data = { t, y, dy, f };
98 lmmin(n_par, par, m_dat, nullptr, &data, lmcurve_evaluate,
104 bool lmfit_exp(int nfit,
108 double parm[], // NOLINT(readability-non-const-parameter)
113 if ((eFitFn < 0) || (eFitFn >= effnNR))
115 fprintf(stderr, "fitfn = %d, should be in the range 0..%d\n",
120 double chisq, ochisq;
124 lm_control_struct control;
125 lm_status_struct *status;
126 int nparam = effnNparams(eFitFn);
130 /* Using default control structure for double precision fitting that
131 * comes with the lmfit package (i.e. from the include file).
133 control = lm_control_double;
134 control.verbosity = (bVerbose ? 1 : 0);
135 control.n_maxpri = 0;
136 control.m_maxpri = 0;
144 printf("%4s %10s Parameters\n", "Step", "chi^2");
146 /* Check whether we have to skip some params */
151 p2 = 1 << (nparam-1);
152 bSkipLast = ((p2 & nfix) == p2);
159 while ((nparam > 0) && (bSkipLast));
162 printf("Using %d out of %d parameters\n", nparam, effnNparams(eFitFn));
168 gmx_lmcurve(nparam, parm, nfit, x, y, dy,
169 lmcurves[eFitFn], &control, status);
170 chisq = gmx::square(status->fnorm);
173 printf("status: fnorm = %g, nfev = %d, userbreak = %d\noutcome = %s\n",
174 status->fnorm, status->nfev, status->userbreak,
175 lm_infmsg[status->outcome]);
180 printf("%4d %8g", j, chisq);
181 for (mmm = 0; (mmm < effnNparams(eFitFn)); mmm++)
183 printf(" %8g", parm[mmm]);
188 bCont = (fabs(ochisq - chisq) > fabs(control.ftol*chisq));
190 while (bCont && (j < maxiter));
194 gmx_fatal(FARGS, "This build of GROMACS was not configured with support "
195 "for lmfit, so the requested fitting cannot be performed. "
196 "See the install guide for instructions on how to build "
197 "GROMACS with lmfit supported.");
198 GMX_UNUSED_VALUE(nfit);
201 GMX_UNUSED_VALUE(dy);
202 GMX_UNUSED_VALUE(parm);
203 GMX_UNUSED_VALUE(bVerbose);
204 GMX_UNUSED_VALUE(eFitFn);
205 GMX_UNUSED_VALUE(nfix);