source: src/Thermostats/Berendsen.cpp@ 343c5a

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 343c5a was f73351, checked in by Frederik Heber <heber@…>, 11 years ago

FIX: Some Thermostats would cause NaN on zero temperature.

  • Berendsen and Woodcock divide through OldTemp which may be zero (because initial velocities are zero). This is now checked.
  • also, default thermostat is now NoThermostat.
  • NOTE: We lack an action to choose a thermostat but this needs some refactoring of the thermostats which are still based on being read from file.
  • TESTFIX: Adapted three regression tests where pcp files were stored that contained still reference to old default Berendsen (with parameter 2.5).
  • Property mode set to 100644
File size: 3.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 * Berendsen.cpp
25 *
26 * Created on: Aug 20, 2010
27 * Author: crueger
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 "Berendsen.hpp"
38
39#include "CodePatterns/Log.hpp"
40#include "config.hpp"
41#include "Element/element.hpp"
42#include "Helpers/defs.hpp"
43#include "Parser/PcpParser_helper.hpp"
44#include "Thermostats/ThermoStatContainer.hpp"
45#include "World.hpp"
46
47Berendsen::Berendsen(double _TempFrequency) :
48 TempFrequency(_TempFrequency)
49{}
50
51Berendsen::Berendsen() :
52 TempFrequency(2.5)
53{}
54
55Berendsen::~Berendsen()
56{}
57
58const char *ThermostatTraits<Berendsen>::name = "Berendsen";
59
60std::string ThermostatTraits<Berendsen>::getName(){
61 return ThermostatTraits<Berendsen>::name;
62}
63
64Thermostat *ThermostatTraits<Berendsen>::make(class ConfigFileBuffer * const fb){
65 double TempFrequency;
66 const int verbose = 0;
67 ParseForParameter(verbose,fb,"Thermostat", 0, 2, 1, double_type, &TempFrequency, 1, critical); // read \tau_T
68 return new Berendsen(TempFrequency);
69}
70
71double Berendsen::scaleAtoms(unsigned int step,double ActualTemp,ATOMSET(std::list) atoms){
72 return doScaleAtoms(step,ActualTemp,atoms.begin(),atoms.end());
73}
74
75double Berendsen::scaleAtoms(unsigned int step,double ActualTemp,ATOMSET(std::vector) atoms){
76 return doScaleAtoms(step,ActualTemp,atoms.begin(),atoms.end());
77}
78
79double Berendsen::scaleAtoms(unsigned int step,double ActualTemp,ATOMSET(std::set) atoms){
80 return doScaleAtoms(step,ActualTemp,atoms.begin(),atoms.end());
81}
82
83template <class ForwardIterator>
84double Berendsen::doScaleAtoms(unsigned int step,double ActualTemp,ForwardIterator begin, ForwardIterator end){
85 LOG(2, "Applying Berendsen-VanGunsteren thermostat...");
86 double ekin=0;
87 if (fabs(ActualTemp) > MYEPSILON) {
88 double ScaleTempFactor = getContainer().TargetTemp/ActualTemp;
89 for(ForwardIterator iter=begin;iter!=end;++iter){
90 Vector U = (*iter)->getAtomicVelocityAtStep(step);
91 if ((*iter)->getFixedIon() == 0) { // even FixedIon moves, only not by other's forces
92 U *= sqrt(1+(World::getInstance().getConfig()->Deltat/TempFrequency)*(ScaleTempFactor-1));
93 ekin += 0.5*(*iter)->getType()->getMass() * U.NormSquared();
94 }
95 (*iter)->setAtomicVelocityAtStep(step, U);
96 }
97 }
98 return ekin;
99}
100
101std::string Berendsen::name(){
102 return ThermostatTraits<Berendsen>::name;
103}
104
105std::string Berendsen::writeParams(){
106 stringstream sstr;
107 sstr << TempFrequency;
108 return sstr.str();
109}
Note: See TracBrowser for help on using the repository browser.