source: src/molecule_geometry.cpp@ 1259df

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 1259df was 8cc22f, checked in by Frederik Heber <heber@…>, 10 years ago

Changed how trajectories are stored, not as vecor but as map.

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