Mantid
Loading...
Searching...
No Matches
Workspaces.cpp
Go to the documentation of this file.
1// Mantid Repository : https://github.com/mantidproject/mantid
2//
3// Copyright © 2018 ISIS Rutherford Appleton Laboratory UKRI,
4// NScD Oak Ridge National Laboratory, European Spallation Source,
5// Institut Laue - Langevin & CSNS, Institute of High Energy Physics, CAS
6// SPDX - License - Identifier: GPL - 3.0 +
7// This code was originally translated from Fortran code on
8// https://ccpforge.cse.rl.ac.uk/gf/project/ral_nlls June 2016
10
12
15 : first_call(0), iter(0), normF0(), normJF0(), normF(), normJF(), normJFold(), normJF_Newton(), Delta(), normd(),
16 use_second_derivatives(false), hybrid_count(0), hybrid_tol(1.0), tr_p(7) {}
17
22void NLLS_workspace::initialize(int n, int m, const nlls_options &options) {
23
24 tr_nu = options.radius_increase;
25
26 J.allocate(m, n);
27 f.allocate(m);
29 hf.allocate(n, n);
30 d.allocate(n);
31 g.allocate(n);
33 y.allocate(n);
34 y.zero();
36 y_sharp.zero();
37
38 if (!options.exact_second_derivatives) {
41 Sks.allocate(n);
43 }
44
45 if (options.output_progress_vectors) {
46 resvec.allocate(options.maxit + 1);
47 gradvec.allocate(options.maxit + 1);
48 }
49
50 if (options.calculate_svd_J) {
51 largest_sv.allocate(options.maxit + 1);
52 smallest_sv.allocate(options.maxit + 1);
53 }
54
55 if (options.model == 3) {
57 }
58}
59
60} // namespace Mantid::CurveFitting::NLLS
void allocate(const int iFrom, const int iTo, const int jFrom, const int jTo)
Resize the matrix.
void allocate(int firstIndex, int lastIndex)
Resize the vector.
NLLS_workspace()
Constructor of the workspace.
Definition: Workspaces.cpp:14
void initialize(int n, int m, const nlls_options &options)
Initialize the workspace.
Definition: Workspaces.cpp:22
int model
specify the model used.
Definition: Workspaces.h:46
int maxit
the maximum number of iterations performed
Definition: Workspaces.h:38
bool output_progress_vectors
Shall we output progess vectors at termination of the routine?
Definition: Workspaces.h:152
double radius_increase
on very successful iterations, the trust-region radius will be increased by the factor ....
Definition: Workspaces.h:99
bool exact_second_derivatives
shall we use explicit second derivatives, or approximate using a secant method
Definition: Workspaces.h:113