source: mulcan_sim.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.4 KB
Line 
1! **************************************************************
2!
3! This file contains the subroutines: mulcan_sim,muca_weight2
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_sim
13!
14! PURPOSE: PERFORM A MULTICANONICAL SIMULATION
15! REQUIRES AS INPUT THE MULTICANONICAL PARAMETER AS CALCULATED
16! BY THE SUBROUTINE mulcan_par
17!
18! CALLS: addang, contacts,energy,metropolis
19!
20 include 'INCL.H'
21
22! external rand
23 external muca_weight2
24
25 logical restart
26
27 parameter(restart=.false.)
28 parameter(kmin=-12,kmax=20,ebin=1.0d0)
29 parameter(nsweep=100000,nequi=100)
30 Parameter(nsave=1000,nmes=10)
31!
32! restart: .true. = restart of simulation
33! .false. = start of simulation with random configuration
34! kmin,kmax: Range of multicanonical parameter
35! ebin: bin size for multicanonical parameter
36! nequi: Number of sweeps for equilibrisation
37! nsweep: Number of sweeps for simulation run
38! nsave: Number of sweeps after which actual configuration is saved
39! for re-starts
40! nmes: Number of sweeps between measurments
41!
42
43 dimension xhist(kmin:kmax),ihist(kmin:kmax)
44 common/muca2/b(kmin:kmax),alpha(kmin:kmax)
45
46
47! FILE with last conformation (for re-starts)
48 open(8,file='EXAMPLES/start.d')
49! File with contact map of reference configuration
50 open(9,file='EXAMPLES/enkefa.ref')
51! File with multicanonical parameter
52 open(10,file='EXAMPLES/muca.d')
53! Result file: Time series of certain quantities
54 open(11, file='EXAMPLES/time.d')
55
56
57 do j=kmin,kmax
58 ihist(j)= 0
59 end do
60
61 nresi=irsml2(1)-irsml1(1) + 1
62! nresi: Number of residues
63
64! READ REFERENCE CONTACT MAP
65 nci = 0
66 do i=1,nresi
67 read(9,*) (iref(i,j) , j=1,nresi)
68 end do
69 do i=1,nresi
70 do j=nresi,i+3,-1
71 if(iref(i,j).eq.1) nci = nci + 1
72 end do
73 end do
74 write (logString, *) 'Number of contacts in reference conformation:',nci
75
76! READ IN FIELDS WITH MULTICANONICAL PARAMETER
77 Do j=kmin,kmax
78 read(10,*) i,b(i),alpha(i)
79 end do
80!
81
82 if(restart) then
83 read(8,*) nswm, eol_old
84 read(8,*) (xhist(j), j=kmin,kmax)
85 do i=1,nvr
86 read(8,*) j, x
87 iv = idvr(j)
88 vlvr(iv) = x
89 end do
90 write (logString, *) 'Last iteration, energy:',nswm,eol_old
91 else
92! _________________________________ random start
93 do i=1,nvr
94 iv=idvr(i) ! provides index of non-fixed variable
95 dv=axvr(iv)*(grnd()-0.5)
96 vr=addang(pi,dv)
97 vlvr(iv)=vr
98 enddo
99 end if
100!
101 eol = energy()
102 write (logString, '(e12.5,/)') eol
103 call contacts(nhy,nhx,dham)
104 write (logString, *) 'Number of contacts in start configuration:',nhy
105 write (logString, *) 'Number of native contacts in start configuration:',
106 & nhx
107 do i=1,nresi
108 write (logString, '(62I1)') (ijcont(i,j), j=1,nresi)
109 end do
110 write (logString, *)
111!
112
113
114 if(.not.restart) then
115! =====================Equilibrization by Metropolis
116 do nsw=1,nequi
117 call metropolis(eol,acz,muca_weight2)
118 end do
119 do i=kmin,kmax
120 ihist(i) = 0
121 end do
122 nswm = 1
123 end if
124
125!======================Simulation
126 acz = 0.0d0
127! LOOP OVER SWEEPS
128 do nsw=nswm,nsweep
129!
130! METROPOLIS UPDATE
131 call metropolis(eol,acz,muca_weight2)
132 muold = min(kmax,max(kmin,int(eol/ebin+sign(0.5d0,eol))))
133 ihist(muold) = ihist(muold) + 1
134!
135! SAVE ACTUAL CONFORMATIONS FOR RE-STARTS:
136 if(mod(nsw,nsave).eq.0) then
137 rewind 8
138 write(8,*) nswm, eol
139 write(8,*) (xhist(j), j=kmin,kmax)
140 do i=1,nvr
141 iv = idvr(i)
142 write(8,*) i,vlvr(iv)
143 end do
144 end if
145! Measurements after NMES sweeps
146 if(mod(nsw,nmes).eq.0) then
147! Take a histogram of energy
148 do i=kmin,kmax
149 xhist(i) = xhist(i) + ihist(i)
150 ihist(i) = 0
151 end do
152! Calculate contacts in actual configuartion and compare with reference
153! configuration
154! call contacts(nhx,nhy,dham)
155! nhx : Number of contcats in actual conformation
156! nhy : Number of contacts which are identical in actual and reference
157! configuration
158! dham: Hamming distance between actual and reference configuration
159!
160 write(11,'(i7,f12.2,2i8,f12.4)') nsw,eol,nhx,nhy,dham
161 end if
162 end do
163! END OF SIMULATION
164
165! FINAL OUTPUT:
166 acz = acz/dble(nsw*nvr)
167 write (logString, *) 'last energy',eol
168 write (logString, *) 'aczeptance rate:',acz
169
170! WRITE DOWN (UN-REWEIGHTED) HISTOGRAM OF MULTICANONICAL SIMULATION
171 do i=kmin,kmax
172 if(xhist(i).gt.0.0d0) then
173 write (logString, *) i,xhist(i)
174 end if
175 end do
176! =====================
177 close(8)
178 close(9)
179 close(10)
180 close(11)
181
182 return
183 end
184
185! ************************************************************
186 real*8 function muca_weight2(x)
187
188 implicit real*8 (a-h,o-z)
189 implicit integer*4 (i-n)
190
191 Parameter(kmin=-12,kmax=20,ebin=1.0d0)
192
193 common/muca2/b(kmin:kmax),alpha(kmin:kmax)
194
195 muold = min(kmax,max(kmin,int(x/ebin+sign(0.5d0,x))))
196 muca_weight2 = b(muold)*x + alpha(muold)
197
198 return
199
200 end
201
202
Note: See TracBrowser for help on using the repository browser.