8b3870e958dda3a06bd374188df77849fdf28f26
[alexxy/gromacs.git] / src / gromacs / correlationfunctions / gmx_lmcurve.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 extern t_lmcurve lmcurves[effnNR+1];
65
66 typedef struct {
67     const double* t;
68     const double* y;
69     const double* dy;
70     double (*f)(const double t, const double* par);
71 } lmcurve_data_struct;
72
73 //! Callback function used by lmmin
74 static void lmcurve_evaluate(
75         const double* par, const int m_dat, const void* data, double* fvec,
76         int* info)
77 {
78     lmcurve_data_struct* D = (lmcurve_data_struct*)data;
79     int                  i;
80     for (i = 0; i < m_dat; i++)
81     {
82         double dy = D->dy[i];
83         if (dy == 0)
84         {
85             dy = 1;
86         }
87         fvec[i] = (D->y[i] - D->f(D->t[i], par))/dy;
88     }
89     *info = 0;
90 }
91
92 //! Calls lmmin with the given data, with callback function \c f.
93 static void gmx_lmcurve(
94         const int n_par, double* par, const int m_dat,
95         const double* t, const double* y, const double *dy,
96         double (*f)(double t, const double* par),
97         const lm_control_struct* control, lm_status_struct* status)
98 {
99     lmcurve_data_struct data = { t, y, dy, f };
100
101     lmmin(n_par, par, m_dat, nullptr, (const void*)&data, lmcurve_evaluate,
102           control, status);
103 }
104
105 #endif
106
107 bool lmfit_exp(int          nfit,
108                const double x[],
109                const double y[],
110                const double dy[],
111                double       parm[],
112                bool         bVerbose,
113                int          eFitFn,
114                int          nfix)
115 {
116     if ((eFitFn < 0) || (eFitFn >= effnNR))
117     {
118         fprintf(stderr, "fitfn = %d, should be in the range 0..%d\n",
119                 eFitFn, effnNR-1);
120         return false;
121     }
122 #if HAVE_LMFIT
123     double             chisq, ochisq;
124     gmx_bool           bCont;
125     int                j;
126     int                maxiter = 100;
127     lm_control_struct  control;
128     lm_status_struct  *status;
129     int                nparam = effnNparams(eFitFn);
130     int                p2;
131     gmx_bool           bSkipLast;
132
133     /* Using default control structure for double precision fitting that
134      * comes with the lmfit package (i.e. from the include file).
135      */
136     control            = lm_control_double;
137     control.verbosity  = (bVerbose ? 1 : 0);
138     control.n_maxpri   = 0;
139     control.m_maxpri   = 0;
140
141     snew(status, 1);
142     /* Initial params */
143     chisq  = 1e12;
144     j      = 0;
145     if (bVerbose)
146     {
147         printf("%4s  %10s  Parameters\n", "Step", "chi^2");
148     }
149     /* Check whether we have to skip some params */
150     if (nfix > 0)
151     {
152         do
153         {
154             p2        = 1 << (nparam-1);
155             bSkipLast = ((p2 & nfix) == p2);
156             if (bSkipLast)
157             {
158                 nparam--;
159                 nfix -= p2;
160             }
161         }
162         while ((nparam > 0) && (bSkipLast));
163         if (bVerbose)
164         {
165             printf("Using %d out of %d parameters\n", nparam, effnNparams(eFitFn));
166         }
167     }
168     do
169     {
170         ochisq = chisq;
171         gmx_lmcurve(nparam, parm, nfit, x, y, dy,
172                     lmcurves[eFitFn], &control, status);
173         chisq = gmx::square(status->fnorm);
174         if (bVerbose)
175         {
176             printf("status: fnorm = %g, nfev = %d, userbreak = %d\noutcome = %s\n",
177                    status->fnorm, status->nfev, status->userbreak,
178                    lm_infmsg[status->outcome]);
179         }
180         if (bVerbose)
181         {
182             int mmm;
183             printf("%4d  %8g", j, chisq);
184             for (mmm = 0; (mmm < effnNparams(eFitFn)); mmm++)
185             {
186                 printf("  %8g", parm[mmm]);
187             }
188             printf("\n");
189         }
190         j++;
191         bCont = (fabs(ochisq - chisq) > fabs(control.ftol*chisq));
192     }
193     while (bCont && (j < maxiter));
194
195     sfree(status);
196 #else
197     gmx_fatal(FARGS, "This build of GROMACS was not configured with support "
198               "for lmfit, so the requested fitting cannot be performed. "
199               "See the install guide for instructions on how to build "
200               "GROMACS with lmfit supported.");
201     GMX_UNUSED_VALUE(nfit);
202     GMX_UNUSED_VALUE(x);
203     GMX_UNUSED_VALUE(y);
204     GMX_UNUSED_VALUE(dy);
205     GMX_UNUSED_VALUE(parm);
206     GMX_UNUSED_VALUE(bVerbose);
207     GMX_UNUSED_VALUE(eFitFn);
208     GMX_UNUSED_VALUE(nfix);
209 #endif
210     return true;
211 }