source: getmol.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: 13.5 KB
Line 
1! **************************
2! **************************************************************
3!
4! This file contains the subroutines: getmol,redres
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 getmol(nml)
15
16! ...................................................................
17! PURPOSE: assemble data for molecule 'nml' according to
18! its sequence using residue library 'reslib'
19!
20! ! Molecules must be assembled in sequential order (1 -> ntlml)
21! (or number of atoms & variables must remain the same)
22!
23! INPUT: irsml1(nml),irsml2(nml),seq(irsml1()...irsml2())
24! nml>1: irsml1(nml-1),iatrs2(irsml2(nml-1))
25! ivrrs1(irsml2(nml-1)),nvrrs(irsml2(nml-1))
26!
27! OUTPUT: molecule - ivrml1,nvrml
28! residues - iatrs1,ixatrs,iatrs2,ivrrs1,nvrrs
29! atoms - nmat,ityat,cgat,blat,baat,csbaat,snbaat,
30! toat,cstoat,sntoat
31! bonds - nbdat,iowat,iyowat,ibdat(1-mxbd,),iybdat(1-mxbd,)
32! ! 1st atom of 'nml': iowat indicates 1st bond
33! to a FOLLOWING atom (not previous) !
34! variables - ityvr,iclvr,iatvr,nmvr
35!
36! CALLS: iopfil,redres,iendst
37! ...................................................................
38
39 include 'INCL.H'
40! arguments
41 integer nml
42! functions
43 integer iopfil, iendst
44
45 integer i, iat, iow, ifirs, ibd, ilars, j, mvr, nbd, n, nbd1, nj
46 integer nh, nrs, ntlat, ntlvr, nat, nxt
47 double precision ba, t, to
48 character res*4
49
50 if (iopfil(lunlib,reslib,'old','formatted').le.izero) then
51 write (logString, '(a,/,a,i3,2a)')
52 & ' getmol> ERROR opening library of residues:',
53 & ' LUN=',lunlib,' FILE=',reslib(1:iendst(reslib))
54 stop
55 endif
56
57 if (nml.eq.1) then
58 ntlat=0
59 ntlvr=0
60 else
61 i=irsml2(nml-1) ! last res. of previous mol.
62 ntlat=iatrs2(i)
63 ntlvr=ivrrs1(i)+nvrrs(i)-1
64 endif
65 ivrml1(nml)=ntlvr + 1
66
67 ilars=irsml2(nml)
68 ifirs=irsml1(nml)
69
70 do nrs=ifirs,ilars ! Residues in molecule
71
72 res = seq(nrs)
73 call tolost(res) ! ensure lower case for residue name
74
75 if (res(:3).eq.'nme'.and.nrs.ne.ilars) then
76 write (logString, '(3a)') ' getmol> residue >',res,
77 & '< allowed at C-terminus only !'
78 close(lunlib)
79 stop
80 elseif (res(:3).eq.'ace'.and.nrs.ne.ifirs) then
81 write (logString, '(3a)') ' getmol> residue >',res,
82 & '< allowed at N-terminus only !'
83 close(lunlib)
84 stop
85 endif
86
87 call redres(res,nat,nxt,nvr)
88
89 if ((nat+ntlat).gt.mxat) then
90 write (logString, '(a,i5)') ' getmol> number of atoms > ',mxat
91 close(lunlib)
92 stop
93 endif
94 if ((nvr+ntlvr).gt.mxvr) then
95 write (logString, '(a,i5)') ' getmol> number of variables > '
96 & ,mxvr
97 close(lunlib)
98 stop
99 endif
100
101 rewind lunlib
102
103! ___________________________________________________________ Atoms
104 do i=1,nat
105 n=i+ntlat
106 nmat(n)=nmath(i)
107 ityat(n)=ityath(i)
108 cgat(n)=cgath(i)
109 blat(n)=blath(i)
110 ba=baath(i)
111 baat(n)=ba
112 csbaat(n)=cos(ba)
113 snbaat(n)=sin(ba)
114 to=toath(i)
115 toat(n)=to
116 cstoat(n)=cos(to)
117 sntoat(n)=sin(to)
118! ______________________________ bonds to previous & following atoms
119 iow=iowath(i)
120 if (iow.eq.0) then ! 1st atom of residue
121 if (nrs.eq.ifirs) then ! 1st atom of 'nml'
122
123 iowat(n)=ibdath(1,i)+ntlat
124 iyowat(n)=iybdath(1,i)
125
126 do j = 1,mxbd-1
127 ibdath(j,i)=ibdath(j+1,i)
128 iybdath(j,i)=iybdath(j+1,i)
129 enddo
130 ibdath(mxbd,i)=0
131 iybdath(mxbd,i)=1
132
133 nbdath(i)=nbdath(i)-1
134 else ! connected with prev. res.
135 nh=ixatrs(nrs-1) ! atom to 'next' residue of prev.res.
136 iowat(n)=nh
137 iyowat(n)=1 !!! only single bonds assumed !!!
138
139! ___________________________ correct atom to 'next' res.
140 nbd=nbdat(nh)
141 if (nbd.eq.mxbd) then
142 write (logString, '(a,i2,a,i4,2a,i4,a)')
143 & ' getmol> need ',(mxbd+2),
144 & 'th bond to connect residues ',
145 & nrs-1,seq(nrs-1),' and ',nrs,seq(nrs)
146 close(lunlib)
147 stop
148 else ! correct atom to 'next' res.
149! _______________________________!! dihedrals for atoms bound to 'nh'
150! are assumed to be phase angles !!
151 do j=1,nbd
152
153 nj=ibdat(j,nh)
154 t=toat(nj)
155
156 if (t.eq.0.0) then
157 write (logString, '(3a,/,2a)')
158 & ' getmol> DIHEDRAL for atom ',nmat(nj),
159 & ' should be PHASE angle with respect to atom ',
160 & nmat(n),' & therefore must be not 0.0 !!'
161 close(lunlib)
162 stop
163 endif
164
165 t=t+to
166 if (abs(t).gt.pi) t=t-sign(pi2,t)
167 toat(nj)=t
168 cstoat(nj)=cos(t)
169 sntoat(nj)=sin(t)
170
171 enddo
172
173 nbd1=nbd+1
174 ibdat(nbd1,nh)=n
175 iybdat(nbd1,nh)=1 ! (only single bonds !)
176 nbdat(nh)=nbd1
177
178 endif
179
180 endif
181
182 else ! connected within res.
183
184 iowat(n)=ntlat+iow
185 iyowat(n)=iyowath(i)
186
187 endif
188
189 nbdat(n)=nbdath(i)
190
191 do j=1,mxbd
192
193 ibd=ibdath(j,i)
194
195 if (ibd.ne.0) then
196 ibdat(j,n)=ibd+ntlat
197 iybdat(j,n)=iybdath(j,i)
198 else
199 ibdat(j,n)=0
200 iybdat(j,n)=1
201 endif
202
203 enddo
204
205 enddo ! ... atoms
206
207! ________________________________________________________ Variables
208 ivrrs1(nrs)=ntlvr+1
209 mvr=0
210
211 do i=1,nvr
212
213 if (nrs.eq.ifirs) then
214
215 iat=iatvrh(i)
216! ____________________________________ Exclude all variables for 1st atom
217! & torsion for atoms bound to it
218 if ( iat.eq.1.or.
219 & (iowath(iat).eq.1.and.ityvrh(i).eq.3)) goto 1
220
221 endif
222
223 mvr=mvr+1
224 ntlvr=ntlvr+1
225 ityvr(ntlvr)=ityvrh(i)
226 iclvr(ntlvr)=iclvrh(i)
227 iatvr(ntlvr)=iatvrh(i)+ntlat
228 nmvr(ntlvr)=nmvrh(i)
229
230 1 enddo ! ... Variables
231
232 nvrrs(nrs)=mvr
233
234 iatrs1(nrs)=ntlat+1 ! first backbone atom of res.
235 ixatrs(nrs)=ntlat+nxt ! last backbone atom
236 ntlat=ntlat+nat
237 iatrs2(nrs)=ntlat ! last atom of res.
238
239 enddo ! ... residues
240
241 close(lunlib)
242
243! _______________________________ Variables
244 if (nml.eq.1) then
245 nvrml(nml)=ntlvr
246 else
247 nvrml(nml)=ntlvr-ivrml1(nml) + 1
248 endif
249
250 return
251 end
252! **************************************
253 subroutine redres(res,nat,nxt,nvrr)
254
255! .......................................................
256! PURPOSE: read atom data for residue 'res' from library
257! (file 'lunlib' 'reslib' opened in routine calling
258! this one)
259!
260! OUTPUT: nat - number of atoms in residue
261! nxt - atom which may bind to following residue
262! nvrr - number of variables in residue
263! for atoms - nmath,blath,baath(rad),toath(rad),
264! ityath,iyowath,iowath (INSIDE residue,
265! =0 if 1st atom)
266! for variables - ityvrh (1=bl/2=ba/3=to),iclvrh,iatvrh,nmvrh
267!
268! LIBRARY: residue-lines:
269! '#', res, nat, nxt; Format: a1,a4,2i4
270! atom-lines:
271! nmat,3{"fix" =' ', clvr,nmvr, blat/baat(deg)/toat(deg)},
272! cgat, ityat, iowat,ibdat1,ibdat2,ibdat3;
273! Format: a4, 3(1x,i2,a1,a3,f9.3), f7.4, i4,4i4
274!
275! CALLS: iendst,tolost
276!
277! .......................................................
278
279 include 'INCL.H'
280! arguments
281 integer nat, nxt, nvrr
282
283! functions
284 integer iendst
285
286 integer icl, ibd, ib3, i, ib1, ib2, ic, iow, iexcp, ity, jow, lg
287 integer nln, j
288 dimension icl(3),ibd(mxbd)
289 character blnk,fix(3),nm(3)*3,res*4,resl*4,line*132
290 data blnk/' '/
291
292 double precision ba, to
293
294 nln=0
295 do i=1,mxbd
296 ibd(i)=0
297 enddo
298
299 resl=res
300
301 call tolost(resl) ! ensure lower case for residue name
302
303! ________________________________ find residue 'resl'
304 1 line=blnk
305 nln=nln+1
306
307 read (lunlib,'(a)',end=2,err=3) line
308 lg=iendst(line)
309
310 if (lg.ge.13.and.line(1:1).eq.'#'.and.line(2:5).eq.resl) then
311! _____________________________________________ read atom data for 'resl'
312 read (line(6:13),'(2i4)',err=3) nat,nxt
313
314 if (nat.gt.mxath) then
315 write (logString, '(a,i5)') ' redres> number of atoms > ',mxath
316 close(lunlib)
317 stop
318 endif
319
320 nvrr=0
321 do i=1,nat
322
323 nln=nln+1
324
325 read (lunlib,'(a4,3(1x,i2,a1,a3,d9.3),d7.4,i4,4i4)',
326 & end=3,err=3)
327 & nmath(i),icl(1),fix(1),nm(1),blath(i),icl(2),fix(2),nm(2),ba,
328 & icl(3),fix(3),nm(3),to,cgath(i),ity,iow,(ibd(j),j=1,mxbd)
329
330 if (ity.le.0.or.ity.gt.mxtyat) goto 6
331 ityath(i)=ity
332
333 jow=abs(iow)
334
335 if (res(:3).eq.'ace'.and.i.eq.1) then ! exception from following check
336 iexcp = 1
337 else
338 iexcp = 0
339 endif
340
341 if (iexcp.eq.0.and.i.le.jow) then
342 if (i.eq.jow) then
343 write(logString,'(5a)') ' redres> atom ',nmath(i),' of ',
344 & resl,' cannot preceed itself '
345 else
346 write (logString, '(5a,i4)')
347 & ' redres> atom ',nmath(i),' of ',
348 & resl,' should be placed AFTER atom #',jow
349 endif
350 goto 5
351 endif
352
353 iowath(i)=jow
354 iyowath(i)=sign(1,iow)
355! ____________________________________ check order & find number of bonds
356! (bonds closing ring must be last !)
357 ib1=abs(ibd(1))
358 ib2=abs(ibd(2))
359 ib3=abs(ibd(3))
360
361 if (ib1.eq.i.or.ib2.eq.i.or.ib3.eq.i) goto 4
362
363 if (ib1.eq.0) then ! no bond to following
364 if (ib2.ne.0.or.ib3.ne.0) goto 4
365 nbdath(i)=0
366 else
367 if (ib1.eq.jow) goto 4
368 if (ib2.eq.0) then
369 if (ib3.ne.0) goto 4
370 nbdath(i)=1
371 else
372 if ( ib2.eq.jow.or.ib2.eq.ib1.or.
373 & (ib2.gt.i.and.ib2.lt.ib1) ) goto 4
374 if (ib3.eq.0) then
375 nbdath(i)=2
376 else
377 if (ib3.eq.jow.or.ib3.eq.ib1.or.ib3.eq.ib2.or.
378 & (ib3.gt.i.and.(ib3.lt.ib1.or.ib3.lt.ib2)) ) goto 4
379 nbdath(i)=3
380 endif
381 endif
382 endif
383
384 do j=1,mxbd
385 ibdath(j,i)=abs(ibd(j))
386 iybdath(j,i)=sign(1,ibd(j))
387 enddo
388
389 baath(i)=ba*cdr ! convert angles into 'radians'
390 toath(i)=to*cdr
391
392! ______________________________ internal degrees of freedom
393 do j=1,3
394 if (fix(j).ne.blnk) then
395 nvrr=nvrr+1
396
397 if (nvrr.gt.mxvrh) then
398 write (logString, '(a,i5)')
399 & ' redres> number of variables > ',
400 & mxvrh
401 close(lunlib)
402 stop
403 endif
404
405 ic=icl(j)
406
407 if ( ic.le.0
408 & .or.(j.eq.3.and.ic.gt.mxtyto) ! dihedral
409 & .or.(j.eq.2.and.ic.gt.mxtyba) ! bond angle
410 & .or.(j.eq.1.and.ic.gt.mxtybl) ) goto 7 ! b. length
411
412 ityvrh(nvrr)=j
413 iclvrh(nvrr)=ic
414 iatvrh(nvrr)=i
415 nmvrh(nvrr)=nm(j)
416
417 endif
418
419 enddo
420 enddo ! ... atoms
421
422 return
423 endif
424
425 goto 1
426
427! ____________________________________________________________ ERRORS
428 2 write (logString, '(4a)') ' redres> residue >',resl,
429 & '< NOT FOUND in ',
430 &reslib(1:iendst(reslib))
431 close(lunlib)
432 stop
433
434 3 write (logString, '(a,i4,2a)') ' redres> ERROR reading line No. ',
435 & nln,
436 & ' in ',reslib(1:iendst(reslib))
437 close(lunlib)
438 stop
439
440 4 write (logString, '(4a)')
441 & ' redres> Incorrect order of bonds for atom ',
442 & nmath(i),' of ',resl
443
444 5 write (logString, '(8x,2a)') '... must correct ',
445 & reslib(1:iendst(reslib))
446 close(lunlib)
447 stop
448
449 6 write (logString, '(a,i2,4a)') ' redres> unknown type :',ity,
450 & ': for atom ',nmath(i),' in residue ',resl
451 close(lunlib)
452 stop
453
454 7 write (logString, '(a,i2,4a)') ' redres> unknown class :',ic,
455 & ': for variable ',nm(j),' in residue ',resl
456 close(lunlib)
457 stop
458
459 end
460
Note: See TracBrowser for help on using the repository browser.