Changeset c0c650 for src/Thermostats
- Timestamp:
- Aug 25, 2010, 9:24:25 AM (15 years ago)
- Branches:
- 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
- Children:
- 3e4162
- Parents:
- 0b882a
- Location:
- src/Thermostats
- Files:
-
- 13 edited
Legend:
- Unmodified
- Added
- Removed
-
src/Thermostats/Berendsen.cpp
r0b882a rc0c650 25 25 {} 26 26 27 Thermostat *ThermostatTraits<class Berendsen>::make(class ConfigFileBuffer * const fb){ 28 double TempFrequency; 29 const int verbose = 0; 30 ParseForParameter(verbose,fb,"Thermostat", 0, 2, 1, double_type, &TempFrequency, 1, critical); // read \tau_T 31 return new class Berendsen(TempFrequency); 32 } 27 33 28 34 double Berendsen::scaleAtoms(unsigned int step,double ActualTemp,ATOMSET(std::list) atoms){ … … 57 63 } 58 64 59 std::string Berendsen::writeParam (){65 std::string Berendsen::writeParams(){ 60 66 stringstream sstr; 61 67 sstr << TempFrequency; -
src/Thermostats/Berendsen.hpp
r0b882a rc0c650 23 23 virtual std::string name(); 24 24 25 virtual std::string writeParam ();25 virtual std::string writeParams(); 26 26 27 27 private: … … 36 36 { 37 37 ThermostatTraits(); 38 virtual Thermostat *make(class ConfigFileBuffer * const fb); 38 39 const char* name; 39 40 }; -
src/Thermostats/GaussianThermostat.cpp
r0b882a rc0c650 28 28 name("Gaussian") 29 29 {} 30 31 Thermostat *ThermostatTraits<GaussianThermostat>::make(class ConfigFileBuffer * const fb){ 32 int ScaleTempStep; 33 const int verbose = 0; 34 ParseForParameter(verbose,fb,"Thermostat", 0, 2, 1, int_type, &ScaleTempStep, 1, critical); // read collision rate 35 return new class GaussianThermostat(ScaleTempStep); 36 } 30 37 31 38 double GaussianThermostat::scaleAtoms(unsigned int step,double ActualTemp,ATOMSET(std::list) atoms){ -
src/Thermostats/GaussianThermostat.hpp
r0b882a rc0c650 41 41 { 42 42 ThermostatTraits(); 43 virtual Thermostat *make(class ConfigFileBuffer * const fb); 43 44 const char* name; 44 45 }; -
src/Thermostats/Langevin.cpp
r0b882a rc0c650 13 13 #include "ThermoStatContainer.hpp" 14 14 15 Langevin::Langevin() 15 Langevin::Langevin(double _TempFrequency,double _alpha) : 16 TempFrequency(_TempFrequency), 17 alpha(_alpha) 16 18 { 17 19 gsl_rng_env_setup(); … … 28 30 name("Langevin") 29 31 {} 32 33 Thermostat *ThermostatTraits<class Langevin>::make(class ConfigFileBuffer * const fb){ 34 double TempFrequency; 35 double alpha; 36 const int verbose = 0; 37 ParseForParameter(verbose,fb,"Thermostat", 0, 2, 1, double_type, &TempFrequency, 1, critical); // read gamma 38 if (ParseForParameter(verbose,fb,"Thermostat", 0, 3, 1, double_type, &alpha, 1, optional)) { 39 DoLog(2) && (Log() << Verbose(2) << "Extended Stochastic Thermostat detected with interpolation coefficient " << alpha << "." << endl); 40 } else { 41 alpha = 1.; 42 } 43 return new class Langevin(TempFrequency,alpha); 44 } 30 45 31 46 double Langevin::scaleAtoms(unsigned int step,double ActualTemp,ATOMSET(std::list) atoms){ … … 51 66 // throw a dice to determine whether it gets hit by a heat bath particle 52 67 if (((((rand()/(double)RAND_MAX))*TempFrequency) < 1.)) { 53 DoLog(3) && (Log() << Verbose(3) << "Particle " << (**iter) << " was hit (sigma " << sigma << "): " << sqrt(U[0]*U[0]+U[1]*U[1]+U[2]*U[2]) << " -> ");68 DoLog(3) && (Log() << Verbose(3) << "Particle " << (**iter) << " was hit (sigma " << sigma << "): " << U.Norm() << " -> "); 54 69 // pick three random numbers from a Boltzmann distribution around the desired temperature T for each momenta axis 55 70 for (int d=0; d<NDIM; d++) { 56 71 U[d] = gsl_ran_gaussian (r, sigma); 57 72 } 58 DoLog(2) && (Log() << Verbose(2) << sqrt(U[0]*U[0]+U[1]*U[1]+U[2]*U[2]) << endl);73 DoLog(2) && (Log() << Verbose(2) << U.Norm() << endl); 59 74 } 60 75 ekin += 0.5*(*iter)->getType()->mass * U.NormSquared(); -
src/Thermostats/Langevin.hpp
r0b882a rc0c650 14 14 { 15 15 public: 16 Langevin( );16 Langevin(double,double); 17 17 virtual ~Langevin(); 18 18 … … 39 39 { 40 40 ThermostatTraits(); 41 virtual Thermostat *make(class ConfigFileBuffer * const fb); 41 42 const char* name; 42 43 }; -
src/Thermostats/NoThermostat.cpp
r0b882a rc0c650 17 17 name("None") 18 18 {} 19 20 Thermostat *ThermostatTraits<class NoThermostat>::make(class ConfigFileBuffer * const fb){ 21 return new class NoThermostat(); 22 } 19 23 20 24 double NoThermostat::scaleAtoms(unsigned int step,double ActualTemp,ATOMSET(std::list) atoms){ -
src/Thermostats/NoThermostat.hpp
r0b882a rc0c650 29 29 { 30 30 ThermostatTraits(); 31 virtual Thermostat *make(class ConfigFileBuffer * const fb); 31 32 const char* name; 32 33 }; -
src/Thermostats/NoseHoover.cpp
r0b882a rc0c650 16 16 17 17 NoseHoover::NoseHoover(double _HooverMass) : 18 HooverMass(_HooverMass) 18 HooverMass(_HooverMass), 19 alpha(0) 19 20 {} 20 21 … … 26 27 {} 27 28 29 Thermostat *ThermostatTraits<class NoseHoover>::make(class ConfigFileBuffer * const fb){ 30 double HooverMass; 31 const int verbose = 0; 32 ParseForParameter(verbose,fb,"Thermostat", 0, 2, 1, double_type, &HooverMass, 1, critical); // read Hoovermass 33 return new class NoseHoover(HooverMass); 34 } 28 35 29 36 double NoseHoover::scaleAtoms(unsigned int step,double ActualTemp,ATOMSET(std::list) atoms){ … … 44 51 init(step,begin,end); 45 52 delta_alpha = (delta_alpha - (3.*count+1.) * getContainer().TargetTemp)/(HooverMass*Units2Electronmass); 46 getContainer().alpha += delta_alpha*World::getInstance().getConfig()->Deltat;47 DoLog(3) && (Log() << Verbose(3) << "alpha = " << delta_alpha << " * " << World::getInstance().getConfig()->Deltat << " = " << getContainer().alpha << "." << endl);53 alpha += delta_alpha*World::getInstance().getConfig()->Deltat; 54 DoLog(3) && (Log() << Verbose(3) << "alpha = " << delta_alpha << " * " << World::getInstance().getConfig()->Deltat << " = " << alpha << "." << endl); 48 55 double ekin =0; 49 56 for(ForwardIterator iter=begin;iter!=end;++iter){ 50 57 Vector &U = (*iter)->Trajectory.U.at(step); 51 58 if ((*iter)->FixedIon == 0) { // even FixedIon moves, only not by other's forces 52 U += World::getInstance().getConfig()->Deltat/(*iter)->getType()->mass * ( getContainer().alpha * (U * (*iter)->getType()->mass));59 U += World::getInstance().getConfig()->Deltat/(*iter)->getType()->mass * (alpha * (U * (*iter)->getType()->mass)); 53 60 ekin += (0.5*(*iter)->getType()->mass) * U.NormSquared(); 54 61 } -
src/Thermostats/NoseHoover.hpp
r0b882a rc0c650 32 32 double HooverMass; 33 33 double delta_alpha; 34 double alpha; 34 35 int count; 35 36 }; … … 39 40 { 40 41 ThermostatTraits(); 42 virtual Thermostat *make(class ConfigFileBuffer * const fb); 41 43 const char* name; 42 44 }; -
src/Thermostats/Thermostat.hpp
r0b882a rc0c650 46 46 struct ThermostatTraits<Thermostat>{ 47 47 ThermostatTraits(); 48 virtual Thermostat *make(class ConfigFileBuffer * const fb)=0; 48 49 const char* name; 49 50 }; -
src/Thermostats/Woodcock.cpp
r0b882a rc0c650 26 26 name("Woodcock") 27 27 {} 28 29 Thermostat *ThermostatTraits<class Woodcock>::make(class ConfigFileBuffer * const fb){ 30 int ScaleTempStep; 31 const int verbose = 0; 32 ParseForParameter(verbose,fb,"Thermostat", 0, 2, 1, int_type, &ScaleTempStep, 1, critical); // read scaling frequency 33 return new class Woodcock(ScaleTempStep); 34 } 28 35 29 36 double Woodcock::scaleAtoms(unsigned int step,double ActualTemp,ATOMSET(std::list) atoms){ -
src/Thermostats/Woodcock.hpp
r0b882a rc0c650 36 36 { 37 37 ThermostatTraits(); 38 virtual Thermostat *make(class ConfigFileBuffer * const fb); 38 39 const char* name; 39 40 };
Note:
See TracChangeset
for help on using the changeset viewer.