1 | ! **************************************************************
|
---|
2 | !
|
---|
3 | ! This file contains the subroutines: init_energy,setpar
|
---|
4 | ! This file contains a BLOCK DATA statement
|
---|
5 | !
|
---|
6 | ! Copyright 2003-2005 Frank Eisenmenger, U.H.E. Hansmann,
|
---|
7 | ! Shura Hayryan, Chin-Ku
|
---|
8 | ! Copyright 2007 Frank Eisenmenger, U.H.E. Hansmann,
|
---|
9 | ! Jan H. Meinke, Sandipan Mohanty
|
---|
10 | !
|
---|
11 | ! **************************************************************
|
---|
12 |
|
---|
13 |
|
---|
14 | subroutine init_energy(libdir)
|
---|
15 |
|
---|
16 | ! ----------------------------------------------
|
---|
17 | ! PURPOSE: initialize energy parameters
|
---|
18 | ! 0 => ECEPP2 or ECEPP3 depending on the value of sh2
|
---|
19 | ! 1 => FLEX
|
---|
20 | ! 2 => Lund force field
|
---|
21 | ! 3 => ECEPP with Abagyan corrections
|
---|
22 | !
|
---|
23 | !
|
---|
24 | ! CALLS: setpar, tessel,iendst
|
---|
25 | !
|
---|
26 | ! contains: BLOCK DATA
|
---|
27 | ! ----------------------------------------------
|
---|
28 |
|
---|
29 | include 'INCL.H'
|
---|
30 |
|
---|
31 | character libdir*(*),tesfil*80
|
---|
32 |
|
---|
33 | if (ientyp.eq.1) then
|
---|
34 | flex = .true.
|
---|
35 | else
|
---|
36 | flex = .false.
|
---|
37 | end if
|
---|
38 |
|
---|
39 | lunlib=10
|
---|
40 | ll=iendst(libdir)
|
---|
41 |
|
---|
42 | if (flex) then ! Flex
|
---|
43 |
|
---|
44 | reslib=libdir(1:ll)//'lib.flx'
|
---|
45 | lunchg=12
|
---|
46 | chgfil=libdir(1:ll)//'charges'
|
---|
47 |
|
---|
48 | else
|
---|
49 |
|
---|
50 | if (sh2) then ! Scheraga (ECEPP/2)
|
---|
51 | reslib=libdir(1:ll)//'lib.sh2'
|
---|
52 | else ! Scheraga (ECEPP/3)
|
---|
53 | reslib=libdir(1:ll)//'lib.sh3'
|
---|
54 | endif
|
---|
55 |
|
---|
56 | endif
|
---|
57 |
|
---|
58 | if (ientyp .eq. 2) then
|
---|
59 | reslib=libdir(1:ll)//'lib.lun'
|
---|
60 | endif
|
---|
61 |
|
---|
62 | call setpar() ! Initialize force field parameters
|
---|
63 |
|
---|
64 |
|
---|
65 | !----Initialize solvation part if necessary
|
---|
66 | write (*,*) 'init_energy: itysol = ',itysol
|
---|
67 | write(*,*) 'init_energy: esol_scaling = ',isolscl
|
---|
68 |
|
---|
69 | its = iabs(itysol)
|
---|
70 |
|
---|
71 | if (its.gt.0.and.its.le.mxtysol) then
|
---|
72 |
|
---|
73 | ll=iendst(libdir)
|
---|
74 | tesfil = libdir(1:ll)//'tes.dat'
|
---|
75 | open(unit=20,file=tesfil,status='old',err=10)
|
---|
76 | call tessel()
|
---|
77 | close(20)
|
---|
78 |
|
---|
79 | else
|
---|
80 |
|
---|
81 | if (itysol.ne.0) then
|
---|
82 | write(*,'(a)') ' init_energy> undefined solvent type !'
|
---|
83 | stop
|
---|
84 | endif
|
---|
85 |
|
---|
86 | endif
|
---|
87 |
|
---|
88 | ! ___________________________ initialise COMMON 'con_r'
|
---|
89 | idloa=ichar('a')
|
---|
90 | idloz=ichar('z')
|
---|
91 | idupa=ichar('A')
|
---|
92 | idupz=ichar('Z')
|
---|
93 |
|
---|
94 | return
|
---|
95 |
|
---|
96 | 10 write (*, '(a)') 'Cannot open the file with surface points'
|
---|
97 | stop
|
---|
98 |
|
---|
99 | end
|
---|
100 |
|
---|
101 | ! *********************
|
---|
102 | subroutine setpar
|
---|
103 |
|
---|
104 | ! __________________________________________________________
|
---|
105 | ! PURPOSE: initialize parameter set for empirical potentials
|
---|
106 | ! depending on variable 'flex'
|
---|
107 | !
|
---|
108 | ! CALLS: None
|
---|
109 | ! __________________________________________________________
|
---|
110 |
|
---|
111 | include 'INCL.H'
|
---|
112 |
|
---|
113 | dimension hbc(mxhbdo,mxhbac),hba(mxhbdo,mxhbac)
|
---|
114 | logical do(mxtyat),ac(mxtyat)
|
---|
115 |
|
---|
116 |
|
---|
117 | tesgrd = .false. ! numerical check of analytical gradients
|
---|
118 |
|
---|
119 | ! ______________________________________ Lennard-Jones parameters
|
---|
120 | if (flex) then
|
---|
121 |
|
---|
122 | plt=plt_f
|
---|
123 | slp=slp_f
|
---|
124 | cohb=cohb_f
|
---|
125 |
|
---|
126 | do i=1,mxtyat
|
---|
127 | do j=1,i
|
---|
128 | cij(i,j)=c_f(i,j)
|
---|
129 | cij(j,i)=cij(i,j)
|
---|
130 | aij(i,j)=a_f(i,j)*1000
|
---|
131 | aij(j,i)=aij(i,j)
|
---|
132 | enddo
|
---|
133 | do j=1,mxtyat
|
---|
134 | a14(i,j)=aij(i,j)
|
---|
135 | enddo
|
---|
136 |
|
---|
137 | do(i)=do_f(i)
|
---|
138 | ac(i)=ac_f(i)
|
---|
139 | enddo
|
---|
140 |
|
---|
141 | do i=1,mxhbdo
|
---|
142 | do j=1,mxhbac
|
---|
143 | hbc(i,j)=chb_f(i,j)
|
---|
144 | hba(i,j)=ahb_f(i,j)
|
---|
145 | enddo
|
---|
146 | enddo
|
---|
147 |
|
---|
148 | do i=1,mxtyto
|
---|
149 | e0to(i)=e0to_f(i)/2.
|
---|
150 | sgto(i)=sgto_f(i)
|
---|
151 | rnto(i)=rnto_f(i)
|
---|
152 | enddo
|
---|
153 |
|
---|
154 | else ! ---------------------- Scheraga:
|
---|
155 |
|
---|
156 | conv=conv/eps_s
|
---|
157 |
|
---|
158 | do i=1,mxtyat
|
---|
159 | atpl(i)=atpl(i)/100
|
---|
160 | efel(i)=efel(i)/100
|
---|
161 | emin(i)=emin(i)/1000
|
---|
162 | rmin(i)=rmin(i)/100
|
---|
163 | enddo
|
---|
164 |
|
---|
165 | do i=1,mxtyat
|
---|
166 | ri=rmin(i)
|
---|
167 | ai=atpl(i)
|
---|
168 | aei=sqrt(ai/efel(i))
|
---|
169 | !c aic=ai/ehm !! ICM
|
---|
170 | !c do j=i,mxtyat !! -"-
|
---|
171 | aic=ai*ehm !! comment for ICM:
|
---|
172 | cij(i,i)=aic*ai/(aei+aei) !! -"-
|
---|
173 | a=.5*cij(i,i)*ri**6
|
---|
174 | aij(i,i)=a !!
|
---|
175 | a14(i,i)=.5*a !!
|
---|
176 | do j=i+1,mxtyat !!
|
---|
177 | aj=atpl(j)
|
---|
178 | ! _______ Constant for 6-12 attractive term (Slater-Kirkwood formula)
|
---|
179 | c=aic*aj/(aei+sqrt(aj/efel(j)))
|
---|
180 | cij(i,j)=c
|
---|
181 | cij(j,i)=c
|
---|
182 | ! ____________________________ repulsive term (form. 3 & 6 of ref 2)
|
---|
183 | rij=.5*(ri+rmin(j))
|
---|
184 | a=.5*c*rij**6
|
---|
185 | aij(i,j)=a
|
---|
186 | aij(j,i)=a
|
---|
187 | a=.5*a
|
---|
188 | a14(i,j)=a
|
---|
189 | a14(j,i)=a
|
---|
190 | enddo
|
---|
191 | do(i)=do_s(i)
|
---|
192 | ac(i)=ac_s(i)
|
---|
193 | enddo
|
---|
194 |
|
---|
195 | ! +++++++++++++++++++++++++++++++++
|
---|
196 | cij(1,1)=45.5d0
|
---|
197 | aij(1,1)=14090.0d0
|
---|
198 | cij(2,2)=45.5d0
|
---|
199 | aij(2,2)=14380.0d0
|
---|
200 | cij(3,3)=45.5d0
|
---|
201 | aij(3,3)=8420.0d0
|
---|
202 | cij(4,4)=45.5d0
|
---|
203 | aij(4,4)=8420.0d0
|
---|
204 | cij(5,5)=45.5d0
|
---|
205 | aij(5,5)=11680.0d0
|
---|
206 | cij(6,6)=45.5d0
|
---|
207 | aij(6,6)=14380.0d0
|
---|
208 | cij(7,7)=370.5d0
|
---|
209 | aij(7,7)=906100.0d0
|
---|
210 | cij(8,8)=766.6d0
|
---|
211 | aij(8,8)=1049000.0d0
|
---|
212 | cij(9,9)=509.5d0
|
---|
213 | aij(9,9)=653600.0d0
|
---|
214 | cij(10,10)=217.2d0
|
---|
215 | aij(10,10)=125600.0d0
|
---|
216 | cij(11,11)=369.0d0
|
---|
217 | aij(11,11)=170200.0d0
|
---|
218 | cij(12,12)=217.2d0
|
---|
219 | aij(12,12)=125600.0d0
|
---|
220 | cij(13,13)=401.3d0
|
---|
221 | aij(13,13)=375200.0d0
|
---|
222 | cij(14,14)=401.3d0
|
---|
223 | aij(14,14)=375200.0d0
|
---|
224 | cij(15,15)=401.3d0
|
---|
225 | aij(15,15)=375200.0d0
|
---|
226 | cij(16,16)=2274.4d0
|
---|
227 | aij(16,16)=5809000.0d0
|
---|
228 | cij(17,17)=45.5d0
|
---|
229 | aij(17,17)=5340.0d0
|
---|
230 | cij(18,18)=370.5d0
|
---|
231 | aij(18,18)=909000.0d0
|
---|
232 | ! +++++++++++++++++++++++++++++++++
|
---|
233 | do i=1,mxtyat
|
---|
234 | a14(i,i)=.5*aij(i,i)
|
---|
235 | enddo
|
---|
236 | ! +++++++++++++++++++++++++++++++++
|
---|
237 |
|
---|
238 | do i=1,mxtyat
|
---|
239 | ! write( *, '(18f14.6)' ) ( a14(i,j), j = 1, mxtyat )
|
---|
240 | enddo
|
---|
241 |
|
---|
242 | do i=1,mxhbdo
|
---|
243 | do j=1,mxhbac
|
---|
244 | hbc(i,j)=chb_s(i,j)
|
---|
245 | hba(i,j)=ahb_s(i,j)
|
---|
246 | enddo
|
---|
247 | enddo
|
---|
248 |
|
---|
249 | do i=1,mxtyto
|
---|
250 | e0to(i)=e0to_s(i)/2.
|
---|
251 | sgto(i)=sgto_s(i)
|
---|
252 | rnto(i)=rnto_s(i)
|
---|
253 | enddo
|
---|
254 |
|
---|
255 |
|
---|
256 | endif
|
---|
257 | ! -------------------------------------------- Hydrogen Bond Parameters
|
---|
258 | do i=1,mxtyat
|
---|
259 | do j=1,mxtyat
|
---|
260 | ihbty(i,j)=0
|
---|
261 | chb(i,j)=0.
|
---|
262 | ahb(i,j)=0.
|
---|
263 | enddo
|
---|
264 | enddo
|
---|
265 |
|
---|
266 | iac=0
|
---|
267 | ido=0
|
---|
268 | do i=1,mxtyat
|
---|
269 | if (do(i)) then
|
---|
270 | ido=ido+1
|
---|
271 | jac=0
|
---|
272 | do j=1,i-1
|
---|
273 | if (ac(j)) then
|
---|
274 | jac=jac+1
|
---|
275 | ihbty(i,j)=1
|
---|
276 | ihbty(j,i)=-1
|
---|
277 | chb(i,j)=hbc(ido,jac)
|
---|
278 | chb(j,i)=chb(i,j)
|
---|
279 | ahb(i,j)=hba(ido,jac)
|
---|
280 | ahb(j,i)=ahb(i,j)
|
---|
281 | endif
|
---|
282 | enddo
|
---|
283 | elseif (ac(i)) then
|
---|
284 | iac=iac+1
|
---|
285 | jdo=0
|
---|
286 | do j=1,i-1
|
---|
287 | if (do(j)) then
|
---|
288 | jdo=jdo+1
|
---|
289 | ihbty(i,j)=-1
|
---|
290 | ihbty(j,i)=1
|
---|
291 | chb(i,j)=hbc(jdo,iac)
|
---|
292 | chb(j,i)=chb(i,j)
|
---|
293 | ahb(i,j)=hba(jdo,iac)
|
---|
294 | ahb(j,i)=ahb(i,j)
|
---|
295 | endif
|
---|
296 | enddo
|
---|
297 | endif
|
---|
298 | enddo
|
---|
299 |
|
---|
300 | do i=1,mxtyto ! eases calculation of torsional derivatives
|
---|
301 | esnto(i)=e0to(i)*sgto(i)*rnto(i)
|
---|
302 | enddo
|
---|
303 |
|
---|
304 | return
|
---|
305 | end
|
---|
306 | ! **************
|
---|
307 | BLOCK DATA
|
---|
308 |
|
---|
309 | include 'INCL.H'
|
---|
310 |
|
---|
311 | ! Atom types ------------------------------------------------------------
|
---|
312 | ! Original types -Scheraga: -Flex:
|
---|
313 | ! H 1 - with aliphatic carbon 1 12
|
---|
314 | ! 2 - with aromatic carbon 3 13
|
---|
315 | ! 3 - with non-sp3 types of nitrogen 2 1
|
---|
316 | ! 4 - with sp3-hybr. nitrogen 2 2
|
---|
317 | ! 5 - with oxygen 4 1
|
---|
318 | ! 6 - with sulfur 3(was 5)1
|
---|
319 | ! C 7 - sp3-hybr. carbon 6,9 3
|
---|
320 | ! 8 - sp2-carbon (carbonyl,carboxyl,carboxylate) 7,11 4
|
---|
321 | ! 9 - aromatic carbon 8,10 4
|
---|
322 | ! O 10 - hydroxyl, ester oxygen (inc. water) 18,19 8
|
---|
323 | ! 11 - carbonyl oxygen 17 9
|
---|
324 | ! 12 - carboxylate oxygen 18,19 10
|
---|
325 | ! N 13 - aliph. nitrogen with 0/1 hydrogen & charged N 13-15 6
|
---|
326 | ! 14 - nitrogen with two hydrogens 13-15 5
|
---|
327 | ! 15 - all other nitrogens (+ sp2-hybrid. in heteroc.) 13-15 7
|
---|
328 | ! S 16 - any sulfur 20,21 18,19
|
---|
329 | ! H 17 - H-delta of Pro, Hyp of ECEPP/3 dataset 5(new) -
|
---|
330 | ! C 18 - C-delta of Pro, Hyp of ECEPP/3 dataset 12(new) -
|
---|
331 |
|
---|
332 | ! Classes for torsional potential ---------------------------------------
|
---|
333 | !
|
---|
334 | ! 1 : 'Omega' = C'(pept.)-N(pept.) [Cpept-Npept]
|
---|
335 | ! 2 : 'Phi' = N(pept.)-C(sp3) [C4-Npept]
|
---|
336 | ! 3 : 'Psi' = C(sp3)-C'(pept.) [C4-Cpept]
|
---|
337 | ! 4 : 'Chi1' = C(sp3)-C(sp3) [C4-C4]
|
---|
338 | ! 5 : C(sp3)-OH (Hydroxyl) [C4-OH]
|
---|
339 | ! 6 : C(sp3)-NH2 [C4-NH2]
|
---|
340 | ! 7 : C(sp3)-NH3+ [C4-NH3+]
|
---|
341 | ! 8 : C(sp3)-NH-(guanidyl) [C4-NHX]
|
---|
342 | ! 9 : C(sp3)-COOH(carboxyl) [C4-COO]
|
---|
343 | ! 10 : C(sp3)-COO-(carboxylate) [C4-COO]
|
---|
344 | ! 11 : C(sp3)-CO(sp2 of amide) [C4-Cpept]
|
---|
345 | ! 12 : C(sp3)-C(aromatic ring) [C4-C3]
|
---|
346 | ! 13 : C(sp3)-S [C4-SC4]
|
---|
347 | ! 14 : C(sp3)-SH [C4-SH]
|
---|
348 | ! 15 : C(aromatic ring)-OH [C3-OH]
|
---|
349 | ! ________________________________________________ "rigid" torsions:
|
---|
350 | ! 16 : C(carboxyl)-OH [C3-OH]
|
---|
351 | ! 17 : -NH-C(sp2 of guanidyl) [C3-NHX]
|
---|
352 | ! 18 : -C(sp3)-NH2 (guanidyl) [not in Flex]
|
---|
353 | ! 19 : -C(sp3)-NH2 (amide) [Cpept-Npept]
|
---|
354 |
|
---|
355 | data conv/332.d0/ ! to convert electrost. energy into [kcal/mole]
|
---|
356 |
|
---|
357 | ! ------------------------- ECEPP/3 potential --------------------------------
|
---|
358 | ! 1) Momany F.A McGuire R.F Burgess A.W Scheraga H.A J Phys Chem v79 2361-2381
|
---|
359 | ! 1975
|
---|
360 | ! 2) Nemethy G Pottle M.S Scheraga H.A, J Phys Chem v87 1883-1887 1983
|
---|
361 | ! 3) Sippl M.J Nemethy G Scheraga H.A J Phys Chem v88 6231-6233 1984
|
---|
362 | ! 4) Nemethy G Gibson K.D Palmer K.A Yoon C.N Paterlini G Zagari A Rumsey S
|
---|
363 | ! Scheraga H.A J Phys Chem v96 6472-6484 1992
|
---|
364 | ! ----------------------------------------------------------------------------
|
---|
365 |
|
---|
366 | data eps_s/2.d0/ ! Distance-INdependent diel. constant
|
---|
367 | ! data eps_s/6.d0/ ! Distance-INdependent diel. constant
|
---|
368 | data plt/78.d0/, slp/0.3d0/ ! Parameters for Epsilon(R)
|
---|
369 |
|
---|
370 | data ehm /362.55d0/ ! Angstrom**2/3 * kcal / mol ! from KONF90
|
---|
371 | ! data ehm /362.09561409d0/ ! Angstrom**2/3 * kcal / mol
|
---|
372 | ! From:
|
---|
373 | ! 1.5
|
---|
374 | ! * elementary charge = 4.80325 *e+2 Angstrom**3/2 * g**1/2 * s**(-1)
|
---|
375 | ! * Planck's constant/2*Pi = 1.0545887 *e-34 Joule * s
|
---|
376 | ! * Avogadro's number = 6.022045 *e+23 mol**(-1)
|
---|
377 | ! / sqrt (mass of electron) = sqrt (9.109534 *e-28 g )
|
---|
378 | ! / thermal equivalent = 4.1868 *e+3 Joule * kcal**(-1)
|
---|
379 | ! data ehm /362.36d0/ ! calculated using Tab II in ref. 2
|
---|
380 | ! data 1/ehm /2.757670d-3/ ! 3*sqrt(m)/(2*e*h) taken from ICM
|
---|
381 |
|
---|
382 | ! ---------------------- atomic polarizabilties (*100,[Angstrom**3])
|
---|
383 | ! 1 2 3 4 5 6 7 8 9 10 11 12
|
---|
384 | data atpl/42.,42.,42.,42.,42.,42.,93.,151.,115.,59.,84.,59.,
|
---|
385 | ! 13 14 15 16 17 18
|
---|
386 | & 93.,93.,93.,220.,42.,93./
|
---|
387 | ! ---------------------- effective numbers of electrons (*100,ref. 2)
|
---|
388 | ! 1 2 3 4 5 6 7 8 9 10 11 12
|
---|
389 | data efel/85.,85.,85.,85.,85.,85.,520.,520.,520.,700.,700.,700.,
|
---|
390 | ! 13 14 15 16 17 18
|
---|
391 | & 610.,610.,610.,1480.,85.,520./
|
---|
392 | ! ------------------------- min. pairwise 6-12 energy (*1000,[kcal/mol])
|
---|
393 | ! 1 2 3 4 5 6 7 8 9 10 11 12
|
---|
394 | data emin/37.,36.,61.,61.,44.,36.,38.,140.,99.,94.,200.,94.,
|
---|
395 | ! 13 14 15 16 17 18
|
---|
396 | & 107.,107.,107.,223.,99.,38./
|
---|
397 | ! ---------------------------- opt. pairwise distance (*100,[Angstrom])
|
---|
398 | ! 1 2 3 4 5 6 7 8 9 10 11
|
---|
399 | data rmin/292.,293.,268.,268.,283.,293.,412.,374.,370.,324.,312.,
|
---|
400 | ! 12 13 14 15 16 17 18
|
---|
401 | & 324.,351.,351.,351.,415.,248.,412./
|
---|
402 | ! ---------------------------------------------- Hydrogen-bond donors
|
---|
403 | ! 1 2 3 4 5 6
|
---|
404 | data do_s/.false.,.false.,.true.,.true.,.true.,.false.,
|
---|
405 | ! 7 8 9 10 11 12
|
---|
406 | & .false.,.false.,.false.,.false.,.false.,.false.,
|
---|
407 | ! 13 14 15 16 17 18
|
---|
408 | & .false.,.false.,.false.,.false.,.false.,.false./
|
---|
409 | ! -------------------------------------------- Hydrogen-bond acceptors
|
---|
410 | ! 1 2 3 4 5 6
|
---|
411 | data ac_s/.false.,.false.,.false.,.false.,.false.,.false.,
|
---|
412 | ! 7 8 9 10 11 12
|
---|
413 | & .false.,.false.,.false.,.true.,.true.,.true.,
|
---|
414 | ! 13 14 15 16 17 18
|
---|
415 | & .true.,.true.,.true.,.false.,.false.,.false./
|
---|
416 | ! & .false.,.true.,.true.,.false.,.false.,.false./ !! ICM
|
---|
417 | ! --------------------------------- HB-parameters (/1000,attraction)
|
---|
418 | data chb_s/2624.,2624.,4610.,.0, ! given as:
|
---|
419 | & 4014.,4014.,5783.,.0, ! (ac_typ x do_typ)
|
---|
420 | & 2624.,2624.,4610.,.0, ! to be used:
|
---|
421 | & 8244.,8244.,8244.,.0, ! (DO_typ x AC_typ)
|
---|
422 | & 8244.,8244.,8244.,.0, ! i.e.:
|
---|
423 | & 8244.,8244.,8244.,.0/ ! ( 3-5 x 10-15 )
|
---|
424 | ! --------------------------------- HB-parameters (/1000,repulsion)
|
---|
425 | data ahb_s/ 5890., 5890.,11220.,.0,
|
---|
426 | & 12040.,12040.,16583.,.0, ! 13344 -> 16583 = Ref. 3
|
---|
427 | & 5890., 5890.,11220.,.0,
|
---|
428 | & 32897.,32897.,32897.,.0,
|
---|
429 | & 32897.,32897.,32897.,.0,
|
---|
430 | & 32897.,32897.,32897.,.0/
|
---|
431 |
|
---|
432 | ! 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15
|
---|
433 | data e0to_s /20.,0.,0.,2.7,.6,1.8,1.8,0.,0.,0.,0.,0.,2.,1.5,3.5
|
---|
434 | ! 16 17 18 19
|
---|
435 | & ,8.,18.,20.,15./
|
---|
436 | data sgto_s /-1.,0.,0., 1.,1., 1., 1.,0.,0.,0.,0.,0.,1., 1.,-1.
|
---|
437 | & ,-1.,-1.,-1.,-1./
|
---|
438 | data rnto_s / 2.,0.,0., 3.,3., 3., 3.,0.,0.,0.,0.,0.,3., 3., 2.
|
---|
439 | & ,2., 2., 2., 2./
|
---|
440 |
|
---|
441 | ! ---------------------------- Flex potential ----------------------------------
|
---|
442 | ! Lavery R Sklenar H Zakrzewska K Pullman B J Biomol Struct Dyn v3 989-1014 1986
|
---|
443 | ! VdW-parameters from: Zhurkin V.B Poltiev V.I Florent'ev V.L Molekulyarnaya
|
---|
444 | ! Biologiya v14 116 1980
|
---|
445 | ! ------------------------------------------------------------------------------
|
---|
446 |
|
---|
447 | data plt_f/78.d0/, slp_f/0.16d0/ ! Parameters for Epsilon(R)
|
---|
448 | data cohb_f/6.d0/ ! Cut-off distance betw. H- & acceptor atom for HB
|
---|
449 | data c_f/ ! ----------- Lennard-Jones C6-parameters (attraction)
|
---|
450 | &40.,40.,40.,40.,40.,40.,100.,126.,126.,86.,121.,105.,126.,105.,
|
---|
451 | &146.,213.1,3*0.,40.,40.,40.,40.,40.,100.,126.,126.,86.,121.,105.,
|
---|
452 | &126.,105.,146.,213.1,4*0.,40.,40.,40.,40.,100.,126.,126.,86.,121.,
|
---|
453 | &105.,126.,105.,146.,213.1,5*0.,40.,40.,40.,100.,126.,126.,86.,121.
|
---|
454 | &,105.,126.,105.,146.,213.1,6*0.,40.,40.,100.,126.,126.,86.,121.,
|
---|
455 | &105.,126.,105.,146.,213.1,7*0.,40.,100.,126.,126.,86.,121.,105.,
|
---|
456 | &126.,105.,146.,213.1,8*0.,250.,316.,316.,217.,305.,264.,316.,264.,
|
---|
457 | &367.,489.,9*0.,400.,400.,274.,385.,334.,400.,334.,464.,537.4,10*0.
|
---|
458 | &,400.,274.,385.,334.,400.,334.,464.,537.4,11*0.,200.,283.,245.,
|
---|
459 | &278.,233.,330.,424.,12*0.,400.,347.,391.,327.,465.,583.,13*0.,300.
|
---|
460 | &,339.,284.,403.,530.,14*0.,400.,334.,467.,556.5,15*0.,280.,391.,
|
---|
461 | &484.,16*0.,550.,673.4,17*0.,246.,38*0./
|
---|
462 | data a_f/ ! ---- Lennard-Jones A12-parameters (/1000,repulsion)
|
---|
463 | &7.74,7.74,7.74,7.74,7.74,7.74,70.6,81.6,81.6,31.3,42.2,36.6,71.4,
|
---|
464 | &62.,78.3,189.3,3*0.,7.74,7.74,7.74,7.74,7.74,70.6,61.7,61.7,31.3,
|
---|
465 | &17.8,15.4,53.7,62.,58.7,189.3,4*0.,7.74,7.74,7.74,7.74,70.6,81.6,
|
---|
466 | &81.6,31.3,42.2,36.6,71.4,62.,78.3,189.3,5*0.,7.74,7.74,7.74,70.6,
|
---|
467 | &61.7,61.7,31.3,17.8,15.4,53.7,62.,58.7,189.3,6*0.,7.74,7.74,70.6,
|
---|
468 | &81.6,81.6,31.3,42.2,36.6,71.4,62.,78.3,189.3,7*0.,7.74,70.6,81.6,
|
---|
469 | &81.6,31.3,42.2,36.6,71.4,62.,78.3,189.3,8*0.,512.,601.,601.,256.,
|
---|
470 | &349.,302.,538.,464.,598.,1196.,9*0.,704.,704.,298.,406.,351.,630.,
|
---|
471 | &544.,699.,1203.5,10*0.,704.,298.,406.,351.,630.,544.,699.,1203.5,
|
---|
472 | &11*0.,129.,176.,153.,269.,233.,303.,561.8,12*0.,240.,208.,366.,
|
---|
473 | &317.,413.,772.5,13*0.,180.,317.,274.,358.,702.3,14*0.,565.,488.,
|
---|
474 | &629.,1105.8,15*0.,421.,544.,976.8,16*0.,705.,1259.5,17*0.,503.3,
|
---|
475 | &38*0./
|
---|
476 | ! ---------------------------------------------- Hydrogen-bond donors
|
---|
477 | ! 1 2 3 4 5 6
|
---|
478 | data do_f/.false.,.false.,.true.,.true.,.true.,.true.,
|
---|
479 | ! 7 8 9 10 11 12
|
---|
480 | & .false.,.false.,.false.,.false.,.false.,.false.,
|
---|
481 | ! 13 14 15 16 17 18
|
---|
482 | & .false.,.false.,.false.,.false.,.false.,.false./
|
---|
483 | ! -------------------------------------------- Hydrogen-bond acceptors
|
---|
484 | ! 1 2 3 4 5 6
|
---|
485 | data ac_f/.false.,.false.,.false.,.false.,.false.,.false.,
|
---|
486 | ! 7 8 9 10 11 12
|
---|
487 | & .false.,.false.,.false.,.true.,.true.,.true.,
|
---|
488 | ! 13 14 15 16 17 18
|
---|
489 | & .false.,.false.,.true.,.true.,.false.,.false./
|
---|
490 | ! --------------------------------- HB-parameters (/1000,attraction)
|
---|
491 | data chb_f/180.,180.,160. ,226.8, ! given as (ac_typ x do_typ)
|
---|
492 | & 175.,175.,150. ,305.8, ! to be used as:
|
---|
493 | & 100.,100., 85. , 85. , ! (Do_typ x Ac_typ)
|
---|
494 | & 165.,165.,150. ,305.8, ! i.e.:
|
---|
495 | & 720.,720.,643.1,845.6, ! ( 3-6 x 10-12,15,16 )
|
---|
496 | & 0., 0., 0. , 0. /
|
---|
497 | ! --------------------------------- HB-parameters (/1000,repulsion)
|
---|
498 | data ahb_f/ 6600., 6600., 6200., 12855.,
|
---|
499 | & 6400., 6400., 7100., 29216.,
|
---|
500 | & 2400., 2400., 2400., 2000.,
|
---|
501 | & 6200., 6200., 7100., 29216.,
|
---|
502 | & 185000.,185000.,172301.,308235.,
|
---|
503 | & 0., 0., 0., 0./
|
---|
504 |
|
---|
505 | ! 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15
|
---|
506 | data e0to_f /20.,2.,.8,2.6,.6,2.1,1.8,3.2,.5,.5,.8,1.2,1.3,1.,6.2
|
---|
507 | ! 16 17 18 19
|
---|
508 | & ,6.2, 8.,0.,20./
|
---|
509 | data sgto_f /-1.,1.,-1., 1.,1., 1., 1.,1.,1.,1.,-1.,1., 1.,1.,-1.
|
---|
510 | & ,-1.,-1.,0.,-1./
|
---|
511 | data rnto_f / 2.,3., 3., 3.,3., 6., 3.,3.,6.,6., 3.,6., 3.,3., 2.
|
---|
512 | & ,2., 2.,0.,2./
|
---|
513 |
|
---|
514 | data (rsnmcd(i),i=1,nrsty)/ ! Names for all amino acid residue types
|
---|
515 | & 'ala ','arg ','arg+','asn ','asp ','asp-','cys ','cyss','gln ',
|
---|
516 | & 'glu ','glu-','gly ','his ','hise','hisd','his+','hyp ','hypu',
|
---|
517 | & 'ile ','leu ','lys ','lys+','met ','phe ','cpro','pro ','cpru',
|
---|
518 | & 'prou','pron','pro+','ser ','thr ','trp ','tyr ','val ' /
|
---|
519 |
|
---|
520 | data (onltcd(i),i=1,nrsty)/ ! One-letter codes for amino acid types
|
---|
521 | & 'A', 'R', 'R', 'N', 'D', 'D', 'C', 'C', 'Q',
|
---|
522 | & 'E', 'E', 'G', 'H', 'H', 'H', 'H', 'P', 'P',
|
---|
523 | & 'I', 'L', 'K', 'K', 'M', 'F', 'P', 'P', 'P',
|
---|
524 | & 'P', 'P', 'P', 'S', 'T', 'W', 'Y', 'V' /
|
---|
525 |
|
---|
526 | ! The vdW radii (in Angstr.) for the atomic groups and
|
---|
527 | ! coefficients for their solvation free energy (kcal/molxA**2)
|
---|
528 |
|
---|
529 | ! Method:
|
---|
530 |
|
---|
531 | ! itysol=1 : OONS --> T.Ooi, et al,
|
---|
532 | ! Proc. Natl. Acad. Sci. USA 8 (1987) 3086-3090.
|
---|
533 | ! itysol=2 : JRF --> J.Vila, et al,
|
---|
534 | ! PROTEINS: Struct Funct Genet 10(1991) 199-218.
|
---|
535 | ! itysol=3 : WE92 --> L.Wesson, D.Eisenberg,
|
---|
536 | ! Protein Science 1 (1992) 227-235.
|
---|
537 | ! itysol=4 : SCH1 --> D.Eisenberg, et al,
|
---|
538 | ! Chem Scrip 29A (1989) 217-221.
|
---|
539 | ! itysol=5 : SCH2 --> A.H.Juffer, et al,
|
---|
540 | ! Proteine Science 4 (1995) 2499-2509.
|
---|
541 | ! itysol=6 : SCH3 --> L.Wesson, D.Eisenberg,
|
---|
542 | ! Protein Science 1 (1992) 227-235.
|
---|
543 | ! itysol=7 : SCH4 --> C.A. Schiffer, et al,
|
---|
544 | ! Mol. Simul. 10(1993) 121-149.
|
---|
545 | ! itysol=8 : EM86 --> D.Eisenberg, A.D. Mclachlan,
|
---|
546 | ! Nature 319 (1986) 199-203.
|
---|
547 | ! itysol=9 : BM --> B. Freyberg, et al,
|
---|
548 | ! J. Mol. Biol. 233 (1993) 275-292.
|
---|
549 |
|
---|
550 | ! ATOM
|
---|
551 | ! TYPE OONS JRF WE92 SCH1 SCH2 SCH3 SCH4 EM86 BM
|
---|
552 |
|
---|
553 | ! 1 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.000 0.000
|
---|
554 | ! 2 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.000 0.000
|
---|
555 | ! 3 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.000 0.000
|
---|
556 | ! 4 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.000 0.000
|
---|
557 | ! 5 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.000 0.000
|
---|
558 | ! 6 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.000 0.000
|
---|
559 | ! 7 0.0080 0.2160 0.0120 0.0180 0.0130 0.0040 0.0325 0.016 1.000
|
---|
560 | ! 8 0.4270 -0.7320 0.0120 0.0180 0.0130 0.0040 0.0325 0.016 1.000
|
---|
561 | ! 9 -0.0080 -0.6780 0.0120 0.0180 0.0130 0.0040 0.0325 0.016 1.000
|
---|
562 | !10 -0.1720 -0.9100 -0.1160 -0.0090 -0.0070 -0.1130 -0.0175 -0.006 0.000
|
---|
563 | !11 -0.0380 -0.2620 -0.1750 -0.0090 -0.0070 -0.1660 -0.2800 -0.006 0.000
|
---|
564 | !12 -0.0380 -0.9100 -0.1750 -0.0370 -0.1120 -0.1660 -0.2800 -0.024 0.000
|
---|
565 | !13 -0.1320 -0.3120 -0.1860 -0.0380 -0.0870 -0.1690 -0.2175 -0.05 0.000
|
---|
566 | !14 -0.1320 -0.3120 -0.1160 -0.0090 -0.0070 -0.1130 -0.0175 -0.006 0.000
|
---|
567 | !15 -0.1320 -0.3120 -0.1160 -0.0090 -0.0070 -0.1130 -0.0175 -0.006 0.000
|
---|
568 | !16 -0.0210 -0.2810 -0.0180 0.0050 -0.0036 -0.0170 -0.0090 0.021 0.000
|
---|
569 | !17 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.000 0.000
|
---|
570 | !18 0.0080 0.2160 0.0120 0.0180 0.0130 0.0040 0.0325 0.016 1.000
|
---|
571 |
|
---|
572 | data coef_sl/54*0.0,
|
---|
573 | & 0.008, 0.216, 0.012, 0.018, 0.013, 0.004, 0.0325, 0.016,1.000,
|
---|
574 | & 0.427,-0.732, 0.012, 0.018, 0.013, 0.004, 0.0325, 0.016,1.000,
|
---|
575 | & -.008,-0.678, 0.012, 0.018, 0.013, 0.004, 0.0325, 0.016,1.000,
|
---|
576 | & -.172,-0.910,-0.116,-0.009,-0.007,-0.113,-0.0175,-0.006,0.000,
|
---|
577 | & -.038,-0.262,-0.116,-0.009,-0.007,-0.113,-0.0175,-0.006,0.000,
|
---|
578 | & -.038,-0.910,-0.175,-0.037,-0.112,-0.166,-0.2800,-0.024,0.000,
|
---|
579 | & -.132,-0.312,-0.186,-0.038,-0.087,-0.169,-0.2175,-0.05, 0.000,
|
---|
580 | & -.132,-0.312,-0.116,-0.009,-0.007,-0.113,-0.0175,-0.006,0.000,
|
---|
581 | & -.132,-0.312,-0.116,-0.009,-0.007,-0.113,-0.0175,-0.006,0.000,
|
---|
582 | & -.021,-0.281,-0.018, 0.005,-.0036,-0.017,-0.0090, 0.021,0.000,
|
---|
583 | & 9*0.0,
|
---|
584 | & 0.008, 0.216, 0.012, 0.018, 0.013, 0.004, 0.0325, 0.016,1.000/
|
---|
585 |
|
---|
586 | data rad_vdw/54*0.,
|
---|
587 | & 2*2.0,6*1.9,2.0,
|
---|
588 | & 2*1.55,6*1.9,1.5,
|
---|
589 | & 2*1.75,6*1.9,1.85,
|
---|
590 | & 9*1.4,
|
---|
591 | & 9*1.4,
|
---|
592 | & 9*1.4,
|
---|
593 | & 2*1.55,6*1.7,1.5,
|
---|
594 | & 2*1.55,6*1.7,1.5,
|
---|
595 | & 2*1.55,6*1.7,1.5,
|
---|
596 | & 2*2.0,6*1.8,1.85,
|
---|
597 | & 9*0.,
|
---|
598 | & 2*2.0,6*1.9,2.0/
|
---|
599 |
|
---|
600 | end
|
---|
601 |
|
---|