3/3 of old-style casting
[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 typedef struct {
65     const double* t;
66     const double* y;
67     const double* dy;
68     double (*f)(const double t, const double* par);
69 } lmcurve_data_struct;
70
71 //! Callback function used by lmmin
72 static void lmcurve_evaluate(
73         const double* par, const int m_dat, const void* data, double* fvec,
74         int* info)
75 {
76     const lmcurve_data_struct* D = reinterpret_cast<const lmcurve_data_struct*>(data);
77     for (int i = 0; i < m_dat; i++)
78     {
79         double dy = D->dy[i];
80         if (dy == 0)
81         {
82             dy = 1;
83         }
84         fvec[i] = (D->y[i] - D->f(D->t[i], par))/dy;
85     }
86     *info = 0;
87 }
88
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)
95 {
96     lmcurve_data_struct data = { t, y, dy, f };
97
98     lmmin(n_par, par, m_dat, nullptr, &data, lmcurve_evaluate,
99           control, status);
100 }
101
102 #endif
103
104 bool lmfit_exp(int          nfit,
105                const double x[],
106                const double y[],
107                const double dy[],
108                double       parm[], // NOLINT(readability-non-const-parameter)
109                bool         bVerbose,
110                int          eFitFn,
111                int          nfix)
112 {
113     if ((eFitFn < 0) || (eFitFn >= effnNR))
114     {
115         fprintf(stderr, "fitfn = %d, should be in the range 0..%d\n",
116                 eFitFn, effnNR-1);
117         return false;
118     }
119 #if HAVE_LMFIT
120     double             chisq, ochisq;
121     gmx_bool           bCont;
122     int                j;
123     int                maxiter = 100;
124     lm_control_struct  control;
125     lm_status_struct  *status;
126     int                nparam = effnNparams(eFitFn);
127     int                p2;
128     gmx_bool           bSkipLast;
129
130     /* Using default control structure for double precision fitting that
131      * comes with the lmfit package (i.e. from the include file).
132      */
133     control            = lm_control_double;
134     control.verbosity  = (bVerbose ? 1 : 0);
135     control.n_maxpri   = 0;
136     control.m_maxpri   = 0;
137
138     snew(status, 1);
139     /* Initial params */
140     chisq  = 1e12;
141     j      = 0;
142     if (bVerbose)
143     {
144         printf("%4s  %10s  Parameters\n", "Step", "chi^2");
145     }
146     /* Check whether we have to skip some params */
147     if (nfix > 0)
148     {
149         do
150         {
151             p2        = 1 << (nparam-1);
152             bSkipLast = ((p2 & nfix) == p2);
153             if (bSkipLast)
154             {
155                 nparam--;
156                 nfix -= p2;
157             }
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,
169                     lmcurves[eFitFn], &control, status);
170         chisq = gmx::square(status->fnorm);
171         if (bVerbose)
172         {
173             printf("status: fnorm = %g, nfev = %d, userbreak = %d\noutcome = %s\n",
174                    status->fnorm, status->nfev, status->userbreak,
175                    lm_infmsg[status->outcome]);
176         }
177         if (bVerbose)
178         {
179             int mmm;
180             printf("%4d  %8g", j, chisq);
181             for (mmm = 0; (mmm < effnNparams(eFitFn)); mmm++)
182             {
183                 printf("  %8g", parm[mmm]);
184             }
185             printf("\n");
186         }
187         j++;
188         bCont = (fabs(ochisq - chisq) > fabs(control.ftol*chisq));
189     }
190     while (bCont && (j < maxiter));
191
192     sfree(status);
193 #else
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);
199     GMX_UNUSED_VALUE(x);
200     GMX_UNUSED_VALUE(y);
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);
206 #endif
207     return true;
208 }