MPQC  2.3.1
coor.h
1 //
2 // coor.h
3 //
4 // Copyright (C) 1996 Limit Point Systems, Inc.
5 //
6 // Author: Curtis Janssen <cljanss@limitpt.com>
7 // Maintainer: LPS
8 //
9 // This file is part of the SC Toolkit.
10 //
11 // The SC Toolkit is free software; you can redistribute it and/or modify
12 // it under the terms of the GNU Library General Public License as published by
13 // the Free Software Foundation; either version 2, or (at your option)
14 // any later version.
15 //
16 // The SC Toolkit is distributed in the hope that it will be useful,
17 // but WITHOUT ANY WARRANTY; without even the implied warranty of
18 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 // GNU Library General Public License for more details.
20 //
21 // You should have received a copy of the GNU Library General Public License
22 // along with the SC Toolkit; see the file COPYING.LIB. If not, write to
23 // the Free Software Foundation, 675 Mass Ave, Cambridge, MA 02139, USA.
24 //
25 // The U.S. Government is granted a limited license as per AL 91-7.
26 //
27 
28 #ifndef _chemistry_molecule_coor_h
29 #define _chemistry_molecule_coor_h
30 
31 #ifdef __GNUC__
32 #pragma interface
33 #endif
34 
35 #include <iostream>
36 #include <vector>
37 
38 #include <math/scmat/matrix.h>
39 #include <math/optimize/transform.h>
40 #include <chemistry/molecule/molecule.h>
41 
42 namespace sc {
43 
46 class IntCoor: public SavableState {
47  protected:
48  // conversion factors from radians, bohr to the preferred units
49  static double bohr_conv;
50  static double radian_conv;
51  char *label_;
52  double value_;
53  public:
54  IntCoor(StateIn&);
55  IntCoor(const IntCoor&);
58  IntCoor(const char* label = 0);
75  IntCoor(const Ref<KeyVal>&);
76 
77  virtual ~IntCoor();
79 
81  virtual const char* label() const;
83  virtual double value() const;
85  virtual void set_value(double);
87  virtual double preferred_value() const;
89  virtual const char* ctype() const = 0;
91  virtual void print(std::ostream & o=ExEnv::out0()) const;
92  virtual void print_details(const Ref<Molecule> &, std::ostream& =ExEnv::out0()) const;
95  virtual double force_constant(Ref<Molecule>&) = 0;
97  virtual void update_value(const Ref<Molecule>&) = 0;
99  virtual void bmat(const Ref<Molecule>&,RefSCVector&bmat,double coef=1.0) = 0;
103  virtual int equivalent(Ref<IntCoor>&) = 0;
104 };
105 
120 class SumIntCoor: public IntCoor {
121  private:
122  std::vector<double> coef_;
123  std::vector<Ref<IntCoor> > coor_;
124  public:
128  SumIntCoor(const char *);
139  SumIntCoor(const Ref<KeyVal>&);
140 
141  ~SumIntCoor();
142  void save_data_state(StateOut&);
143 
145  int n();
148  void add(Ref<IntCoor>&,double coef);
150  void normalize();
151 
152  // IntCoor overrides
154  double preferred_value() const;
156  const char* ctype() const;
158  void print_details(const Ref<Molecule> &, std::ostream& =ExEnv::out0()) const;
160  double force_constant(Ref<Molecule>&);
162  void update_value(const Ref<Molecule>&);
164  void bmat(const Ref<Molecule>&,RefSCVector&bmat,double coef = 1.0);
166  int equivalent(Ref<IntCoor>&);
167 };
168 
189 class SetIntCoor: public SavableState {
190  private:
191  std::vector<Ref<IntCoor> > coor_;
192  public:
193  SetIntCoor();
205  SetIntCoor(const Ref<KeyVal>&);
206 
207  virtual ~SetIntCoor();
208  void save_data_state(StateOut&);
209 
211  void add(const Ref<IntCoor>&);
213  void add(const Ref<SetIntCoor>&);
215  void pop();
217  void clear();
219  int n() const;
221  Ref<IntCoor> coor(int i) const;
223  virtual void fd_bmat(const Ref<Molecule>&,RefSCMatrix&);
225  virtual void bmat(const Ref<Molecule>&, RefSCMatrix&);
231  virtual void print_details(const Ref<Molecule> &,std::ostream& =ExEnv::out0()) const;
233  virtual void update_values(const Ref<Molecule>&);
235  virtual void values_to_vector(const RefSCVector&);
236 };
237 
238 
239 // ////////////////////////////////////////////////////////////////////////
240 
241 class BitArrayLTri;
242 
246 {
247  protected:
248  Ref<Molecule> molecule_;
249 
250  int linear_bends_;
251  int linear_lbends_;
252  int linear_tors_;
253  int linear_stors_;
254  int nextra_bonds_;
255  int *extra_bonds_;
256  double linear_bend_thres_;
257  double linear_tors_thres_;
258  double radius_scale_factor_;
259 
260  void init_constants();
261 
262  double cos_ijk(Molecule& m, int i, int j, int k);
263  int hterminal(Molecule& m, BitArrayLTri& bonds, int i);
264  int nearest_contact(int i, Molecule& m);
265 
266  void add_bonds(const Ref<SetIntCoor>& list, BitArrayLTri& bonds, Molecule& m);
267  void add_bends(const Ref<SetIntCoor>& list, BitArrayLTri& bonds, Molecule& m);
268  void add_tors(const Ref<SetIntCoor>& list, BitArrayLTri& bonds, Molecule& m);
269  void add_out(const Ref<SetIntCoor>& list, BitArrayLTri& bonds, Molecule& m);
270  public:
274  IntCoorGen(const Ref<Molecule>&, int nextra=0, int *extra=0);
316  IntCoorGen(const Ref<KeyVal>&);
318 
319  ~IntCoorGen();
320 
322  void save_data_state(StateOut&);
323 
325  virtual void generate(const Ref<SetIntCoor>&);
326 
328  virtual void print(std::ostream& out=ExEnv::out0()) const;
329 };
330 
331 
332 // ////////////////////////////////////////////////////////////////////////
333 
334 
339 {
340  protected:
341  Ref<Molecule> molecule_;
342  RefSCDimension dnatom3_; // the number of atoms x 3
343  Ref<SCMatrixKit> matrixkit_; // used to construct matrices
344 
345  int debug_;
346  public:
365  MolecularCoor(const Ref<KeyVal>&);
366 
367  virtual ~MolecularCoor();
368 
369  void save_data_state(StateOut&);
370 
373  RefSCDimension dim_natom3() { return dnatom3_; }
374 
376  Ref<Molecule> molecule() const { return molecule_; }
377 
379  virtual void print(std::ostream& =ExEnv::out0()) const = 0;
380  virtual void print_simples(std::ostream& =ExEnv::out0()) const = 0;
381 
385  virtual RefSCDimension dim() = 0;
386 
390  int to_cartesian(const RefSCVector&internal);
391  virtual int to_cartesian(const Ref<Molecule>&mol,
392  const RefSCVector&internal) = 0;
393 
397  virtual int to_internal(RefSCVector&internal) = 0;
398 
403  virtual int to_cartesian(RefSCVector&cartesian,RefSCVector&internal) = 0;
404 
409  virtual int to_internal(RefSCVector&internal,RefSCVector&cartesian) = 0;
410 
414  virtual int to_cartesian(RefSymmSCMatrix&cartesian,
415  RefSymmSCMatrix&internal) =0;
416 
420  virtual int to_internal(RefSymmSCMatrix&internal,
421  RefSymmSCMatrix&cartesian) = 0;
422 
425  virtual void guess_hessian(RefSymmSCMatrix&hessian) = 0;
426 
430 
432  virtual int nconstrained();
433 
438 
439  Ref<SCMatrixKit> matrixkit() const { return matrixkit_; }
440 };
441 
442 
446 {
447  protected:
448  Ref<IntCoorGen> generator_;
449 
450  void form_K_matrix(RefSCDimension& dredundant,
451  RefSCDimension& dfixed,
452  RefSCMatrix& K,
453  int*& is_totally_symmetric);
454 
455  RefSCDimension dim_; // corresponds to the number of variable coordinates
456  RefSCDimension dvc_; // the number of variable + constant coordinates
457 
458  Ref<SetIntCoor> variable_; // the variable internal coordinates
459  Ref<SetIntCoor> constant_; // the constant internal coordinates
460 
461  Ref<SetIntCoor> fixed_;
462  Ref<SetIntCoor> watched_;
463  Ref<IntCoor> followed_;
464 
465  // these are all of the basic coordinates
466  Ref<SetIntCoor> bonds_;
467  Ref<SetIntCoor> bends_;
468  Ref<SetIntCoor> tors_;
469  Ref<SetIntCoor> outs_;
470  // these are provided by the user or generated coordinates that
471  // could not be assigned to any of the above catagories
472  Ref<SetIntCoor> extras_;
473 
474  Ref<SetIntCoor> all_;
475 
476  // Useful relationships
477  // variable_->n() + constant_->n() = 3N-6(5)
478  // symm_->n() + asymm_->n() = 3N-6(5)
479 
480  int update_bmat_; // if 1 recompute the b matrix during to_cartesian
481  int only_totally_symmetric_; // only coors with tot. symm comp. are varied
482  double symmetry_tolerance_; // tol used to find coors with tot. sym. comp.
483  double simple_tolerance_; // tol used to see if a simple is included
484  double coordinate_tolerance_; // tol used to see if a coor is included
485  double cartesian_tolerance_; // tol used in intco->cart transformation
486  double scale_bonds_; // scale factor for bonds
487  double scale_bends_; // scale factor for bends
488  double scale_tors_; // scale factor for tors
489  double scale_outs_; // scale factor for outs
490 
491  int nextra_bonds_;
492  int* extra_bonds_;
493 
494  int given_fixed_values_; // if true make molecule have given fixed values
495 
496  int decouple_bonds_;
497  int decouple_bends_;
498 
499  int max_update_steps_;
500  double max_update_disp_;
501 
505  virtual void init();
508  virtual void new_coords();
510  virtual void read_keyval(const Ref<KeyVal>&);
511 
512  // control whether or not to print coordinates when they are formed
513  int form_print_simples_;
514  int form_print_variable_;
515  int form_print_constant_;
516  int form_print_molecule_;
517  public:
625 
626  virtual ~IntMolecularCoor();
627  void save_data_state(StateOut&);
628 
631  virtual void form_coordinates(int keep_variable=0) =0;
632 
635  virtual int all_to_cartesian(const Ref<Molecule> &,RefSCVector&internal);
638  virtual int all_to_internal(const Ref<Molecule> &,RefSCVector&internal);
639 
642  virtual RefSCDimension dim();
643  virtual int to_cartesian(const Ref<Molecule> &,const RefSCVector&internal);
644  virtual int to_internal(RefSCVector&internal);
645  virtual int to_cartesian(RefSCVector&cartesian,RefSCVector&internal);
646  virtual int to_internal(RefSCVector&internal,RefSCVector&cartesian);
647  virtual int to_cartesian(RefSymmSCMatrix&cart,RefSymmSCMatrix&internal);
648  virtual int to_internal(RefSymmSCMatrix&internal,RefSymmSCMatrix&cart);
649  virtual void print(std::ostream& =ExEnv::out0()) const;
650  virtual void print_simples(std::ostream& =ExEnv::out0()) const;
651  virtual void print_variable(std::ostream& =ExEnv::out0()) const;
652  virtual void print_constant(std::ostream& =ExEnv::out0()) const;
653  int nconstrained();
654 };
655 
656 // ///////////////////////////////////////////////////////////////////////
657 
668 {
669  protected:
670  // true if coordinates should be changed during optimization
671  int change_coordinates_;
672  // true if hessian should be transformed too (should always be true)
673  int transform_hessian_;
674  // max value for the condition number if coordinates can be changed
675  double max_kappa2_;
676 
677  void init();
678  public:
700 
701  virtual ~SymmMolecularCoor();
702  void save_data_state(StateOut&);
703 
706  void form_coordinates(int keep_variable=0);
707 
709  void guess_hessian(RefSymmSCMatrix&hessian);
712 
716 
717  void print(std::ostream& =ExEnv::out0()) const;
718 };
719 
720 // ///////////////////////////////////////////////////////////////////////
721 
725 {
726 
727  public:
732 
733  virtual ~RedundMolecularCoor();
734  void save_data_state(StateOut&);
735 
738  void form_coordinates(int keep_variable=0);
740  void guess_hessian(RefSymmSCMatrix&hessian);
743 };
744 
745 // ///////////////////////////////////////////////////////////////////////
746 
752 {
753  private:
754  protected:
755  RefSCDimension dim_; // the number of atoms x 3
756 
758  virtual void init();
759  public:
764 
765  virtual ~CartMolecularCoor();
766 
767  void save_data_state(StateOut&);
768 
770  virtual RefSCDimension dim();
771  virtual int to_cartesian(const Ref<Molecule>&,const RefSCVector&internal);
772  virtual int to_internal(RefSCVector&internal);
773  virtual int to_cartesian(RefSCVector&cartesian,RefSCVector&internal);
774  virtual int to_internal(RefSCVector&internal,RefSCVector&cartesian);
775  virtual int to_cartesian(RefSymmSCMatrix&cart,RefSymmSCMatrix&internal);
776  virtual int to_internal(RefSymmSCMatrix&internal,RefSymmSCMatrix&cart);
777  virtual void print(std::ostream& =ExEnv::out0()) const;
778  virtual void print_simples(std::ostream& =ExEnv::out0()) const;
779  void guess_hessian(RefSymmSCMatrix&hessian);
781 };
782 
783 }
784 
785 #endif
786 
787 // Local Variables:
788 // mode: c++
789 // c-file-style: "CLJ"
790 // End:
sc::MolecularCoor::save_data_state
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
sc::IntMolecularCoor::all_to_internal
virtual int all_to_internal(const Ref< Molecule > &, RefSCVector &internal)
Like to_internal(), except all internal coordinates are considered, not just the variable ones.
sc::CartMolecularCoor::init
virtual void init()
Initializes the dimensions.
sc::CartMolecularCoor::inverse_hessian
RefSymmSCMatrix inverse_hessian(RefSymmSCMatrix &)
Given an Hessian, return the inverse of that hessian.
sc::SetIntCoor
The SetIntCoor class describes a set of internal coordinates.
Definition: coor.h:189
sc::IntCoor::update_value
virtual void update_value(const Ref< Molecule > &)=0
Recalculate the value of the coordinate.
sc::SymmMolecularCoor::inverse_hessian
RefSymmSCMatrix inverse_hessian(RefSymmSCMatrix &)
Invert the hessian.
sc::SetIntCoor::clear
void clear()
Removes all coordinates from the set.
sc::SetIntCoor::bmat
virtual void bmat(const Ref< Molecule > &, RefSCMatrix &)
Compute the B matrix the old-fashioned way.
sc::RedundMolecularCoor
The RedundMolecularCoor class provides a redundant set of simple internal coordinates.
Definition: coor.h:724
sc::SumIntCoor::force_constant
double force_constant(Ref< Molecule > &)
Returns the weighted sum of the individual force constants.
sc::BitArrayLTri
Definition: bitarray.h:68
sc::IntCoor::ctype
virtual const char * ctype() const =0
Returns a string representation of the type of coordinate this is.
sc::RefSymmSCMatrix
The RefSymmSCMatrix class is a smart pointer to an SCSymmSCMatrix specialization.
Definition: matrix.h:261
sc::Molecule
The Molecule class contains information about molecules.
Definition: molecule.h:127
sc::IntMolecularCoor::to_internal
virtual int to_internal(RefSCVector &internal)
Fill in the vector `‘internal’' with the current internal coordinates.
sc::CartMolecularCoor::guess_hessian
void guess_hessian(RefSymmSCMatrix &hessian)
Calculate an approximate hessian and place the result in `‘hessian’'.
sc::SumIntCoor::add
void add(Ref< IntCoor > &, double coef)
Add a coordinate to the linear combination.
sc::RedundMolecularCoor::inverse_hessian
RefSymmSCMatrix inverse_hessian(RefSymmSCMatrix &)
Invert the hessian.
sc::SetIntCoor::n
int n() const
Returns the number of coordinates in the set.
sc::RefSCMatrix
The RefSCMatrix class is a smart pointer to an SCMatrix specialization.
Definition: matrix.h:135
sc::MolecularCoor::to_cartesian
int to_cartesian(const RefSCVector &internal)
Given a set of displaced internal coordinates, update the cartesian coordinates of the Molecule conta...
sc::SetIntCoor::update_values
virtual void update_values(const Ref< Molecule > &)
Recalculate the values of the internal coordinates in the set.
sc::SetIntCoor::pop
void pop()
Removes the last coordinate from this set.
sc::Ref
A template class that maintains references counts.
Definition: ref.h:332
sc::SumIntCoor
SumIntCoor is used to construct linear combinations of internal coordinates.
Definition: coor.h:120
sc::SumIntCoor::ctype
const char * ctype() const
Always returns `‘SUM’'.
sc::MolecularCoor::dim
virtual RefSCDimension dim()=0
Returns a smart reference to an SCDimension equal to the number of coordinates (be they Cartesian,...
sc::MolecularCoor::guess_hessian
virtual void guess_hessian(RefSymmSCMatrix &hessian)=0
Calculate an approximate hessian and place the result in `‘hessian’'.
sc::SymmMolecularCoor::save_data_state
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
sc::SymmMolecularCoor::print
void print(std::ostream &=ExEnv::out0()) const
Print the coordinate.
sc::SymmMolecularCoor::guess_hessian
void guess_hessian(RefSymmSCMatrix &hessian)
Form the approximate hessian.
sc::RedundMolecularCoor::save_data_state
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
sc::SymmMolecularCoor::init
void init()
This is called by the constructors of classes derived from IntMolecularCoor.
sc::IntMolecularCoor::nconstrained
int nconstrained()
Returns the number of constrained coordinates.
sc::SetIntCoor::save_data_state
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
sc::CartMolecularCoor::save_data_state
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
sc::IntMolecularCoor::new_coords
virtual void new_coords()
Allocates memory for the SetIntCoor's used to store the simple and internal coordinates.
sc::IntCoor::preferred_value
virtual double preferred_value() const
Returns the value of the coordinate in more familiar units.
sc::StateIn
Restores objects that derive from SavableState.
Definition: statein.h:70
sc::IntCoorGen::save_data_state
void save_data_state(StateOut &)
Standard member.
sc::RefSCDimension
The RefSCDimension class is a smart pointer to an SCDimension specialization.
Definition: dim.h:156
sc::MolecularCoor::inverse_hessian
virtual RefSymmSCMatrix inverse_hessian(RefSymmSCMatrix &)=0
Given an Hessian, return the inverse of that hessian.
sc::SumIntCoor::print_details
void print_details(const Ref< Molecule > &, std::ostream &=ExEnv::out0()) const
Print the individual coordinates in the sum with their coefficients.
sc::SetIntCoor::print_details
virtual void print_details(const Ref< Molecule > &, std::ostream &=ExEnv::out0()) const
Print the coordinates in the set.
sc::IntCoor::label
virtual const char * label() const
Returns the string containing the label for the internal coordinate.
sc::SumIntCoor::preferred_value
double preferred_value() const
Returns the value of the coordinate in a.u. and radians.
sc::IntCoor::print
virtual void print(std::ostream &o=ExEnv::out0()) const
Print information about the coordinate.
sc::RedundMolecularCoor::form_coordinates
void form_coordinates(int keep_variable=0)
Actually form the variable and constant internal coordinates from the simple internal coordinates.
sc::IntCoor::save_data_state
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
sc::IntCoorGen
IntCoorGen generates a set of simple internal coordinates for a molecule.
Definition: coor.h:245
sc::CartMolecularCoor::dim
virtual RefSCDimension dim()
These implement the virtual functions inherited from MolecularCoor.
sc::RedundMolecularCoor::guess_hessian
void guess_hessian(RefSymmSCMatrix &hessian)
Form the approximate hessian.
sc::SumIntCoor::n
int n()
Returns the number of coordinates in this linear combination.
sc::MolecularCoor::to_internal
virtual int to_internal(RefSCVector &internal)=0
Fill in the vector `‘internal’' with the current internal coordinates.
sc::IntMolecularCoor
The IntMolecularCoor abstract class describes a molecule's coordinates in terms of internal coordinat...
Definition: coor.h:445
sc::IntMolecularCoor::dim
virtual RefSCDimension dim()
These implement the virtual functions inherited from MolecularCoor.
sc::IntMolecularCoor::read_keyval
virtual void read_keyval(const Ref< KeyVal > &)
Reads the KeyVal input.
sc::IntCoor::force_constant
virtual double force_constant(Ref< Molecule > &)=0
Returns the value of the force constant associated with this coordinate.
sc::SetIntCoor::coor
Ref< IntCoor > coor(int i) const
Returns a reference to the i'th coordinate in the set.
sc::IntMolecularCoor::all_to_cartesian
virtual int all_to_cartesian(const Ref< Molecule > &, RefSCVector &internal)
Like to_cartesians(), except all internal coordinates are considered, not just the variable ones.
sc::IntCoor::bmat
virtual void bmat(const Ref< Molecule > &, RefSCVector &bmat, double coef=1.0)=0
Fill in a row the the B matrix.
sc::MolecularCoor::nconstrained
virtual int nconstrained()
Returns the number of constrained coordinates.
sc::RefSCVector
The RefSCVector class is a smart pointer to an SCVector specialization.
Definition: matrix.h:55
sc::SymmMolecularCoor::form_coordinates
void form_coordinates(int keep_variable=0)
Actually form the variable and constant internal coordinates from simple internal coordinates.
sc::StateOut
Serializes objects that derive from SavableState.
Definition: stateout.h:61
sc::SetIntCoor::guess_hessian
virtual void guess_hessian(Ref< Molecule > &, RefSymmSCMatrix &)
Create an approximate Hessian for this set of coordinates.
sc::MolecularCoor::molecule
Ref< Molecule > molecule() const
Returns the molecule.
Definition: coor.h:376
sc::IntMolecularCoor::save_data_state
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
sc::SumIntCoor::save_data_state
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
sc::IntCoorGen::print
virtual void print(std::ostream &out=ExEnv::out0()) const
Print out information about this.
sc::MolecularCoor
The MolecularCoor abstract class describes the coordinate system used to describe a molecule.
Definition: coor.h:338
sc::CartMolecularCoor::to_internal
virtual int to_internal(RefSCVector &internal)
Fill in the vector `‘internal’' with the current internal coordinates.
sc::SumIntCoor::normalize
void normalize()
This function normalizes all the coefficients.
sc::IntCoorGen::IntCoorGen
IntCoorGen(const Ref< Molecule > &, int nextra=0, int *extra=0)
Create an IntCoorGen given a Molecule and, optionally, extra bonds.
sc::ExEnv::out0
static std::ostream & out0()
Return an ostream that writes from node 0.
sc::IntCoor
The IntCoor abstract class describes an internal coordinate of a molecule.
Definition: coor.h:46
sc::MolecularCoor::dim_natom3
RefSCDimension dim_natom3()
Returns a smart reference to an SCDimension equal to the number of atoms in the molecule times 3.
Definition: coor.h:373
sc::MolecularCoor::print
virtual void print(std::ostream &=ExEnv::out0()) const =0
Print the coordinate.
sc::MolecularCoor::change_coordinates
virtual Ref< NonlinearTransform > change_coordinates()
When this is called, MoleculeCoor may select a new internal coordinate system and return a transform ...
sc::SavableState
Base class for objects that can save/restore state.
Definition: state.h:46
sc::IntMolecularCoor::print
virtual void print(std::ostream &=ExEnv::out0()) const
Print the coordinate.
sc::SetIntCoor::add
void add(const Ref< IntCoor > &)
Adds an internal coordinate to the set.
sc::SetIntCoor::values_to_vector
virtual void values_to_vector(const RefSCVector &)
Copy the values of the internal coordinates to a vector.
sc::IntCoor::equivalent
virtual int equivalent(Ref< IntCoor > &)=0
Test to see if this internal coordinate is equivalent to that one.
sc::IntCoorGen::generate
virtual void generate(const Ref< SetIntCoor > &)
This generates a set of internal coordinates.
sc::SymmMolecularCoor
The SymmMolecularCoor class derives from IntMolecularCoor.
Definition: coor.h:667
sc::CartMolecularCoor
The CartMolecularCoor class implements Cartesian coordinates in a way suitable for use in geometry op...
Definition: coor.h:751
sc::IntMolecularCoor::init
virtual void init()
This is called by the constructors of classes derived from IntMolecularCoor.
sc::SumIntCoor::equivalent
int equivalent(Ref< IntCoor > &)
Always returns 0.
sc::SumIntCoor::update_value
void update_value(const Ref< Molecule > &)
Recalculate the value of the coordinate.
sc::SymmMolecularCoor::change_coordinates
Ref< NonlinearTransform > change_coordinates()
This overrides MoleculeCoor's change_coordinates and might transform to a new set of coordinates.
sc::CartMolecularCoor::print
virtual void print(std::ostream &=ExEnv::out0()) const
Print the coordinate.
sc::IntCoor::value
virtual double value() const
Returns the value of the coordinate in atomic units or radians.
sc::IntCoor::set_value
virtual void set_value(double)
Sets the value of the coordinate in atomic units or radians.
sc::SumIntCoor::bmat
void bmat(const Ref< Molecule > &, RefSCVector &bmat, double coef=1.0)
Fill in a row the the B matrix.
sc::IntMolecularCoor::form_coordinates
virtual void form_coordinates(int keep_variable=0)=0
Actually form the variable and constant internal coordinates from the simple internal coordinates.
sc::SetIntCoor::fd_bmat
virtual void fd_bmat(const Ref< Molecule > &, RefSCMatrix &)
Compute the B matrix by finite displacements.

Generated at Sun Jan 26 2020 23:33:03 for MPQC 2.3.1 using the documentation package Doxygen 1.8.16.