source: src/molecule_geometry.cpp@ 2d701e

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 2d701e was 2d701e, checked in by Frederik Heber <heber@…>, 10 years ago

FIX: molecule::Translate() used in CenterInBox(), CenterEdge(), CenterPeriodic(), and CenterAtVector().

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