Rosetta 3.5
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ReplicaExchangeMC.cc
Go to the documentation of this file.
1 // -*- mode:c++;tab-width:2;indent-tabs-mode:t;show-trailing-whitespace:t;rm-trailing-spaces:t -*-
2 // vi: set ts=2 noet:
3 //
4 // (c) Copyright Rosetta Commons Member Institutions.
5 // (c) This file is part of the Rosetta software suite and is made available under license.
6 // (c) The Rosetta software is developed by the contributing members of the Rosetta Commons.
7 // (c) For more information, see http://www.rosettacommons.org. Questions about this can be
8 // (c) addressed to University of Washington UW TechTransfer, email: license@u.washington.edu.
9 
10 /// @file protocols/moves/ReplicaExchangeMC.cc
11 /// @brief implementing a Recplica Exchange Monte Carlo Mover
12 /// @author Yuan Liu (wendao@u.washington.edu
13 
14 #ifdef USEMPI
15 #include <mpi.h>
16 #endif //USEMPI
17 
18 #include <cmath>
19 #include <numeric/random/random.hh>
20 #include <basic/Tracer.hh>
22 #include <utility/exit.hh>
23 
24 #include <utility/vector1.hh>
25 
26 
27 static basic::Tracer TR("protocols.ReplicaExchangeMC");
28 static numeric::random::RandomGenerator re_RG(1801); // <- Magic number, do not change it!!!
29 
30 namespace protocols {
31 namespace moves {
32 
34  Pose const & init_pose, // PoseCOP init_pose,
35  ScoreFunction const & scorefxn, // ScoreFunctionCOP scorefxn,
36  utility::vector1<core::Real> const &tlist,
37  core::Size nint
38 ):MonteCarlo(init_pose, scorefxn, 0.0),
39 rank_(0),
40 size_(1),
41 nreplica_frequency_(nint),
42 //noutput_(200),
43 ntrials_(0),
44 Tlist_(tlist),
45 last_energylist(NULL),
46 T_tag(NULL),
47 T_rev(NULL),
48 T_ndx(0)
49 {
50  init();
51 }
52 
54  ScoreFunction const & scorefxn, // ScoreFunctionCOP scorefxn,
55  utility::vector1<core::Real> const &tlist,
56  core::Size nint
57 ):MonteCarlo(scorefxn, 0.0),
58 rank_(0),
59 size_(1),
60 nreplica_frequency_(nint),
61 //noutput_(200),
62 ntrials_(0),
63 Tlist_(tlist),
64 last_energylist(NULL),
65 T_tag(NULL),
66 T_rev(NULL),
67 T_ndx(0)
68 {
69  init();
70 }
71 
73 {
74 #ifdef USEMPI
75  //init mpi parameter
76  MPI_Comm_rank( MPI_COMM_WORLD, &rank_ );
77  MPI_Comm_size( MPI_COMM_WORLD, &size_ );
78 #endif
79 
80  //make sure the tlist match the nproc
81  if( size_ != (int)Tlist_.size() ){
82  std::cerr << "size " << size_ << " and Tlist-size: " << Tlist_.size() << " don't match up!" << std::endl;
83  }
84  runtime_assert(size_ == static_cast<int>(Tlist_.size()));
85 
86  //setup init list
87  last_energylist = new double[size_];
88  T_tag = new int[size_];
89  T_rev = new int[size_];
90 
91  for (int i=0; i<size_; i++)
92  {
93  T_tag[i] = i;
94  T_rev[i] = i;
95  }
96 
97  //set temperature
98  T_ndx = rank_;
100 
101  TR << "I am proc:" << rank_ << ", of all " << size_ << " procs!";
102  TR << " T=" << Tlist_[T_tag[T_ndx]+1] << " ninterval=" << nreplica_frequency_ << std::endl;
103 
104  if (rank_>0) return;
105 
106  //setup exchange_schedule
108  for (int i=0; i<size_-1; i+=2)
109  {
110  std::pair<int, int> elem(i, i+1);
111  list.push_back(elem);
112  }
113  exchange_schedule.push_back(list);
114  list.clear();
115  for (int i=1; i<size_-1; i+=2)
116  {
117  std::pair<int, int> elem(i, i+1);
118  list.push_back(elem);
119  }
120  exchange_schedule.push_back(list);
121 
122  //debug
123  TR << "Building exchange schedule 1" << std::endl;
124  for (Size i=1; i<=exchange_schedule[1].size(); i++)
125  {
126  TR << exchange_schedule[1][i].first << "<-->" << exchange_schedule[1][i].second << std::endl;
127  }
128  TR << "Building exchange schedule 2" << std::endl;
129  for (Size i=1; i<=exchange_schedule[2].size(); i++)
130  {
131  TR << exchange_schedule[2][i].first << "<-->" << exchange_schedule[2][i].second << std::endl;
132  }
133 }
134 
136 {
137  if (last_energylist!=NULL) delete [] last_energylist;
138  if (T_tag!=NULL) delete [] T_tag;
139  if (T_rev!=NULL) delete [] T_rev;
140 }
141 
143 {
144  static int flag=-1;
145  int nlist=(3+flag)/2; // 1 or 2, switch
146 
147  TR << "Building " << nlist << std::endl;
148  for(Size i=1; i<=exchange_schedule[nlist].size(); i++)
149  {
150  Size node1=T_rev[exchange_schedule[nlist][i].first];
151  Size node2=T_rev[exchange_schedule[nlist][i].second];
152 
153  //TR << node1 << "<==>" << node2 << std::endl;
154  //TR << "proc" << node1 << ": e=" << elist[node1] <<" T=" << Tlist_[T_tag[node1]+1] << std::endl;
155  //TR << "proc" << node2 << ": e=" << elist[node2] <<" T=" << Tlist_[T_tag[node2]+1] << std::endl;
156 
157  Real delta=(1.0/Tlist_[T_tag[node1]+1]-1.0/Tlist_[T_tag[node2]+1])*(elist[node2]-elist[node1]);
158  Real probability = 0;
159  //if (delta>0) probability = 1.0;
160  //else probability = std::exp( std::max(-40.0,delta) );
161  if (delta<0) probability = 1.0;
162  else probability = std::exp( std::max(-40.0, -delta) );
163  TR << "Try:" << Tlist_[exchange_schedule[nlist][i].first+1]
164  << " <==> " << Tlist_[exchange_schedule[nlist][i].second+1]
165  << std::endl;
166  //TR << "Delta=" << delta << " Prob=" << probability << std::endl;
167 
168  if (re_RG.uniform()<probability)
169  {
170  TR << "Switch:" << Tlist_[exchange_schedule[nlist][i].first+1] << "<==>" << Tlist_[exchange_schedule[nlist][i].second+1] << std::endl;
171  //switch
172  Size tmp=T_tag[node1];
173  T_tag[node1]=T_tag[node2];
174  T_tag[node2]=tmp;
175  tmp=T_rev[exchange_schedule[nlist][i].first];
176  T_rev[exchange_schedule[nlist][i].first]=T_rev[exchange_schedule[nlist][i].second];
177  T_rev[exchange_schedule[nlist][i].second]=tmp;
178  }
179  }
180  TR << "Done" << std::endl;
181 
182  //switch
183  flag = -flag;
184 }
185 
186 bool
188  Pose & pose,//PoseOP pose,
189  std::string const & move_type,
190  core::Real const proposal_density_ratio)
191 {
192  ntrials_++;
193 
194 #ifdef USEMPI
195  MPI_Barrier(MPI_COMM_WORLD);
196  if (ntrials_ % (nreplica_frequency_/2) == 0)
197  {
198  double last_energy_ = last_accepted_score();
199 
200  //get infomation
201  MPI_Gather(&last_energy_, 1, MPI_DOUBLE, last_energylist, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD);
202  //change the T_tag and T_rev at node0
204  //public the new T_tag
205  MPI_Scatter(T_tag, 1, MPI_INT, &T_ndx, 1, MPI_INT, 0, MPI_COMM_WORLD);
207  }
208 #endif
209 
210  //if (ntrials_ % noutput_ == 0 ) {
211  //this line should be moved to another condition after debug, and controld by another option
212  //TR << "proc=" << rank_ << " step=" << ntrials_ << " T=" << temperature() << " E=" << last_accepted_score() << std::endl;
213  //}
214 
215  return Parent::boltzmann( pose, move_type, proposal_density_ratio );
216 }
217 
218 } // moves
219 } // prot
220 
221