8 #ifndef INCLUDE_MOLASSEMBLER_DG_EIGEN_REFINEMENT_PROBLEM_H
9 #define INCLUDE_MOLASSEMBLER_DG_EIGEN_REFINEMENT_PROBLEM_H
11 #include <Eigen/Dense>
16 namespace Molassembler {
17 namespace DistanceGeometry {
27 template<
unsigned dimensionality,
typename FloatType,
bool SIMD>
30 dimensionality == 3 || dimensionality == 4,
31 "EigenRefinementProblem is only suitable for three or four dimensions"
38 using VectorType = Eigen::Matrix<FloatType, Eigen::Dynamic, 1>;
54 void distanceTerm(
unsigned ,
unsigned ,
double ) {}
55 void chiralTerm(
const ChiralConstraint::SiteSequence& ,
double ) {}
56 void dihedralTerm(
const DihedralConstraint::SiteSequence& ,
double ) {}
57 void fourthDimensionTerm(
unsigned ,
double ) {}
67 static std::string
name() {
68 std::string
name =
"Eigen";
70 name +=
"RefinementProblem<dimensionality=";
71 name += std::to_string(dimensionality);
74 if(std::is_same<FloatType, float>::value) {
81 name +=
", SIMD=true>";
83 name +=
", SIMD=false>";
97 return positions.template segment<dimensionality>(dimensionality * index);
108 return positions.template segment<3>(dimensionality * index);
118 const ChiralConstraint::AtomListType& atomList
120 if(atomList.size() == 1) {
126 sum += positions.template segment<dimensionality>(dimensionality * i);
129 return sum / atomList.size();
139 const ChiralConstraint::AtomListType& atomList
141 if(atomList.size() == 1) {
147 sum += positions.template segment<3>(dimensionality * i);
150 return sum / atomList.size();
180 mutable double proportionChiralConstraintsCorrectSign = 0.0;
186 const Eigen::MatrixXd& squaredBounds,
187 std::vector<ChiralConstraint> passChiralConstraints,
188 std::vector<DihedralConstraint> passDihedralConstraints
192 const unsigned N = squaredBounds.cols();
193 const unsigned strictlyUpperTriangularElements = N * (N - 1) / 2;
198 for(
unsigned linearIndex = 0, i = 0; i < N; ++i) {
199 for(
unsigned j = i + 1; j < N; ++j, ++linearIndex) {
209 for(
unsigned i = 0; i < C; ++i) {
219 for(
unsigned i = 0; i < D; ++i) {
236 Eigen::Ref<VectorType> gradient
249 Eigen::Ref<VectorType> gradient
261 Eigen::Ref<VectorType> gradient
271 Eigen::Ref<VectorType> gradient
286 assert(parameters.size() == gradient.size());
287 assert(parameters.size() % dimensionality == 0);
304 unsigned nonZeroChiralConstraints = 0;
305 unsigned incorrectNonZeroChiralConstraints = 0;
314 if(!constraint.targetVolumeIsZero()) {
315 nonZeroChiralConstraints += 1;
319 const double volume = (
330 (volume < 0 && constraint.lower > 0)
331 || (volume > 0 && constraint.lower < 0)
333 incorrectNonZeroChiralConstraints += 1;
338 if(nonZeroChiralConstraints == 0) {
339 proportionChiralConstraintsCorrectSign = 1.0;
341 proportionChiralConstraintsCorrectSign =
static_cast<double>(
342 nonZeroChiralConstraints - incorrectNonZeroChiralConstraints
343 ) / nonZeroChiralConstraints;
346 return proportionChiralConstraintsCorrectSign;
349 template<
typename Visitor>
350 auto visitTerms(
const VectorType& positions, Visitor&& visitor)
const {
351 FloatType dummyValue = 0;
353 dummyGradient.resize(positions.size());
354 dummyGradient.setZero();
380 template<
typename Visitor>
386 const unsigned N = positions.size() / dimensionality;
389 for(
unsigned i = 0; i < N; ++i) {
390 for(
unsigned j = i + 1; j < N; ++j) {
391 const double ijDistance = (
392 positions.template segment<dimensionality>(dimensionality * j)
393 - positions.template segment<dimensionality>(dimensionality * i)
397 ijDistance - bounds.
upperBound(i, j) > visitor.deviationThreshold
398 || bounds.
lowerBound(i, j) - ijDistance > visitor.deviationThreshold
400 visitor.distanceOverThreshold(i, j, ijDistance);
401 if(visitor.earlyExit) {
402 return visitor.value;
410 const double volume = (
424 volume - constraint.upper > visitor.deviationThreshold
425 || constraint.lower - volume > visitor.deviationThreshold
427 visitor.chiralOverThreshold(constraint, volume);
428 if(visitor.earlyExit) {
429 return visitor.value;
448 const double dihedral = std::atan2(
449 a.cross(b).dot(-g.normalized()),
453 const double constraintSumHalved = (constraint.lower + constraint.upper) / 2;
456 if(dihedral < constraintSumHalved - M_PI) {
457 phi = dihedral + 2 * M_PI;
458 }
else if(dihedral > constraintSumHalved + M_PI) {
459 phi = dihedral - 2 * M_PI;
464 const double term = std::fabs(phi - constraintSumHalved) - (constraint.upper - constraint.lower) / 2;
466 if(term > visitor.deviationThreshold) {
467 visitor.dihedralOverThreshold(constraint, dihedral, term);
468 if(visitor.earlyExit) {
469 return visitor.value;
474 return visitor.value;
483 template<
class Visitor,
bool dependent = SIMD, std::enable_if_t<!dependent,
int>...>
487 Eigen::Ref<VectorType> gradient,
490 assert(positions.size() == gradient.size());
491 const unsigned N = positions.size() / dimensionality;
493 for(
unsigned linearIndex = 0, i = 0; i < N - 1; ++i) {
494 for(
unsigned j = i + 1; j < N; ++j, ++linearIndex) {
497 assert(lowerBoundSquared <= upperBoundSquared);
501 positions.template segment<dimensionality>(dimensionality * i)
502 - positions.template segment<dimensionality>(dimensionality * j)
505 const FloatType squareDistance = positionDifference.squaredNorm();
508 const FloatType upperTerm = squareDistance / upperBoundSquared - 1;
511 const FloatType value = upperTerm * upperTerm;
513 visitor.distanceTerm(i, j, value);
517 gradient.template segment<dimensionality>(dimensionality * i) += f;
518 gradient.template segment<dimensionality>(dimensionality * j) -= f;
521 const FloatType quotient = lowerBoundSquared + squareDistance;
522 const FloatType lowerTerm = 2 * lowerBoundSquared / quotient - 1;
525 const FloatType value = lowerTerm * lowerTerm;
527 visitor.distanceTerm(i, j, value);
537 gradient.template segment<dimensionality>(dimensionality * i) -= g;
538 gradient.template segment<dimensionality>(dimensionality * j) += g;
540 visitor.distanceTerm(i, j, 0.0);
550 template<
class Visitor,
bool dependent = SIMD, std::enable_if_t<dependent,
int>...>
554 Eigen::Ref<VectorType> gradient,
557 assert(positions.size() == gradient.size());
592 const unsigned N = positions.size() / dimensionality;
593 const unsigned differencesCount = N * (N - 1) / 2;
598 for(
unsigned i = 0; i < N - 1; ++i) {
599 auto iPosition = positions.template segment<dimensionality>(dimensionality * i);
600 for(
unsigned j = i + 1; j < N; ++j) {
601 positionDifferences.col(offset) = (
603 - positions.template segment<dimensionality>(dimensionality * j)
611 const VectorType squareDistances = positionDifferences.colwise().squaredNorm();
620 for(
unsigned i = 0; i < N - 1; ++i) {
621 for(
unsigned j = i + 1; j < N; ++j) {
623 positions.template segment<dimensionality>(dimensionality * i)
624 - positions.template segment<dimensionality>(dimensionality * j)
627 assert(positionDifferences.col(offset) == positionDifference);
629 const FloatType squareDistance = positionDifference.squaredNorm();
631 assert(squareDistance == squareDistances(offset));
636 assert(upperTerm == upperTerms(offset));
645 for(
unsigned iOffset = 0, i = 0; i < N - 1; ++i) {
646 const unsigned crossTerms = N - i - 1;
647 for(
unsigned jOffset = 0, j = i + 1; j < N; ++j, ++jOffset) {
648 const FloatType upperTerm = upperTerms(iOffset + jOffset);
650 error += upperTerm * upperTerm;
657 4 * positionDifferences.col(iOffset + jOffset) * upperTerm
662 gradient.template segment<dimensionality>(dimensionality * i) += f;
663 gradient.template segment<dimensionality>(dimensionality * j) -= f;
668 const FloatType quotient = lowerBoundSquared + squareDistances(iOffset + jOffset);
669 const FloatType lowerTerm = 2 * lowerBoundSquared / quotient - 1;
673 error += lowerTerm * lowerTerm;
677 8 * lowerBoundSquared * positionDifferences.col(iOffset + jOffset) * lowerTerm / (
683 gradient.template segment<dimensionality>(dimensionality * i) -= g;
684 gradient.template segment<dimensionality>(dimensionality * j) += g;
689 iOffset += crossTerms;
696 template<
typename Visitor>
700 Eigen::Ref<VectorType> gradient,
703 unsigned nonZeroChiralConstraints = 0;
704 unsigned incorrectNonZeroChiralConstraints = 0;
716 const FloatType volume = (alphaMinusDelta).dot(
717 betaMinusDelta.cross(gammaMinusDelta)
720 if(!constraint.targetVolumeIsZero()) {
721 ++nonZeroChiralConstraints;
723 ( volume < 0 && constraint.lower > 0)
724 || (volume > 0 && constraint.lower < 0)
726 ++incorrectNonZeroChiralConstraints;
731 const FloatType upperTerm =
static_cast<FloatType
>(constraint.weight) * (volume - static_cast<FloatType>(constraint.upper));
733 const FloatType value = upperTerm * upperTerm;
735 visitor.chiralTerm(constraint.sites, value);
739 const FloatType lowerTerm =
static_cast<FloatType
>(constraint.weight) * (static_cast<FloatType>(constraint.lower) - volume);
741 const FloatType value = lowerTerm * lowerTerm;
743 visitor.chiralTerm(constraint.sites, value);
746 const FloatType factor = 2 * (
747 std::max(FloatType {0}, upperTerm)
748 - std::max(FloatType {0}, lowerTerm)
760 factor * betaMinusDelta.cross(gammaMinusDelta)
761 / constraint.sites[0].size()
765 factor * gammaMinusDelta.cross(alphaMinusDelta)
766 / constraint.sites[1].size()
770 factor * alphaMinusDelta.cross(betaMinusDelta)
771 / constraint.sites[2].size()
775 factor * betaMinusGamma.cross(alphaMinusGamma)
776 / constraint.sites[3].size()
779 for(
const AtomIndex alphaI : constraint.sites[0]) {
780 gradient.template segment<3>(dimensionality * alphaI) += iContribution;
783 for(
const AtomIndex betaI : constraint.sites[1]) {
784 gradient.template segment<3>(dimensionality * betaI) += jContribution;
787 for(
const AtomIndex gammaI : constraint.sites[2]) {
788 gradient.template segment<3>(dimensionality * gammaI) += kContribution;
791 for(
const AtomIndex deltaI : constraint.sites[3]) {
792 gradient.template segment<3>(dimensionality * deltaI) += lContribution;
795 visitor.chiralTerm(constraint.sites, 0.0);
800 if(nonZeroChiralConstraints == 0) {
801 proportionChiralConstraintsCorrectSign = 1;
803 proportionChiralConstraintsCorrectSign =
static_cast<double>(
804 nonZeroChiralConstraints - incorrectNonZeroChiralConstraints
805 ) / nonZeroChiralConstraints;
812 template<
class Visitor>
816 Eigen::Ref<VectorType> gradient,
819 assert(positions.size() == gradient.size());
820 constexpr FloatType reductionFactor = 1.0 / 10;
836 FloatType phi = std::atan2(
837 a.cross(b).dot(-g.normalized()),
841 const FloatType constraintSumHalved = (
static_cast<FloatType
>(constraint.lower) + static_cast<FloatType>(constraint.upper)) / 2;
843 constexpr FloatType pi {M_PI};
845 if(phi < constraintSumHalved - pi) {
847 }
else if(phi > constraintSumHalved + pi) {
851 const FloatType w_phi = phi - constraintSumHalved;
853 FloatType h_phi = std::fabs(w_phi) - (
static_cast<FloatType
>(constraint.upper) - static_cast<FloatType>(constraint.lower)) / 2;
859 visitor.dihedralTerm(constraint.sites, 0.0);
864 const FloatType errorContribution = h_phi * h_phi * reductionFactor;
865 error += errorContribution;
866 visitor.dihedralTerm(constraint.sites, errorContribution);
869 h_phi *=
static_cast<int>(0 < w_phi) - static_cast<int>(w_phi < 0);
872 h_phi *= 2 * reductionFactor;
875 const FloatType gLength = g.norm();
877 throw std::runtime_error(
"Encountered zero-length g vector in dihedral errors");
880 const FloatType aLengthSq = a.squaredNorm();
882 throw std::runtime_error(
"Encountered zero-length a vector in dihedral gradient contributions");
885 const FloatType bLengthSq = b.squaredNorm();
887 throw std::runtime_error(
"Encountered zero-length b vector in dihedral gradient contributions");
890 const FloatType fDotG = f.dot(g);
891 const FloatType gDotH = g.dot(h);
894 -(h_phi / constraint.sites[0].size()) * gLength * a
898 (h_phi / constraint.sites[1].size()) * (
899 (gLength + fDotG / gLength) * a
900 - (gDotH / gLength) * b
905 (h_phi / constraint.sites[2].size()) * (
906 (gDotH / gLength - gLength) * b
907 - (fDotG / gLength) * a
912 (h_phi / constraint.sites[3].size()) * gLength * b
916 !iContribution.array().isFinite().all()
917 || !jContribution.array().isFinite().all()
918 || !kContribution.array().isFinite().all()
919 || !lContribution.array().isFinite().all()
921 throw std::runtime_error(
"Encountered non-finite dihedral gradient contributions");
925 for(
const AtomIndex alphaConstitutingIndex : constraint.sites[0]) {
926 gradient.template segment<3>(dimensionality * alphaConstitutingIndex) += iContribution;
929 for(
const AtomIndex betaConstitutingIndex: constraint.sites[1]) {
930 gradient.template segment<3>(dimensionality * betaConstitutingIndex) += jContribution;
933 for(
const AtomIndex gammaConstitutingIndex : constraint.sites[2]) {
934 gradient.template segment<3>(dimensionality * gammaConstitutingIndex) += kContribution;
937 for(
const AtomIndex deltaConstitutingIndex : constraint.sites[3]) {
938 gradient.template segment<3>(dimensionality * deltaConstitutingIndex) += lContribution;
947 template<
class Visitor = DefaultTermVisitor>
951 Eigen::Ref<VectorType> gradient,
952 Visitor&& visitor = {}
954 if(dimensionality != 4) {
959 const unsigned N = positions.size() / dimensionality;
961 for(
unsigned i = 0; i < N; ++i) {
962 const FloatType fourthDimensionalValue = positions(4 * i + 3);
963 const FloatType value = fourthDimensionalValue * fourthDimensionalValue;
965 visitor.fourthDimensionTerm(i, value);
966 gradient(4 * i + 3) += 2 * fourthDimensionalValue;
969 for(
unsigned i = 0; i < N; ++i) {
970 visitor.fourthDimensionTerm(i, 0.0);
981 template<
typename RefinementType>
985 template<
unsigned,
typename,
bool>
class BaseClass,
986 unsigned dimensionality,
987 typename FloatingPointType,
989 >
static auto templateArgumentHelper(BaseClass<dimensionality, FloatingPointType, SIMD> )
991 std::integral_constant<unsigned, dimensionality>,
993 std::integral_constant<bool, SIMD>
999 using TemplateArgumentsTuple = decltype(templateArgumentHelper(std::declval<RefinementType>()));
1000 using DimensionalityConstant = std::tuple_element_t<0, TemplateArgumentsTuple>;
1001 using FloatingPointType = std::tuple_element_t<1, TemplateArgumentsTuple>;
1002 using SimdConstant = std::tuple_element_t<2, TemplateArgumentsTuple>;
double upper
Upper bound on signed volume.
Definition: DistanceGeometry.h:41
void chiralContributionsImpl(const VectorType &positions, FloatType &error, Eigen::Ref< VectorType > gradient, Visitor &&visitor) const
Adds chiral error and gradient contributions.
Definition: EigenRefinement.h:697
Definition: DistanceBoundsMatrix.h:28
bool dihedralTerms
Whether to enable dihedral terms.
Definition: EigenRefinement.h:175
bool compressFourthDimension
Whether to compress the fourth dimension.
Definition: EigenRefinement.h:173
static ThreeDimensionalVector getAveragePosition3D(const VectorType &positions, const ChiralConstraint::AtomListType &atomList)
Calculate the average three-dimensional position of several atoms.
Definition: EigenRefinement.h:137
void dihedralContributions(const VectorType &positions, FloatType &error, Eigen::Ref< VectorType > gradient) const
Adds dihedral error and gradient contributions.
Definition: EigenRefinement.h:258
Data struct representing a dihedral constraint.
Definition: DistanceGeometry.h:61
Eigen-based refinement error function.
Definition: EigenRefinement.h:28
VectorType chiralLowerConstraints
Chiral lower constraints, in sequence of chiralConstraints.
Definition: EigenRefinement.h:163
double calculateProportionChiralConstraintsCorrectSign(const VectorType &positions) const
Calculates the number of chiral constraints with correct sign.
Definition: EigenRefinement.h:303
std::vector< DihedralConstraint > dihedralConstraints
List of dihedral constraints.
Definition: EigenRefinement.h:171
VectorType chiralUpperConstraints
Chiral upper constraints, in sequence of chiralConstraints.
Definition: EigenRefinement.h:161
static double & lowerBound(Eigen::Ref< Eigen::MatrixXd > matrix, const AtomIndex i, const AtomIndex j)
Uses Floyd's algorithm to smooth the matrix.
Definition: DistanceBoundsMatrix.h:35
double lower
Lower bound on signed volume.
Definition: DistanceGeometry.h:39
double lower
Lower bound on dihedral angle.
Definition: DistanceGeometry.h:68
void fourthDimensionContributions(const VectorType &positions, FloatType &error, Eigen::Ref< VectorType > gradient) const
Adds pairwise distance error and gradient contributions.
Definition: EigenRefinement.h:268
void distanceContributionsImpl(const VectorType &positions, FloatType &error, Eigen::Ref< VectorType > gradient, Visitor &&) const
SIMD implementation of distance contributions.
Definition: EigenRefinement.h:551
void distanceContributionsImpl(const VectorType &positions, FloatType &error, Eigen::Ref< VectorType > gradient, Visitor &&visitor) const
Adds distance error and gradient contributions (non-SIMD)
Definition: EigenRefinement.h:484
Eigen::Matrix< FloatType, 3, 1 > ThreeDimensionalVector
Three dimensional vector.
Definition: EigenRefinement.h:41
static FullDimensionalVector getPosition(const VectorType &positions, const AtomIndex index)
Copy a full dimensional position of an atom.
Definition: EigenRefinement.h:93
VectorType lowerDistanceBoundsSquared
Lower distance bounds squared, linearized in i < j.
Definition: EigenRefinement.h:159
std::size_t AtomIndex
Unsigned integer atom index type. Used to refer to particular atoms.
Definition: Types.h:51
static double & upperBound(Eigen::Ref< Eigen::MatrixXd > matrix, const AtomIndex i, const AtomIndex j)
Uses Floyd's algorithm to smooth the matrix.
Definition: DistanceBoundsMatrix.h:43
Eigen::Matrix< FloatType, 3, Eigen::Dynamic > ThreeDimensionalMatrixType
Three-row dynamic column matrix.
Definition: EigenRefinement.h:46
auto visitUnfulfilledConstraints(const DistanceBoundsMatrix &bounds, const VectorType &positions, Visitor &&visitor) const
Visit all unfulfilled constraints.
Definition: EigenRefinement.h:381
Eigen::Matrix< FloatType, dimensionality, Eigen::Dynamic > FullDimensionalMatrixType
Three or four-row dynamic column matrix.
Definition: EigenRefinement.h:48
void distanceContributions(const VectorType &positions, FloatType &error, Eigen::Ref< VectorType > gradient) const
Adds pairwise distance error and gradient contributions.
Definition: EigenRefinement.h:233
void chiralContributions(const VectorType &positions, FloatType &error, Eigen::Ref< VectorType > gradient) const
Adds chiral error and gradient contributions.
Definition: EigenRefinement.h:246
Eigen::Matrix< FloatType, Eigen::Dynamic, 1 > VectorType
Vector layout of positions.
Definition: EigenRefinement.h:38
VectorType upperDistanceBoundsSquared
Upper distance bounds squared, linearized in i < j.
Definition: EigenRefinement.h:157
This is for when you have a fully qualified Refinement typename but want to get the template argument...
Definition: EigenRefinement.h:982
Visitor that can be passed to visit all terms.
Definition: EigenRefinement.h:53
VectorType dihedralConstraintDiffsHalved
Dihedral bounds upper minus lower, halved, in sequence.
Definition: EigenRefinement.h:167
void operator()(const VectorType ¶meters, FloatType &value, Eigen::Ref< VectorType > gradient) const
Calculates the error value and gradient for all contributions.
Definition: EigenRefinement.h:285
FloatType FloatingPointType
Template argument specifying floating-point type.
Definition: EigenRefinement.h:50
VectorType dihedralConstraintSumsHalved
Dihedral bounds' averages, in sequence of dihedralConstraints.
Definition: EigenRefinement.h:165
std::vector< ChiralConstraint > chiralConstraints
List of chiral constraints.
Definition: EigenRefinement.h:169
static FullDimensionalVector getAveragePosition(const VectorType &positions, const ChiralConstraint::AtomListType &atomList)
Calculate the average full-dimensional position of several atoms.
Definition: EigenRefinement.h:116
constexpr T sum(const ArrayType< T, size > &array)
Sum up all elements of an array-like class.
Definition: Containers.h:111
static std::string name()
Get a string representation of the type name.
Definition: EigenRefinement.h:67
Eigen::Matrix< FloatType, dimensionality, 1 > FullDimensionalVector
Three or four dimensional vector (depending on dimensionality)
Definition: EigenRefinement.h:43
Data struct representing a chiral constraint.
Definition: DistanceGeometry.h:32
static ThreeDimensionalVector getPosition3D(const VectorType &positions, const AtomIndex index)
Copy a three-dimensional position of an atom.
Definition: EigenRefinement.h:104
void dihedralContributionsImpl(const VectorType &positions, FloatType &error, Eigen::Ref< VectorType > gradient, Visitor &&visitor) const
Adds dihedral error and gradient contributions.
Definition: EigenRefinement.h:813
void fourthDimensionContributionsImpl(const VectorType &positions, FloatType &error, Eigen::Ref< VectorType > gradient, Visitor &&visitor={}) const
Adds fourth dimension error and gradient contributions.
Definition: EigenRefinement.h:948
Class storing atom-pairwise distance bounds.
double upper
Upper bound on dihedral angle.
Definition: DistanceGeometry.h:70