Skip to content

Commit

Permalink
scaling: trying to estimate b_star before actual optimization
Browse files Browse the repository at this point in the history
it's not finished, but it doesn't seem to improve the results
  • Loading branch information
wojdyr committed May 25, 2024
1 parent 1658ac3 commit f3c57ea
Show file tree
Hide file tree
Showing 4 changed files with 46 additions and 6 deletions.
9 changes: 6 additions & 3 deletions include/gemmi/levmar.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -25,9 +25,7 @@ namespace gemmi {
/// a is n x n matrix (in vector)
/// b is vector of length n,
/// This function returns vector x[] in b[], and 1-matrix in a[].
inline void jordan_solve(std::vector<double>& a, std::vector<double>& b) {
assert(a.size() == b.size() * b.size());
int n = (int) b.size();
inline void jordan_solve(double* a, double* b, int n) {
for (int i = 0; i < n; i++) {
// looking for a pivot element
int maxnr = -1;
Expand Down Expand Up @@ -70,6 +68,11 @@ inline void jordan_solve(std::vector<double>& a, std::vector<double>& b) {
}
}

inline void jordan_solve(std::vector<double>& a, std::vector<double>& b) {
assert(a.size() == b.size() * b.size());
jordan_solve(a.data(), b.data(), (int)b.size());
}


#ifdef GEMMI_DEBUG_LEVMAR
inline void debug_print(const std::string& name, std::vector<double> &a) {
Expand Down
41 changes: 38 additions & 3 deletions include/gemmi/scaling.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -217,7 +217,7 @@ struct Scaling {
return p.fcmol + (Real)get_solvent_scale(p.stol2) * p.fmask;
}

// quick linear fit (ignoring sigma) to get initial parameters
// quick linear fit (ignoring sigma) to get initial k_overall and isotropic B
void fit_isotropic_b_approximately() {
double sx = 0, sy = 0, sxx = 0, sxy = 0;
int n = 0;
Expand All @@ -242,6 +242,7 @@ struct Scaling {
set_b_overall({b_iso, b_iso, b_iso, 0, 0, 0});
}

// least-squares fitting of k_overall only
double lsq_k_overall() const {
double sxx = 0, sxy = 0;
for (const Point& p : points) {
Expand All @@ -255,9 +256,43 @@ struct Scaling {
return sxx != 0. ? sxy / sxx : 1.;
}

void fit_parameters() {
// For testing only, don't use it.
// Estimates anisotropic b_star using other parameters (incl. isotropic B),
// following P. Afonine et al, doi:10.1107/S0907444913000462 sec. 2.1.
// The symmetry constraints are not implemented - don't use it!
void fit_b_star_approximately() {
double b_iso = 1/3. * get_b_overall().trace();
//size_t nc = constraint_matrix.size();
double M[36] = {};
double b[6] = {};
//std::vector<double> Vc(nc);
for (const Point& p : points) {
double fcalc = std::abs(get_fcalc(p));
// the factor 1 / 2 pi^2 will be taken into account later
double Z = std::log(p.fobs / (k_overall * fcalc)) + b_iso * p.stol2;
Vec3 h(p.hkl);
double V[6] = {h.x * h.x, h.y * h.y, h.z * h.z,
2 * h.x * h.y, 2 * h.x * h.z, 2 * h.y * h.z};
//for (size_t i = 0; i < nc; ++i)
// Vc[i] = vec6_dot(constraint_matrix[i], V);
for (size_t i = 0; i < 6; ++i) {
for (size_t j = 0; j < 6; ++j)
M[6*i+j] += V[i] * V[j];
b[i] -= Z * V[i];
}
}
jordan_solve(M, b, 6);
double b_star_iso = 1/3. * b_star.trace();
SMat33<double> u_star{b[0], b[1], b[2], b[3], b[4], b[5]};
// u_to_b() / (2 pi^2) = 8 pi^2 / 2 pi^2 = 4
b_star = u_star.scaled(4.0).added_kI(b_star_iso);
//auto e = get_b_overall().elements_pdb();
//printf("fitted B = {%g %g %g %g %g %g}\n", e[0], e[1], e[2], e[3], e[4], e[5]);
}

double fit_parameters() {
LevMar levmar;
levmar.fit(*this);
return levmar.fit(*this);
}


Expand Down
1 change: 1 addition & 0 deletions prog/sfcalc.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -301,6 +301,7 @@ void process_with_fft(const gemmi::Structure& st,
scaling.k_overall = scaling.lsq_k_overall();
//fprintf(stderr, "initial k_ov=%g\n", scaling.k_overall);
}
//scaling.fit_b_star_approximately();
scaling.fit_parameters();
if (scaling.use_solvent)
fprintf(stderr, "Bulk solvent parameters: k_sol=%g B_sol=%g\n",
Expand Down
1 change: 1 addition & 0 deletions python/scaling.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -23,6 +23,7 @@ void add_scaling(py::module& m) {
.def("prepare_points", &Scaling::prepare_points,
py::arg("calc"), py::arg("obs"), py::arg("mask")=static_cast<FPhiData*>(nullptr))
.def("fit_isotropic_b_approximately", &Scaling::fit_isotropic_b_approximately)
.def("fit_b_star_approximately", &Scaling::fit_b_star_approximately)
.def("fit_parameters", &Scaling::fit_parameters)
.def("get_overall_scale_factor", &Scaling::get_overall_scale_factor, py::arg("hkl"))
.def("get_overall_scale_factor", [](const Scaling& self, py::array_t<int> hkl) {
Expand Down

0 comments on commit f3c57ea

Please sign in to comment.