00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "Minuit2/Numerical2PGradientCalculator.h"
00011 #include "Minuit2/InitialGradientCalculator.h"
00012 #include "Minuit2/MnFcn.h"
00013 #include "Minuit2/MnUserTransformation.h"
00014 #include "Minuit2/MnMachinePrecision.h"
00015 #include "Minuit2/MinimumParameters.h"
00016 #include "Minuit2/FunctionGradient.h"
00017 #include "Minuit2/MnStrategy.h"
00018
00019
00020
00021 #if defined(DEBUG) || defined(WARNINGMSG)
00022 #include "Minuit2/MnPrint.h"
00023 #ifdef _OPENMP
00024 #include <omp.h>
00025 #include <iomanip>
00026 #ifdef DEBUG
00027 #define DEBUG_MP
00028 #endif
00029 #endif
00030 #endif
00031
00032 #include <math.h>
00033
00034 #include "Minuit2/MPIProcess.h"
00035
00036 namespace ROOT {
00037
00038 namespace Minuit2 {
00039
00040
00041 FunctionGradient Numerical2PGradientCalculator::operator()(const MinimumParameters& par) const {
00042
00043
00044 InitialGradientCalculator gc(fFcn, fTransformation, fStrategy);
00045 FunctionGradient gra = gc(par);
00046
00047 return (*this)(par, gra);
00048 }
00049
00050
00051
00052 FunctionGradient Numerical2PGradientCalculator::operator()(const std::vector<double>& params) const {
00053
00054
00055 int npar = params.size();
00056
00057 MnAlgebraicVector par(npar);
00058 for (int i = 0; i < npar; ++i) {
00059 par(i) = params[i];
00060 }
00061
00062 double fval = Fcn()(par);
00063
00064 MinimumParameters minpars = MinimumParameters(par, fval);
00065
00066 return (*this)(minpars);
00067
00068 }
00069
00070
00071
00072 FunctionGradient Numerical2PGradientCalculator::operator()(const MinimumParameters& par, const FunctionGradient& Gradient) const {
00073
00074
00075
00076
00077
00078
00079
00080 assert(par.IsValid());
00081
00082
00083 double fcnmin = par.Fval();
00084
00085
00086 double eps2 = Precision().Eps2();
00087 double eps = Precision().Eps();
00088
00089 double dfmin = 8.*eps2*(fabs(fcnmin)+Fcn().Up());
00090 double vrysml = 8.*eps*eps;
00091
00092
00093
00094
00095
00096 unsigned int n = (par.Vec()).size();
00097 unsigned int ncycle = Ncycle();
00098
00099 MnAlgebraicVector grd = Gradient.Grad();
00100 MnAlgebraicVector g2 = Gradient.G2();
00101 MnAlgebraicVector gstep = Gradient.Gstep();
00102
00103 #ifndef _OPENMP
00104 MPIProcess mpiproc(n,0);
00105 #endif
00106
00107 #ifdef DEBUG
00108 std::cout << "Calculating Gradient at x = " << par.Vec() << std::endl;
00109 #endif
00110
00111 #ifndef _OPENMP
00112
00113 MnAlgebraicVector x = par.Vec();
00114
00115 unsigned int startElementIndex = mpiproc.StartElementIndex();
00116 unsigned int endElementIndex = mpiproc.EndElementIndex();
00117
00118 for(unsigned int i = startElementIndex; i < endElementIndex; i++) {
00119
00120 #else
00121
00122
00123
00124 #pragma omp parallel
00125 #pragma omp for
00126
00127
00128 for(int i = 0; i < int(n); i++) {
00129
00130 #endif
00131
00132 #ifdef DEBUG_MP
00133 int ith = omp_get_thread_num();
00134
00135 #endif
00136
00137 #ifdef _OPENMP
00138
00139 MnAlgebraicVector x = par.Vec();
00140 #endif
00141
00142 double xtf = x(i);
00143 double epspri = eps2 + fabs(grd(i)*eps2);
00144 double stepb4 = 0.;
00145 for(unsigned int j = 0; j < ncycle; j++) {
00146 double optstp = sqrt(dfmin/(fabs(g2(i))+epspri));
00147 double step = std::max(optstp, fabs(0.1*gstep(i)));
00148
00149 if(Trafo().Parameter(Trafo().ExtOfInt(i)).HasLimits()) {
00150 if(step > 0.5) step = 0.5;
00151 }
00152 double stpmax = 10.*fabs(gstep(i));
00153 if(step > stpmax) step = stpmax;
00154
00155 double stpmin = std::max(vrysml, 8.*fabs(eps2*x(i)));
00156 if(step < stpmin) step = stpmin;
00157
00158
00159 if(fabs((step-stepb4)/step) < StepTolerance()) {
00160
00161
00162
00163 break;
00164 }
00165 gstep(i) = step;
00166 stepb4 = step;
00167
00168
00169
00170
00171
00172 x(i) = xtf + step;
00173 double fs1 = Fcn()(x);
00174 x(i) = xtf - step;
00175 double fs2 = Fcn()(x);
00176 x(i) = xtf;
00177
00178 double grdb4 = grd(i);
00179 grd(i) = 0.5*(fs1 - fs2)/step;
00180 g2(i) = (fs1 + fs2 - 2.*fcnmin)/step/step;
00181
00182 if(fabs(grdb4-grd(i))/(fabs(grd(i))+dfmin/step) < GradTolerance()) {
00183
00184
00185
00186
00187 break;
00188 }
00189 }
00190
00191
00192 #ifdef DEBUG_MP
00193 #pragma omp critical
00194 {
00195 std::cout << "Gradient for thread " << ith << " " << i << " " << std::setprecision(15) << grd(i) << " " << g2(i) << std::endl;
00196 }
00197 #endif
00198
00199
00200
00201
00202 }
00203
00204 #ifdef DEBUG
00205 std::cout << "Gradient = " << grd << std::endl;
00206 #endif
00207
00208
00209
00210
00211 #ifndef _OPENMP
00212 mpiproc.SyncVector(grd);
00213 mpiproc.SyncVector(g2);
00214 mpiproc.SyncVector(gstep);
00215 #endif
00216
00217 return FunctionGradient(grd, g2, gstep);
00218 }
00219
00220 const MnMachinePrecision& Numerical2PGradientCalculator::Precision() const {
00221
00222 return fTransformation.Precision();
00223 }
00224
00225 unsigned int Numerical2PGradientCalculator::Ncycle() const {
00226
00227 return Strategy().GradientNCycles();
00228 }
00229
00230 double Numerical2PGradientCalculator::StepTolerance() const {
00231
00232 return Strategy().GradientStepTolerance();
00233 }
00234
00235 double Numerical2PGradientCalculator::GradTolerance() const {
00236
00237 return Strategy().GradientTolerance();
00238 }
00239
00240
00241 }
00242
00243 }