1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857
|
/*
Project: Adun
Copyright (C) 2007 Michael Johnston & Jordi Villa-Freixa
Author: Michael Johnston
This application is free software; you can redistribute it and/or
modify it under the terms of the GNU General Public
License as published by the Free Software Foundation; either
version 2 of the License, or (at your option) any later version.
This application is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
Library General Public License for more details.
You should have received a copy of the GNU General Public
License along with this library; if not, write to the Free
Software Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111 USA.
*/
#include "AdunKernel/AdunMinimiser.h"
#include <gsl/gsl_blas.h>
static NSArray* allowedAlgorithms;
/**
Methods that wrap calls to framework object so
pointers to the methods can be used with the GSL minimisation scheme
*/
@interface AdMinimiser (GSLWrappers)
/**
Calculates the energy of the systems in the system collection
*/
- (double) energyOfConfiguration: (const gsl_vector*) configuration usingParameters: (void*) params;
/**
Updates the configuration of the systems in system collection using \e configuration.
The sections of \e configuration corresponding to each system are given by rangeForSystem:()
Calculates the gradient (negative of the force) acting on the elements of the systems
and writes it to \e gradient
*/
- (void) gradientOfConfiguration: (const gsl_vector*) configuration
usingParameters: (void*) params
gradientBuffer: (gsl_vector*) gradient;
/**
Updates the configuration of the systems in system collection using \e configuration.
The sections of \e configuration corresponding to each system are given by rangeForSystem:().
Calculates the gradient (negative of the force) acting on the elements of the systems and the energy.
Stores the energy in \e energyBuffer and the gradient in \e gradient
*/
- (void) energyAndGradientOfConfiguration: (const gsl_vector*) configuration
usingParameters: (void*) params
energyBuffer: (double*) energyBuffer
gradientBuffer: (gsl_vector*) gradient;
/**
Initialises the minimistiation function structure used by the gsl functions.
*/
- (void) initGSLMinimisationFunctionStruct;
@end
@implementation AdMinimiser
+ (void) initialize
{
allowedAlgorithms = [NSArray arrayWithObjects:
@"SteepestDescent",
@"FletcherReevesCG",
@"PolakRibiereCG",
@"BFGS",
nil];
[allowedAlgorithms retain];
}
- (void) _storeConfigurationInVector: (gsl_vector*) vector
{
int i, j, k;
AdMatrix* coordinates;
NSRange range;
NSEnumerator* systemEnum;
NSValue* pointer;
id system;
systemEnum = [systemRanges keyEnumerator];
while((pointer = [systemEnum nextObject]))
{
system = (id)[pointer pointerValue];
range = [[systemRanges objectForKey: pointer] rangeValue];
coordinates = [system coordinates];
k = range.location*3;
for(i=0; i<coordinates->no_rows; i++)
for(j=0; j<3; j++, k++)
gsl_vector_set(vector, k, coordinates->matrix[i][j]);
}
}
- (void) emptyPool
{
[pool release];
pool = [[NSAutoreleasePool alloc] init];
AdLogMemoryUsage();
}
- (void) endProduction
{
endSimulation = YES;
}
- (id) init
{
return [self initWithForceFields: nil];
}
- (id) initWithForceFields: (AdForceFieldCollection*) aForceFieldCollection
{
return [self initWithSystems: nil
forceFields: aForceFieldCollection];
}
- (id) initWithSystems: (AdSystemCollection*) aSystemCollection
forceFields: (AdForceFieldCollection*) aForceFieldCollection
{
return [self initWithSystems: aSystemCollection
forceFields: aForceFieldCollection
absoluteTolerance: 0.005
numberOfSteps: 1E2];
}
- (id) initWithSystems: (AdSystemCollection*) aSystemCollection
forceFields: (AdForceFieldCollection*) aForceFieldCollection
absoluteTolerance: (double) absTol
numberOfSteps: (unsigned int) steps
{
return [self initWithSystems: aSystemCollection
forceFields: aForceFieldCollection
absoluteTolerance: absTol
numberOfSteps: steps
algorithm: @"SteepestDescent"
stepSize: 0.1
tolerance: 0.1];
}
- (id) initWithSystems: (AdSystemCollection*) aSystemCollection
forceFields: (AdForceFieldCollection*) aForceFieldCollection
absoluteTolerance: (double) absTol
numberOfSteps: (unsigned int) steps
algorithm: (NSString*) algorithm
stepSize: (double) stepsize
tolerance: (double) tol
{
return [self initWithSystems: aSystemCollection
forceFields: aForceFieldCollection
absoluteTolerance: absTol
numberOfSteps: steps
algorithm: algorithm
stepSize: stepsize
tolerance: tol
checkFPErrorInterval: 1000];
}
- (id) initWithSystems: (AdSystemCollection*) aSystemCollection
forceFields: (AdForceFieldCollection*) aForceFieldCollection
absoluteTolerance: (double) absTol
numberOfSteps: (unsigned int) steps
algorithm: (NSString*) aString
stepSize: (double) step
tolerance: (double) tol
checkFPErrorInterval: (unsigned int) interval
{
if((self = [super init]))
{
endSimulation = NO;
checkFPErrorInterval = interval;
currentStep = 0;
absoluteTolerance = fabs(absTol);
tolerance = fabs(tol);
stepSize = step;
numberOfSteps = steps;
systemRanges = [NSMutableDictionary new];
currentEnergy = gradientNorm = 0;
converged = NO;
constraints = [NSMutableDictionary new];
NSDebugLLog(@"AdSimulator", @"The maximium number of iterations is %d", numberOfSteps);
if(aString == nil)
aString = @"SteepestDescent";
[self setAlgorithm: aString];
[self setForceFields: aForceFieldCollection];
[self setSystems: aSystemCollection];
timer = [AdMainLoopTimer mainLoopTimer];
}
return self;
}
- (void) dealloc
{
gsl_matrix* matrix;
NSEnumerator* keyEnum;
NSValue* key, *value;
//Free any constraint matrices
keyEnum = [constraints keyEnumerator];
while((key = [keyEnum nextObject]))
{
value = [constraints objectForKey: key];
matrix = [value pointerValue];
gsl_matrix_free(matrix);
}
[constraints release];
[systemCollection release];
[systemRanges release];
[forceFieldCollection release];
[algorithmName release];
[fullSystems release];
[super dealloc];
}
- (BOOL) production: (NSError**) error
{
int status;
struct tms start;
struct tms end;
gsl_multimin_fdfminimizer *minimiser;
gsl_vector* initialConfiguration;
converged = NO;
gradientNorm = 0;
//Allocate and set the minimiser
minimiser = gsl_multimin_fdfminimizer_alloc(algorithmType, numberOfElements*3);
[self initGSLMinimisationFunctionStruct];
if(minimiser == NULL)
{
NSWarnLog(@"Error while initialising minimiser");
return NO;
}
initialConfiguration = gsl_vector_calloc(numberOfElements*3);
[self _storeConfigurationInVector: initialConfiguration];
currentEnergy = [self energyOfConfiguration: initialConfiguration
usingParameters: NULL];
NSDebugLLog(@"AdMinimiser",
@"Initial energy %lf\n", currentEnergy);
GSPrintf(stdout,
@"\nBeginning minimisation using %@ algorithm - Max steps %d\n",
algorithmName, numberOfSteps);
GSPrintf(stdout,
@"Initial energy %lf\n", currentEnergy);
gsl_multimin_fdfminimizer_set(minimiser,
&minimistationFunctions,
initialConfiguration,
stepSize,
tolerance);
//Send a notification
[[NSNotificationCenter defaultCenter]
postNotificationName: @"AdConfigurationGeneratorWillBeginProductionNotification"
object: self];
//Setup timers
[timer sendMessage: @selector(emptyPool)
toObject: self
interval: 100
name: @"Autorelease"];
[timer sendMessage: @selector(checkFloatingPointErrors)
toObject: self
interval: checkFPErrorInterval
name: @"FloatingPointErrors"];
//Minimsation
times(&start);
NSDebugLLog(@"AdMinimiser", @"Beginning minimisation");
pool = [[NSAutoreleasePool alloc] init];
for(currentStep=0; currentStep < numberOfSteps; currentStep++)
{
status = gsl_multimin_fdfminimizer_iterate(minimiser);
currentEnergy = gsl_multimin_fdfminimizer_minimum(minimiser);
gradientNorm = gsl_blas_dnrm2(gsl_multimin_fdfminimizer_gradient(minimiser));
[timer increment];
NSDebugLLog(@"SimulationLoop",
@"Finished minimisation - step %d",
currentStep);
NSDebugLLog(@"SimulationLoop",
@"Current energy %lf. Gradient norm %lf",
currentEnergy,
gradientNorm);
if(status)
{
GSPrintf(stdout, @"Error during iteration\n");
GSPrintf(stdout, @"Description - %s\n", gsl_strerror(status));
NSWarnLog(@"Error during iteration");
NSWarnLog(@"Description - %s", gsl_strerror(status));
break;
}
status = gsl_multimin_test_gradient(
gsl_multimin_fdfminimizer_gradient(minimiser),
absoluteTolerance);
if(status == GSL_SUCCESS)
{
GSPrintf(stdout, @"Reached covergence\n");
converged = YES;
break;
}
if(endSimulation)
{
GSPrintf(stdout, @"Exiting on user request\n");
break;
}
}
times(&end);
GSPrintf(stdout, @"Minimistation ended - step %d\n",
currentStep);
if(currentStep == 0)
currentStep = 1;
GSPrintf(stdout,
@"Final energy %lf. Gradient norm %lf. Target accuracy %lf\n\n",
currentEnergy,
gradientNorm,
absoluteTolerance);
NSDebugLLog(@"SimulationLoop",
@"Minimistation complete");
AdLogTimingInformation(&start, &end, currentStep);
fflush(stdout);
[timer removeMessageWithName: @"FloatingPointErrors"];
[timer removeMessageWithName: @"Autorelease"];
gsl_multimin_fdfminimizer_free(minimiser);
gsl_vector_free(initialConfiguration);
minimiser = NULL;
initialConfiguration = NULL;
[pool release];
return YES;
}
- (BOOL) restartFrom: (int) step error: (NSError**) error
{
NSWarnLog(@"AdMinimiser instances do not support restarting");
return NO;
}
/*
* Object Accessors
*/
- (AdSystemCollection*) systems
{
return [[systemCollection retain] autorelease];
}
- (void) setSystems: (AdSystemCollection*) aCollection
{
int start, end;
NSRange range;
NSEnumerator* systemEnum;
id system;
//We need to assign ranges to the systems
[systemCollection release];
systemCollection = [aCollection retain];
[fullSystems release];
fullSystems = [systemCollection fullSystems];
[fullSystems retain];
[systemRanges removeAllObjects];
numberOfElements= start = end = 0;
systemEnum = [[systemCollection fullSystems] objectEnumerator];
while((system = [systemEnum nextObject]))
{
numberOfElements += [system numberOfElements];
start = end;
end = start + [system numberOfElements];
range = NSMakeRange(start, [system numberOfElements]);
[systemRanges setObject: [NSValue valueWithRange: range]
forKey: [NSValue valueWithPointer: system]];
}
NSDebugLLog(@"AdMinimiser",
@"There are %d elements in total", numberOfElements);
NSDebugLLog(@"AdMinimiser",
@"Systems and ranges %@", systemRanges);
}
- (AdForceFieldCollection*) forceFields
{
return [[forceFieldCollection retain] autorelease];
}
- (void) setForceFields: (AdForceFieldCollection*) object
{
[forceFieldCollection release];
forceFieldCollection = [object retain];
}
- (void) setNumberOfSteps: (unsigned int) anInt
{
numberOfSteps = anInt;
}
- (unsigned int) numberOfSteps
{
return numberOfSteps;
}
- (double) stepSize
{
return stepSize;
}
- (void) setStepSize: (double) aDouble
{
stepSize = aDouble;
}
- (unsigned int) currentStep
{
return currentStep;
}
- (void) setAlgorithm: (NSString*) aString
{
if(![allowedAlgorithms containsObject: aString])
[NSException raise: NSInvalidArgumentException
format: @"Unknown minimistation algorithm %@", aString];
[algorithmName release];
algorithmName = [aString retain];
if([aString isEqual: @"FletcherReevesCG"])
algorithmType = gsl_multimin_fdfminimizer_conjugate_fr;
else if([aString isEqual: @"PolakRibiereCG"])
algorithmType = gsl_multimin_fdfminimizer_conjugate_pr;
else if([aString isEqual: @"SteepestDescent"])
algorithmType = gsl_multimin_fdfminimizer_steepest_descent;
else if([aString isEqual: @"BFGS"])
algorithmType = gsl_multimin_fdfminimizer_vector_bfgs;
}
- (NSString*) algorithm
{
return [[algorithmName retain] autorelease];
}
- (double) tolerance
{
return tolerance;
}
- (void) setTolerance: (double) aDouble
{
tolerance = fabs(aDouble);
}
- (double) absoluteTolerance
{
return absoluteTolerance;
}
- (void) setAbsoluteTolerance: (double) aDouble;
{
absoluteTolerance = fabs(aDouble);
}
- (double) currentEnergy
{
return currentEnergy;
}
- (BOOL) converged
{
return converged;
}
- (double) gradientNorm
{
return gradientNorm;
}
- (void) setConstraints: (AdDataMatrix*) aMatrix forSystem: (AdSystem*) system
{
gsl_matrix* constraintMatrix;
//Check system exists
if(![[systemCollection fullSystems] containsObject: system])
[NSException raise: NSInvalidArgumentException
format: @"System %@ is not one of the full systems known by this minimiser object",
system];
//Check they are the right dimension
if((int)[aMatrix numberOfRows] != 3*numberOfElements)
[NSException raise: NSInvalidArgumentException
format: @"Constraint matrix has incorrect number of rows (%@, Expected - %d)"
,[aMatrix numberOfRows], numberOfElements*3];
//Remove any previous constraints set for system (if any)
[self removeConstraintsForSystem: system];
constraintMatrix = AdGSLMatrixFromAdDataMatrix(aMatrix);
[constraints setObject: [NSValue valueWithPointer: constraintMatrix]
forKey: [NSValue valueWithPointer: system]];
}
- (AdDataMatrix*) constraintsForSystem: (AdSystem*) system
{
gsl_matrix* constraintMatrix;
NSValue* value;
value = [constraints objectForKey: [NSValue valueWithPointer: system]];
//Check if any constraints were set
if(value == nil)
return nil;
constraintMatrix = [value pointerValue];
return [AdDataMatrix matrixFromGSLMatrix: constraintMatrix];
}
- (void) removeConstraintsForSystem: (AdSystem*) system
{
gsl_matrix* constraintMatrix;
NSValue* value;
value = [constraints objectForKey: [NSValue valueWithPointer: system]];
//Check if any constraints were set
if(value != nil)
{
[constraints removeObjectForKey: [NSValue valueWithPointer: system]];
constraintMatrix = [value pointerValue];
gsl_matrix_free(constraintMatrix);
}
}
@end
//The GSL functions wrapping the methods in the category below.
double energyFunction(const gsl_vector* configuration, void* params);
void gradientFunction(const gsl_vector* configuration, void* params, gsl_vector* gradient);
void energyAndGradientFunction(const gsl_vector* configuration,
void* params,
double* energyBuffer,
gsl_vector* gradient);
double energyFunction(const gsl_vector* configuration, void* params)
{
return [(AdMinimiser*)params energyOfConfiguration: configuration
usingParameters: params];
}
void gradientFunction(const gsl_vector* configuration, void* params, gsl_vector* gradient)
{
[(AdMinimiser*)params gradientOfConfiguration: configuration
usingParameters: params
gradientBuffer: gradient];
}
void energyAndGradientFunction(const gsl_vector* configuration,
void* params, double*
energyBuffer,
gsl_vector* gradient)
{
[(AdMinimiser*)params energyAndGradientOfConfiguration: configuration
usingParameters: params
energyBuffer: energyBuffer
gradientBuffer: gradient];
}
@implementation AdMinimiser (GSLWrappers)
- (void) _updateSystemsWithConfiguration: (const gsl_vector*) configuration
{
int i,j,k;
double movement = 0;
AdMatrix* coordinates;
NSRange range;
NSEnumerator* systemEnum;
id system;
//Set the new configuration
systemEnum = [[systemCollection fullSystems] objectEnumerator];
while((system = [systemEnum nextObject]))
{
coordinates = [system coordinates];
[system object: self willBeginWritingToMatrix: coordinates];
range = [[systemRanges objectForKey: [NSValue valueWithPointer: system]]
rangeValue];
for(i=0, k = range.location*3; i < coordinates->no_rows; i++)
for(j = 0; j < 3; j++, k++)
{ movement += pow(coordinates->matrix[i][j] - gsl_vector_get(configuration, k), 2);
coordinates->matrix[i][j] = gsl_vector_get(configuration, k);
}
[system object: self didFinishWritingToMatrix: coordinates];
}
}
/**
Removes the projection along a set of constraint vectors
from \e vector.
The constraint vectors are retreived from the contraints ivar which is a dictionary.
If the dictionary contains constraints for \e system this method does nothing.
*/
- (void) _applyConstraintsToForceVector: (gsl_vector*) vector forSystem: (AdSystem*) system
{
unsigned int i;
double projection;
gsl_vector* constraintVector;
gsl_matrix* constraintMatrix;
gsl_vector_view view;
NSValue* value;
value = [constraints objectForKey: [NSValue valueWithPointer: system]];
if(value != nil)
{
constraintMatrix = [value pointerValue];
//A vector structure to hold the vectors
constraintVector = gsl_vector_alloc(vector->size);
//Iterate over the constraint vectors
for(i=0; i<constraintMatrix->size2; i++)
{
//Copy the constraint vector into the vector struct
view = gsl_matrix_column(constraintMatrix, i);
gsl_vector_memcpy(constraintVector, &view.vector);
//Find the projection
gsl_blas_ddot(vector, constraintVector, &projection);
//Remove the projection along the constraint vector
gsl_vector_scale(constraintVector, -projection);
gsl_vector_sub(vector, constraintVector);
}
gsl_vector_free(constraintVector);
}
}
- (double) energyOfConfiguration: (const gsl_vector*) configuration usingParameters: (void*) params
{
double energy = 0;
NSEnumerator * forceFieldsEnum, *systemEnum;
NSArray* forceFields;
id forceField, system;
NSMutableSet* forceFieldSet;
[self _updateSystemsWithConfiguration: configuration];
forceFieldSet = [NSMutableSet set];
[forceFieldCollection evaluateEnergies];
systemEnum = [[systemCollection fullSystems] objectEnumerator];
while((system = [systemEnum nextObject]))
{
forceFields = [forceFieldCollection forceFieldsForSystem: system
activityFlag: AdActiveForceFields];
[forceFieldSet addObjectsFromArray: forceFields];
}
forceFieldsEnum = [forceFieldSet objectEnumerator];
while((forceField = [forceFieldsEnum nextObject]))
energy += [forceField totalEnergy];
return energy;
}
- (void) gradientOfConfiguration: (const gsl_vector*) configuration
usingParameters: (void*) params
gradientBuffer: (gsl_vector*) gradient
{
double energy = 0;
[self energyAndGradientOfConfiguration: configuration
usingParameters: params
energyBuffer: &energy
gradientBuffer: gradient];
}
/**
Convienience function called from energyAndGradientOfConfiguration:usingParameters:energyBuffer:gradientBuffer:
It adds the gradient information in vector, which corresponds to the atoms in \e aSystem
to the \e gradient vector.
It applies any constraints to the force-vector first.
See the main method for more.
*/
- (void) _updateGradient: (gsl_vector*) gradient
ofSystem: (AdSystem*) aSystem
withForces: (gsl_vector*) vector
{
unsigned int i,k;
double value;
NSRange gradientRange;
[self _applyConstraintsToForceVector: vector
forSystem: aSystem];
//Now update the gradient vector with the force vector
//
//systemRange will tell us the range assigned to the
//atoms from interactingSystem.
//For example, there are 10000 total atoms, and the
//ones from 560 to 6780 are part of this system
//The elements in the gradient vector corresponding to
//this system are then 560*3 - 6780*3. (3DOF's per atoms)
gradientRange = [[systemRanges objectForKey:
[NSValue valueWithPointer: aSystem]]
rangeValue];
for(i=0, k = gradientRange.location*3; i<vector->size; i++, k++)
{
value = gsl_vector_get(gradient, k);
//Force is the negative of the gradient
value -= gsl_vector_get(vector, i);
gsl_vector_set(gradient, k, value);
}
}
- (void) energyAndGradientOfConfiguration: (const gsl_vector*) configuration
usingParameters: (void*) params
energyBuffer: (double*) energyBuffer
gradientBuffer: (gsl_vector*) gradient
{
unsigned int i;
int j, k;
double energy;
AdMatrix* forces;
gsl_vector* myGradient, *forceVector;
NSRange interactionRange;
NSEnumerator * forceFieldsEnum, *systemEnum, *interactingSystemsEnum;
NSArray* forceFields;
NSMutableSet* forceFieldSet;
id forceField, system, interactingSystem;
/*
* The main points here
* 1) Adun gets the forces on a system by system basis
* 2) They are returned as matrices
* So if there are three systems there will be three matrices.
* 3) The gsl functions require one vector containing the total
* gradient on all atoms - this is the gradient parameter
* passed above
*
* This means we have to do two things
* 1) Change the matrices into vectors
* 2) For each system update the correct elements in the gradientVector
* To do this we use the systemRange ivar which tells us which elements
* in the gradient vector correspond to the atoms of each system.
*/
//We cant be sure whats in the gradient buffer already.
//Since we need to add directly to the gradient elements
//we work on our own vector and copy the results to the buffer
//at the end
myGradient = gsl_vector_calloc(numberOfElements*3);
[self _updateSystemsWithConfiguration: configuration];
forceFieldSet = [NSMutableSet set];
[forceFieldCollection evaluateForces];
systemEnum = [[systemCollection fullSystems] objectEnumerator];
while((system = [systemEnum nextObject]))
{
forceFields = [forceFieldCollection forceFieldsForSystem: system
activityFlag: AdActiveForceFields];
[forceFieldSet addObjectsFromArray: forceFields];
}
energy = 0;
forceFieldsEnum = [forceFieldSet objectEnumerator];
while((forceField = [forceFieldsEnum nextObject]))
{
forces = [forceField forces];
system = [forceField system];
//We have to differentiate between forces from
//a normal system and forces from an interaction system
if([system isKindOfClass: [AdInteractionSystem class]])
{
//Only use the parts of the interaction system forces that
//correspond to full systems in the system collection
//This usually would be both but there could be e.g. a dipole
//system whose elements can't be moved.
interactingSystemsEnum = [[(AdInteractionSystem*)system systems] objectEnumerator];
while((interactingSystem = [interactingSystemsEnum nextObject]))
{
if([fullSystems containsObject: interactingSystem])
{
//Get the range for this system in the joint
//interaction system matrices
interactionRange = [system rangeForSystem: interactingSystem];
//Create a vector from the part of the force-matrix
//corresponding to interactingSystem
forceVector = gsl_vector_alloc(interactionRange.length*3);
for(i=interactionRange.location, k = 0; i < NSMaxRange(interactionRange); i++)
for(j=0;j<3; j++, k++)
gsl_vector_set(forceVector, k, forces->matrix[i][j]);
//This will also apply any constraints etc.
[self _updateGradient: myGradient
ofSystem: interactingSystem
withForces: forceVector];
gsl_vector_free(forceVector);
}
}
}
else
{
//Create a force-vector for the system
forceVector = gsl_vector_alloc(forces->no_rows*3);
for(i=0, k = 0; i < (unsigned int)forces->no_rows; i++)
for(j=0;j<3; j++, k++)
gsl_vector_set(forceVector, k, forces->matrix[i][j]);
[self _updateGradient: myGradient
ofSystem: system
withForces: forceVector];
gsl_vector_free(forceVector);
}
energy += [forceField totalEnergy];
}
gsl_vector_memcpy(gradient, myGradient);
gsl_vector_free(myGradient);
*energyBuffer = energy;
}
- (void) initGSLMinimisationFunctionStruct
{
minimistationFunctions.f = (double (*)(const gsl_vector*, void*))energyFunction;
minimistationFunctions.df = (void (*)(const gsl_vector*, void*, gsl_vector*))gradientFunction;
minimistationFunctions.fdf = (void (*)(const gsl_vector*, void*, double*, gsl_vector*))energyAndGradientFunction;
minimistationFunctions.n = numberOfElements*3;
//Set the param to be ourself to the functions can call the above methods
minimistationFunctions.params = (void*)self;
}
@end
|