source: src/molecule_geometry.cpp@ 72d108

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 72d108 was 1f91f4, checked in by Frederik Heber <heber@…>, 14 years ago

factored functionality in PrincipalAxisSystemAction() and RotateToPrincipalAxisSystemAction() into class molecule:

  • Property mode set to 100644
File size: 17.3 KB
Line 
1/*
2 * Project: MoleCuilder
3 * Description: creates and alters molecular systems
4 * Copyright (C) 2010 University of Bonn. All rights reserved.
5 * Please see the LICENSE file or "Copyright notice" in builder.cpp for details.
6 */
7
8/*
9 * molecule_geometry.cpp
10 *
11 * Created on: Oct 5, 2009
12 * Author: heber
13 */
14
15// include config.h
16#ifdef HAVE_CONFIG_H
17#include <config.h>
18#endif
19
20#include "CodePatterns/MemDebug.hpp"
21
22#include "Helpers/helpers.hpp"
23#include "CodePatterns/Log.hpp"
24#include "CodePatterns/Verbose.hpp"
25#include "LinearAlgebra/Line.hpp"
26#include "LinearAlgebra/RealSpaceMatrix.hpp"
27#include "LinearAlgebra/Plane.hpp"
28
29#include "atom.hpp"
30#include "bond.hpp"
31#include "config.hpp"
32#include "element.hpp"
33#include "LinearAlgebra/leastsquaremin.hpp"
34#include "molecule.hpp"
35#include "World.hpp"
36
37#include "Box.hpp"
38
39#include <boost/foreach.hpp>
40
41#include <gsl/gsl_eigen.h>
42#include <gsl/gsl_multimin.h>
43
44
45/************************************* Functions for class molecule *********************************/
46
47
48/** Centers the molecule in the box whose lengths are defined by vector \a *BoxLengths.
49 * \param *out output stream for debugging
50 */
51bool molecule::CenterInBox()
52{
53 bool status = true;
54 const Vector *Center = DetermineCenterOfAll();
55 const Vector *CenterBox = DetermineCenterOfBox();
56 Box &domain = World::getInstance().getDomain();
57
58 // go through all atoms
59 BOOST_FOREACH(atom* iter, atoms){
60 *iter -= *Center;
61 *iter -= *CenterBox;
62 }
63 atoms.transformNodes(boost::bind(&Box::WrapPeriodically,domain,_1));
64
65 delete(Center);
66 delete(CenterBox);
67 return status;
68};
69
70
71/** Bounds the molecule in the box whose lengths are defined by vector \a *BoxLengths.
72 * \param *out output stream for debugging
73 */
74bool molecule::BoundInBox()
75{
76 bool status = true;
77 Box &domain = World::getInstance().getDomain();
78
79 // go through all atoms
80 atoms.transformNodes(boost::bind(&Box::WrapPeriodically,domain,_1));
81
82 return status;
83};
84
85/** Centers the edge of the atoms at (0,0,0).
86 * \param *out output stream for debugging
87 * \param *max coordinates of other edge, specifying box dimensions.
88 */
89void molecule::CenterEdge(Vector *max)
90{
91 Vector *min = new Vector;
92
93// Log() << Verbose(3) << "Begin of CenterEdge." << endl;
94 molecule::const_iterator iter = begin(); // start at first in list
95 if (iter != end()) { //list not empty?
96 for (int i=NDIM;i--;) {
97 max->at(i) = (*iter)->at(i);
98 min->at(i) = (*iter)->at(i);
99 }
100 for (; iter != end(); ++iter) {// continue with second if present
101 //(*iter)->Output(1,1,out);
102 for (int i=NDIM;i--;) {
103 max->at(i) = (max->at(i) < (*iter)->at(i)) ? (*iter)->at(i) : max->at(i);
104 min->at(i) = (min->at(i) > (*iter)->at(i)) ? (*iter)->at(i) : min->at(i);
105 }
106 }
107// Log() << Verbose(4) << "Maximum is ";
108// max->Output(out);
109// Log() << Verbose(0) << ", Minimum is ";
110// min->Output(out);
111// Log() << Verbose(0) << endl;
112 min->Scale(-1.);
113 (*max) += (*min);
114 Translate(min);
115 }
116 delete(min);
117// Log() << Verbose(3) << "End of CenterEdge." << endl;
118};
119
120/** Centers the center of the atoms at (0,0,0).
121 * \param *out output stream for debugging
122 * \param *center return vector for translation vector
123 */
124void molecule::CenterOrigin()
125{
126 int Num = 0;
127 molecule::const_iterator iter = begin(); // start at first in list
128 Vector Center;
129
130 Center.Zero();
131 if (iter != end()) { //list not empty?
132 for (; iter != end(); ++iter) { // continue with second if present
133 Num++;
134 Center += (*iter)->getPosition();
135 }
136 Center.Scale(-1./(double)Num); // divide through total number (and sign for direction)
137 Translate(&Center);
138 }
139};
140
141/** Returns vector pointing to center of all atoms.
142 * \return pointer to center of all vector
143 */
144Vector * molecule::DetermineCenterOfAll() const
145{
146 molecule::const_iterator iter = begin(); // start at first in list
147 Vector *a = new Vector();
148 double Num = 0;
149
150 a->Zero();
151
152 if (iter != end()) { //list not empty?
153 for (; iter != end(); ++iter) { // continue with second if present
154 Num++;
155 (*a) += (*iter)->getPosition();
156 }
157 a->Scale(1./(double)Num); // divide through total mass (and sign for direction)
158 }
159 return a;
160};
161
162/** Returns vector pointing to center of the domain.
163 * \return pointer to center of the domain
164 */
165Vector * molecule::DetermineCenterOfBox() const
166{
167 Vector *a = new Vector(0.5,0.5,0.5);
168 const RealSpaceMatrix &M = World::getInstance().getDomain().getM();
169 (*a) *= M;
170 return a;
171};
172
173/** Returns vector pointing to center of gravity.
174 * \param *out output stream for debugging
175 * \return pointer to center of gravity vector
176 */
177Vector * molecule::DetermineCenterOfGravity() const
178{
179 molecule::const_iterator iter = begin(); // start at first in list
180 Vector *a = new Vector();
181 Vector tmp;
182 double Num = 0;
183
184 a->Zero();
185
186 if (iter != end()) { //list not empty?
187 for (; iter != end(); ++iter) { // continue with second if present
188 Num += (*iter)->getType()->getMass();
189 tmp = (*iter)->getType()->getMass() * (*iter)->getPosition();
190 (*a) += tmp;
191 }
192 a->Scale(1./Num); // divide through total mass
193 }
194// Log() << Verbose(1) << "Resulting center of gravity: ";
195// a->Output(out);
196// Log() << Verbose(0) << endl;
197 return a;
198};
199
200/** Centers the center of gravity of the atoms at (0,0,0).
201 * \param *out output stream for debugging
202 * \param *center return vector for translation vector
203 */
204void molecule::CenterPeriodic()
205{
206 Vector NewCenter;
207 DeterminePeriodicCenter(NewCenter);
208 // go through all atoms
209 BOOST_FOREACH(atom* iter, atoms){
210 *iter -= NewCenter;
211 }
212};
213
214
215/** Centers the center of gravity of the atoms at (0,0,0).
216 * \param *out output stream for debugging
217 * \param *center return vector for translation vector
218 */
219void molecule::CenterAtVector(Vector *newcenter)
220{
221 // go through all atoms
222 BOOST_FOREACH(atom* iter, atoms){
223 *iter -= *newcenter;
224 }
225};
226
227/** Calculate the inertia tensor of a the molecule.
228 *
229 * @return inertia tensor
230 */
231RealSpaceMatrix molecule::getInertiaTensor() const
232{
233 RealSpaceMatrix InertiaTensor;
234 Vector *CenterOfGravity = DetermineCenterOfGravity();
235
236 // reset inertia tensor
237 InertiaTensor.setZero();
238
239 // sum up inertia tensor
240 for (molecule::const_iterator iter = begin(); iter != end(); ++iter) {
241 Vector x = (*iter)->getPosition();
242 x -= *CenterOfGravity;
243 const double mass = (*iter)->getType()->getMass();
244 InertiaTensor.at(0,0) += mass*(x[1]*x[1] + x[2]*x[2]);
245 InertiaTensor.at(0,1) += mass*(-x[0]*x[1]);
246 InertiaTensor.at(0,2) += mass*(-x[0]*x[2]);
247 InertiaTensor.at(1,0) += mass*(-x[1]*x[0]);
248 InertiaTensor.at(1,1) += mass*(x[0]*x[0] + x[2]*x[2]);
249 InertiaTensor.at(1,2) += mass*(-x[1]*x[2]);
250 InertiaTensor.at(2,0) += mass*(-x[2]*x[0]);
251 InertiaTensor.at(2,1) += mass*(-x[2]*x[1]);
252 InertiaTensor.at(2,2) += mass*(x[0]*x[0] + x[1]*x[1]);
253 }
254 // print InertiaTensor
255 DoLog(0) && (Log() << Verbose(0) << "The inertia tensor of molecule "
256 << getName() << " is:"
257 << InertiaTensor << endl);
258
259 delete CenterOfGravity;
260 return InertiaTensor;
261}
262
263/** Rotates the molecule in such a way that biggest principal axis corresponds
264 * to given \a Axis.
265 *
266 * @param Axis Axis to align with biggest principal axis
267 */
268void molecule::RotateToPrincipalAxisSystem(Vector &Axis)
269{
270 Vector *CenterOfGravity = DetermineCenterOfGravity();
271 RealSpaceMatrix InertiaTensor = getInertiaTensor();
272
273 // diagonalize to determine principal axis system
274 Vector Eigenvalues = InertiaTensor.transformToEigenbasis();
275
276 for(int i=0;i<NDIM;i++)
277 DoLog(0) && (Log() << Verbose(0) << "eigenvalue = " << Eigenvalues[i] << ", eigenvector = " << InertiaTensor.column(i) << endl);
278
279 DoLog(0) && (Log() << Verbose(0) << "Transforming to PAS ... ");
280
281 // obtain first column, eigenvector to biggest eigenvalue
282 Vector BiggestEigenvector(InertiaTensor.column(Eigenvalues.SmallestComponent()));
283 Vector DesiredAxis(Axis);
284
285 // Creation Line that is the rotation axis
286 DesiredAxis.VectorProduct(BiggestEigenvector);
287 Line RotationAxis(Vector(0.,0.,0.), DesiredAxis);
288
289 // determine angle
290 const double alpha = BiggestEigenvector.Angle(Axis);
291
292 DoLog(0) && (Log() << Verbose(0) << "Rotation angle is " << alpha << endl);
293
294 // and rotate
295 for (molecule::iterator iter = begin(); iter != end(); ++iter) {
296 *(*iter) -= *CenterOfGravity;
297 (*iter)->setPosition(RotationAxis.rotateVector((*iter)->getPosition(), alpha));
298 *(*iter) += *CenterOfGravity;
299 }
300 DoLog(0) && (Log() << Verbose(0) << "done." << endl);
301
302 delete CenterOfGravity;
303}
304
305/** Scales all atoms by \a *factor.
306 * \param *factor pointer to scaling factor
307 *
308 * TODO: Is this realy what is meant, i.e.
309 * x=(x[0]*factor[0],x[1]*factor[1],x[2]*factor[2]) (current impl)
310 * or rather
311 * x=(**factor) * x (as suggested by comment)
312 */
313void molecule::Scale(const double ** const factor)
314{
315 for (molecule::const_iterator iter = begin(); iter != end(); ++iter) {
316 for (int j=0;j<MDSteps;j++)
317 (*iter)->Trajectory.R.at(j).ScaleAll(*factor);
318 (*iter)->ScaleAll(*factor);
319 }
320};
321
322/** Translate all atoms by given vector.
323 * \param trans[] translation vector.
324 */
325void molecule::Translate(const Vector *trans)
326{
327 for (molecule::const_iterator iter = begin(); iter != end(); ++iter) {
328 for (int j=0;j<MDSteps;j++)
329 (*iter)->Trajectory.R.at(j) += (*trans);
330 *(*iter) += (*trans);
331 }
332};
333
334/** Translate the molecule periodically in the box.
335 * \param trans[] translation vector.
336 * TODO treatment of trajetories missing
337 */
338void molecule::TranslatePeriodically(const Vector *trans)
339{
340 Box &domain = World::getInstance().getDomain();
341
342 // go through all atoms
343 BOOST_FOREACH(atom* iter, atoms){
344 *iter += *trans;
345 }
346 atoms.transformNodes(boost::bind(&Box::WrapPeriodically,domain,_1));
347
348};
349
350
351/** Mirrors all atoms against a given plane.
352 * \param n[] normal vector of mirror plane.
353 */
354void molecule::Mirror(const Vector *n)
355{
356 OBSERVE;
357 Plane p(*n,0);
358 atoms.transformNodes(boost::bind(&Plane::mirrorVector,p,_1));
359};
360
361/** Determines center of molecule (yet not considering atom masses).
362 * \param center reference to return vector
363 */
364void molecule::DeterminePeriodicCenter(Vector &center)
365{
366 const RealSpaceMatrix &matrix = World::getInstance().getDomain().getM();
367 const RealSpaceMatrix &inversematrix = World::getInstance().getDomain().getM();
368 double tmp;
369 bool flag;
370 Vector Testvector, Translationvector;
371 Vector Center;
372
373 do {
374 Center.Zero();
375 flag = true;
376 for (molecule::const_iterator iter = begin(); iter != end(); ++iter) {
377#ifdef ADDHYDROGEN
378 if ((*iter)->getType()->getAtomicNumber() != 1) {
379#endif
380 Testvector = inversematrix * (*iter)->getPosition();
381 Translationvector.Zero();
382 for (BondList::const_iterator Runner = (*iter)->ListOfBonds.begin(); Runner != (*iter)->ListOfBonds.end(); (++Runner)) {
383 if ((*iter)->nr < (*Runner)->GetOtherAtom((*iter))->nr) // otherwise we shift one to, the other fro and gain nothing
384 for (int j=0;j<NDIM;j++) {
385 tmp = (*iter)->at(j) - (*Runner)->GetOtherAtom(*iter)->at(j);
386 if ((fabs(tmp)) > BondDistance) {
387 flag = false;
388 DoLog(0) && (Log() << Verbose(0) << "Hit: atom " << (*iter)->getName() << " in bond " << *(*Runner) << " has to be shifted due to " << tmp << "." << endl);
389 if (tmp > 0)
390 Translationvector[j] -= 1.;
391 else
392 Translationvector[j] += 1.;
393 }
394 }
395 }
396 Testvector += Translationvector;
397 Testvector *= matrix;
398 Center += Testvector;
399 Log() << Verbose(1) << "vector is: " << Testvector << endl;
400#ifdef ADDHYDROGEN
401 // now also change all hydrogens
402 for (BondList::const_iterator Runner = (*iter)->ListOfBonds.begin(); Runner != (*iter)->ListOfBonds.end(); (++Runner)) {
403 if ((*Runner)->GetOtherAtom((*iter))->getType()->getAtomicNumber() == 1) {
404 Testvector = inversematrix * (*Runner)->GetOtherAtom((*iter))->getPosition();
405 Testvector += Translationvector;
406 Testvector *= matrix;
407 Center += Testvector;
408 Log() << Verbose(1) << "Hydrogen vector is: " << Testvector << endl;
409 }
410 }
411 }
412#endif
413 }
414 } while (!flag);
415
416 Center.Scale(1./static_cast<double>(getAtomCount()));
417 CenterAtVector(&Center);
418};
419
420/** Align all atoms in such a manner that given vector \a *n is along z axis.
421 * \param n[] alignment vector.
422 */
423void molecule::Align(Vector *n)
424{
425 double alpha, tmp;
426 Vector z_axis;
427 z_axis[0] = 0.;
428 z_axis[1] = 0.;
429 z_axis[2] = 1.;
430
431 // rotate on z-x plane
432 DoLog(0) && (Log() << Verbose(0) << "Begin of Aligning all atoms." << endl);
433 alpha = atan(-n->at(0)/n->at(2));
434 DoLog(1) && (Log() << Verbose(1) << "Z-X-angle: " << alpha << " ... ");
435 for (molecule::const_iterator iter = begin(); iter != end(); ++iter) {
436 tmp = (*iter)->at(0);
437 (*iter)->set(0, cos(alpha) * tmp + sin(alpha) * (*iter)->at(2));
438 (*iter)->set(2, -sin(alpha) * tmp + cos(alpha) * (*iter)->at(2));
439 for (int j=0;j<MDSteps;j++) {
440 tmp = (*iter)->Trajectory.R.at(j)[0];
441 (*iter)->Trajectory.R.at(j)[0] = cos(alpha) * tmp + sin(alpha) * (*iter)->Trajectory.R.at(j)[2];
442 (*iter)->Trajectory.R.at(j)[2] = -sin(alpha) * tmp + cos(alpha) * (*iter)->Trajectory.R.at(j)[2];
443 }
444 }
445 // rotate n vector
446 tmp = n->at(0);
447 n->at(0) = cos(alpha) * tmp + sin(alpha) * n->at(2);
448 n->at(2) = -sin(alpha) * tmp + cos(alpha) * n->at(2);
449 DoLog(1) && (Log() << Verbose(1) << "alignment vector after first rotation: " << n << endl);
450
451 // rotate on z-y plane
452 alpha = atan(-n->at(1)/n->at(2));
453 DoLog(1) && (Log() << Verbose(1) << "Z-Y-angle: " << alpha << " ... ");
454 for (molecule::const_iterator iter = begin(); iter != end(); ++iter) {
455 tmp = (*iter)->at(1);
456 (*iter)->set(1, cos(alpha) * tmp + sin(alpha) * (*iter)->at(2));
457 (*iter)->set(2, -sin(alpha) * tmp + cos(alpha) * (*iter)->at(2));
458 for (int j=0;j<MDSteps;j++) {
459 tmp = (*iter)->Trajectory.R.at(j)[1];
460 (*iter)->Trajectory.R.at(j)[1] = cos(alpha) * tmp + sin(alpha) * (*iter)->Trajectory.R.at(j)[2];
461 (*iter)->Trajectory.R.at(j)[2] = -sin(alpha) * tmp + cos(alpha) * (*iter)->Trajectory.R.at(j)[2];
462 }
463 }
464 // rotate n vector (for consistency check)
465 tmp = n->at(1);
466 n->at(1) = cos(alpha) * tmp + sin(alpha) * n->at(2);
467 n->at(2) = -sin(alpha) * tmp + cos(alpha) * n->at(2);
468
469
470 DoLog(1) && (Log() << Verbose(1) << "alignment vector after second rotation: " << n << endl);
471 DoLog(0) && (Log() << Verbose(0) << "End of Aligning all atoms." << endl);
472};
473
474
475/** Calculates sum over least square distance to line hidden in \a *x.
476 * \param *x offset and direction vector
477 * \param *params pointer to lsq_params structure
478 * \return \f$ sum_i^N | y_i - (a + t_i b)|^2\f$
479 */
480double LeastSquareDistance (const gsl_vector * x, void * params)
481{
482 double res = 0, t;
483 Vector a,b,c,d;
484 struct lsq_params *par = (struct lsq_params *)params;
485
486 // initialize vectors
487 a[0] = gsl_vector_get(x,0);
488 a[1] = gsl_vector_get(x,1);
489 a[2] = gsl_vector_get(x,2);
490 b[0] = gsl_vector_get(x,3);
491 b[1] = gsl_vector_get(x,4);
492 b[2] = gsl_vector_get(x,5);
493 // go through all atoms
494 for (molecule::const_iterator iter = par->mol->begin(); iter != par->mol->end(); ++iter) {
495 if ((*iter)->getType() == ((struct lsq_params *)params)->type) { // for specific type
496 c = (*iter)->getPosition() - a;
497 t = c.ScalarProduct(b); // get direction parameter
498 d = t*b; // and create vector
499 c -= d; // ... yielding distance vector
500 res += d.ScalarProduct(d); // add squared distance
501 }
502 }
503 return res;
504};
505
506/** By minimizing the least square distance gains alignment vector.
507 * \bug this is not yet working properly it seems
508 */
509void molecule::GetAlignvector(struct lsq_params * par) const
510{
511 int np = 6;
512
513 const gsl_multimin_fminimizer_type *T =
514 gsl_multimin_fminimizer_nmsimplex;
515 gsl_multimin_fminimizer *s = NULL;
516 gsl_vector *ss;
517 gsl_multimin_function minex_func;
518
519 size_t iter = 0, i;
520 int status;
521 double size;
522
523 /* Initial vertex size vector */
524 ss = gsl_vector_alloc (np);
525
526 /* Set all step sizes to 1 */
527 gsl_vector_set_all (ss, 1.0);
528
529 /* Starting point */
530 par->x = gsl_vector_alloc (np);
531 par->mol = this;
532
533 gsl_vector_set (par->x, 0, 0.0); // offset
534 gsl_vector_set (par->x, 1, 0.0);
535 gsl_vector_set (par->x, 2, 0.0);
536 gsl_vector_set (par->x, 3, 0.0); // direction
537 gsl_vector_set (par->x, 4, 0.0);
538 gsl_vector_set (par->x, 5, 1.0);
539
540 /* Initialize method and iterate */
541 minex_func.f = &LeastSquareDistance;
542 minex_func.n = np;
543 minex_func.params = (void *)par;
544
545 s = gsl_multimin_fminimizer_alloc (T, np);
546 gsl_multimin_fminimizer_set (s, &minex_func, par->x, ss);
547
548 do
549 {
550 iter++;
551 status = gsl_multimin_fminimizer_iterate(s);
552
553 if (status)
554 break;
555
556 size = gsl_multimin_fminimizer_size (s);
557 status = gsl_multimin_test_size (size, 1e-2);
558
559 if (status == GSL_SUCCESS)
560 {
561 printf ("converged to minimum at\n");
562 }
563
564 printf ("%5d ", (int)iter);
565 for (i = 0; i < (size_t)np; i++)
566 {
567 printf ("%10.3e ", gsl_vector_get (s->x, i));
568 }
569 printf ("f() = %7.3f size = %.3f\n", s->fval, size);
570 }
571 while (status == GSL_CONTINUE && iter < 100);
572
573 for (i=0;i<(size_t)np;i++)
574 gsl_vector_set(par->x, i, gsl_vector_get(s->x, i));
575 //gsl_vector_free(par->x);
576 gsl_vector_free(ss);
577 gsl_multimin_fminimizer_free (s);
578};
Note: See TracBrowser for help on using the repository browser.