source: src/vector_ops.cpp@ 436f04

Action_Thermostats Add_AtomRandomPerturbation Add_FitFragmentPartialChargesAction Add_RotateAroundBondAction Add_SelectAtomByNameAction Added_ParseSaveFragmentResults AddingActions_SaveParseParticleParameters Adding_Graph_to_ChangeBondActions Adding_MD_integration_tests Adding_ParticleName_to_Atom Adding_StructOpt_integration_tests AtomFragments Automaking_mpqc_open AutomationFragmentation_failures Candidate_v1.5.4 Candidate_v1.6.0 Candidate_v1.6.1 ChangeBugEmailaddress ChangingTestPorts ChemicalSpaceEvaluator CombiningParticlePotentialParsing Combining_Subpackages Debian_Package_split Debian_package_split_molecuildergui_only Disabling_MemDebug Docu_Python_wait EmpiricalPotential_contain_HomologyGraph EmpiricalPotential_contain_HomologyGraph_documentation Enable_parallel_make_install Enhance_userguide Enhanced_StructuralOptimization Enhanced_StructuralOptimization_continued Example_ManyWaysToTranslateAtom Exclude_Hydrogens_annealWithBondGraph FitPartialCharges_GlobalError Fix_BoundInBox_CenterInBox_MoleculeActions Fix_ChargeSampling_PBC Fix_ChronosMutex Fix_FitPartialCharges Fix_FitPotential_needs_atomicnumbers Fix_ForceAnnealing Fix_IndependentFragmentGrids Fix_ParseParticles Fix_ParseParticles_split_forward_backward_Actions Fix_PopActions Fix_QtFragmentList_sorted_selection Fix_Restrictedkeyset_FragmentMolecule Fix_StatusMsg Fix_StepWorldTime_single_argument Fix_Verbose_Codepatterns Fix_fitting_potentials Fixes ForceAnnealing_goodresults ForceAnnealing_oldresults ForceAnnealing_tocheck ForceAnnealing_with_BondGraph ForceAnnealing_with_BondGraph_continued ForceAnnealing_with_BondGraph_continued_betteresults ForceAnnealing_with_BondGraph_contraction-expansion FragmentAction_writes_AtomFragments FragmentMolecule_checks_bonddegrees GeometryObjects Gui_Fixes Gui_displays_atomic_force_velocity ImplicitCharges IndependentFragmentGrids IndependentFragmentGrids_IndividualZeroInstances IndependentFragmentGrids_IntegrationTest IndependentFragmentGrids_Sole_NN_Calculation JobMarket_RobustOnKillsSegFaults JobMarket_StableWorkerPool JobMarket_unresolvable_hostname_fix MoreRobust_FragmentAutomation ODR_violation_mpqc_open PartialCharges_OrthogonalSummation PdbParser_setsAtomName PythonUI_with_named_parameters QtGui_reactivate_TimeChanged_changes Recreated_GuiChecks Rewrite_FitPartialCharges RotateToPrincipalAxisSystem_UndoRedo SaturateAtoms_findBestMatching SaturateAtoms_singleDegree StoppableMakroAction Subpackage_CodePatterns Subpackage_JobMarket Subpackage_LinearAlgebra Subpackage_levmar Subpackage_mpqc_open Subpackage_vmg Switchable_LogView ThirdParty_MPQC_rebuilt_buildsystem TrajectoryDependenant_MaxOrder TremoloParser_IncreasedPrecision TremoloParser_MultipleTimesteps TremoloParser_setsAtomName Ubuntu_1604_changes stable
Last change on this file since 436f04 was 112b09, checked in by Tillmann Crueger <crueger@…>, 15 years ago

Added #include "Helpers/MemDebug.hpp" to all .cpp files

  • Property mode set to 100644
File size: 2.7 KB
Line 
1/*
2 * vector_ops.cpp
3 *
4 * Created on: Apr 1, 2010
5 * Author: crueger
6 */
7
8#include "Helpers/MemDebug.hpp"
9
10#include "vector.hpp"
11#include "Plane.hpp"
12#include "log.hpp"
13#include "verbose.hpp"
14#include "gslmatrix.hpp"
15#include "leastsquaremin.hpp"
16#include "info.hpp"
17#include "Helpers/fast_functions.hpp"
18#include "Exceptions/LinearDependenceException.hpp"
19#include "Exceptions/SkewException.hpp"
20
21#include <gsl/gsl_linalg.h>
22#include <gsl/gsl_matrix.h>
23#include <gsl/gsl_permutation.h>
24#include <gsl/gsl_vector.h>
25
26/**
27 * !@file
28 * These files defines several common operation on vectors that should not
29 * become part of the main vector class, because they are either to complex
30 * or need methods from other subsystems that should not be moved to
31 * the LinAlg-Subsystem
32 */
33
34/** Creates a new vector as the one with least square distance to a given set of \a vectors.
35 * \param *vectors set of vectors
36 * \param num number of vectors
37 * \return true if success, false if failed due to linear dependency
38 */
39bool LSQdistance(Vector &res,const Vector **vectors, int num)
40{
41 int j;
42
43 for (j=0;j<num;j++) {
44 Log() << Verbose(1) << j << "th atom's vector: " << vectors[j] << endl;
45 }
46
47 int np = 3;
48 struct LSQ_params par;
49
50 const gsl_multimin_fminimizer_type *T =
51 gsl_multimin_fminimizer_nmsimplex;
52 gsl_multimin_fminimizer *s = NULL;
53 gsl_vector *ss, *y;
54 gsl_multimin_function minex_func;
55
56 size_t iter = 0, i;
57 int status;
58 double size;
59
60 /* Initial vertex size vector */
61 ss = gsl_vector_alloc (np);
62 y = gsl_vector_alloc (np);
63
64 /* Set all step sizes to 1 */
65 gsl_vector_set_all (ss, 1.0);
66
67 /* Starting point */
68 par.vectors = vectors;
69 par.num = num;
70
71 for (i=NDIM;i--;)
72 gsl_vector_set(y, i, (vectors[0]->at(i) - vectors[1]->at(i))/2.);
73
74 /* Initialize method and iterate */
75 minex_func.f = &LSQ;
76 minex_func.n = np;
77 minex_func.params = (void *)&par;
78
79 s = gsl_multimin_fminimizer_alloc (T, np);
80 gsl_multimin_fminimizer_set (s, &minex_func, y, ss);
81
82 do
83 {
84 iter++;
85 status = gsl_multimin_fminimizer_iterate(s);
86
87 if (status)
88 break;
89
90 size = gsl_multimin_fminimizer_size (s);
91 status = gsl_multimin_test_size (size, 1e-2);
92
93 if (status == GSL_SUCCESS)
94 {
95 printf ("converged to minimum at\n");
96 }
97
98 printf ("%5d ", (int)iter);
99 for (i = 0; i < (size_t)np; i++)
100 {
101 printf ("%10.3e ", gsl_vector_get (s->x, i));
102 }
103 printf ("f() = %7.3f size = %.3f\n", s->fval, size);
104 }
105 while (status == GSL_CONTINUE && iter < 100);
106
107 for (i=(size_t)np;i--;)
108 res[i] = gsl_vector_get(s->x, i);
109 gsl_vector_free(y);
110 gsl_vector_free(ss);
111 gsl_multimin_fminimizer_free (s);
112
113 return true;
114};
Note: See TracBrowser for help on using the repository browser.