lsst.meas.modelfit  13.0-11-g2fa83af+13
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Pages
optimizer.cc
Go to the documentation of this file.
1 // -*- lsst-c++ -*-
2 /*
3  * LSST Data Management System
4  * Copyright 2008-2013 LSST Corporation.
5  *
6  * This product includes software developed by the
7  * LSST Project (http://www.lsst.org/).
8  *
9  * This program is free software: you can redistribute it and/or modify
10  * it under the terms of the GNU General Public License as published by
11  * the Free Software Foundation, either version 3 of the License, or
12  * (at your option) any later version.
13  *
14  * This program 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
17  * GNU General Public License for more details.
18  *
19  * You should have received a copy of the LSST License Statement and
20  * the GNU General Public License along with this program. If not,
21  * see <http://www.lsstcorp.org/LegalNotices/>.
22  */
23 
24 #include "pybind11/pybind11.h"
25 
26 #include "numpy/arrayobject.h"
27 #include "ndarray/pybind11.h"
28 
29 #include "lsst/pex/config/python.h"
30 
35 #include "lsst/afw/table/Catalog.h"
36 
37 namespace py = pybind11;
38 using namespace pybind11::literals;
39 
40 namespace lsst {
41 namespace meas {
42 namespace modelfit {
43 namespace {
44 
45 using PyOptimizerObjective = py::class_<OptimizerObjective, std::shared_ptr<OptimizerObjective>>;
46 using PyOptimizerControl = py::class_<OptimizerControl, std::shared_ptr<OptimizerControl>>;
47 using PyOptimizerHistoryRecorder =
48  py::class_<OptimizerHistoryRecorder, std::shared_ptr<OptimizerHistoryRecorder>>;
49 using PyOptimizer = py::class_<Optimizer, std::shared_ptr<Optimizer>>;
50 
51 static PyOptimizerObjective declareOptimizerObjective(py::module &mod) {
52  PyOptimizerObjective cls(mod, "OptimizerObjective");
53  // Class is abstract, so no constructor.
54  cls.def_readonly("dataSize", &OptimizerObjective::dataSize);
55  cls.def_readonly("parameterSize", &OptimizerObjective::parameterSize);
56  cls.def_static("makeFromLikelihood", &OptimizerObjective::makeFromLikelihood, "likelihood"_a,
57  "prior"_a = nullptr);
58  // class is abstract and not subclassable in Python, so we don't wrap the ctor
59  cls.def("fillObjectiveValueGrid", &OptimizerObjective::fillObjectiveValueGrid, "parameters"_a,
60  "output"_a);
61  cls.def("computeResiduals", &OptimizerObjective::computeResiduals, "parameters"_a, "residuals"_a);
62  cls.def("differentiateResiduals", &OptimizerObjective::differentiateResiduals, "parameters"_a,
63  "derivatives"_a);
64  cls.def("hasPrior", &OptimizerObjective::hasPrior);
65  cls.def("computePrior", &OptimizerObjective::computePrior, "parameters"_a);
66  cls.def("differentiatePrior", &OptimizerObjective::differentiatePrior, "parameters"_a, "gradient"_a,
67  "hessian"_a);
68  return cls;
69 }
70 
71 static PyOptimizerControl declareOptimizerControl(py::module &mod) {
72  PyOptimizerControl cls(mod, "OptimizerControl");
73  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, noSR1Term);
74  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, skipSR1UpdateThreshold);
75  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, minTrustRadiusThreshold);
76  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, gradientThreshold);
77  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, numDiffRelStep);
78  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, numDiffAbsStep);
79  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, numDiffTrustRadiusStep);
80  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, stepAcceptThreshold);
81  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, trustRegionInitialSize);
82  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, trustRegionGrowReductionRatio);
83  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, trustRegionGrowStepFraction);
84  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, trustRegionGrowFactor);
85  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, trustRegionShrinkReductionRatio);
86  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, trustRegionShrinkFactor);
87  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, trustRegionSolverTolerance);
88  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, maxInnerIterations);
89  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, maxOuterIterations);
90  LSST_DECLARE_CONTROL_FIELD(cls, OptimizerControl, doSaveIterations);
91  cls.def(py::init<>());
92  return cls;
93 }
94 
95 static PyOptimizerHistoryRecorder declareOptimizerHistoryRecorder(py::module &mod) {
96  PyOptimizerHistoryRecorder cls(mod, "OptimizerHistoryRecorder");
97  cls.def(py::init<afw::table::Schema &, std::shared_ptr<Model>, bool>(), "schema"_a, "model"_a,
98  "doRecordDerivatives"_a);
99  cls.def(py::init<afw::table::Schema const &>(), "schema"_a);
100  cls.def("apply", &OptimizerHistoryRecorder::apply, "outerIterCount"_a, "innerIterCount"_a, "history"_a,
101  "optimizer"_a);
102  cls.def("unpackDerivatives",
103  (void (OptimizerHistoryRecorder::*)(ndarray::Array<Scalar const, 1, 1> const &,
104  ndarray::Array<Scalar, 1, 1> const &,
105  ndarray::Array<Scalar, 2, 2> const &) const) &
106  OptimizerHistoryRecorder::unpackDerivatives,
107  "nested"_a, "gradient"_a, "hessian"_a);
108  cls.def("unpackDerivatives", (void (OptimizerHistoryRecorder::*)(
109  afw::table::BaseRecord const &, ndarray::Array<Scalar, 1, 1> const &,
110  ndarray::Array<Scalar, 2, 2> const &) const) &
111  OptimizerHistoryRecorder::unpackDerivatives,
112  "record"_a, "gradient"_a, "hessian"_a);
113  // Other unpackDerivatives overloads do the same thing but with Eigen types,
114  // which makes them redundant in Python where it's all just NumPy.
115  cls.def("fillObjectiveModelGrid", &OptimizerHistoryRecorder::fillObjectiveModelGrid, "record"_a,
116  "parameters"_a, "output"_a);
117  cls.def_readonly("outer", &OptimizerHistoryRecorder::outer);
118  cls.def_readonly("inner", &OptimizerHistoryRecorder::inner);
119  cls.def_readonly("state", &OptimizerHistoryRecorder::state);
120  cls.def_readonly("objective", &OptimizerHistoryRecorder::objective);
121  cls.def_readonly("prior", &OptimizerHistoryRecorder::prior);
122  cls.def_readonly("trust", &OptimizerHistoryRecorder::trust);
123  cls.def_readonly("parameters", &OptimizerHistoryRecorder::parameters);
124  cls.def_readonly("derivatives", &OptimizerHistoryRecorder::derivatives);
125  return cls;
126 }
127 
128 static PyOptimizer declareOptimizer(py::module &mod) {
129  PyOptimizer cls(mod, "Optimizer");
130  // StateFlags enum is used as bitflag, so we wrap values as int class attributes.
131  cls.attr("CONVERGED_GRADZERO") = py::cast(int(Optimizer::CONVERGED_GRADZERO));
132  cls.attr("CONVERGED_TR_SMALL") = py::cast(int(Optimizer::CONVERGED_TR_SMALL));
133  cls.attr("CONVERGED") = py::cast(int(Optimizer::CONVERGED));
134  cls.attr("FAILED_MAX_INNER_ITERATIONS") = py::cast(int(Optimizer::FAILED_MAX_INNER_ITERATIONS));
135  cls.attr("FAILED_MAX_OUTER_ITERATIONS") = py::cast(int(Optimizer::FAILED_MAX_OUTER_ITERATIONS));
136  cls.attr("FAILED_MAX_ITERATIONS") = py::cast(int(Optimizer::FAILED_MAX_ITERATIONS));
137  cls.attr("FAILED_EXCEPTION") = py::cast(int(Optimizer::FAILED_EXCEPTION));
138  cls.attr("FAILED_NAN") = py::cast(int(Optimizer::FAILED_NAN));
139  cls.attr("FAILED") = py::cast(int(Optimizer::FAILED));
140  cls.attr("STATUS_STEP_REJECTED") = py::cast(int(Optimizer::STATUS_STEP_REJECTED));
141  cls.attr("STATUS_STEP_ACCEPTED") = py::cast(int(Optimizer::STATUS_STEP_ACCEPTED));
142  cls.attr("STATUS_STEP") = py::cast(int(Optimizer::STATUS_STEP));
143  cls.attr("STATUS_TR_UNCHANGED") = py::cast(int(Optimizer::STATUS_TR_UNCHANGED));
144  cls.attr("STATUS_TR_DECREASED") = py::cast(int(Optimizer::STATUS_TR_DECREASED));
145  cls.attr("STATUS_TR_INCREASED") = py::cast(int(Optimizer::STATUS_TR_INCREASED));
146  cls.attr("STATUS_TR") = py::cast(int(Optimizer::STATUS_TR));
147  cls.attr("STATUS") = py::cast(int(Optimizer::STATUS));
148  cls.def(py::init<std::shared_ptr<Optimizer::Objective const>, ndarray::Array<Scalar const, 1, 1> const &,
150  "objective"_a, "parameters"_a, "ctrl"_a);
151  cls.def("getObjective", &Optimizer::getObjective);
152  cls.def("getControl", &Optimizer::getControl, py::return_value_policy::copy);
153  cls.def("step", (bool (Optimizer::*)()) & Optimizer::step);
154  cls.def("step", (bool (Optimizer::*)(Optimizer::HistoryRecorder const &, afw::table::BaseCatalog &)) &
155  Optimizer::step,
156  "recorder"_a, "history"_a);
157  cls.def("run", (int (Optimizer::*)()) & Optimizer::run);
158  cls.def("run", (int (Optimizer::*)(Optimizer::HistoryRecorder const &, afw::table::BaseCatalog &)) &
159  Optimizer::run,
160  "recorder"_a, "history"_a);
161  cls.def("getState", &Optimizer::getState);
162  cls.def("getObjectiveValue", &Optimizer::getObjectiveValue);
163  cls.def("getParameters", &Optimizer::getParameters);
164  cls.def("getResiduals", &Optimizer::getResiduals);
165  cls.def("getGradient", &Optimizer::getGradient);
166  cls.def("getHessian", &Optimizer::getHessian);
167  cls.def("removeSR1Term", &Optimizer::removeSR1Term);
168  return cls;
169 }
170 
171 PYBIND11_PLUGIN(optimizer) {
172  py::module::import("lsst.meas.modelfit.model");
173  py::module::import("lsst.meas.modelfit.likelihood");
174  py::module::import("lsst.meas.modelfit.priors");
175 
176  py::module mod("optimizer");
177 
178  if (_import_array() < 0) {
179  PyErr_SetString(PyExc_ImportError, "numpy.core.multiarray failed to import");
180  return nullptr;
181  }
182 
183  auto clsObjective = declareOptimizerObjective(mod);
184  auto clsControl = declareOptimizerControl(mod);
185  auto clsHistoryRecorder = declareOptimizerHistoryRecorder(mod);
186  auto cls = declareOptimizer(mod);
187  cls.attr("Objective") = clsObjective;
188  cls.attr("Control") = clsControl;
189  cls.attr("HistoryRecorder") = clsHistoryRecorder;
190 
191  mod.def("solveTrustRegion", &solveTrustRegion, "x"_a, "F"_a, "g"_a, "r"_a, "tolerance"_a);
192 
193  return mod.ptr();
194 }
195 }
196 }
197 }
198 } // namespace lsst::meas::modelfit::anonymous
void solveTrustRegion(ndarray::Array< Scalar, 1, 1 > const &x, ndarray::Array< Scalar const, 2, 1 > const &F, ndarray::Array< Scalar const, 1, 1 > const &g, double r, double tolerance)
Solve a symmetric quadratic matrix equation with a ball constraint.