source: mulcan_par.f

Last change on this file was 38d77eb, checked in by baerbaer <baerbaer@…>, 14 years ago

Redirected standard out to logString.

SMMP produced a lot of log messages. This became an issue when run in massively
parallel environments. I replaced all writes to standard out to a write to logString.
The next step is to pass this string to a function that writes the messages to a log
file according to their importance and the chosen log level.

git-svn-id: svn+ssh://svn.berlios.de/svnroot/repos/smmp/trunk@34 26dc1dd8-5c4e-0410-9ffe-d298b4865968

  • Property mode set to 100644
File size: 5.3 KB
Line 
1! **************************************************************
2!
3! This file contains the subroutines: mulcan_par, muca_weight
4!
5! Copyright 2003-2005 Frank Eisenmenger, U.H.E. Hansmann,
6! Shura Hayryan, Chin-Ku
7! Copyright 2007 Frank Eisenmenger, U.H.E. Hansmann,
8! Jan H. Meinke, Sandipan Mohanty
9!
10! **************************************************************
11
12 subroutine mulcan_par
13!
14! PURPOSE: CALCULATION OF ESTIMATORS FOR MULTICANONICAL WEIGHTS
15! USING A METHOD DESCRIBED FIRST IN: Bernd Berg,
16! {\it J.~Stat.~Phys.} {\bf 82}, 331~(1996). THE METHOD
17! IS STABLE, BUT MAY NOT ALWAYS BE THE FASTEST WAY OF
18! GETTING PARAMETERS
19!
20! CALLS: addang,energy,metropolis
21!
22
23 include 'INCL.H'
24
25! external rand
26 external muca_weight
27
28 logical l_iter
29
30 parameter(l_iter=.false.)
31 parameter(kmin=-12,kmax=20,ebin=1.0d0)
32 parameter(xmin = kmin,xmax=kmax)
33 parameter(nsweep=100000,nup=5000)
34 parameter(temp=1000.0)
35!
36! l_iter: .true. = for iteration of multicanonical parameters
37! .false. = starts new calculation of multicanonica parameters
38! kmin,kmax: Range of multicanonical parameter
39! ebin: bin size for multicanonical parameter
40! nup: Number of sweeps between updates of multicanonical parameters
41! nsweep: Number of sweeps for simulation run
42! temp: temperature for initial canonical run
43!
44
45 dimension xhist(kmin:kmax),g1(kmin:kmax),ent(kmin:kmax)
46 common/muca/b(kmin:kmax),alpha(kmin:kmax),beta
47 common/jhist/ihist(kmin:kmax)
48
49
50! Files with multicanonical parameter
51 open(11,file='EXAMPLES/mpar_full.d')
52 open(12,file='EXAMPLES/muca.d')
53
54
55 beta=1.0/ ( temp * 1.98773d-3 )
56 do j=kmin,kmax
57 ihist(j)= 0
58 end do
59
60
61 if(l_iter) then
62! READ IN FIELDS WITH MULTICANONICAL PARAMETER
63 Do j=kmin,kmax
64 read(11,*) i,xhist(i),g1(i),b(i),alpha(i),ent(i)
65 end do
66 else
67 do j=kmin,kmax
68 b(j) = 0.0d0
69 alpha(j) = 0.0d0
70 g1(j) = 0.0d0
71 xhist(j) = 0.0d0
72 ent(j) = 0.0d0
73 end do
74 end if
75!
76! _________________________________ random start
77 do i=1,nvr
78 iv=idvr(i) ! provides index of non-fixed variable
79 dv=axvr(iv)*(grnd()-0.5)
80 vr=addang(pi,dv)
81 vlvr(iv)=vr
82 enddo
83!
84 eol = energy()
85 write (logString, '(a,e12.5,/)') 'Energy of start configuration: ',eol
86 write (logString, *)
87
88 call outpdb(1, 'start.pdb')
89!
90
91
92!======================Simulation
93 acz = 0.0d0
94! LOOP OVER SWEEPS
95 do nsw=1,nsweep
96!
97! METROPOLIS UPDATE
98 call metropolis(eol,acz,muca_weight)
99
100 muold = int(min(xmax,max(xmin,eol/ebin+sign(0.5d0,eol))))
101 ihist(muold) = ihist(muold) + 1
102 write (logString, *) nsw, eol, acz
103!
104! ITERATE MULTICANONICAL WEIGHTS EVERY NUP SWEEPS
105 if(mod(nsw,nup).eq.0) then
106 xhist(kmax) = xhist(kmax) + ihist(kmax)
107 xhist(kmax-1) = xhist(kmax-1) + ihist(kmax-1)
108 rewind 11
109 do i=kmax-2,kmin,-1
110 xhist(i) = xhist(i) + ihist(i)
111 if((ihist(i).gt.1).and.(ihist(i+1).gt.1)) then
112 g0 = float(ihist(i+1))*float(ihist(i))
113 g0 = g0/float(ihist(i)+ihist(i+1))
114 g1(i) = g1(i) + g0
115 g0 = g0/g1(i)
116 bi_1 = log(float(ihist(i+1)))
117 bi = log(float(ihist(i)))
118 b(i) = b(i) + g0*(bi_1 - bi)/ebin
119 else
120 if(xhist(i).lt.2.0) then
121 b(i) = max(0.0d0,b(i+1))
122 end if
123 end if
124 end do
125 b(kmax) = b(kmax-2)
126 b(kmax-1) = b(kmax-2)
127 alpha(kmax) = 0
128 ent(kmax) = (beta+b(kmax))*float(kmax)*ebin + alpha(kmax)
129 write(11,'(i5,5f15.4)') kmax,xhist(kmax),g1(kmax),
130 & b(kmax),alpha(kmax),ent(kmax)
131 ihist(kmax) = 0
132 do i=kmax-1,kmin,-1
133 ihist(i) = 0
134 alpha(i) = alpha(i+1) + (b(i+1) - b(i))*(i+0.5d0)*ebin
135 ent(i) = (beta+b(i))*float(i)*ebin + alpha(i)
136 write(11,'(i5,5f15.4)') i,xhist(i),g1(i),b(i),alpha(i),ent(i)
137 end do
138 end if
139!
140 end do
141! END OF SIMULATION
142
143! FINAL OUTPUT:
144 acz = acz/dble(nsw*nvr)
145 write (logString, *) 'last energy',eol
146 write (logString, *) 'aczeptance rate:',acz
147!
148! OUTPUT OF FINAL HISTOGRAM
149!
150 write (logString, *) 'Histogram:'
151 do i=kmin,kmax
152 if(xhist(i).gt.0.0d0) then
153 write (logString, *) i,xhist(i)
154 end if
155 end do
156!
157! WRITE DOWN MULTICANONICAL PARAMETER
158 alpha(kmax)=0.0d0
159 b(kmax) = b(kmax)+beta
160 do i= kmax-1,kmin,-1
161 b(i) = b(i) + beta
162 alpha(i) = alpha(i+1) + (b(i+1) - b(i))*(i+0.5d0)*ebin
163 end do
164 rewind 12
165 do i=kmin,kmax
166 write(12,*) i,b(i),alpha(i)
167 end do
168
169 close(11)
170 close(12)
171!
172! =====================
173!
174 return
175 end
176
177!**********************************************************************
178
179 real*8 function muca_weight(x)
180
181 implicit real*8 (a-h,o-z)
182 implicit integer*8 (i-n)
183
184 parameter(kmin=-12,kmax=20,ebin=1.0d0)
185
186 common/muca/b(kmin:kmax),alpha(kmin:kmax),beta
187
188 xmin = kmin
189 xmax = kmax
190 muold = int(min(xmax,max(xmin,x/ebin+sign(0.5d0,x))))
191 muca_weight = (beta+b(muold))*x + alpha(muold)
192
193 return
194
195 end
196
Note: See TracBrowser for help on using the repository browser.