| 1 | /*
 | 
|---|
| 2 |  * Langevin.cpp
 | 
|---|
| 3 |  *
 | 
|---|
| 4 |  *  Created on: Aug 20, 2010
 | 
|---|
| 5 |  *      Author: crueger
 | 
|---|
| 6 |  */
 | 
|---|
| 7 | 
 | 
|---|
| 8 | // include config.h
 | 
|---|
| 9 | #ifdef HAVE_CONFIG_H
 | 
|---|
| 10 | #include <config.h>
 | 
|---|
| 11 | #endif
 | 
|---|
| 12 | 
 | 
|---|
| 13 | #include "Helpers/MemDebug.hpp"
 | 
|---|
| 14 | 
 | 
|---|
| 15 | #include "Langevin.hpp"
 | 
|---|
| 16 | #include "element.hpp"
 | 
|---|
| 17 | #include "config.hpp"
 | 
|---|
| 18 | #include "Helpers/Verbose.hpp"
 | 
|---|
| 19 | #include "Helpers/Log.hpp"
 | 
|---|
| 20 | #include "ThermoStatContainer.hpp"
 | 
|---|
| 21 | 
 | 
|---|
| 22 | Langevin::Langevin(double _TempFrequency,double _alpha) :
 | 
|---|
| 23 |   TempFrequency(_TempFrequency),
 | 
|---|
| 24 |   alpha(_alpha)
 | 
|---|
| 25 | {
 | 
|---|
| 26 |   gsl_rng_env_setup();
 | 
|---|
| 27 |   T = gsl_rng_default;
 | 
|---|
| 28 |   r = gsl_rng_alloc (T);
 | 
|---|
| 29 | }
 | 
|---|
| 30 | 
 | 
|---|
| 31 | Langevin::Langevin() :
 | 
|---|
| 32 |   TempFrequency(2.5),
 | 
|---|
| 33 |   alpha(0.)
 | 
|---|
| 34 | {}
 | 
|---|
| 35 | 
 | 
|---|
| 36 | Langevin::~Langevin()
 | 
|---|
| 37 | {
 | 
|---|
| 38 |   gsl_rng_free (r);
 | 
|---|
| 39 | }
 | 
|---|
| 40 | 
 | 
|---|
| 41 | const char *ThermostatTraits<Langevin>::name = "Langevin";
 | 
|---|
| 42 | 
 | 
|---|
| 43 | std::string ThermostatTraits<Langevin>::getName(){
 | 
|---|
| 44 |   return ThermostatTraits<Langevin>::name;
 | 
|---|
| 45 | }
 | 
|---|
| 46 | 
 | 
|---|
| 47 | Thermostat *ThermostatTraits<Langevin>::make(class ConfigFileBuffer * const fb){
 | 
|---|
| 48 |   double TempFrequency;
 | 
|---|
| 49 |   double alpha;
 | 
|---|
| 50 |   const int verbose = 0;
 | 
|---|
| 51 |   ParseForParameter(verbose,fb,"Thermostat", 0, 2, 1, double_type, &TempFrequency, 1, critical); // read gamma
 | 
|---|
| 52 |   if (ParseForParameter(verbose,fb,"Thermostat", 0, 3, 1, double_type, &alpha, 1, optional)) {
 | 
|---|
| 53 |     DoLog(2) && (Log() << Verbose(2) << "Extended Stochastic Thermostat detected with interpolation coefficient " << alpha << "." << endl);
 | 
|---|
| 54 |   } else {
 | 
|---|
| 55 |     alpha = 1.;
 | 
|---|
| 56 |   }
 | 
|---|
| 57 |   return new Langevin(TempFrequency,alpha);
 | 
|---|
| 58 | }
 | 
|---|
| 59 | 
 | 
|---|
| 60 | double Langevin::scaleAtoms(unsigned int step,double ActualTemp,ATOMSET(std::list) atoms){
 | 
|---|
| 61 |   return doScaleAtoms(step,ActualTemp,atoms.begin(),atoms.end());
 | 
|---|
| 62 | }
 | 
|---|
| 63 | 
 | 
|---|
| 64 | double Langevin::scaleAtoms(unsigned int step,double ActualTemp,ATOMSET(std::vector) atoms){
 | 
|---|
| 65 |   return doScaleAtoms(step,ActualTemp,atoms.begin(),atoms.end());
 | 
|---|
| 66 | }
 | 
|---|
| 67 | 
 | 
|---|
| 68 | double Langevin::scaleAtoms(unsigned int step,double ActualTemp,ATOMSET(std::set) atoms){
 | 
|---|
| 69 |   return doScaleAtoms(step,ActualTemp,atoms.begin(),atoms.end());
 | 
|---|
| 70 | }
 | 
|---|
| 71 | 
 | 
|---|
| 72 | template <class ForwardIterator>
 | 
|---|
| 73 | double Langevin::doScaleAtoms(unsigned int step,double ActualTemp,ForwardIterator begin, ForwardIterator end){
 | 
|---|
| 74 |   DoLog(2) && (Log() << Verbose(2) <<  "Applying Langevin thermostat..." << endl);
 | 
|---|
| 75 |   double ekin=0;
 | 
|---|
| 76 |   for(ForwardIterator iter=begin;iter!=end;++iter){
 | 
|---|
| 77 |     double sigma  = sqrt(getContainer().TargetTemp/(*iter)->getType()->getMass()); // sigma = (k_b T)/m (Hartree/atomicmass = atomiclength/atomictime)
 | 
|---|
| 78 |     Vector &U = (*iter)->Trajectory.U.at(step);
 | 
|---|
| 79 |     if ((*iter)->FixedIon == 0) { // even FixedIon moves, only not by other's forces
 | 
|---|
| 80 |       // throw a dice to determine whether it gets hit by a heat bath particle
 | 
|---|
| 81 |       if (((((rand()/(double)RAND_MAX))*TempFrequency) < 1.)) {
 | 
|---|
| 82 |         DoLog(3) && (Log() << Verbose(3) << "Particle " << (**iter) << " was hit (sigma " << sigma << "): " << U.Norm() << " -> ");
 | 
|---|
| 83 |         // pick three random numbers from a Boltzmann distribution around the desired temperature T for each momenta axis
 | 
|---|
| 84 |         for (int d=0; d<NDIM; d++) {
 | 
|---|
| 85 |           U[d] = gsl_ran_gaussian (r, sigma);
 | 
|---|
| 86 |         }
 | 
|---|
| 87 |         DoLog(2) && (Log() << Verbose(2) << U.Norm() << endl);
 | 
|---|
| 88 |       }
 | 
|---|
| 89 |       ekin += 0.5*(*iter)->getType()->getMass() * U.NormSquared();
 | 
|---|
| 90 |     }
 | 
|---|
| 91 |   }
 | 
|---|
| 92 |   return ekin;
 | 
|---|
| 93 | }
 | 
|---|
| 94 | 
 | 
|---|
| 95 | std::string Langevin::name(){
 | 
|---|
| 96 |   return ThermostatTraits<Langevin>::name;
 | 
|---|
| 97 | }
 | 
|---|
| 98 | 
 | 
|---|
| 99 | std::string Langevin::writeParams(){
 | 
|---|
| 100 |   stringstream sstr;
 | 
|---|
| 101 |   sstr << TempFrequency << "\t" << alpha;
 | 
|---|
| 102 |   return sstr.str();
 | 
|---|
| 103 | }
 | 
|---|