source: pysmmp/ParallelTempering.py

Last change on this file was 225475c, checked in by Jan Meinke <j.meinke@…>, 14 years ago

Move pysmmp in its own subdirectory.

  • Property mode set to 100644
File size: 4.7 KB
Line 
1#!/usr/bin/env python
2
3#import smmp
4import sys, copy, random, math, numpy as np
5
6from mpi4py import MPI
7#from algorithms import CanonicalMonteCarlo
8
9class ParallelTempering(object):
10 """Parallel tempering is a generalized Monte Carlo algorithm. A set of
11 replica is simulated at different temperatures. Each replica performs a
12 canoncial Monte Carlo simulation. After a given number of sweeps, an exchange
13 is attempted between neighboring replicas.
14 """
15
16 def __init__(self, T, myUniverse):
17 """Initialize ParallelTempering with a set of tempertures T.
18 """
19 self.kB = 1.98773e-3
20
21 self._T = copy.deepcopy(T)
22 self._numberOfReplica = len(T)
23
24 # Calculate the number of processors available per replica
25 numberOfProcessors = MPI.COMM_WORLD.Get_size()
26 numberOfProcessorsPerReplica = numberOfProcessors // self._numberOfReplica
27 numberOfProcessors = numberOfProcessorsPerReplica * self._numberOfReplica
28 self._myReplica = MPI.COMM_WORLD.Get_rank() // numberOfProcessorsPerReplica
29
30 # Create the replica internal communicators for the energy function
31 groupWorld = MPI.COMM_WORLD.Get_group()
32 group = [MPI.Group.Range_incl(groupWorld, [[i * numberOfProcessorsPerReplica, (i + 1) * numberOfProcessorsPerReplica - 1, 1], ]) for i in xrange(self._numberOfReplica)]
33 comm = [MPI.COMM_WORLD.Create(g) for g in group]
34
35 # Create the communicator needed to exchange information about replica
36 ptGroup = MPI.Group.Incl(groupWorld, range(0, numberOfProcessors, numberOfProcessorsPerReplica))
37 self._ptComm = MPI.COMM_WORLD.Create(ptGroup)
38
39 self._localMaster = self._ptComm != MPI.COMM_NULL
40 if self._localMaster:
41 self._myPTRank = self._ptComm.Get_rank()
42
43 self._oddCycle = False
44
45
46 def exchange(self):
47 """Collect temperatures and energies from all replicas, determine neighbor and do replica exchange move.
48
49 The exchange step is the heart of parallel tempering. The probability of exchange between two
50 temperatures i and j is given by P_{ij} = \Delta \beta \Delta E, where \beta = 1/k_B T.
51 """
52 if self._localMaster:
53# energies = ptComm.allgather(self._myUniverse.energy())
54 myE = [250, 255, 75, 125, 120, 185][self._myReplica]
55 myT = [420, 370, 250, 275, 300, 330][self._myReplica]
56 energies = self._ptComm.allgather(myE)
57# temperatures = ptComm.allgather(self._myUniverse.temperature())
58 temperatures = self._ptComm.allgather(myT)
59
60 temperatureRank = zip(temperatures, np.arange(self._numberOfReplica))
61 temperatureRank.sort()
62
63 if (self._myPTRank % 2):
64 if self._oddCycle:
65 rank0 = self._myPTRank
66 else:
67 rank0 = self._myPTRank - 1
68
69 rank1 = (rank0 + 1) % self._ptComm.Get_size()
70
71 self._oddCycle = not self._oddCycle
72
73 E0 = energies[temperatureRank[rank0][1]]
74 E1 = energies[temperatureRank[rank1][1]]
75 T0 = self._T[rank0]
76 T1 = self._T[rank1]
77 beta0 = 1.0 / (T0 * self.kB)
78 beta1 = 1.0 / (T1 * self.kB)
79
80 P = math.exp((beta0 - beta1) * (E0 - E1))
81 print self._myPTRank, P
82
83 if P >= random.random(): # accept exchange
84 rank = self._ptComm.allgather(np.array([temperatureRank[rank1][1], temperatureRank[rank0][1]]))
85 else:
86 rank = self._ptComm.allgather(np.array([temperatureRank[rank0][1], temperatureRank[rank1][1]]))
87 else:
88 random.random() # Do one random number to stay in sync with the other processes for the same replica.
89 rank = self._ptComm.allgather(np.array([-1, -1]))
90
91 rank = np.array(rank[1::2]).flatten()
92 print(rank)
93 rankTemperature = zip(rank, self._T)
94 rankTemperature.sort()
95 print self._myPTRank, rankTemperature[self._myReplica][1]
96# self._myUniverse.setTemperature(rankTemperature[self._myReplica][1])
97 else:
98 random.random() # Do one random number to stay in sync with the master process for the same replica.
99
100
101
102
103if __name__ == "__main__":
104
105 n = int(sys.argv[1])
106 T = [250, 275, 300, 330, 370, 420]
107
108 myPT = ParallelTempering(T, [])
109 myPT.exchange()
110 myPT.exchange()
111
112
113
114
115
116
117
Note: See TracBrowser for help on using the repository browser.