Main Page   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Related Pages  

colvarvalue.h

Go to the documentation of this file.
00001 // -*- c++ -*-
00002 
00003 // This file is part of the Collective Variables module (Colvars).
00004 // The original version of Colvars and its updates are located at:
00005 // https://github.com/Colvars/colvars
00006 // Please update all Colvars source files before making any changes.
00007 // If you wish to distribute your changes, please submit them to the
00008 // Colvars repository at GitHub.
00009 
00010 #ifndef COLVARVALUE_H
00011 #define COLVARVALUE_H
00012 
00013 #include "colvarmodule.h"
00014 #include "colvartypes.h"
00015 
00016 
00039 
00040 
00041 class colvarvalue {
00042 
00043 public:
00044 
00050   enum Type {
00052     type_notset,
00054     type_scalar,
00056     type_3vector,
00058     type_unit3vector,
00060     type_unit3vectorderiv,
00062     type_quaternion,
00064     type_quaternionderiv,
00066     type_vector,
00068     type_all
00069   };
00070 
00072   Type value_type;
00073 
00075   cvm::real real_value;
00076 
00078   cvm::rvector rvector_value;
00079 
00081   cvm::quaternion quaternion_value;
00082 
00084   cvm::vector1d<cvm::real> vector1d_value;
00085 
00088   std::vector<Type> elem_types;
00089 
00092   std::vector<int> elem_indices;
00093 
00096   std::vector<int> elem_sizes;
00097 
00099   static inline bool type_checking()
00100   {
00101     return true;
00102   }
00103 
00105   static std::string const type_desc(Type t);
00106 
00108   static std::string const type_keyword(Type t);
00109 
00111   static size_t num_df(Type t);
00112 
00114   static size_t num_dimensions(Type t);
00115 
00117   size_t size() const;
00118 
00121   colvarvalue();
00122 
00124   colvarvalue(Type const &vti);
00125 
00127   colvarvalue(cvm::real const &x);
00128 
00132   colvarvalue(cvm::rvector const &v, Type vti = type_3vector);
00133 
00135   colvarvalue(cvm::quaternion const &q, Type vti = type_quaternion);
00136 
00138   colvarvalue(cvm::vector1d<cvm::real> const &v, Type vti = type_vector);
00139 
00141   colvarvalue(colvarvalue const &x);
00142 
00143 
00145   void reset();
00146 
00151   void apply_constraints();
00152 
00154   inline Type type() const
00155   {
00156     return value_type;
00157   }
00158 
00160   void type(Type const &vti);
00161 
00163   void type(colvarvalue const &x);
00164 
00167   void is_derivative();
00168 
00170   cvm::real norm2() const;
00171 
00173   inline cvm::real norm() const
00174   {
00175     return cvm::sqrt(this->norm2());
00176   }
00177 
00179   cvm::real sum() const;
00180 
00182   colvarvalue ones() const;
00183 
00185   cvm::real dist2(colvarvalue const &x2) const;
00186 
00188   colvarvalue dist2_grad(colvarvalue const &x2) const;
00189 
00192   static colvarvalue const interpolate(colvarvalue const &x1,
00193                                        colvarvalue const &x2,
00194                                        cvm::real const lambda = 0.5);
00195 
00197   colvarvalue & operator = (colvarvalue const &x);
00198 
00199   void operator += (colvarvalue const &x);
00200   void operator -= (colvarvalue const &x);
00201   void operator *= (cvm::real const &a);
00202   void operator /= (cvm::real const &a);
00203 
00204   // Binary operators (return values)
00205   friend colvarvalue operator + (colvarvalue const &x1, colvarvalue const &x2);
00206   friend colvarvalue operator - (colvarvalue const &x1, colvarvalue const &x2);
00207   friend colvarvalue operator * (colvarvalue const &x, cvm::real const &a);
00208   friend colvarvalue operator * (cvm::real const &a, colvarvalue const &x);
00209   friend colvarvalue operator / (colvarvalue const &x, cvm::real const &a);
00210 
00212   friend cvm::real operator * (colvarvalue const &x1, colvarvalue const &x2);
00213 
00214   // Cast to scalar
00215   inline operator cvm::real() const
00216   {
00217     if (value_type != type_scalar) {
00218       cvm::error("Error: trying to use a variable of type \""+
00219                  type_desc(value_type)+"\" as one of type \""+
00220                  type_desc(type_scalar)+"\".\n");
00221     }
00222     return real_value;
00223   }
00224 
00225   // Cast to 3-vector
00226   inline operator cvm::rvector() const
00227   {
00228     if ((value_type != type_3vector) &&
00229         (value_type != type_unit3vector) &&
00230         (value_type != type_unit3vectorderiv)) {
00231       cvm::error("Error: trying to use a variable of type \""+
00232                  type_desc(value_type)+"\" as one of type \""+
00233                  type_desc(type_3vector)+"\".\n");
00234     }
00235     return rvector_value;
00236   }
00237 
00238   // Cast to quaternion
00239   inline operator cvm::quaternion() const
00240   {
00241     if ((value_type != type_quaternion) &&
00242         (value_type != type_quaternionderiv)) {
00243       cvm::error("Error: trying to use a variable of type \""+
00244                  type_desc(value_type)+"\" as one of type \""+
00245                  type_desc(type_quaternion)+"\".\n");
00246     }
00247     return quaternion_value;
00248   }
00249 
00250   // Create a n-dimensional vector from one of the basic types, or return the existing vector
00251   cvm::vector1d<cvm::real> const as_vector() const;
00252 
00253 
00255   inline bool is_scalar() const
00256   {
00257     return (value_type == type_scalar);
00258   }
00259 
00260 
00264   void add_elem(colvarvalue const &x);
00265 
00267   colvarvalue const get_elem(int const i_begin, int const i_end, Type const vt) const;
00268 
00270   colvarvalue const get_elem(int const icv) const;
00271 
00274   void set_elem(int const icv, colvarvalue const &x);
00275 
00277   void set_elem(int const i_begin, int const i_end, colvarvalue const &x);
00278 
00280   void set_random();
00281 
00283   void set_ones(cvm::real assigned_value = 1.0);
00284 
00286   cvm::real operator [] (int const i) const;
00287 
00289   cvm::real & operator [] (int const i);
00290 
00292   static int check_types(colvarvalue const &x1, colvarvalue const &x2);
00293 
00295   static int check_types_assign(Type const &vt1, Type const &vt2);
00296 
00298   void undef_op() const;
00299 
00300 
00302   friend std::ostream & operator << (std::ostream &os, colvarvalue const &q);
00303 
00305   friend std::istream & operator >> (std::istream &is, colvarvalue &q);
00306 
00310   size_t output_width(size_t const &real_width) const;
00311 
00313   std::string to_simple_string() const;
00314 
00316   int from_simple_string(std::string const &s);
00317 
00318 
00319   // optimized routines for operations on arrays of colvar values;
00320   // xv and result are assumed to have the same number of elements
00321 
00324   static void inner_opt(colvarvalue const                        &x,
00325                         std::vector<colvarvalue>::iterator       &xv,
00326                         std::vector<colvarvalue>::iterator const &xv_end,
00327                         std::vector<cvm::real>::iterator         &result);
00328 
00331   static void inner_opt(colvarvalue const                        &x,
00332                         std::list<colvarvalue>::iterator         &xv,
00333                         std::list<colvarvalue>::iterator const   &xv_end,
00334                         std::vector<cvm::real>::iterator         &result);
00335 
00339   static void p2leg_opt(colvarvalue const                        &x,
00340                         std::vector<colvarvalue>::iterator       &xv,
00341                         std::vector<colvarvalue>::iterator const &xv_end,
00342                         std::vector<cvm::real>::iterator         &result);
00343 
00346   static void p2leg_opt(colvarvalue const                        &x,
00347                         std::list<colvarvalue>::iterator         &xv,
00348                         std::list<colvarvalue>::iterator const   &xv_end,
00349                         std::vector<cvm::real>::iterator         &result);
00350 
00351 };
00352 
00353 
00354 inline size_t colvarvalue::size() const
00355 {
00356   switch (value_type) {
00357   case colvarvalue::type_notset:
00358   default:
00359     return 0; break;
00360   case colvarvalue::type_scalar:
00361     return 1; break;
00362   case colvarvalue::type_3vector:
00363   case colvarvalue::type_unit3vector:
00364   case colvarvalue::type_unit3vectorderiv:
00365     return 3; break;
00366   case colvarvalue::type_quaternion:
00367   case colvarvalue::type_quaternionderiv:
00368     return 4; break;
00369   case colvarvalue::type_vector:
00370     return vector1d_value.size(); break;
00371   }
00372 }
00373 
00374 
00375 inline cvm::real colvarvalue::operator [] (int const i) const
00376 {
00377   switch (value_type) {
00378   case colvarvalue::type_notset:
00379   default:
00380     cvm::error("Error: trying to access a colvar value "
00381                "that is not initialized.\n", COLVARS_BUG_ERROR);
00382     return 0.0; break;
00383   case colvarvalue::type_scalar:
00384     return real_value; break;
00385   case colvarvalue::type_3vector:
00386   case colvarvalue::type_unit3vector:
00387   case colvarvalue::type_unit3vectorderiv:
00388     return rvector_value[i]; break;
00389   case colvarvalue::type_quaternion:
00390   case colvarvalue::type_quaternionderiv:
00391     return quaternion_value[i]; break;
00392   case colvarvalue::type_vector:
00393     return vector1d_value[i]; break;
00394   }
00395 }
00396 
00397 
00398 inline cvm::real & colvarvalue::operator [] (int const i)
00399 {
00400   switch (value_type) {
00401   case colvarvalue::type_notset:
00402   default:
00403     cvm::error("Error: trying to access a colvar value "
00404                "that is not initialized.\n", COLVARS_BUG_ERROR);
00405     return real_value; break;
00406   case colvarvalue::type_scalar:
00407     return real_value; break;
00408   case colvarvalue::type_3vector:
00409   case colvarvalue::type_unit3vector:
00410   case colvarvalue::type_unit3vectorderiv:
00411     return rvector_value[i]; break;
00412   case colvarvalue::type_quaternion:
00413   case colvarvalue::type_quaternionderiv:
00414     return quaternion_value[i]; break;
00415   case colvarvalue::type_vector:
00416     return vector1d_value[i]; break;
00417   }
00418 }
00419 
00420 
00421 inline int colvarvalue::check_types(colvarvalue const &x1,
00422                                     colvarvalue const &x2)
00423 {
00424   if (!colvarvalue::type_checking()) {
00425     return COLVARS_OK;
00426   }
00427 
00428   if (x1.type() != x2.type()) {
00429     if (((x1.type() == type_unit3vector) &&
00430          (x2.type() == type_unit3vectorderiv)) ||
00431         ((x2.type() == type_unit3vector) &&
00432          (x1.type() == type_unit3vectorderiv)) ||
00433         ((x1.type() == type_quaternion) &&
00434          (x2.type() == type_quaternionderiv)) ||
00435         ((x2.type() == type_quaternion) &&
00436          (x1.type() == type_quaternionderiv))) {
00437       return COLVARS_OK;
00438     } else {
00439       cvm::error("Trying to perform an operation between two colvar "
00440                  "values with different types, \""+
00441                  colvarvalue::type_desc(x1.type())+
00442                  "\" and \""+
00443                  colvarvalue::type_desc(x2.type())+
00444                  "\".\n");
00445       return COLVARS_ERROR;
00446     }
00447   }
00448 
00449   if (x1.type() == type_vector) {
00450     if (x1.vector1d_value.size() != x2.vector1d_value.size()) {
00451       cvm::error("Trying to perform an operation between two vector colvar "
00452                  "values with different sizes, "+
00453                  cvm::to_str(x1.vector1d_value.size())+
00454                  " and "+
00455                  cvm::to_str(x2.vector1d_value.size())+
00456                  ".\n");
00457       return COLVARS_ERROR;
00458     }
00459   }
00460   return COLVARS_OK;
00461 }
00462 
00463 
00464 inline int colvarvalue::check_types_assign(colvarvalue::Type const &vt1,
00465                                            colvarvalue::Type const &vt2)
00466 {
00467   if (!colvarvalue::type_checking()) {
00468     return COLVARS_OK;
00469   }
00470 
00471   if (vt1 != type_notset) {
00472     if (((vt1 == type_unit3vector) &&
00473          (vt2 == type_unit3vectorderiv)) ||
00474         ((vt2 == type_unit3vector) &&
00475          (vt1 == type_unit3vectorderiv)) ||
00476         ((vt1 == type_quaternion) &&
00477          (vt2 == type_quaternionderiv)) ||
00478         ((vt2 == type_quaternion) &&
00479          (vt1 == type_quaternionderiv))) {
00480       return COLVARS_OK;
00481     } else {
00482       if (vt1 != vt2) {
00483         cvm::error("Trying to assign a colvar value with type \""+
00484                    type_desc(vt2)+"\" to one with type \""+
00485                    type_desc(vt1)+"\".\n");
00486         return COLVARS_ERROR;
00487       }
00488     }
00489   }
00490   return COLVARS_OK;
00491 }
00492 
00493 
00494 inline colvarvalue & colvarvalue::operator = (colvarvalue const &x)
00495 {
00496   check_types_assign(this->type(), x.type());
00497   value_type = x.type();
00498 
00499   switch (this->type()) {
00500   case colvarvalue::type_scalar:
00501     this->real_value = x.real_value;
00502     break;
00503   case colvarvalue::type_3vector:
00504   case colvarvalue::type_unit3vector:
00505   case colvarvalue::type_unit3vectorderiv:
00506     this->rvector_value = x.rvector_value;
00507     break;
00508   case colvarvalue::type_quaternion:
00509   case colvarvalue::type_quaternionderiv:
00510     this->quaternion_value = x.quaternion_value;
00511     break;
00512   case colvarvalue::type_vector:
00513     vector1d_value = x.vector1d_value;
00514     elem_types = x.elem_types;
00515     elem_indices = x.elem_indices;
00516     elem_sizes = x.elem_sizes;
00517     break;
00518   case colvarvalue::type_notset:
00519   default:
00520     undef_op();
00521     break;
00522   }
00523   return *this;
00524 }
00525 
00526 
00527 inline void colvarvalue::operator += (colvarvalue const &x)
00528 {
00529   colvarvalue::check_types(*this, x);
00530 
00531   switch (this->type()) {
00532   case colvarvalue::type_scalar:
00533     this->real_value += x.real_value;
00534     break;
00535   case colvarvalue::type_3vector:
00536   case colvarvalue::type_unit3vector:
00537   case colvarvalue::type_unit3vectorderiv:
00538     this->rvector_value += x.rvector_value;
00539     break;
00540   case colvarvalue::type_quaternion:
00541   case colvarvalue::type_quaternionderiv:
00542     this->quaternion_value += x.quaternion_value;
00543     break;
00544   case colvarvalue::type_vector:
00545     this->vector1d_value += x.vector1d_value;
00546     break;
00547   case colvarvalue::type_notset:
00548   default:
00549     undef_op();
00550   }
00551 }
00552 
00553 
00554 inline void colvarvalue::operator -= (colvarvalue const &x)
00555 {
00556   colvarvalue::check_types(*this, x);
00557 
00558   switch (value_type) {
00559   case colvarvalue::type_scalar:
00560     real_value -= x.real_value;
00561     break;
00562   case colvarvalue::type_3vector:
00563   case colvarvalue::type_unit3vector:
00564   case colvarvalue::type_unit3vectorderiv:
00565     rvector_value -= x.rvector_value;
00566     break;
00567   case colvarvalue::type_quaternion:
00568   case colvarvalue::type_quaternionderiv:
00569     quaternion_value -= x.quaternion_value;
00570     break;
00571   case colvarvalue::type_vector:
00572     this->vector1d_value -= x.vector1d_value;
00573     break;
00574   case colvarvalue::type_notset:
00575   default:
00576     undef_op();
00577   }
00578 }
00579 
00580 
00581 inline void colvarvalue::operator *= (cvm::real const &a)
00582 {
00583   switch (value_type) {
00584   case colvarvalue::type_scalar:
00585     real_value *= a;
00586     break;
00587   case colvarvalue::type_3vector:
00588   case colvarvalue::type_unit3vectorderiv:
00589     rvector_value *= a;
00590     break;
00591   case colvarvalue::type_quaternion:
00592   case colvarvalue::type_quaternionderiv:
00593     quaternion_value *= a;
00594     break;
00595   case colvarvalue::type_vector:
00596     this->vector1d_value *= a;
00597     break;
00598   case colvarvalue::type_notset:
00599   default:
00600     undef_op();
00601   }
00602 }
00603 
00604 
00605 inline void colvarvalue::operator /= (cvm::real const &a)
00606 {
00607   switch (value_type) {
00608   case colvarvalue::type_scalar:
00609     real_value /= a; break;
00610   case colvarvalue::type_3vector:
00611   case colvarvalue::type_unit3vector:
00612   case colvarvalue::type_unit3vectorderiv:
00613     rvector_value /= a; break;
00614   case colvarvalue::type_quaternion:
00615   case colvarvalue::type_quaternionderiv:
00616     quaternion_value /= a; break;
00617   case colvarvalue::type_vector:
00618     this->vector1d_value /= a;
00619     break;
00620   case colvarvalue::type_notset:
00621   default:
00622     undef_op();
00623   }
00624 }
00625 
00626 
00627 inline cvm::vector1d<cvm::real> const colvarvalue::as_vector() const
00628 {
00629   switch (value_type) {
00630   case colvarvalue::type_scalar:
00631     {
00632       cvm::vector1d<cvm::real> v(1);
00633       v[0] = real_value;
00634       return v;
00635     }
00636   case colvarvalue::type_3vector:
00637   case colvarvalue::type_unit3vector:
00638   case colvarvalue::type_unit3vectorderiv:
00639     return rvector_value.as_vector();
00640   case colvarvalue::type_quaternion:
00641   case colvarvalue::type_quaternionderiv:
00642     return quaternion_value.as_vector();
00643   case colvarvalue::type_vector:
00644     return vector1d_value;
00645   case colvarvalue::type_notset:
00646   default:
00647     return cvm::vector1d<cvm::real>(0);
00648   }
00649 }
00650 
00651 
00652 inline cvm::real colvarvalue::norm2() const
00653 {
00654   switch (value_type) {
00655   case colvarvalue::type_scalar:
00656     return (this->real_value)*(this->real_value);
00657   case colvarvalue::type_3vector:
00658   case colvarvalue::type_unit3vector:
00659   case colvarvalue::type_unit3vectorderiv:
00660     return (this->rvector_value).norm2();
00661   case colvarvalue::type_quaternion:
00662   case colvarvalue::type_quaternionderiv:
00663     return (this->quaternion_value).norm2();
00664   case colvarvalue::type_vector:
00665     if (elem_types.size() > 0) {
00666       // if we have information about non-scalar types, use it
00667       cvm::real result = 0.0;
00668       size_t i;
00669       for (i = 0; i < elem_types.size(); i++) {
00670         result += (this->get_elem(i)).norm2();
00671       }
00672       return result;
00673     } else {
00674       return vector1d_value.norm2();
00675     }
00676     break;
00677   case colvarvalue::type_notset:
00678   default:
00679     return 0.0;
00680   }
00681 }
00682 
00683 
00684 inline cvm::real colvarvalue::sum() const
00685 {
00686   switch (value_type) {
00687   case colvarvalue::type_scalar:
00688     return (this->real_value);
00689   case colvarvalue::type_3vector:
00690   case colvarvalue::type_unit3vector:
00691   case colvarvalue::type_unit3vectorderiv:
00692     return (this->rvector_value).x + (this->rvector_value).y +
00693       (this->rvector_value).z;
00694   case colvarvalue::type_quaternion:
00695   case colvarvalue::type_quaternionderiv:
00696     return (this->quaternion_value).q0 + (this->quaternion_value).q1 +
00697       (this->quaternion_value).q2 + (this->quaternion_value).q3;
00698   case colvarvalue::type_vector:
00699     return (this->vector1d_value).sum();
00700   case colvarvalue::type_notset:
00701   default:
00702     return 0.0;
00703   }
00704 }
00705 
00706 
00707 inline cvm::real colvarvalue::dist2(colvarvalue const &x2) const
00708 {
00709   colvarvalue::check_types(*this, x2);
00710 
00711   switch (this->type()) {
00712   case colvarvalue::type_scalar:
00713     return (this->real_value - x2.real_value)*(this->real_value - x2.real_value);
00714   case colvarvalue::type_3vector:
00715     return (this->rvector_value - x2.rvector_value).norm2();
00716   case colvarvalue::type_unit3vector:
00717   case colvarvalue::type_unit3vectorderiv:
00718     // angle between (*this) and x2 is the distance
00719     return cvm::acos(this->rvector_value * x2.rvector_value) * cvm::acos(this->rvector_value * x2.rvector_value);
00720   case colvarvalue::type_quaternion:
00721   case colvarvalue::type_quaternionderiv:
00722     // angle between (*this) and x2 is the distance, the quaternion
00723     // object has it implemented internally
00724     return this->quaternion_value.dist2(x2.quaternion_value);
00725   case colvarvalue::type_vector:
00726     return (this->vector1d_value - x2.vector1d_value).norm2();
00727   case colvarvalue::type_notset:
00728   default:
00729     this->undef_op();
00730     return 0.0;
00731   };
00732 }
00733 
00734 
00735 #endif

Generated on Mon Apr 22 04:26:23 2024 for VMD (current) by doxygen1.2.14 written by Dimitri van Heesch, © 1997-2002