1 /*
2 * This file is part of the GROMACS molecular simulation package.
3 *
4 * Copyright (c) 2016,2018,2019, 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.
8 *
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.
13 *
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.
18 *
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.
23 *
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.
31 *
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.
34 */
35 /*! \internal
36 * \file
37 * \brief
38 * Defines a driver routine for lmfit, and a callback for it to use.
39 *
40 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
41 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \ingroup module_correlationfunctions
43 */
44 #include "gmxpre.h"
45
46 #include "gmx_lmcurve.h"
47
48 #include "config.h"
49
50 #include <cmath>
51
52 #if HAVE_LMFIT
53 # include <lmmin.h>
54 # include <lmstruct.h>
55 #endif
56
57 #include "gromacs/correlationfunctions/expfit.h"
58 #include "gromacs/math/functions.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/smalloc.h"
61
62 #if HAVE_LMFIT
63
64 typedef struct
65 {
66 const double* t;
67 const double* y;
68 const double* dy;
69 double (*f)(const double t, const double* par);
70 } lmcurve_data_struct;
71
72 //! Callback function used by lmmin
lmcurve_evaluate(const double * par,const int m_dat,const void * data,double * fvec,int * info)73 static void lmcurve_evaluate(const double* par, const int m_dat, const void* data, double* fvec, int* info)
74 {
75 const lmcurve_data_struct* D = reinterpret_cast<const lmcurve_data_struct*>(data);
76 for (int i = 0; i < m_dat; i++)
77 {
78 double dy = D->dy[i];
79 if (dy == 0)
80 {
81 dy = 1;
82 }
83 fvec[i] = (D->y[i] - D->f(D->t[i], par)) / dy;
84 }
85 *info = 0;
86 }
87
88 //! Calls lmmin with the given data, with callback function \c f.
gmx_lmcurve(const int n_par,double * par,const int m_dat,const double * t,const double * y,const double * dy,double (* f)(double t,const double * par),const lm_control_struct * control,lm_status_struct * status)89 static void gmx_lmcurve(const int n_par,
90 double* par,
91 const int m_dat,
92 const double* t,
93 const double* y,
94 const double* dy,
95 double (*f)(double t, const double* par),
96 const lm_control_struct* control,
97 lm_status_struct* status)
98 {
99 lmcurve_data_struct data = { t, y, dy, f };
100
101 lmmin(n_par, par, m_dat, nullptr, &data, lmcurve_evaluate, control, status);
102 }
103
104 #endif
105
lmfit_exp(int nfit,const double x[],const double y[],const double dy[],double parm[],bool bVerbose,int eFitFn,int nfix)106 bool lmfit_exp(int nfit,
107 const double x[],
108 const double y[],
109 const double dy[],
110 double parm[], // NOLINT(readability-non-const-parameter)
111 bool bVerbose,
112 int eFitFn,
113 int nfix)
114 {
115 if ((eFitFn < 0) || (eFitFn >= effnNR))
116 {
117 fprintf(stderr, "fitfn = %d, should be in the range 0..%d\n", eFitFn, effnNR - 1);
118 return false;
119 }
120 #if HAVE_LMFIT
121 double chisq, ochisq;
122 gmx_bool bCont;
123 int j;
124 int maxiter = 100;
125 lm_control_struct control;
126 lm_status_struct* status;
127 int nparam = effnNparams(eFitFn);
128 int p2;
129 gmx_bool bSkipLast;
130
131 /* Using default control structure for double precision fitting that
132 * comes with the lmfit package (i.e. from the include file).
133 */
134 control = lm_control_double;
135 control.verbosity = (bVerbose ? 1 : 0);
136 control.n_maxpri = 0;
137 control.m_maxpri = 0;
138
139 snew(status, 1);
140 /* Initial params */
141 chisq = 1e12;
142 j = 0;
143 if (bVerbose)
144 {
145 printf("%4s %10s Parameters\n", "Step", "chi^2");
146 }
147 /* Check whether we have to skip some params */
148 if (nfix > 0)
149 {
150 do
151 {
152 p2 = 1 << (nparam - 1);
153 bSkipLast = ((p2 & nfix) == p2);
154 if (bSkipLast)
155 {
156 nparam--;
157 nfix -= p2;
158 }
159 } while ((nparam > 0) && (bSkipLast));
160 if (bVerbose)
161 {
162 printf("Using %d out of %d parameters\n", nparam, effnNparams(eFitFn));
163 }
164 }
165 do
166 {
167 ochisq = chisq;
168 gmx_lmcurve(nparam, parm, nfit, x, y, dy, lmcurves[eFitFn], &control, status);
169 chisq = gmx::square(status->fnorm);
170 if (bVerbose)
171 {
172 printf("status: fnorm = %g, nfev = %d, userbreak = %d\noutcome = %s\n", status->fnorm,
173 status->nfev, status->userbreak, lm_infmsg[status->outcome]);
174 }
175 if (bVerbose)
176 {
177 int mmm;
178 printf("%4d %8g", j, chisq);
179 for (mmm = 0; (mmm < effnNparams(eFitFn)); mmm++)
180 {
181 printf(" %8g", parm[mmm]);
182 }
183 printf("\n");
184 }
185 j++;
186 bCont = (fabs(ochisq - chisq) > fabs(control.ftol * chisq));
187 } while (bCont && (j < maxiter));
188
189 sfree(status);
190 #else
191 gmx_fatal(FARGS,
192 "This build of GROMACS was not configured with support "
193 "for lmfit, so the requested fitting cannot be performed. "
194 "See the install guide for instructions on how to build "
195 "GROMACS with lmfit supported.");
196 GMX_UNUSED_VALUE(nfit);
197 GMX_UNUSED_VALUE(x);
198 GMX_UNUSED_VALUE(y);
199 GMX_UNUSED_VALUE(dy);
200 GMX_UNUSED_VALUE(parm);
201 GMX_UNUSED_VALUE(bVerbose);
202 GMX_UNUSED_VALUE(eFitFn);
203 GMX_UNUSED_VALUE(nfix);
204 #endif
205 return true;
206 }
207