flexiblesusy is hosted by Hepforge, IPPP Durham
FlexibleSUSY
gsl_multiroot_fsolver.cpp
Go to the documentation of this file.
1// ====================================================================
2// This file is part of FlexibleSUSY.
3//
4// FlexibleSUSY is free software: you can redistribute it and/or modify
5// it under the terms of the GNU General Public License as published
6// by the Free Software Foundation, either version 3 of the License,
7// or (at your option) any later version.
8//
9// FlexibleSUSY is distributed in the hope that it will be useful, but
10// WITHOUT ANY WARRANTY; without even the implied warranty of
11// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12// General Public License for more details.
13//
14// You should have received a copy of the GNU General Public License
15// along with FlexibleSUSY. If not, see
16// <http://www.gnu.org/licenses/>.
17// ====================================================================
18
20#include "logger.hpp"
21#include "error.hpp"
22#include <string>
23
24namespace flexiblesusy {
25
27 const gsl_multiroot_fsolver_type* type, std::size_t dim,
28 gsl_multiroot_function* f, const GSL_vector& start)
29{
30 solver = gsl_multiroot_fsolver_alloc(type, dim);
31
32 if (!solver) {
33 throw OutOfMemoryError(
34 std::string("Cannot allocate gsl_multiroot_fsolver ") +
35 gsl_multiroot_fsolver_name(solver));
36 }
37
38 gsl_multiroot_fsolver_set(solver, f, start.raw());
39}
40
42{
43 gsl_multiroot_fsolver_free(solver);
44}
45
47
49{
50 return gsl_multiroot_fsolver_iterate(solver);
51}
52
53void GSL_multiroot_fsolver::print_state(std::size_t iteration) const
54{
55 VERBOSE_MSG("\t\t\tIteration " << iteration
56 << ": x = " << GSL_vector(solver->x)
57 << ", f(x) = " << GSL_vector(solver->f));
58}
59
60int GSL_multiroot_fsolver::test_residual(double precision) const noexcept
61{
62 return gsl_multiroot_test_residual(solver->f, precision);
63}
64
65} // namespace flexiblesusy
void print_state(std::size_t iteration) const
GSL_multiroot_fsolver(const gsl_multiroot_fsolver_type *type, std::size_t dim, gsl_multiroot_function *f, const GSL_vector &start)
int test_residual(double precision) const noexcept
const gsl_vector * raw() const noexcept
get raw pointer
Definition: gsl_vector.cpp:191
Not enough memory.
Definition: error.hpp:194
#define VERBOSE_MSG(msg)
Definition: logger.hpp:57
std::complex< double > f(double tau) noexcept