1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747
|
// -*- c++ -*-
// This file is part of the Collective Variables module (Colvars).
// The original version of Colvars and its updates are located at:
// https://github.com/Colvars/colvars
// Please update all Colvars source files before making any changes.
// If you wish to distribute your changes, please submit them to the
// Colvars repository at GitHub.
#ifndef COLVARVALUE_H
#define COLVARVALUE_H
#include "colvarmodule.h"
#include "colvartypes.h"
/// \brief Value of a collective variable: this is a metatype which
/// can be set at runtime. By default it is set to be a scalar
/// number, and can be treated as such in all operations (this is
/// done by most \link colvar::cvc \endlink implementations).
///
/// \link colvarvalue \endlink allows \link colvar \endlink to be
/// treat different data types. By default, a \link colvarvalue
/// \endlink variable is a scalar number. To use it as
/// another type, declare and initialize it as
/// `colvarvalue x(colvarvalue::type_xxx)`, use `x.type (colvarvalue::type_xxx)`
/// at a later stage, or if unset,
/// assign the type with `x = y;`, provided y is correctly set.
///
/// All operators (either unary or binary) on a \link
/// colvarvalue \endlink object performs one or more checks on the
/// \link Type \endlink, except when reading from a stream, when there is no way to
/// detect the \link Type \endlink. To use `is >> x;` x \b MUST
/// already have a type correcly set up for properly parsing the
/// stream. No problem of course with the output streams: `os << x;`
///
/// \em Note \em on \em performance: to avoid type checks in a long array of \link
/// colvarvalue \endlink objects, use one of the existing "_opt" functions or implement a new one
class colvarvalue {
public:
/// \brief Possible types of value
///
/// These three cover most possibilities of data type one can
/// devise. If you need to implement a new colvar with a very
/// complex data type, it's better to put an allocatable class here
enum Type {
/// Undefined type
type_notset,
/// Scalar number, implemented as \link colvarmodule::real \endlink (default)
type_scalar,
/// 3-dimensional vector, implemented as \link colvarmodule::rvector \endlink
type_3vector,
/// 3-dimensional unit vector, implemented as \link colvarmodule::rvector \endlink
type_unit3vector,
/// 3-dimensional vector that is a derivative of a unitvector
type_unit3vectorderiv,
/// 4-dimensional unit vector representing a rotation, implemented as \link colvarmodule::quaternion \endlink
type_quaternion,
/// 4-dimensional vector that is a derivative of a quaternion
type_quaternionderiv,
/// vector (arbitrary dimension)
type_vector,
/// Needed to iterate through enum
type_all
};
/// Current type of this colvarvalue object
Type value_type;
/// \brief Real data member
cvm::real real_value;
/// \brief 3-dimensional vector data member
cvm::rvector rvector_value;
/// \brief Quaternion data member
cvm::quaternion quaternion_value;
/// \brief Generic vector data member
cvm::vector1d<cvm::real> vector1d_value;
/// \brief If \link vector1d_value \endlink is a concatenation of colvarvalues,
/// keep track of the individual types
std::vector<Type> elem_types;
/// \brief If \link vector1d_value \endlink is a concatenation of colvarvalues,
/// these mark the initial components of each colvarvalue
std::vector<int> elem_indices;
/// \brief If \link vector1d_value \endlink is a concatenation of colvarvalues,
/// these mark how many components for each colvarvalue
std::vector<int> elem_sizes;
/// \brief Whether or not the type check is enforced
static inline bool type_checking()
{
return true;
}
/// Runtime description of value types
static std::string const type_desc(Type t);
/// User keywords for specifying value types in the configuration
static std::string const type_keyword(Type t);
/// Number of degrees of freedom for each supported type
static size_t num_df(Type t);
/// Number of dimensions for each supported type (used to allocate vector1d_value)
static size_t num_dimensions(Type t);
/// Number of dimensions of this variable
size_t size() const;
/// \brief Default constructor: this class defaults to a scalar
/// number and always behaves like it unless you change its type
inline colvarvalue()
: value_type(type_scalar), real_value(0.0)
{}
/// Constructor from a type specification
inline colvarvalue(Type const &vti)
: value_type(vti), real_value(0.0)
{
reset();
}
/// Copy constructor from real base type
inline colvarvalue(cvm::real const &x)
: value_type(type_scalar), real_value(x)
{}
/// \brief Copy constructor from rvector base type (Note: this sets
/// by default a type \link type_3vector \endlink , if you want a
/// \link type_unit3vector \endlink you must set it explicitly)
inline colvarvalue(cvm::rvector const &v, Type vti = type_3vector)
: value_type(vti), real_value(0.0), rvector_value(v)
{}
/// \brief Copy constructor from quaternion base type
inline colvarvalue(cvm::quaternion const &q, Type vti = type_quaternion)
: value_type(vti), real_value(0.0), quaternion_value(q)
{}
/// Copy constructor from vector1d base type
colvarvalue(cvm::vector1d<cvm::real> const &v, Type vti = type_vector);
/// Copy constructor from another \link colvarvalue \endlink
colvarvalue(colvarvalue const &x);
/// Set to the null value for the data type currently defined
void reset();
/// \brief If the variable has constraints (e.g. unitvector or
/// quaternion), transform it to satisfy them; this function needs
/// to be called only when the \link colvarvalue \endlink
/// is calculated outside of \link colvar::cvc \endlink objects
void apply_constraints();
/// Get the current type
inline Type type() const
{
return value_type;
}
/// Set the type explicitly
void type(Type const &vti);
/// Set the type after another \link colvarvalue \endlink
void type(colvarvalue const &x);
/// Make the type a derivative of the original type
/// (so that its constraints do not apply)
void is_derivative();
/// Square norm of this colvarvalue
cvm::real norm2() const;
/// Norm of this colvarvalue
inline cvm::real norm() const
{
return cvm::sqrt(this->norm2());
}
/// Sum of the components of this colvarvalue (if more than one dimension)
cvm::real sum() const;
/// Return a colvarvalue object of the same type and all components set to 1
colvarvalue ones() const;
/// Square distance between this \link colvarvalue \endlink and another
cvm::real dist2(colvarvalue const &x2) const;
/// Derivative with respect to this \link colvarvalue \endlink of the square distance
colvarvalue dist2_grad(colvarvalue const &x2) const;
/// Return the midpoint between x1 and x2, optionally weighted by lambda
/// (which must be between 0.0 and 1.0)
static colvarvalue const interpolate(colvarvalue const &x1,
colvarvalue const &x2,
cvm::real const lambda = 0.5);
/// Assignment operator (type of x is checked)
colvarvalue & operator = (colvarvalue const &x);
void operator += (colvarvalue const &x);
void operator -= (colvarvalue const &x);
void operator *= (cvm::real const &a);
void operator /= (cvm::real const &a);
// Binary operators (return values)
friend colvarvalue operator + (colvarvalue const &x1, colvarvalue const &x2);
friend colvarvalue operator - (colvarvalue const &x1, colvarvalue const &x2);
friend colvarvalue operator * (colvarvalue const &x, cvm::real const &a);
friend colvarvalue operator * (cvm::real const &a, colvarvalue const &x);
friend colvarvalue operator / (colvarvalue const &x, cvm::real const &a);
/// Inner product
friend cvm::real operator * (colvarvalue const &x1, colvarvalue const &x2);
// Cast to scalar
inline operator cvm::real() const
{
if (value_type != type_scalar) {
cvm::error("Error: trying to use a variable of type \""+
type_desc(value_type)+"\" as one of type \""+
type_desc(type_scalar)+"\".\n");
}
return real_value;
}
// Cast to 3-vector
inline operator cvm::rvector() const
{
if ((value_type != type_3vector) &&
(value_type != type_unit3vector) &&
(value_type != type_unit3vectorderiv)) {
cvm::error("Error: trying to use a variable of type \""+
type_desc(value_type)+"\" as one of type \""+
type_desc(type_3vector)+"\".\n");
}
return rvector_value;
}
// Cast to quaternion
inline operator cvm::quaternion() const
{
if ((value_type != type_quaternion) &&
(value_type != type_quaternionderiv)) {
cvm::error("Error: trying to use a variable of type \""+
type_desc(value_type)+"\" as one of type \""+
type_desc(type_quaternion)+"\".\n");
}
return quaternion_value;
}
// Create a n-dimensional vector from one of the basic types, or return the existing vector
cvm::vector1d<cvm::real> const as_vector() const;
/// Whether this variable is a real number
inline bool is_scalar() const
{
return (value_type == type_scalar);
}
/// Add an element to the vector (requires that type_vector is already set).
/// This is only needed to use this object as a vector of "complex" colvar values.
/// To use it instead as a plain n-dimensional vector, access vector1d_value directly.
void add_elem(colvarvalue const &x);
/// Get a single colvarvalue out of elements of the vector
colvarvalue const get_elem(int const i_begin, int const i_end, Type const vt) const;
/// Get a single colvarvalue out of elements of the vector
colvarvalue const get_elem(int const icv) const;
/// Set elements of the vector from a single colvarvalue (uses the rank of x
/// to compute the length)
void set_elem(int const icv, colvarvalue const &x);
/// Set elements of the vector from a single colvarvalue
void set_elem(int const i_begin, int const i_end, colvarvalue const &x);
/// Make each element a random number in N(0,1)
void set_random();
/// Make each element equal to the given argument
void set_ones(cvm::real assigned_value = 1.0);
/// Get a scalar number out of an element of the vector
cvm::real operator [] (int const i) const;
/// Use an element of the vector as a scalar number
cvm::real & operator [] (int const i);
/// Ensure that the two types are the same within a binary operator
static int check_types(colvarvalue const &x1, colvarvalue const &x2);
/// Ensure that the two types are the same within an assignment, or that the left side is type_notset
static int check_types_assign(Type const &vt1, Type const &vt2);
/// Undefined operation
void undef_op() const;
/// \brief Formatted output operator
friend std::ostream & operator << (std::ostream &os, colvarvalue const &q);
/// \brief Formatted input operator
friend std::istream & operator >> (std::istream &is, colvarvalue &q);
/// Give the number of characters required to output this
/// colvarvalue, given the current type assigned and the number of
/// characters for a real number
size_t output_width(size_t const &real_width) const;
/// Formats value as a script-friendly string (space separated list)
std::string to_simple_string() const;
/// Parses value from a script-friendly string (space separated list)
int from_simple_string(std::string const &s);
// optimized routines for operations on arrays of colvar values;
// xv and result are assumed to have the same number of elements
/// \brief Optimized routine for the inner product of one collective
/// variable with an array
static void inner_opt(colvarvalue const &x,
std::vector<colvarvalue>::iterator &xv,
std::vector<colvarvalue>::iterator const &xv_end,
std::vector<cvm::real>::iterator &result);
/// \brief Optimized routine for the inner product of one collective
/// variable with an array
static void inner_opt(colvarvalue const &x,
std::list<colvarvalue>::iterator &xv,
std::list<colvarvalue>::iterator const &xv_end,
std::vector<cvm::real>::iterator &result);
/// \brief Optimized routine for the second order Legendre
/// polynomial, (3cos^2(w)-1)/2, of one collective variable with an
/// array
static void p2leg_opt(colvarvalue const &x,
std::vector<colvarvalue>::iterator &xv,
std::vector<colvarvalue>::iterator const &xv_end,
std::vector<cvm::real>::iterator &result);
/// \brief Optimized routine for the second order Legendre
/// polynomial of one collective variable with an array
static void p2leg_opt(colvarvalue const &x,
std::list<colvarvalue>::iterator &xv,
std::list<colvarvalue>::iterator const &xv_end,
std::vector<cvm::real>::iterator &result);
};
inline size_t colvarvalue::size() const
{
switch (value_type) {
case colvarvalue::type_notset:
default:
return 0; break;
case colvarvalue::type_scalar:
return 1; break;
case colvarvalue::type_3vector:
case colvarvalue::type_unit3vector:
case colvarvalue::type_unit3vectorderiv:
return 3; break;
case colvarvalue::type_quaternion:
case colvarvalue::type_quaternionderiv:
return 4; break;
case colvarvalue::type_vector:
return vector1d_value.size(); break;
}
}
inline cvm::real colvarvalue::operator [] (int const i) const
{
switch (value_type) {
case colvarvalue::type_notset:
default:
cvm::error("Error: trying to access a colvar value "
"that is not initialized.\n", BUG_ERROR);
return 0.0; break;
case colvarvalue::type_scalar:
return real_value; break;
case colvarvalue::type_3vector:
case colvarvalue::type_unit3vector:
case colvarvalue::type_unit3vectorderiv:
return rvector_value[i]; break;
case colvarvalue::type_quaternion:
case colvarvalue::type_quaternionderiv:
return quaternion_value[i]; break;
case colvarvalue::type_vector:
return vector1d_value[i]; break;
}
}
inline cvm::real & colvarvalue::operator [] (int const i)
{
switch (value_type) {
case colvarvalue::type_notset:
default:
cvm::error("Error: trying to access a colvar value "
"that is not initialized.\n", BUG_ERROR);
return real_value; break;
case colvarvalue::type_scalar:
return real_value; break;
case colvarvalue::type_3vector:
case colvarvalue::type_unit3vector:
case colvarvalue::type_unit3vectorderiv:
return rvector_value[i]; break;
case colvarvalue::type_quaternion:
case colvarvalue::type_quaternionderiv:
return quaternion_value[i]; break;
case colvarvalue::type_vector:
return vector1d_value[i]; break;
}
}
inline int colvarvalue::check_types(colvarvalue const &x1,
colvarvalue const &x2)
{
if (!colvarvalue::type_checking()) {
return COLVARS_OK;
}
if (x1.type() != x2.type()) {
if (((x1.type() == type_unit3vector) &&
(x2.type() == type_unit3vectorderiv)) ||
((x2.type() == type_unit3vector) &&
(x1.type() == type_unit3vectorderiv)) ||
((x1.type() == type_quaternion) &&
(x2.type() == type_quaternionderiv)) ||
((x2.type() == type_quaternion) &&
(x1.type() == type_quaternionderiv))) {
return COLVARS_OK;
} else {
cvm::error("Trying to perform an operation between two colvar "
"values with different types, \""+
colvarvalue::type_desc(x1.type())+
"\" and \""+
colvarvalue::type_desc(x2.type())+
"\".\n");
return COLVARS_ERROR;
}
}
if (x1.type() == type_vector) {
if (x1.vector1d_value.size() != x2.vector1d_value.size()) {
cvm::error("Trying to perform an operation between two vector colvar "
"values with different sizes, "+
cvm::to_str(x1.vector1d_value.size())+
" and "+
cvm::to_str(x2.vector1d_value.size())+
".\n");
return COLVARS_ERROR;
}
}
return COLVARS_OK;
}
inline int colvarvalue::check_types_assign(colvarvalue::Type const &vt1,
colvarvalue::Type const &vt2)
{
if (!colvarvalue::type_checking()) {
return COLVARS_OK;
}
if (vt1 != type_notset) {
if (((vt1 == type_unit3vector) &&
(vt2 == type_unit3vectorderiv)) ||
((vt2 == type_unit3vector) &&
(vt1 == type_unit3vectorderiv)) ||
((vt1 == type_quaternion) &&
(vt2 == type_quaternionderiv)) ||
((vt2 == type_quaternion) &&
(vt1 == type_quaternionderiv))) {
return COLVARS_OK;
} else {
if (vt1 != vt2) {
cvm::error("Trying to assign a colvar value with type \""+
type_desc(vt2)+"\" to one with type \""+
type_desc(vt1)+"\".\n");
return COLVARS_ERROR;
}
}
}
return COLVARS_OK;
}
inline colvarvalue & colvarvalue::operator = (colvarvalue const &x)
{
check_types_assign(this->type(), x.type());
value_type = x.type();
switch (this->type()) {
case colvarvalue::type_scalar:
this->real_value = x.real_value;
break;
case colvarvalue::type_3vector:
case colvarvalue::type_unit3vector:
case colvarvalue::type_unit3vectorderiv:
this->rvector_value = x.rvector_value;
break;
case colvarvalue::type_quaternion:
case colvarvalue::type_quaternionderiv:
this->quaternion_value = x.quaternion_value;
break;
case colvarvalue::type_vector:
vector1d_value = x.vector1d_value;
elem_types = x.elem_types;
elem_indices = x.elem_indices;
elem_sizes = x.elem_sizes;
break;
case colvarvalue::type_notset:
default:
undef_op();
break;
}
return *this;
}
inline void colvarvalue::operator += (colvarvalue const &x)
{
colvarvalue::check_types(*this, x);
switch (this->type()) {
case colvarvalue::type_scalar:
this->real_value += x.real_value;
break;
case colvarvalue::type_3vector:
case colvarvalue::type_unit3vector:
case colvarvalue::type_unit3vectorderiv:
this->rvector_value += x.rvector_value;
break;
case colvarvalue::type_quaternion:
case colvarvalue::type_quaternionderiv:
this->quaternion_value += x.quaternion_value;
break;
case colvarvalue::type_vector:
this->vector1d_value += x.vector1d_value;
break;
case colvarvalue::type_notset:
default:
undef_op();
}
}
inline void colvarvalue::operator -= (colvarvalue const &x)
{
colvarvalue::check_types(*this, x);
switch (value_type) {
case colvarvalue::type_scalar:
real_value -= x.real_value;
break;
case colvarvalue::type_3vector:
case colvarvalue::type_unit3vector:
case colvarvalue::type_unit3vectorderiv:
rvector_value -= x.rvector_value;
break;
case colvarvalue::type_quaternion:
case colvarvalue::type_quaternionderiv:
quaternion_value -= x.quaternion_value;
break;
case colvarvalue::type_vector:
this->vector1d_value -= x.vector1d_value;
break;
case colvarvalue::type_notset:
default:
undef_op();
}
}
inline void colvarvalue::operator *= (cvm::real const &a)
{
switch (value_type) {
case colvarvalue::type_scalar:
real_value *= a;
break;
case colvarvalue::type_3vector:
case colvarvalue::type_unit3vectorderiv:
rvector_value *= a;
break;
case colvarvalue::type_quaternion:
case colvarvalue::type_quaternionderiv:
quaternion_value *= a;
break;
case colvarvalue::type_vector:
this->vector1d_value *= a;
break;
case colvarvalue::type_notset:
default:
undef_op();
}
}
inline void colvarvalue::operator /= (cvm::real const &a)
{
switch (value_type) {
case colvarvalue::type_scalar:
real_value /= a; break;
case colvarvalue::type_3vector:
case colvarvalue::type_unit3vector:
case colvarvalue::type_unit3vectorderiv:
rvector_value /= a; break;
case colvarvalue::type_quaternion:
case colvarvalue::type_quaternionderiv:
quaternion_value /= a; break;
case colvarvalue::type_vector:
this->vector1d_value /= a;
break;
case colvarvalue::type_notset:
default:
undef_op();
}
}
inline cvm::vector1d<cvm::real> const colvarvalue::as_vector() const
{
switch (value_type) {
case colvarvalue::type_scalar:
{
cvm::vector1d<cvm::real> v(1);
v[0] = real_value;
return v;
}
case colvarvalue::type_3vector:
case colvarvalue::type_unit3vector:
case colvarvalue::type_unit3vectorderiv:
return rvector_value.as_vector();
case colvarvalue::type_quaternion:
case colvarvalue::type_quaternionderiv:
return quaternion_value.as_vector();
case colvarvalue::type_vector:
return vector1d_value;
case colvarvalue::type_notset:
default:
return cvm::vector1d<cvm::real>(0);
}
}
inline cvm::real colvarvalue::norm2() const
{
switch (value_type) {
case colvarvalue::type_scalar:
return (this->real_value)*(this->real_value);
case colvarvalue::type_3vector:
case colvarvalue::type_unit3vector:
case colvarvalue::type_unit3vectorderiv:
return (this->rvector_value).norm2();
case colvarvalue::type_quaternion:
case colvarvalue::type_quaternionderiv:
return (this->quaternion_value).norm2();
case colvarvalue::type_vector:
if (elem_types.size() > 0) {
// if we have information about non-scalar types, use it
cvm::real result = 0.0;
size_t i;
for (i = 0; i < elem_types.size(); i++) {
result += (this->get_elem(i)).norm2();
}
return result;
} else {
return vector1d_value.norm2();
}
break;
case colvarvalue::type_notset:
default:
return 0.0;
}
}
inline cvm::real colvarvalue::sum() const
{
switch (value_type) {
case colvarvalue::type_scalar:
return (this->real_value);
case colvarvalue::type_3vector:
case colvarvalue::type_unit3vector:
case colvarvalue::type_unit3vectorderiv:
return (this->rvector_value).x + (this->rvector_value).y +
(this->rvector_value).z;
case colvarvalue::type_quaternion:
case colvarvalue::type_quaternionderiv:
return (this->quaternion_value).q0 + (this->quaternion_value).q1 +
(this->quaternion_value).q2 + (this->quaternion_value).q3;
case colvarvalue::type_vector:
return (this->vector1d_value).sum();
case colvarvalue::type_notset:
default:
return 0.0;
}
}
inline cvm::real colvarvalue::dist2(colvarvalue const &x2) const
{
colvarvalue::check_types(*this, x2);
switch (this->type()) {
case colvarvalue::type_scalar:
return (this->real_value - x2.real_value)*(this->real_value - x2.real_value);
case colvarvalue::type_3vector:
return (this->rvector_value - x2.rvector_value).norm2();
case colvarvalue::type_unit3vector:
case colvarvalue::type_unit3vectorderiv:
// angle between (*this) and x2 is the distance
return cvm::acos(this->rvector_value * x2.rvector_value) * cvm::acos(this->rvector_value * x2.rvector_value);
case colvarvalue::type_quaternion:
case colvarvalue::type_quaternionderiv:
// angle between (*this) and x2 is the distance, the quaternion
// object has it implemented internally
return this->quaternion_value.dist2(x2.quaternion_value);
case colvarvalue::type_vector:
return (this->vector1d_value - x2.vector1d_value).norm2();
case colvarvalue::type_notset:
default:
this->undef_op();
return 0.0;
};
}
#endif
|