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