MPICH2 and energy fixes in WHAM-M and Cluster-M.
[unres.git] / source / wham / src-M / energy_p_new.F
1       subroutine etotal(energia,fact)
2       implicit real*8 (a-h,o-z)
3       include 'DIMENSIONS'
4       include 'DIMENSIONS.ZSCOPT'
5
6 #ifndef ISNAN
7       external proc_proc
8 #endif
9 #ifdef WINPGI
10 cMS$ATTRIBUTES C ::  proc_proc
11 #endif
12
13       include 'COMMON.IOUNITS'
14       double precision energia(0:max_ene),energia1(0:max_ene+1)
15 #ifdef MPL
16       include 'COMMON.INFO'
17       external d_vadd
18       integer ready
19 #endif
20       include 'COMMON.FFIELD'
21       include 'COMMON.DERIV'
22       include 'COMMON.INTERACT'
23       include 'COMMON.SBRIDGE'
24       include 'COMMON.CHAIN'
25       double precision fact(6)
26 cd      write(iout, '(a,i2)')'Calling etotal ipot=',ipot
27 cd    print *,'nnt=',nnt,' nct=',nct
28 C
29 C Compute the side-chain and electrostatic interaction energy
30 C
31       goto (101,102,103,104,105) ipot
32 C Lennard-Jones potential.
33   101 call elj(evdw,evdw_t)
34 cd    print '(a)','Exit ELJ'
35       goto 106
36 C Lennard-Jones-Kihara potential (shifted).
37   102 call eljk(evdw,evdw_t)
38       goto 106
39 C Berne-Pechukas potential (dilated LJ, angular dependence).
40   103 call ebp(evdw,evdw_t)
41       goto 106
42 C Gay-Berne potential (shifted LJ, angular dependence).
43   104 call egb(evdw,evdw_t)
44       goto 106
45 C Gay-Berne-Vorobjev potential (shifted LJ, angular dependence).
46   105 call egbv(evdw,evdw_t)
47 C
48 C Calculate electrostatic (H-bonding) energy of the main chain.
49 C
50   106 call eelec(ees,evdw1,eel_loc,eello_turn3,eello_turn4)
51 C
52 C Calculate excluded-volume interaction energy between peptide groups
53 C and side chains.
54 C
55       call escp(evdw2,evdw2_14)
56 c
57 c Calculate the bond-stretching energy
58 c
59       call ebond(estr)
60 c      write (iout,*) "estr",estr
61
62 C Calculate the disulfide-bridge and other energy and the contributions
63 C from other distance constraints.
64 cd    print *,'Calling EHPB'
65       call edis(ehpb)
66 cd    print *,'EHPB exitted succesfully.'
67 C
68 C Calculate the virtual-bond-angle energy.
69 C
70       call ebend(ebe)
71 cd    print *,'Bend energy finished.'
72 C
73 C Calculate the SC local energy.
74 C
75       call esc(escloc)
76 cd    print *,'SCLOC energy finished.'
77 C
78 C Calculate the virtual-bond torsional energy.
79 C
80 cd    print *,'nterm=',nterm
81       call etor(etors,edihcnstr,fact(1))
82 C
83 C 6/23/01 Calculate double-torsional energy
84 C
85       call etor_d(etors_d,fact(2))
86 C
87 C 21/5/07 Calculate local sicdechain correlation energy
88 C
89       call eback_sc_corr(esccor)
90
91 C 12/1/95 Multi-body terms
92 C
93       n_corr=0
94       n_corr1=0
95       if (wcorr4.gt.0.0d0 .or. wcorr5.gt.0.0d0 .or. wcorr6.gt.0.0d0 
96      &    .or. wturn6.gt.0.0d0) then
97 c         print *,"calling multibody_eello"
98          call multibody_eello(ecorr,ecorr5,ecorr6,eturn6,n_corr,n_corr1)
99 c         write (*,*) 'n_corr=',n_corr,' n_corr1=',n_corr1
100 c         print *,ecorr,ecorr5,ecorr6,eturn6
101       endif
102       if (wcorr4.eq.0.0d0 .and. wcorr.gt.0.0d0) then
103          call multibody_hb(ecorr,ecorr5,ecorr6,n_corr,n_corr1)
104       endif
105 c      write (iout,*) "ft(6)",fact(6)," evdw",evdw," evdw_t",evdw_t
106 #ifdef SPLITELE
107       etot=wsc*(evdw+fact(6)*evdw_t)+wscp*evdw2+welec*fact(1)*ees
108      & +wvdwpp*evdw1
109      & +wang*ebe+wtor*fact(1)*etors+wscloc*escloc
110      & +wstrain*ehpb+nss*ebr+wcorr*fact(3)*ecorr+wcorr5*fact(4)*ecorr5
111      & +wcorr6*fact(5)*ecorr6+wturn4*fact(3)*eello_turn4
112      & +wturn3*fact(2)*eello_turn3+wturn6*fact(5)*eturn6
113      & +wel_loc*fact(2)*eel_loc+edihcnstr+wtor_d*fact(2)*etors_d
114      & +wbond*estr+wsccor*fact(1)*esccor
115 #else
116       etot=wsc*(evdw+fact(6)*evdw_t)+wscp*evdw2
117      & +welec*fact(1)*(ees+evdw1)
118      & +wang*ebe+wtor*fact(1)*etors+wscloc*escloc
119      & +wstrain*ehpb+nss*ebr+wcorr*fact(3)*ecorr+wcorr5*fact(4)*ecorr5
120      & +wcorr6*fact(5)*ecorr6+wturn4*fact(3)*eello_turn4
121      & +wturn3*fact(2)*eello_turn3+wturn6*fact(5)*eturn6
122      & +wel_loc*fact(2)*eel_loc+edihcnstr+wtor_d*fact(2)*etors_d
123      & +wbond*estr+wsccor*fact(1)*esccor
124 #endif
125       energia(0)=etot
126       energia(1)=evdw
127 #ifdef SCP14
128       energia(2)=evdw2-evdw2_14
129       energia(17)=evdw2_14
130 #else
131       energia(2)=evdw2
132       energia(17)=0.0d0
133 #endif
134 #ifdef SPLITELE
135       energia(3)=ees
136       energia(16)=evdw1
137 #else
138       energia(3)=ees+evdw1
139       energia(16)=0.0d0
140 #endif
141       energia(4)=ecorr
142       energia(5)=ecorr5
143       energia(6)=ecorr6
144       energia(7)=eel_loc
145       energia(8)=eello_turn3
146       energia(9)=eello_turn4
147       energia(10)=eturn6
148       energia(11)=ebe
149       energia(12)=escloc
150       energia(13)=etors
151       energia(14)=etors_d
152       energia(15)=ehpb
153       energia(18)=estr
154       energia(19)=esccor
155       energia(20)=edihcnstr
156       energia(21)=evdw_t
157 c detecting NaNQ
158 #ifdef ISNAN
159 #ifdef AIX
160       if (isnan(etot).ne.0) energia(0)=1.0d+99
161 #else
162       if (isnan(etot)) energia(0)=1.0d+99
163 #endif
164 #else
165       i=0
166 #ifdef WINPGI
167       idumm=proc_proc(etot,i)
168 #else
169       call proc_proc(etot,i)
170 #endif
171       if(i.eq.1)energia(0)=1.0d+99
172 #endif
173 #ifdef MPL
174 c     endif
175 #endif
176       if (calc_grad) then
177 C
178 C Sum up the components of the Cartesian gradient.
179 C
180 #ifdef SPLITELE
181       do i=1,nct
182         do j=1,3
183           gradc(j,i,icg)=wsc*gvdwc(j,i)+wscp*gvdwc_scp(j,i)+
184      &                welec*fact(1)*gelc(j,i)+wvdwpp*gvdwpp(j,i)+
185      &                wbond*gradb(j,i)+
186      &                wstrain*ghpbc(j,i)+
187      &                wcorr*fact(3)*gradcorr(j,i)+
188      &                wel_loc*fact(2)*gel_loc(j,i)+
189      &                wturn3*fact(2)*gcorr3_turn(j,i)+
190      &                wturn4*fact(3)*gcorr4_turn(j,i)+
191      &                wcorr5*fact(4)*gradcorr5(j,i)+
192      &                wcorr6*fact(5)*gradcorr6(j,i)+
193      &                wturn6*fact(5)*gcorr6_turn(j,i)+
194      &                wsccor*fact(2)*gsccorc(j,i)
195           gradx(j,i,icg)=wsc*gvdwx(j,i)+wscp*gradx_scp(j,i)+
196      &                  wbond*gradbx(j,i)+
197      &                  wstrain*ghpbx(j,i)+wcorr*gradxorr(j,i)+
198      &                  wsccor*fact(2)*gsccorx(j,i)
199         enddo
200 #else
201       do i=1,nct
202         do j=1,3
203           gradc(j,i,icg)=wsc*gvdwc(j,i)+wscp*gvdwc_scp(j,i)+
204      &                welec*fact(1)*gelc(j,i)+wstrain*ghpbc(j,i)+
205      &                wbond*gradb(j,i)+
206      &                wcorr*fact(3)*gradcorr(j,i)+
207      &                wel_loc*fact(2)*gel_loc(j,i)+
208      &                wturn3*fact(2)*gcorr3_turn(j,i)+
209      &                wturn4*fact(3)*gcorr4_turn(j,i)+
210      &                wcorr5*fact(4)*gradcorr5(j,i)+
211      &                wcorr6*fact(5)*gradcorr6(j,i)+
212      &                wturn6*fact(5)*gcorr6_turn(j,i)+
213      &                wsccor*fact(2)*gsccorc(j,i)
214           gradx(j,i,icg)=wsc*gvdwx(j,i)+wscp*gradx_scp(j,i)+
215      &                  wbond*gradbx(j,i)+
216      &                  wstrain*ghpbx(j,i)+wcorr*gradxorr(j,i)+
217      &                  wsccor*fact(1)*gsccorx(j,i)
218         enddo
219 #endif
220       enddo
221
222
223       do i=1,nres-3
224         gloc(i,icg)=gloc(i,icg)+wcorr*fact(3)*gcorr_loc(i)
225      &   +wcorr5*fact(4)*g_corr5_loc(i)
226      &   +wcorr6*fact(5)*g_corr6_loc(i)
227      &   +wturn4*fact(3)*gel_loc_turn4(i)
228      &   +wturn3*fact(2)*gel_loc_turn3(i)
229      &   +wturn6*fact(5)*gel_loc_turn6(i)
230      &   +wel_loc*fact(2)*gel_loc_loc(i)
231 c     &   +wsccor*fact(1)*gsccor_loc(i)
232 c BYLA ROZNICA Z CLUSTER< OSTATNIA LINIA DODANA
233       enddo
234       endif
235       return
236       end
237 C------------------------------------------------------------------------
238       subroutine enerprint(energia,fact)
239       implicit real*8 (a-h,o-z)
240       include 'DIMENSIONS'
241       include 'DIMENSIONS.ZSCOPT'
242       include 'COMMON.IOUNITS'
243       include 'COMMON.FFIELD'
244       include 'COMMON.SBRIDGE'
245       double precision energia(0:max_ene),fact(6)
246       etot=energia(0)
247       evdw=energia(1)+fact(6)*energia(21)
248 #ifdef SCP14
249       evdw2=energia(2)+energia(17)
250 #else
251       evdw2=energia(2)
252 #endif
253       ees=energia(3)
254 #ifdef SPLITELE
255       evdw1=energia(16)
256 #endif
257       ecorr=energia(4)
258       ecorr5=energia(5)
259       ecorr6=energia(6)
260       eel_loc=energia(7)
261       eello_turn3=energia(8)
262       eello_turn4=energia(9)
263       eello_turn6=energia(10)
264       ebe=energia(11)
265       escloc=energia(12)
266       etors=energia(13)
267       etors_d=energia(14)
268       ehpb=energia(15)
269       esccor=energia(19)
270       edihcnstr=energia(20)
271       estr=energia(18)
272 #ifdef SPLITELE
273       write (iout,10) evdw,wsc,evdw2,wscp,ees,welec*fact(1),evdw1,
274      &  wvdwpp,
275      &  estr,wbond,ebe,wang,escloc,wscloc,etors,wtor*fact(1),
276      &  etors_d,wtor_d*fact(2),ehpb,wstrain,
277      &  ecorr,wcorr*fact(3),ecorr5,wcorr5*fact(4),ecorr6,wcorr6*fact(5),
278      &  eel_loc,wel_loc*fact(2),eello_turn3,wturn3*fact(2),
279      &  eello_turn4,wturn4*fact(3),eello_turn6,wturn6*fact(5),
280      &  esccor,wsccor*fact(1),edihcnstr,ebr*nss,etot
281    10 format (/'Virtual-chain energies:'//
282      & 'EVDW=  ',1pE16.6,' WEIGHT=',1pD16.6,' (SC-SC)'/
283      & 'EVDW2= ',1pE16.6,' WEIGHT=',1pD16.6,' (SC-p)'/
284      & 'EES=   ',1pE16.6,' WEIGHT=',1pD16.6,' (p-p elec)'/
285      & 'EVDWPP=',1pE16.6,' WEIGHT=',1pD16.6,' (p-p VDW)'/
286      & 'ESTR=  ',1pE16.6,' WEIGHT=',1pD16.6,' (stretching)'/
287      & 'EBE=   ',1pE16.6,' WEIGHT=',1pD16.6,' (bending)'/
288      & 'ESC=   ',1pE16.6,' WEIGHT=',1pD16.6,' (SC local)'/
289      & 'ETORS= ',1pE16.6,' WEIGHT=',1pD16.6,' (torsional)'/
290      & 'ETORSD=',1pE16.6,' WEIGHT=',1pD16.6,' (double torsional)'/
291      & 'EHBP=  ',1pE16.6,' WEIGHT=',1pD16.6,
292      & ' (SS bridges & dist. cnstr.)'/
293      & 'ECORR4=',1pE16.6,' WEIGHT=',1pD16.6,' (multi-body)'/
294      & 'ECORR5=',1pE16.6,' WEIGHT=',1pD16.6,' (multi-body)'/
295      & 'ECORR6=',1pE16.6,' WEIGHT=',1pD16.6,' (multi-body)'/
296      & 'EELLO= ',1pE16.6,' WEIGHT=',1pD16.6,' (electrostatic-local)'/
297      & 'ETURN3=',1pE16.6,' WEIGHT=',1pD16.6,' (turns, 3rd order)'/
298      & 'ETURN4=',1pE16.6,' WEIGHT=',1pD16.6,' (turns, 4th order)'/
299      & 'ETURN6=',1pE16.6,' WEIGHT=',1pD16.6,' (turns, 6th order)'/
300      & 'ESCCOR=',1pE16.6,' WEIGHT=',1pD16.6,' (backbone-rotamer corr)'/
301      & 'EDIHC= ',1pE16.6,' (dihedral angle constraints)'/
302      & 'ESS=   ',1pE16.6,' (disulfide-bridge intrinsic energy)'/ 
303      & 'ETOT=  ',1pE16.6,' (total)')
304 #else
305       write (iout,10) evdw,wsc,evdw2,wscp,ees,welec*fact(1),estr,wbond,
306      &  ebe,wang,escloc,wscloc,etors,wtor*fact(1),etors_d,wtor_d*fact2,
307      &  ehpb,wstrain,ecorr,wcorr*fact(3),ecorr5,wcorr5*fact(4),
308      &  ecorr6,wcorr6*fact(5),eel_loc,wel_loc*fact(2),
309      &  eello_turn3,wturn3*fact(2),eello_turn4,wturn4*fact(3),
310      &  eello_turn6,wturn6*fact(5),esccor*fact(1),wsccor,
311      &  edihcnstr,ebr*nss,etot
312    10 format (/'Virtual-chain energies:'//
313      & 'EVDW=  ',1pE16.6,' WEIGHT=',1pD16.6,' (SC-SC)'/
314      & 'EVDW2= ',1pE16.6,' WEIGHT=',1pD16.6,' (SC-p)'/
315      & 'EES=   ',1pE16.6,' WEIGHT=',1pD16.6,' (p-p)'/
316      & 'ESTR=  ',1pE16.6,' WEIGHT=',1pD16.6,' (stretching)'/
317      & 'EBE=   ',1pE16.6,' WEIGHT=',1pD16.6,' (bending)'/
318      & 'ESC=   ',1pE16.6,' WEIGHT=',1pD16.6,' (SC local)'/
319      & 'ETORS= ',1pE16.6,' WEIGHT=',1pD16.6,' (torsional)'/
320      & 'ETORSD=',1pE16.6,' WEIGHT=',1pD16.6,' (double torsional)'/
321      & 'EHBP=  ',1pE16.6,' WEIGHT=',1pD16.6,
322      & ' (SS bridges & dist. cnstr.)'/
323      & 'ECORR4=',1pE16.6,' WEIGHT=',1pD16.6,' (multi-body)'/
324      & 'ECORR5=',1pE16.6,' WEIGHT=',1pD16.6,' (multi-body)'/
325      & 'ECORR6=',1pE16.6,' WEIGHT=',1pD16.6,' (multi-body)'/
326      & 'EELLO= ',1pE16.6,' WEIGHT=',1pD16.6,' (electrostatic-local)'/
327      & 'ETURN3=',1pE16.6,' WEIGHT=',1pD16.6,' (turns, 3rd order)'/
328      & 'ETURN4=',1pE16.6,' WEIGHT=',1pD16.6,' (turns, 4th order)'/
329      & 'ETURN6=',1pE16.6,' WEIGHT=',1pD16.6,' (turns, 6th order)'/
330      & 'ESCCOR=',1pE16.6,' WEIGHT=',1pD16.6,' (backbone-rotamer corr)'/
331      & 'EDIHC= ',1pE16.6,' (dihedral angle constraints)'/
332      & 'ESS=   ',1pE16.6,' (disulfide-bridge intrinsic energy)'/ 
333      & 'ETOT=  ',1pE16.6,' (total)')
334 #endif
335       return
336       end
337 C-----------------------------------------------------------------------
338       subroutine elj(evdw,evdw_t)
339 C
340 C This subroutine calculates the interaction energy of nonbonded side chains
341 C assuming the LJ potential of interaction.
342 C
343       implicit real*8 (a-h,o-z)
344       include 'DIMENSIONS'
345       include 'DIMENSIONS.ZSCOPT'
346       include "DIMENSIONS.COMPAR"
347       parameter (accur=1.0d-10)
348       include 'COMMON.GEO'
349       include 'COMMON.VAR'
350       include 'COMMON.LOCAL'
351       include 'COMMON.CHAIN'
352       include 'COMMON.DERIV'
353       include 'COMMON.INTERACT'
354       include 'COMMON.TORSION'
355       include 'COMMON.ENEPS'
356       include 'COMMON.SBRIDGE'
357       include 'COMMON.NAMES'
358       include 'COMMON.IOUNITS'
359       include 'COMMON.CONTACTS'
360       dimension gg(3)
361       integer icant
362       external icant
363 cd    print *,'Entering ELJ nnt=',nnt,' nct=',nct,' expon=',expon
364 c ROZNICA z cluster
365       do i=1,210
366         do j=1,2
367           eneps_temp(j,i)=0.0d0
368         enddo
369       enddo
370 cROZNICA
371
372       evdw=0.0D0
373       evdw_t=0.0d0
374       do i=iatsc_s,iatsc_e
375         itypi=iabs(itype(i))
376         if (itypi.eq.ntyp1) cycle
377         itypi1=iabs(itype(i+1))
378         xi=c(1,nres+i)
379         yi=c(2,nres+i)
380         zi=c(3,nres+i)
381 C Change 12/1/95
382         num_conti=0
383 C
384 C Calculate SC interaction energy.
385 C
386         do iint=1,nint_gr(i)
387 cd        write (iout,*) 'i=',i,' iint=',iint,' istart=',istart(i,iint),
388 cd   &                  'iend=',iend(i,iint)
389           do j=istart(i,iint),iend(i,iint)
390             itypj=iabs(itype(j))
391             if (itypj.eq.ntyp1) cycle
392             xj=c(1,nres+j)-xi
393             yj=c(2,nres+j)-yi
394             zj=c(3,nres+j)-zi
395 C Change 12/1/95 to calculate four-body interactions
396             rij=xj*xj+yj*yj+zj*zj
397             rrij=1.0D0/rij
398 c           write (iout,*)'i=',i,' j=',j,' itypi=',itypi,' itypj=',itypj
399             eps0ij=eps(itypi,itypj)
400             fac=rrij**expon2
401             e1=fac*fac*aa(itypi,itypj)
402             e2=fac*bb(itypi,itypj)
403             evdwij=e1+e2
404             ij=icant(itypi,itypj)
405 c ROZNICA z cluster
406             eneps_temp(1,ij)=eneps_temp(1,ij)+e1/dabs(eps0ij)
407             eneps_temp(2,ij)=eneps_temp(2,ij)+e2/eps0ij
408 c
409
410 cd          sigm=dabs(aa(itypi,itypj)/bb(itypi,itypj))**(1.0D0/6.0D0)
411 cd          epsi=bb(itypi,itypj)**2/aa(itypi,itypj)
412 cd          write (iout,'(2(a3,i3,2x),6(1pd12.4)/2(3(1pd12.4),5x)/)')
413 cd   &        restyp(itypi),i,restyp(itypj),j,aa(itypi,itypj),
414 cd   &        bb(itypi,itypj),1.0D0/dsqrt(rrij),evdwij,epsi,sigm,
415 cd   &        (c(k,i),k=1,3),(c(k,j),k=1,3)
416             if (bb(itypi,itypj).gt.0.0d0) then
417               evdw=evdw+evdwij
418             else
419               evdw_t=evdw_t+evdwij
420             endif
421             if (calc_grad) then
422
423 C Calculate the components of the gradient in DC and X
424 C
425             fac=-rrij*(e1+evdwij)
426             gg(1)=xj*fac
427             gg(2)=yj*fac
428             gg(3)=zj*fac
429             do k=1,3
430               gvdwx(k,i)=gvdwx(k,i)-gg(k)
431               gvdwx(k,j)=gvdwx(k,j)+gg(k)
432             enddo
433             do k=i,j-1
434               do l=1,3
435                 gvdwc(l,k)=gvdwc(l,k)+gg(l)
436               enddo
437             enddo
438             endif
439 C
440 C 12/1/95, revised on 5/20/97
441 C
442 C Calculate the contact function. The ith column of the array JCONT will 
443 C contain the numbers of atoms that make contacts with the atom I (of numbers
444 C greater than I). The arrays FACONT and GACONT will contain the values of
445 C the contact function and its derivative.
446 C
447 C Uncomment next line, if the correlation interactions include EVDW explicitly.
448 c           if (j.gt.i+1 .and. evdwij.le.0.0D0) then
449 C Uncomment next line, if the correlation interactions are contact function only
450             if (j.gt.i+1.and. eps0ij.gt.0.0D0) then
451               rij=dsqrt(rij)
452               sigij=sigma(itypi,itypj)
453               r0ij=rs0(itypi,itypj)
454 C
455 C Check whether the SC's are not too far to make a contact.
456 C
457               rcut=1.5d0*r0ij
458               call gcont(rij,rcut,1.0d0,0.2d0*rcut,fcont,fprimcont)
459 C Add a new contact, if the SC's are close enough, but not too close (r<sigma).
460 C
461               if (fcont.gt.0.0D0) then
462 C If the SC-SC distance if close to sigma, apply spline.
463 cAdam           call gcont(-rij,-1.03d0*sigij,2.0d0*sigij,1.0d0,
464 cAdam &             fcont1,fprimcont1)
465 cAdam           fcont1=1.0d0-fcont1
466 cAdam           if (fcont1.gt.0.0d0) then
467 cAdam             fprimcont=fprimcont*fcont1+fcont*fprimcont1
468 cAdam             fcont=fcont*fcont1
469 cAdam           endif
470 C Uncomment following 4 lines to have the geometric average of the epsilon0's
471 cga             eps0ij=1.0d0/dsqrt(eps0ij)
472 cga             do k=1,3
473 cga               gg(k)=gg(k)*eps0ij
474 cga             enddo
475 cga             eps0ij=-evdwij*eps0ij
476 C Uncomment for AL's type of SC correlation interactions.
477 cadam           eps0ij=-evdwij
478                 num_conti=num_conti+1
479                 jcont(num_conti,i)=j
480                 facont(num_conti,i)=fcont*eps0ij
481                 fprimcont=eps0ij*fprimcont/rij
482                 fcont=expon*fcont
483 cAdam           gacont(1,num_conti,i)=-fprimcont*xj+fcont*gg(1)
484 cAdam           gacont(2,num_conti,i)=-fprimcont*yj+fcont*gg(2)
485 cAdam           gacont(3,num_conti,i)=-fprimcont*zj+fcont*gg(3)
486 C Uncomment following 3 lines for Skolnick's type of SC correlation.
487                 gacont(1,num_conti,i)=-fprimcont*xj
488                 gacont(2,num_conti,i)=-fprimcont*yj
489                 gacont(3,num_conti,i)=-fprimcont*zj
490 cd              write (iout,'(2i5,2f10.5)') i,j,rij,facont(num_conti,i)
491 cd              write (iout,'(2i3,3f10.5)') 
492 cd   &           i,j,(gacont(kk,num_conti,i),kk=1,3)
493               endif
494             endif
495           enddo      ! j
496         enddo        ! iint
497 C Change 12/1/95
498         num_cont(i)=num_conti
499       enddo          ! i
500       if (calc_grad) then
501       do i=1,nct
502         do j=1,3
503           gvdwc(j,i)=expon*gvdwc(j,i)
504           gvdwx(j,i)=expon*gvdwx(j,i)
505         enddo
506       enddo
507       endif
508 C******************************************************************************
509 C
510 C                              N O T E !!!
511 C
512 C To save time, the factor of EXPON has been extracted from ALL components
513 C of GVDWC and GRADX. Remember to multiply them by this factor before further 
514 C use!
515 C
516 C******************************************************************************
517       return
518       end
519 C-----------------------------------------------------------------------------
520       subroutine eljk(evdw,evdw_t)
521 C
522 C This subroutine calculates the interaction energy of nonbonded side chains
523 C assuming the LJK potential of interaction.
524 C
525       implicit real*8 (a-h,o-z)
526       include 'DIMENSIONS'
527       include 'DIMENSIONS.ZSCOPT'
528       include "DIMENSIONS.COMPAR"
529       include 'COMMON.GEO'
530       include 'COMMON.VAR'
531       include 'COMMON.LOCAL'
532       include 'COMMON.CHAIN'
533       include 'COMMON.DERIV'
534       include 'COMMON.INTERACT'
535       include 'COMMON.ENEPS'
536       include 'COMMON.IOUNITS'
537       include 'COMMON.NAMES'
538       dimension gg(3)
539       logical scheck
540       integer icant
541       external icant
542 c     print *,'Entering ELJK nnt=',nnt,' nct=',nct,' expon=',expon
543       do i=1,210
544         do j=1,2
545           eneps_temp(j,i)=0.0d0
546         enddo
547       enddo
548       evdw=0.0D0
549       evdw_t=0.0d0
550       do i=iatsc_s,iatsc_e
551         itypi=iabs(itype(i))
552         if (itypi.eq.ntyp1) cycle
553         itypi1=iabs(itype(i+1))
554         xi=c(1,nres+i)
555         yi=c(2,nres+i)
556         zi=c(3,nres+i)
557 C
558 C Calculate SC interaction energy.
559 C
560         do iint=1,nint_gr(i)
561           do j=istart(i,iint),iend(i,iint)
562             itypj=iabs(itype(j))
563             if (itypj.eq.ntyp1) cycle
564             xj=c(1,nres+j)-xi
565             yj=c(2,nres+j)-yi
566             zj=c(3,nres+j)-zi
567             rrij=1.0D0/(xj*xj+yj*yj+zj*zj)
568             fac_augm=rrij**expon
569             e_augm=augm(itypi,itypj)*fac_augm
570             r_inv_ij=dsqrt(rrij)
571             rij=1.0D0/r_inv_ij 
572             r_shift_inv=1.0D0/(rij+r0(itypi,itypj)-sigma(itypi,itypj))
573             fac=r_shift_inv**expon
574             e1=fac*fac*aa(itypi,itypj)
575             e2=fac*bb(itypi,itypj)
576             evdwij=e_augm+e1+e2
577             ij=icant(itypi,itypj)
578             eneps_temp(1,ij)=eneps_temp(1,ij)+(e1+a_augm)
579      &        /dabs(eps(itypi,itypj))
580             eneps_temp(2,ij)=eneps_temp(2,ij)+e2/eps(itypi,itypj)
581 cd          sigm=dabs(aa(itypi,itypj)/bb(itypi,itypj))**(1.0D0/6.0D0)
582 cd          epsi=bb(itypi,itypj)**2/aa(itypi,itypj)
583 cd          write (iout,'(2(a3,i3,2x),8(1pd12.4)/2(3(1pd12.4),5x)/)')
584 cd   &        restyp(itypi),i,restyp(itypj),j,aa(itypi,itypj),
585 cd   &        bb(itypi,itypj),augm(itypi,itypj),epsi,sigm,
586 cd   &        sigma(itypi,itypj),1.0D0/dsqrt(rrij),evdwij,
587 cd   &        (c(k,i),k=1,3),(c(k,j),k=1,3)
588             if (bb(itypi,itypj).gt.0.0d0) then
589               evdw=evdw+evdwij
590             else 
591               evdw_t=evdw_t+evdwij
592             endif
593             if (calc_grad) then
594
595 C Calculate the components of the gradient in DC and X
596 C
597             fac=-2.0D0*rrij*e_augm-r_inv_ij*r_shift_inv*(e1+e1+e2)
598             gg(1)=xj*fac
599             gg(2)=yj*fac
600             gg(3)=zj*fac
601             do k=1,3
602               gvdwx(k,i)=gvdwx(k,i)-gg(k)
603               gvdwx(k,j)=gvdwx(k,j)+gg(k)
604             enddo
605             do k=i,j-1
606               do l=1,3
607                 gvdwc(l,k)=gvdwc(l,k)+gg(l)
608               enddo
609             enddo
610             endif
611           enddo      ! j
612         enddo        ! iint
613       enddo          ! i
614       if (calc_grad) then
615       do i=1,nct
616         do j=1,3
617           gvdwc(j,i)=expon*gvdwc(j,i)
618           gvdwx(j,i)=expon*gvdwx(j,i)
619         enddo
620       enddo
621       endif
622       return
623       end
624 C-----------------------------------------------------------------------------
625       subroutine ebp(evdw,evdw_t)
626 C
627 C This subroutine calculates the interaction energy of nonbonded side chains
628 C assuming the Berne-Pechukas potential of interaction.
629 C
630       implicit real*8 (a-h,o-z)
631       include 'DIMENSIONS'
632       include 'DIMENSIONS.ZSCOPT'
633       include "DIMENSIONS.COMPAR"
634       include 'COMMON.GEO'
635       include 'COMMON.VAR'
636       include 'COMMON.LOCAL'
637       include 'COMMON.CHAIN'
638       include 'COMMON.DERIV'
639       include 'COMMON.NAMES'
640       include 'COMMON.INTERACT'
641       include 'COMMON.ENEPS'
642       include 'COMMON.IOUNITS'
643       include 'COMMON.CALC'
644       common /srutu/ icall
645 c     double precision rrsave(maxdim)
646       logical lprn
647       integer icant
648       external icant
649       do i=1,210
650         do j=1,2
651           eneps_temp(j,i)=0.0d0
652         enddo
653       enddo
654       evdw=0.0D0
655       evdw_t=0.0d0
656 c     print *,'Entering EBP nnt=',nnt,' nct=',nct,' expon=',expon
657 c     if (icall.eq.0) then
658 c       lprn=.true.
659 c     else
660         lprn=.false.
661 c     endif
662       ind=0
663       do i=iatsc_s,iatsc_e
664         itypi=iabs(itype(i))
665         if (itypi.eq.ntyp1) cycle
666         itypi1=iabs(itype(i+1))
667         xi=c(1,nres+i)
668         yi=c(2,nres+i)
669         zi=c(3,nres+i)
670         dxi=dc_norm(1,nres+i)
671         dyi=dc_norm(2,nres+i)
672         dzi=dc_norm(3,nres+i)
673         dsci_inv=vbld_inv(i+nres)
674 C
675 C Calculate SC interaction energy.
676 C
677         do iint=1,nint_gr(i)
678           do j=istart(i,iint),iend(i,iint)
679             ind=ind+1
680             itypj=iabs(itype(j))
681             if (itypj.eq.ntyp1) cycle
682             dscj_inv=vbld_inv(j+nres)
683             chi1=chi(itypi,itypj)
684             chi2=chi(itypj,itypi)
685             chi12=chi1*chi2
686             chip1=chip(itypi)
687             chip2=chip(itypj)
688             chip12=chip1*chip2
689             alf1=alp(itypi)
690             alf2=alp(itypj)
691             alf12=0.5D0*(alf1+alf2)
692 C For diagnostics only!!!
693 c           chi1=0.0D0
694 c           chi2=0.0D0
695 c           chi12=0.0D0
696 c           chip1=0.0D0
697 c           chip2=0.0D0
698 c           chip12=0.0D0
699 c           alf1=0.0D0
700 c           alf2=0.0D0
701 c           alf12=0.0D0
702             xj=c(1,nres+j)-xi
703             yj=c(2,nres+j)-yi
704             zj=c(3,nres+j)-zi
705             dxj=dc_norm(1,nres+j)
706             dyj=dc_norm(2,nres+j)
707             dzj=dc_norm(3,nres+j)
708             rrij=1.0D0/(xj*xj+yj*yj+zj*zj)
709 cd          if (icall.eq.0) then
710 cd            rrsave(ind)=rrij
711 cd          else
712 cd            rrij=rrsave(ind)
713 cd          endif
714             rij=dsqrt(rrij)
715 C Calculate the angle-dependent terms of energy & contributions to derivatives.
716             call sc_angular
717 C Calculate whole angle-dependent part of epsilon and contributions
718 C to its derivatives
719             fac=(rrij*sigsq)**expon2
720             e1=fac*fac*aa(itypi,itypj)
721             e2=fac*bb(itypi,itypj)
722             evdwij=eps1*eps2rt*eps3rt*(e1+e2)
723             eps2der=evdwij*eps3rt
724             eps3der=evdwij*eps2rt
725             evdwij=evdwij*eps2rt*eps3rt
726             ij=icant(itypi,itypj)
727             aux=eps1*eps2rt**2*eps3rt**2
728             eneps_temp(1,ij)=eneps_temp(1,ij)+e1*aux
729      &        /dabs(eps(itypi,itypj))
730             eneps_temp(2,ij)=eneps_temp(2,ij)+e2*aux/eps(itypi,itypj)
731             if (bb(itypi,itypj).gt.0.0d0) then
732               evdw=evdw+evdwij
733             else
734               evdw_t=evdw_t+evdwij
735             endif
736             if (calc_grad) then
737             if (lprn) then
738             sigm=dabs(aa(itypi,itypj)/bb(itypi,itypj))**(1.0D0/6.0D0)
739             epsi=bb(itypi,itypj)**2/aa(itypi,itypj)
740             write (iout,'(2(a3,i3,2x),15(0pf7.3))')
741      &        restyp(itypi),i,restyp(itypj),j,
742      &        epsi,sigm,chi1,chi2,chip1,chip2,
743      &        eps1,eps2rt**2,eps3rt**2,1.0D0/dsqrt(sigsq),
744      &        om1,om2,om12,1.0D0/dsqrt(rrij),
745      &        evdwij
746             endif
747 C Calculate gradient components.
748             e1=e1*eps1*eps2rt**2*eps3rt**2
749             fac=-expon*(e1+evdwij)
750             sigder=fac/sigsq
751             fac=rrij*fac
752 C Calculate radial part of the gradient
753             gg(1)=xj*fac
754             gg(2)=yj*fac
755             gg(3)=zj*fac
756 C Calculate the angular part of the gradient and sum add the contributions
757 C to the appropriate components of the Cartesian gradient.
758             call sc_grad
759             endif
760           enddo      ! j
761         enddo        ! iint
762       enddo          ! i
763 c     stop
764       return
765       end
766 C-----------------------------------------------------------------------------
767       subroutine egb(evdw,evdw_t)
768 C
769 C This subroutine calculates the interaction energy of nonbonded side chains
770 C assuming the Gay-Berne potential of interaction.
771 C
772       implicit real*8 (a-h,o-z)
773       include 'DIMENSIONS'
774       include 'DIMENSIONS.ZSCOPT'
775       include "DIMENSIONS.COMPAR"
776       include 'COMMON.GEO'
777       include 'COMMON.VAR'
778       include 'COMMON.LOCAL'
779       include 'COMMON.CHAIN'
780       include 'COMMON.DERIV'
781       include 'COMMON.NAMES'
782       include 'COMMON.INTERACT'
783       include 'COMMON.ENEPS'
784       include 'COMMON.IOUNITS'
785       include 'COMMON.CALC'
786       logical lprn
787       common /srutu/icall
788       integer icant
789       external icant
790       do i=1,210
791         do j=1,2
792           eneps_temp(j,i)=0.0d0
793         enddo
794       enddo
795 c     print *,'Entering EGB nnt=',nnt,' nct=',nct,' expon=',expon
796       evdw=0.0D0
797       evdw_t=0.0d0
798       lprn=.false.
799 c      if (icall.gt.0) lprn=.true.
800       ind=0
801       do i=iatsc_s,iatsc_e
802         itypi=iabs(itype(i))
803         if (itypi.eq.ntyp1) cycle
804         itypi1=iabs(itype(i+1))
805         xi=c(1,nres+i)
806         yi=c(2,nres+i)
807         zi=c(3,nres+i)
808         dxi=dc_norm(1,nres+i)
809         dyi=dc_norm(2,nres+i)
810         dzi=dc_norm(3,nres+i)
811         dsci_inv=vbld_inv(i+nres)
812 C
813 C Calculate SC interaction energy.
814 C
815         do iint=1,nint_gr(i)
816           do j=istart(i,iint),iend(i,iint)
817             ind=ind+1
818             itypj=iabs(itype(j))
819             if (itypj.eq.ntyp1) cycle
820             dscj_inv=vbld_inv(j+nres)
821             sig0ij=sigma(itypi,itypj)
822             chi1=chi(itypi,itypj)
823             chi2=chi(itypj,itypi)
824             chi12=chi1*chi2
825             chip1=chip(itypi)
826             chip2=chip(itypj)
827             chip12=chip1*chip2
828             alf1=alp(itypi)
829             alf2=alp(itypj)
830             alf12=0.5D0*(alf1+alf2)
831 C For diagnostics only!!!
832 c           chi1=0.0D0
833 c           chi2=0.0D0
834 c           chi12=0.0D0
835 c           chip1=0.0D0
836 c           chip2=0.0D0
837 c           chip12=0.0D0
838 c           alf1=0.0D0
839 c           alf2=0.0D0
840 c           alf12=0.0D0
841             xj=c(1,nres+j)-xi
842             yj=c(2,nres+j)-yi
843             zj=c(3,nres+j)-zi
844             dxj=dc_norm(1,nres+j)
845             dyj=dc_norm(2,nres+j)
846             dzj=dc_norm(3,nres+j)
847 c            write (iout,*) i,j,xj,yj,zj
848             rrij=1.0D0/(xj*xj+yj*yj+zj*zj)
849             rij=dsqrt(rrij)
850 C Calculate angle-dependent terms of energy and contributions to their
851 C derivatives.
852             call sc_angular
853             sigsq=1.0D0/sigsq
854             sig=sig0ij*dsqrt(sigsq)
855             rij_shift=1.0D0/rij-sig+sig0ij
856 C I hate to put IF's in the loops, but here don't have another choice!!!!
857             if (rij_shift.le.0.0D0) then
858               evdw=1.0D20
859               return
860             endif
861             sigder=-sig*sigsq
862 c---------------------------------------------------------------
863             rij_shift=1.0D0/rij_shift 
864             fac=rij_shift**expon
865             e1=fac*fac*aa(itypi,itypj)
866             e2=fac*bb(itypi,itypj)
867             evdwij=eps1*eps2rt*eps3rt*(e1+e2)
868             eps2der=evdwij*eps3rt
869             eps3der=evdwij*eps2rt
870             evdwij=evdwij*eps2rt*eps3rt
871             if (bb(itypi,itypj).gt.0) then
872               evdw=evdw+evdwij
873             else
874               evdw_t=evdw_t+evdwij
875             endif
876             ij=icant(itypi,itypj)
877             aux=eps1*eps2rt**2*eps3rt**2
878             eneps_temp(1,ij)=eneps_temp(1,ij)+aux*e1
879      &        /dabs(eps(itypi,itypj))
880             eneps_temp(2,ij)=eneps_temp(2,ij)+aux*e2/eps(itypi,itypj)
881 c            write (iout,*) "i",i," j",j," itypi",itypi," itypj",itypj,
882 c     &         " ij",ij," eneps",aux*e1/dabs(eps(itypi,itypj)),
883 c     &         aux*e2/eps(itypi,itypj)
884 c            if (lprn) then
885             sigm=dabs(aa(itypi,itypj)/bb(itypi,itypj))**(1.0D0/6.0D0)
886             epsi=bb(itypi,itypj)**2/aa(itypi,itypj)
887 #ifdef DEBUG
888             write (iout,'(2(a3,i3,2x),17(0pf7.3))')
889      &        restyp(itypi),i,restyp(itypj),j,
890      &        epsi,sigm,chi1,chi2,chip1,chip2,
891      &        eps1,eps2rt**2,eps3rt**2,sig,sig0ij,
892      &        om1,om2,om12,1.0D0/rij,1.0D0/rij_shift,
893      &        evdwij
894              write (iout,*) "partial sum", evdw, evdw_t
895 #endif
896 c            endif
897             if (calc_grad) then
898 C Calculate gradient components.
899             e1=e1*eps1*eps2rt**2*eps3rt**2
900             fac=-expon*(e1+evdwij)*rij_shift
901             sigder=fac*sigder
902             fac=rij*fac
903 C Calculate the radial part of the gradient
904             gg(1)=xj*fac
905             gg(2)=yj*fac
906             gg(3)=zj*fac
907 C Calculate angular part of the gradient.
908             call sc_grad
909             endif
910           enddo      ! j
911         enddo        ! iint
912       enddo          ! i
913       return
914       end
915 C-----------------------------------------------------------------------------
916       subroutine egbv(evdw,evdw_t)
917 C
918 C This subroutine calculates the interaction energy of nonbonded side chains
919 C assuming the Gay-Berne-Vorobjev potential of interaction.
920 C
921       implicit real*8 (a-h,o-z)
922       include 'DIMENSIONS'
923       include 'DIMENSIONS.ZSCOPT'
924       include "DIMENSIONS.COMPAR"
925       include 'COMMON.GEO'
926       include 'COMMON.VAR'
927       include 'COMMON.LOCAL'
928       include 'COMMON.CHAIN'
929       include 'COMMON.DERIV'
930       include 'COMMON.NAMES'
931       include 'COMMON.INTERACT'
932       include 'COMMON.ENEPS'
933       include 'COMMON.IOUNITS'
934       include 'COMMON.CALC'
935       common /srutu/ icall
936       logical lprn
937       integer icant
938       external icant
939       do i=1,210
940         do j=1,2
941           eneps_temp(j,i)=0.0d0
942         enddo
943       enddo
944       evdw=0.0D0
945       evdw_t=0.0d0
946 c     print *,'Entering EGB nnt=',nnt,' nct=',nct,' expon=',expon
947       evdw=0.0D0
948       lprn=.false.
949 c      if (icall.gt.0) lprn=.true.
950       ind=0
951       do i=iatsc_s,iatsc_e
952         itypi=iabs(itype(i))
953         if (itypi.eq.ntyp1) cycle
954         itypi1=iabs(itype(i+1))
955         xi=c(1,nres+i)
956         yi=c(2,nres+i)
957         zi=c(3,nres+i)
958         dxi=dc_norm(1,nres+i)
959         dyi=dc_norm(2,nres+i)
960         dzi=dc_norm(3,nres+i)
961         dsci_inv=vbld_inv(i+nres)
962 C
963 C Calculate SC interaction energy.
964 C
965         do iint=1,nint_gr(i)
966           do j=istart(i,iint),iend(i,iint)
967             ind=ind+1
968             itypj=iabs(itype(j))
969             if (itypj.eq.ntyp1) cycle
970             dscj_inv=vbld_inv(j+nres)
971             sig0ij=sigma(itypi,itypj)
972             r0ij=r0(itypi,itypj)
973             chi1=chi(itypi,itypj)
974             chi2=chi(itypj,itypi)
975             chi12=chi1*chi2
976             chip1=chip(itypi)
977             chip2=chip(itypj)
978             chip12=chip1*chip2
979             alf1=alp(itypi)
980             alf2=alp(itypj)
981             alf12=0.5D0*(alf1+alf2)
982 C For diagnostics only!!!
983 c           chi1=0.0D0
984 c           chi2=0.0D0
985 c           chi12=0.0D0
986 c           chip1=0.0D0
987 c           chip2=0.0D0
988 c           chip12=0.0D0
989 c           alf1=0.0D0
990 c           alf2=0.0D0
991 c           alf12=0.0D0
992             xj=c(1,nres+j)-xi
993             yj=c(2,nres+j)-yi
994             zj=c(3,nres+j)-zi
995             dxj=dc_norm(1,nres+j)
996             dyj=dc_norm(2,nres+j)
997             dzj=dc_norm(3,nres+j)
998             rrij=1.0D0/(xj*xj+yj*yj+zj*zj)
999             rij=dsqrt(rrij)
1000 C Calculate angle-dependent terms of energy and contributions to their
1001 C derivatives.
1002             call sc_angular
1003             sigsq=1.0D0/sigsq
1004             sig=sig0ij*dsqrt(sigsq)
1005             rij_shift=1.0D0/rij-sig+r0ij
1006 C I hate to put IF's in the loops, but here don't have another choice!!!!
1007             if (rij_shift.le.0.0D0) then
1008               evdw=1.0D20
1009               return
1010             endif
1011             sigder=-sig*sigsq
1012 c---------------------------------------------------------------
1013             rij_shift=1.0D0/rij_shift 
1014             fac=rij_shift**expon
1015             e1=fac*fac*aa(itypi,itypj)
1016             e2=fac*bb(itypi,itypj)
1017             evdwij=eps1*eps2rt*eps3rt*(e1+e2)
1018             eps2der=evdwij*eps3rt
1019             eps3der=evdwij*eps2rt
1020             fac_augm=rrij**expon
1021             e_augm=augm(itypi,itypj)*fac_augm
1022             evdwij=evdwij*eps2rt*eps3rt
1023             if (bb(itypi,itypj).gt.0.0d0) then
1024               evdw=evdw+evdwij+e_augm
1025             else
1026               evdw_t=evdw_t+evdwij+e_augm
1027             endif
1028             ij=icant(itypi,itypj)
1029             aux=eps1*eps2rt**2*eps3rt**2
1030             eneps_temp(1,ij)=eneps_temp(1,ij)+aux*(e1+e_augm)
1031      &        /dabs(eps(itypi,itypj))
1032             eneps_temp(2,ij)=eneps_temp(2,ij)+aux*e2/eps(itypi,itypj)
1033 c            eneps_temp(ij)=eneps_temp(ij)
1034 c     &         +(evdwij+e_augm)/eps(itypi,itypj)
1035 c            if (lprn) then
1036 c            sigm=dabs(aa(itypi,itypj)/bb(itypi,itypj))**(1.0D0/6.0D0)
1037 c            epsi=bb(itypi,itypj)**2/aa(itypi,itypj)
1038 c            write (iout,'(2(a3,i3,2x),17(0pf7.3))')
1039 c     &        restyp(itypi),i,restyp(itypj),j,
1040 c     &        epsi,sigm,sig,(augm(itypi,itypj)/epsi)**(1.0D0/12.0D0),
1041 c     &        chi1,chi2,chip1,chip2,
1042 c     &        eps1,eps2rt**2,eps3rt**2,
1043 c     &        om1,om2,om12,1.0D0/rij,1.0D0/rij_shift,
1044 c     &        evdwij+e_augm
1045 c            endif
1046             if (calc_grad) then
1047 C Calculate gradient components.
1048             e1=e1*eps1*eps2rt**2*eps3rt**2
1049             fac=-expon*(e1+evdwij)*rij_shift
1050             sigder=fac*sigder
1051             fac=rij*fac-2*expon*rrij*e_augm
1052 C Calculate the radial part of the gradient
1053             gg(1)=xj*fac
1054             gg(2)=yj*fac
1055             gg(3)=zj*fac
1056 C Calculate angular part of the gradient.
1057             call sc_grad
1058             endif
1059           enddo      ! j
1060         enddo        ! iint
1061       enddo          ! i
1062       return
1063       end
1064 C-----------------------------------------------------------------------------
1065       subroutine sc_angular
1066 C Calculate eps1,eps2,eps3,sigma, and parts of their derivatives in om1,om2,
1067 C om12. Called by ebp, egb, and egbv.
1068       implicit none
1069       include 'COMMON.CALC'
1070       erij(1)=xj*rij
1071       erij(2)=yj*rij
1072       erij(3)=zj*rij
1073       om1=dxi*erij(1)+dyi*erij(2)+dzi*erij(3)
1074       om2=dxj*erij(1)+dyj*erij(2)+dzj*erij(3)
1075       om12=dxi*dxj+dyi*dyj+dzi*dzj
1076       chiom12=chi12*om12
1077 C Calculate eps1(om12) and its derivative in om12
1078       faceps1=1.0D0-om12*chiom12
1079       faceps1_inv=1.0D0/faceps1
1080       eps1=dsqrt(faceps1_inv)
1081 C Following variable is eps1*deps1/dom12
1082       eps1_om12=faceps1_inv*chiom12
1083 C Calculate sigma(om1,om2,om12) and the derivatives of sigma**2 in om1,om2,
1084 C and om12.
1085       om1om2=om1*om2
1086       chiom1=chi1*om1
1087       chiom2=chi2*om2
1088       facsig=om1*chiom1+om2*chiom2-2.0D0*om1om2*chiom12
1089       sigsq=1.0D0-facsig*faceps1_inv
1090       sigsq_om1=(chiom1-chiom12*om2)*faceps1_inv
1091       sigsq_om2=(chiom2-chiom12*om1)*faceps1_inv
1092       sigsq_om12=-chi12*(om1om2*faceps1-om12*facsig)*faceps1_inv**2
1093 C Calculate eps2 and its derivatives in om1, om2, and om12.
1094       chipom1=chip1*om1
1095       chipom2=chip2*om2
1096       chipom12=chip12*om12
1097       facp=1.0D0-om12*chipom12
1098       facp_inv=1.0D0/facp
1099       facp1=om1*chipom1+om2*chipom2-2.0D0*om1om2*chipom12
1100 C Following variable is the square root of eps2
1101       eps2rt=1.0D0-facp1*facp_inv
1102 C Following three variables are the derivatives of the square root of eps
1103 C in om1, om2, and om12.
1104       eps2rt_om1=-4.0D0*(chipom1-chipom12*om2)*facp_inv
1105       eps2rt_om2=-4.0D0*(chipom2-chipom12*om1)*facp_inv
1106       eps2rt_om12=4.0D0*chip12*(om1om2*facp-om12*facp1)*facp_inv**2 
1107 C Evaluate the "asymmetric" factor in the VDW constant, eps3
1108       eps3rt=1.0D0-alf1*om1+alf2*om2-alf12*om12 
1109 C Calculate whole angle-dependent part of epsilon and contributions
1110 C to its derivatives
1111       return
1112       end
1113 C----------------------------------------------------------------------------
1114       subroutine sc_grad
1115       implicit real*8 (a-h,o-z)
1116       include 'DIMENSIONS'
1117       include 'DIMENSIONS.ZSCOPT'
1118       include 'COMMON.CHAIN'
1119       include 'COMMON.DERIV'
1120       include 'COMMON.CALC'
1121       double precision dcosom1(3),dcosom2(3)
1122       eom1=eps2der*eps2rt_om1-2.0D0*alf1*eps3der+sigder*sigsq_om1
1123       eom2=eps2der*eps2rt_om2+2.0D0*alf2*eps3der+sigder*sigsq_om2
1124       eom12=evdwij*eps1_om12+eps2der*eps2rt_om12
1125      &     -2.0D0*alf12*eps3der+sigder*sigsq_om12
1126       do k=1,3
1127         dcosom1(k)=rij*(dc_norm(k,nres+i)-om1*erij(k))
1128         dcosom2(k)=rij*(dc_norm(k,nres+j)-om2*erij(k))
1129       enddo
1130       do k=1,3
1131         gg(k)=gg(k)+eom1*dcosom1(k)+eom2*dcosom2(k)
1132       enddo 
1133       do k=1,3
1134         gvdwx(k,i)=gvdwx(k,i)-gg(k)
1135      &            +(eom12*(dc_norm(k,nres+j)-om12*dc_norm(k,nres+i))
1136      &            +eom1*(erij(k)-om1*dc_norm(k,nres+i)))*dsci_inv
1137         gvdwx(k,j)=gvdwx(k,j)+gg(k)
1138      &            +(eom12*(dc_norm(k,nres+i)-om12*dc_norm(k,nres+j))
1139      &            +eom2*(erij(k)-om2*dc_norm(k,nres+j)))*dscj_inv
1140       enddo
1141
1142 C Calculate the components of the gradient in DC and X
1143 C
1144       do k=i,j-1
1145         do l=1,3
1146           gvdwc(l,k)=gvdwc(l,k)+gg(l)
1147         enddo
1148       enddo
1149       return
1150       end
1151 c------------------------------------------------------------------------------
1152       subroutine vec_and_deriv
1153       implicit real*8 (a-h,o-z)
1154       include 'DIMENSIONS'
1155       include 'DIMENSIONS.ZSCOPT'
1156       include 'COMMON.IOUNITS'
1157       include 'COMMON.GEO'
1158       include 'COMMON.VAR'
1159       include 'COMMON.LOCAL'
1160       include 'COMMON.CHAIN'
1161       include 'COMMON.VECTORS'
1162       include 'COMMON.DERIV'
1163       include 'COMMON.INTERACT'
1164       dimension uyder(3,3,2),uzder(3,3,2),vbld_inv_temp(2)
1165 C Compute the local reference systems. For reference system (i), the
1166 C X-axis points from CA(i) to CA(i+1), the Y axis is in the 
1167 C CA(i)-CA(i+1)-CA(i+2) plane, and the Z axis is perpendicular to this plane.
1168       do i=1,nres-1
1169 c          if (i.eq.nres-1 .or. itel(i+1).eq.0) then
1170           if (i.eq.nres-1) then
1171 C Case of the last full residue
1172 C Compute the Z-axis
1173             call vecpr(dc_norm(1,i),dc_norm(1,i-1),uz(1,i))
1174             costh=dcos(pi-theta(nres))
1175             fac=1.0d0/dsqrt(1.0d0-costh*costh)
1176             do k=1,3
1177               uz(k,i)=fac*uz(k,i)
1178             enddo
1179             if (calc_grad) then
1180 C Compute the derivatives of uz
1181             uzder(1,1,1)= 0.0d0
1182             uzder(2,1,1)=-dc_norm(3,i-1)
1183             uzder(3,1,1)= dc_norm(2,i-1) 
1184             uzder(1,2,1)= dc_norm(3,i-1)
1185             uzder(2,2,1)= 0.0d0
1186             uzder(3,2,1)=-dc_norm(1,i-1)
1187             uzder(1,3,1)=-dc_norm(2,i-1)
1188             uzder(2,3,1)= dc_norm(1,i-1)
1189             uzder(3,3,1)= 0.0d0
1190             uzder(1,1,2)= 0.0d0
1191             uzder(2,1,2)= dc_norm(3,i)
1192             uzder(3,1,2)=-dc_norm(2,i) 
1193             uzder(1,2,2)=-dc_norm(3,i)
1194             uzder(2,2,2)= 0.0d0
1195             uzder(3,2,2)= dc_norm(1,i)
1196             uzder(1,3,2)= dc_norm(2,i)
1197             uzder(2,3,2)=-dc_norm(1,i)
1198             uzder(3,3,2)= 0.0d0
1199             endif
1200 C Compute the Y-axis
1201             facy=fac
1202             do k=1,3
1203               uy(k,i)=fac*(dc_norm(k,i-1)-costh*dc_norm(k,i))
1204             enddo
1205             if (calc_grad) then
1206 C Compute the derivatives of uy
1207             do j=1,3
1208               do k=1,3
1209                 uyder(k,j,1)=2*dc_norm(k,i-1)*dc_norm(j,i)
1210      &                        -dc_norm(k,i)*dc_norm(j,i-1)
1211                 uyder(k,j,2)=-dc_norm(j,i)*dc_norm(k,i)
1212               enddo
1213               uyder(j,j,1)=uyder(j,j,1)-costh
1214               uyder(j,j,2)=1.0d0+uyder(j,j,2)
1215             enddo
1216             do j=1,2
1217               do k=1,3
1218                 do l=1,3
1219                   uygrad(l,k,j,i)=uyder(l,k,j)
1220                   uzgrad(l,k,j,i)=uzder(l,k,j)
1221                 enddo
1222               enddo
1223             enddo 
1224             call unormderiv(uy(1,i),uyder(1,1,1),facy,uygrad(1,1,1,i))
1225             call unormderiv(uy(1,i),uyder(1,1,2),facy,uygrad(1,1,2,i))
1226             call unormderiv(uz(1,i),uzder(1,1,1),fac,uzgrad(1,1,1,i))
1227             call unormderiv(uz(1,i),uzder(1,1,2),fac,uzgrad(1,1,2,i))
1228             endif
1229           else
1230 C Other residues
1231 C Compute the Z-axis
1232             call vecpr(dc_norm(1,i),dc_norm(1,i+1),uz(1,i))
1233             costh=dcos(pi-theta(i+2))
1234             fac=1.0d0/dsqrt(1.0d0-costh*costh)
1235             do k=1,3
1236               uz(k,i)=fac*uz(k,i)
1237             enddo
1238             if (calc_grad) then
1239 C Compute the derivatives of uz
1240             uzder(1,1,1)= 0.0d0
1241             uzder(2,1,1)=-dc_norm(3,i+1)
1242             uzder(3,1,1)= dc_norm(2,i+1) 
1243             uzder(1,2,1)= dc_norm(3,i+1)
1244             uzder(2,2,1)= 0.0d0
1245             uzder(3,2,1)=-dc_norm(1,i+1)
1246             uzder(1,3,1)=-dc_norm(2,i+1)
1247             uzder(2,3,1)= dc_norm(1,i+1)
1248             uzder(3,3,1)= 0.0d0
1249             uzder(1,1,2)= 0.0d0
1250             uzder(2,1,2)= dc_norm(3,i)
1251             uzder(3,1,2)=-dc_norm(2,i) 
1252             uzder(1,2,2)=-dc_norm(3,i)
1253             uzder(2,2,2)= 0.0d0
1254             uzder(3,2,2)= dc_norm(1,i)
1255             uzder(1,3,2)= dc_norm(2,i)
1256             uzder(2,3,2)=-dc_norm(1,i)
1257             uzder(3,3,2)= 0.0d0
1258             endif
1259 C Compute the Y-axis
1260             facy=fac
1261             do k=1,3
1262               uy(k,i)=facy*(dc_norm(k,i+1)-costh*dc_norm(k,i))
1263             enddo
1264             if (calc_grad) then
1265 C Compute the derivatives of uy
1266             do j=1,3
1267               do k=1,3
1268                 uyder(k,j,1)=2*dc_norm(k,i+1)*dc_norm(j,i)
1269      &                        -dc_norm(k,i)*dc_norm(j,i+1)
1270                 uyder(k,j,2)=-dc_norm(j,i)*dc_norm(k,i)
1271               enddo
1272               uyder(j,j,1)=uyder(j,j,1)-costh
1273               uyder(j,j,2)=1.0d0+uyder(j,j,2)
1274             enddo
1275             do j=1,2
1276               do k=1,3
1277                 do l=1,3
1278                   uygrad(l,k,j,i)=uyder(l,k,j)
1279                   uzgrad(l,k,j,i)=uzder(l,k,j)
1280                 enddo
1281               enddo
1282             enddo 
1283             call unormderiv(uy(1,i),uyder(1,1,1),facy,uygrad(1,1,1,i))
1284             call unormderiv(uy(1,i),uyder(1,1,2),facy,uygrad(1,1,2,i))
1285             call unormderiv(uz(1,i),uzder(1,1,1),fac,uzgrad(1,1,1,i))
1286             call unormderiv(uz(1,i),uzder(1,1,2),fac,uzgrad(1,1,2,i))
1287           endif
1288           endif
1289       enddo
1290       if (calc_grad) then
1291       do i=1,nres-1
1292         vbld_inv_temp(1)=vbld_inv(i+1)
1293         if (i.lt.nres-1) then
1294           vbld_inv_temp(2)=vbld_inv(i+2)
1295         else
1296           vbld_inv_temp(2)=vbld_inv(i)
1297         endif
1298         do j=1,2
1299           do k=1,3
1300             do l=1,3
1301               uygrad(l,k,j,i)=vbld_inv_temp(j)*uygrad(l,k,j,i)
1302               uzgrad(l,k,j,i)=vbld_inv_temp(j)*uzgrad(l,k,j,i)
1303             enddo
1304           enddo
1305         enddo
1306       enddo
1307       endif
1308       return
1309       end
1310 C-----------------------------------------------------------------------------
1311       subroutine vec_and_deriv_test
1312       implicit real*8 (a-h,o-z)
1313       include 'DIMENSIONS'
1314       include 'DIMENSIONS.ZSCOPT'
1315       include 'COMMON.IOUNITS'
1316       include 'COMMON.GEO'
1317       include 'COMMON.VAR'
1318       include 'COMMON.LOCAL'
1319       include 'COMMON.CHAIN'
1320       include 'COMMON.VECTORS'
1321       dimension uyder(3,3,2),uzder(3,3,2)
1322 C Compute the local reference systems. For reference system (i), the
1323 C X-axis points from CA(i) to CA(i+1), the Y axis is in the 
1324 C CA(i)-CA(i+1)-CA(i+2) plane, and the Z axis is perpendicular to this plane.
1325       do i=1,nres-1
1326           if (i.eq.nres-1) then
1327 C Case of the last full residue
1328 C Compute the Z-axis
1329             call vecpr(dc_norm(1,i),dc_norm(1,i-1),uz(1,i))
1330             costh=dcos(pi-theta(nres))
1331             fac=1.0d0/dsqrt(1.0d0-costh*costh)
1332 c            write (iout,*) 'fac',fac,
1333 c     &        1.0d0/dsqrt(scalar(uz(1,i),uz(1,i)))
1334             fac=1.0d0/dsqrt(scalar(uz(1,i),uz(1,i)))
1335             do k=1,3
1336               uz(k,i)=fac*uz(k,i)
1337             enddo
1338 C Compute the derivatives of uz
1339             uzder(1,1,1)= 0.0d0
1340             uzder(2,1,1)=-dc_norm(3,i-1)
1341             uzder(3,1,1)= dc_norm(2,i-1) 
1342             uzder(1,2,1)= dc_norm(3,i-1)
1343             uzder(2,2,1)= 0.0d0
1344             uzder(3,2,1)=-dc_norm(1,i-1)
1345             uzder(1,3,1)=-dc_norm(2,i-1)
1346             uzder(2,3,1)= dc_norm(1,i-1)
1347             uzder(3,3,1)= 0.0d0
1348             uzder(1,1,2)= 0.0d0
1349             uzder(2,1,2)= dc_norm(3,i)
1350             uzder(3,1,2)=-dc_norm(2,i) 
1351             uzder(1,2,2)=-dc_norm(3,i)
1352             uzder(2,2,2)= 0.0d0
1353             uzder(3,2,2)= dc_norm(1,i)
1354             uzder(1,3,2)= dc_norm(2,i)
1355             uzder(2,3,2)=-dc_norm(1,i)
1356             uzder(3,3,2)= 0.0d0
1357 C Compute the Y-axis
1358             do k=1,3
1359               uy(k,i)=fac*(dc_norm(k,i-1)-costh*dc_norm(k,i))
1360             enddo
1361             facy=fac
1362             facy=1.0d0/dsqrt(scalar(dc_norm(1,i),dc_norm(1,i))*
1363      &       (scalar(dc_norm(1,i-1),dc_norm(1,i-1))**2-
1364      &        scalar(dc_norm(1,i),dc_norm(1,i-1))**2))
1365             do k=1,3
1366 c              uy(k,i)=facy*(dc_norm(k,i+1)-costh*dc_norm(k,i))
1367               uy(k,i)=
1368 c     &        facy*(
1369      &        dc_norm(k,i-1)*scalar(dc_norm(1,i),dc_norm(1,i))
1370      &        -scalar(dc_norm(1,i),dc_norm(1,i-1))*dc_norm(k,i)
1371 c     &        )
1372             enddo
1373 c            write (iout,*) 'facy',facy,
1374 c     &       1.0d0/dsqrt(scalar(uy(1,i),uy(1,i)))
1375             facy=1.0d0/dsqrt(scalar(uy(1,i),uy(1,i)))
1376             do k=1,3
1377               uy(k,i)=facy*uy(k,i)
1378             enddo
1379 C Compute the derivatives of uy
1380             do j=1,3
1381               do k=1,3
1382                 uyder(k,j,1)=2*dc_norm(k,i-1)*dc_norm(j,i)
1383      &                        -dc_norm(k,i)*dc_norm(j,i-1)
1384                 uyder(k,j,2)=-dc_norm(j,i)*dc_norm(k,i)
1385               enddo
1386 c              uyder(j,j,1)=uyder(j,j,1)-costh
1387 c              uyder(j,j,2)=1.0d0+uyder(j,j,2)
1388               uyder(j,j,1)=uyder(j,j,1)
1389      &          -scalar(dc_norm(1,i),dc_norm(1,i-1))
1390               uyder(j,j,2)=scalar(dc_norm(1,i),dc_norm(1,i))
1391      &          +uyder(j,j,2)
1392             enddo
1393             do j=1,2
1394               do k=1,3
1395                 do l=1,3
1396                   uygrad(l,k,j,i)=uyder(l,k,j)
1397                   uzgrad(l,k,j,i)=uzder(l,k,j)
1398                 enddo
1399               enddo
1400             enddo 
1401             call unormderiv(uy(1,i),uyder(1,1,1),facy,uygrad(1,1,1,i))
1402             call unormderiv(uy(1,i),uyder(1,1,2),facy,uygrad(1,1,2,i))
1403             call unormderiv(uz(1,i),uzder(1,1,1),fac,uzgrad(1,1,1,i))
1404             call unormderiv(uz(1,i),uzder(1,1,2),fac,uzgrad(1,1,2,i))
1405           else
1406 C Other residues
1407 C Compute the Z-axis
1408             call vecpr(dc_norm(1,i),dc_norm(1,i+1),uz(1,i))
1409             costh=dcos(pi-theta(i+2))
1410             fac=1.0d0/dsqrt(1.0d0-costh*costh)
1411             fac=1.0d0/dsqrt(scalar(uz(1,i),uz(1,i)))
1412             do k=1,3
1413               uz(k,i)=fac*uz(k,i)
1414             enddo
1415 C Compute the derivatives of uz
1416             uzder(1,1,1)= 0.0d0
1417             uzder(2,1,1)=-dc_norm(3,i+1)
1418             uzder(3,1,1)= dc_norm(2,i+1) 
1419             uzder(1,2,1)= dc_norm(3,i+1)
1420             uzder(2,2,1)= 0.0d0
1421             uzder(3,2,1)=-dc_norm(1,i+1)
1422             uzder(1,3,1)=-dc_norm(2,i+1)
1423             uzder(2,3,1)= dc_norm(1,i+1)
1424             uzder(3,3,1)= 0.0d0
1425             uzder(1,1,2)= 0.0d0
1426             uzder(2,1,2)= dc_norm(3,i)
1427             uzder(3,1,2)=-dc_norm(2,i) 
1428             uzder(1,2,2)=-dc_norm(3,i)
1429             uzder(2,2,2)= 0.0d0
1430             uzder(3,2,2)= dc_norm(1,i)
1431             uzder(1,3,2)= dc_norm(2,i)
1432             uzder(2,3,2)=-dc_norm(1,i)
1433             uzder(3,3,2)= 0.0d0
1434 C Compute the Y-axis
1435             facy=fac
1436             facy=1.0d0/dsqrt(scalar(dc_norm(1,i),dc_norm(1,i))*
1437      &       (scalar(dc_norm(1,i+1),dc_norm(1,i+1))**2-
1438      &        scalar(dc_norm(1,i),dc_norm(1,i+1))**2))
1439             do k=1,3
1440 c              uy(k,i)=facy*(dc_norm(k,i+1)-costh*dc_norm(k,i))
1441               uy(k,i)=
1442 c     &        facy*(
1443      &        dc_norm(k,i+1)*scalar(dc_norm(1,i),dc_norm(1,i))
1444      &        -scalar(dc_norm(1,i),dc_norm(1,i+1))*dc_norm(k,i)
1445 c     &        )
1446             enddo
1447 c            write (iout,*) 'facy',facy,
1448 c     &       1.0d0/dsqrt(scalar(uy(1,i),uy(1,i)))
1449             facy=1.0d0/dsqrt(scalar(uy(1,i),uy(1,i)))
1450             do k=1,3
1451               uy(k,i)=facy*uy(k,i)
1452             enddo
1453 C Compute the derivatives of uy
1454             do j=1,3
1455               do k=1,3
1456                 uyder(k,j,1)=2*dc_norm(k,i+1)*dc_norm(j,i)
1457      &                        -dc_norm(k,i)*dc_norm(j,i+1)
1458                 uyder(k,j,2)=-dc_norm(j,i)*dc_norm(k,i)
1459               enddo
1460 c              uyder(j,j,1)=uyder(j,j,1)-costh
1461 c              uyder(j,j,2)=1.0d0+uyder(j,j,2)
1462               uyder(j,j,1)=uyder(j,j,1)
1463      &          -scalar(dc_norm(1,i),dc_norm(1,i+1))
1464               uyder(j,j,2)=scalar(dc_norm(1,i),dc_norm(1,i))
1465      &          +uyder(j,j,2)
1466             enddo
1467             do j=1,2
1468               do k=1,3
1469                 do l=1,3
1470                   uygrad(l,k,j,i)=uyder(l,k,j)
1471                   uzgrad(l,k,j,i)=uzder(l,k,j)
1472                 enddo
1473               enddo
1474             enddo 
1475             call unormderiv(uy(1,i),uyder(1,1,1),facy,uygrad(1,1,1,i))
1476             call unormderiv(uy(1,i),uyder(1,1,2),facy,uygrad(1,1,2,i))
1477             call unormderiv(uz(1,i),uzder(1,1,1),fac,uzgrad(1,1,1,i))
1478             call unormderiv(uz(1,i),uzder(1,1,2),fac,uzgrad(1,1,2,i))
1479           endif
1480       enddo
1481       do i=1,nres-1
1482         do j=1,2
1483           do k=1,3
1484             do l=1,3
1485               uygrad(l,k,j,i)=vblinv*uygrad(l,k,j,i)
1486               uzgrad(l,k,j,i)=vblinv*uzgrad(l,k,j,i)
1487             enddo
1488           enddo
1489         enddo
1490       enddo
1491       return
1492       end
1493 C-----------------------------------------------------------------------------
1494       subroutine check_vecgrad
1495       implicit real*8 (a-h,o-z)
1496       include 'DIMENSIONS'
1497       include 'DIMENSIONS.ZSCOPT'
1498       include 'COMMON.IOUNITS'
1499       include 'COMMON.GEO'
1500       include 'COMMON.VAR'
1501       include 'COMMON.LOCAL'
1502       include 'COMMON.CHAIN'
1503       include 'COMMON.VECTORS'
1504       dimension uygradt(3,3,2,maxres),uzgradt(3,3,2,maxres)
1505       dimension uyt(3,maxres),uzt(3,maxres)
1506       dimension uygradn(3,3,2),uzgradn(3,3,2),erij(3)
1507       double precision delta /1.0d-7/
1508       call vec_and_deriv
1509 cd      do i=1,nres
1510 crc          write(iout,'(2i5,2(3f10.5,5x))') i,1,dc_norm(:,i)
1511 crc          write(iout,'(2i5,2(3f10.5,5x))') i,2,uy(:,i)
1512 crc          write(iout,'(2i5,2(3f10.5,5x)/)')i,3,uz(:,i)
1513 cd          write(iout,'(2i5,2(3f10.5,5x))') i,1,
1514 cd     &     (dc_norm(if90,i),if90=1,3)
1515 cd          write(iout,'(2i5,2(3f10.5,5x))') i,2,(uy(if90,i),if90=1,3)
1516 cd          write(iout,'(2i5,2(3f10.5,5x)/)')i,3,(uz(if90,i),if90=1,3)
1517 cd          write(iout,'(a)')
1518 cd      enddo
1519       do i=1,nres
1520         do j=1,2
1521           do k=1,3
1522             do l=1,3
1523               uygradt(l,k,j,i)=uygrad(l,k,j,i)
1524               uzgradt(l,k,j,i)=uzgrad(l,k,j,i)
1525             enddo
1526           enddo
1527         enddo
1528       enddo
1529       call vec_and_deriv
1530       do i=1,nres
1531         do j=1,3
1532           uyt(j,i)=uy(j,i)
1533           uzt(j,i)=uz(j,i)
1534         enddo
1535       enddo
1536       do i=1,nres
1537 cd        write (iout,*) 'i=',i
1538         do k=1,3
1539           erij(k)=dc_norm(k,i)
1540         enddo
1541         do j=1,3
1542           do k=1,3
1543             dc_norm(k,i)=erij(k)
1544           enddo
1545           dc_norm(j,i)=dc_norm(j,i)+delta
1546 c          fac=dsqrt(scalar(dc_norm(1,i),dc_norm(1,i)))
1547 c          do k=1,3
1548 c            dc_norm(k,i)=dc_norm(k,i)/fac
1549 c          enddo
1550 c          write (iout,*) (dc_norm(k,i),k=1,3)
1551 c          write (iout,*) (erij(k),k=1,3)
1552           call vec_and_deriv
1553           do k=1,3
1554             uygradn(k,j,1)=(uy(k,i)-uyt(k,i))/delta
1555             uygradn(k,j,2)=(uy(k,i-1)-uyt(k,i-1))/delta
1556             uzgradn(k,j,1)=(uz(k,i)-uzt(k,i))/delta
1557             uzgradn(k,j,2)=(uz(k,i-1)-uzt(k,i-1))/delta
1558           enddo 
1559 c          write (iout,'(i5,3f8.5,3x,3f8.5,5x,3f8.5,3x,3f8.5)') 
1560 c     &      j,(uzgradt(k,j,1,i),k=1,3),(uzgradn(k,j,1),k=1,3),
1561 c     &      (uzgradt(k,j,2,i-1),k=1,3),(uzgradn(k,j,2),k=1,3)
1562         enddo
1563         do k=1,3
1564           dc_norm(k,i)=erij(k)
1565         enddo
1566 cd        do k=1,3
1567 cd          write (iout,'(i5,3f8.5,3x,3f8.5,5x,3f8.5,3x,3f8.5)') 
1568 cd     &      k,(uygradt(k,l,1,i),l=1,3),(uygradn(k,l,1),l=1,3),
1569 cd     &      (uygradt(k,l,2,i-1),l=1,3),(uygradn(k,l,2),l=1,3)
1570 cd          write (iout,'(i5,3f8.5,3x,3f8.5,5x,3f8.5,3x,3f8.5)') 
1571 cd     &      k,(uzgradt(k,l,1,i),l=1,3),(uzgradn(k,l,1),l=1,3),
1572 cd     &      (uzgradt(k,l,2,i-1),l=1,3),(uzgradn(k,l,2),l=1,3)
1573 cd          write (iout,'(a)')
1574 cd        enddo
1575       enddo
1576       return
1577       end
1578 C--------------------------------------------------------------------------
1579       subroutine set_matrices
1580       implicit real*8 (a-h,o-z)
1581       include 'DIMENSIONS'
1582       include 'DIMENSIONS.ZSCOPT'
1583       include 'COMMON.IOUNITS'
1584       include 'COMMON.GEO'
1585       include 'COMMON.VAR'
1586       include 'COMMON.LOCAL'
1587       include 'COMMON.CHAIN'
1588       include 'COMMON.DERIV'
1589       include 'COMMON.INTERACT'
1590       include 'COMMON.CONTACTS'
1591       include 'COMMON.TORSION'
1592       include 'COMMON.VECTORS'
1593       include 'COMMON.FFIELD'
1594       double precision auxvec(2),auxmat(2,2)
1595 C
1596 C Compute the virtual-bond-torsional-angle dependent quantities needed
1597 C to calculate the el-loc multibody terms of various order.
1598 C
1599       do i=3,nres+1
1600         if (i .lt. nres+1) then
1601           sin1=dsin(phi(i))
1602           cos1=dcos(phi(i))
1603           sintab(i-2)=sin1
1604           costab(i-2)=cos1
1605           obrot(1,i-2)=cos1
1606           obrot(2,i-2)=sin1
1607           sin2=dsin(2*phi(i))
1608           cos2=dcos(2*phi(i))
1609           sintab2(i-2)=sin2
1610           costab2(i-2)=cos2
1611           obrot2(1,i-2)=cos2
1612           obrot2(2,i-2)=sin2
1613           Ug(1,1,i-2)=-cos1
1614           Ug(1,2,i-2)=-sin1
1615           Ug(2,1,i-2)=-sin1
1616           Ug(2,2,i-2)= cos1
1617           Ug2(1,1,i-2)=-cos2
1618           Ug2(1,2,i-2)=-sin2
1619           Ug2(2,1,i-2)=-sin2
1620           Ug2(2,2,i-2)= cos2
1621         else
1622           costab(i-2)=1.0d0
1623           sintab(i-2)=0.0d0
1624           obrot(1,i-2)=1.0d0
1625           obrot(2,i-2)=0.0d0
1626           obrot2(1,i-2)=0.0d0
1627           obrot2(2,i-2)=0.0d0
1628           Ug(1,1,i-2)=1.0d0
1629           Ug(1,2,i-2)=0.0d0
1630           Ug(2,1,i-2)=0.0d0
1631           Ug(2,2,i-2)=1.0d0
1632           Ug2(1,1,i-2)=0.0d0
1633           Ug2(1,2,i-2)=0.0d0
1634           Ug2(2,1,i-2)=0.0d0
1635           Ug2(2,2,i-2)=0.0d0
1636         endif
1637         if (i .gt. 3 .and. i .lt. nres+1) then
1638           obrot_der(1,i-2)=-sin1
1639           obrot_der(2,i-2)= cos1
1640           Ugder(1,1,i-2)= sin1
1641           Ugder(1,2,i-2)=-cos1
1642           Ugder(2,1,i-2)=-cos1
1643           Ugder(2,2,i-2)=-sin1
1644           dwacos2=cos2+cos2
1645           dwasin2=sin2+sin2
1646           obrot2_der(1,i-2)=-dwasin2
1647           obrot2_der(2,i-2)= dwacos2
1648           Ug2der(1,1,i-2)= dwasin2
1649           Ug2der(1,2,i-2)=-dwacos2
1650           Ug2der(2,1,i-2)=-dwacos2
1651           Ug2der(2,2,i-2)=-dwasin2
1652         else
1653           obrot_der(1,i-2)=0.0d0
1654           obrot_der(2,i-2)=0.0d0
1655           Ugder(1,1,i-2)=0.0d0
1656           Ugder(1,2,i-2)=0.0d0
1657           Ugder(2,1,i-2)=0.0d0
1658           Ugder(2,2,i-2)=0.0d0
1659           obrot2_der(1,i-2)=0.0d0
1660           obrot2_der(2,i-2)=0.0d0
1661           Ug2der(1,1,i-2)=0.0d0
1662           Ug2der(1,2,i-2)=0.0d0
1663           Ug2der(2,1,i-2)=0.0d0
1664           Ug2der(2,2,i-2)=0.0d0
1665         endif
1666         if (i.gt. nnt+2 .and. i.lt.nct+2) then
1667           if (itype(i-2).le.ntyp) then
1668             iti = itortyp(itype(i-2))
1669           else 
1670             iti=ntortyp+1
1671           endif
1672         else
1673           iti=ntortyp+1
1674         endif
1675         if (i.gt. nnt+1 .and. i.lt.nct+1) then
1676           if (itype(i-1).le.ntyp) then
1677             iti1 = itortyp(itype(i-1))
1678           else
1679             iti1=ntortyp+1
1680           endif
1681         else
1682           iti1=ntortyp+1
1683         endif
1684 cd        write (iout,*) '*******i',i,' iti1',iti
1685 cd        write (iout,*) 'b1',b1(:,iti)
1686 cd        write (iout,*) 'b2',b2(:,iti)
1687 cd        write (iout,*) 'Ug',Ug(:,:,i-2)
1688 c        print *,"itilde1 i iti iti1",i,iti,iti1
1689         if (i .gt. iatel_s+2) then
1690           call matvec2(Ug(1,1,i-2),b2(1,iti),Ub2(1,i-2))
1691           call matmat2(EE(1,1,iti),Ug(1,1,i-2),EUg(1,1,i-2))
1692           call matmat2(CC(1,1,iti),Ug(1,1,i-2),CUg(1,1,i-2))
1693           call matmat2(DD(1,1,iti),Ug(1,1,i-2),DUg(1,1,i-2))
1694           call matmat2(Dtilde(1,1,iti),Ug2(1,1,i-2),DtUg2(1,1,i-2))
1695           call matvec2(Ctilde(1,1,iti1),obrot(1,i-2),Ctobr(1,i-2))
1696           call matvec2(Dtilde(1,1,iti),obrot2(1,i-2),Dtobr2(1,i-2))
1697         else
1698           do k=1,2
1699             Ub2(k,i-2)=0.0d0
1700             Ctobr(k,i-2)=0.0d0 
1701             Dtobr2(k,i-2)=0.0d0
1702             do l=1,2
1703               EUg(l,k,i-2)=0.0d0
1704               CUg(l,k,i-2)=0.0d0
1705               DUg(l,k,i-2)=0.0d0
1706               DtUg2(l,k,i-2)=0.0d0
1707             enddo
1708           enddo
1709         endif
1710 c        print *,"itilde2 i iti iti1",i,iti,iti1
1711         call matvec2(Ugder(1,1,i-2),b2(1,iti),Ub2der(1,i-2))
1712         call matmat2(EE(1,1,iti),Ugder(1,1,i-2),EUgder(1,1,i-2))
1713         call matmat2(CC(1,1,iti1),Ugder(1,1,i-2),CUgder(1,1,i-2))
1714         call matmat2(DD(1,1,iti),Ugder(1,1,i-2),DUgder(1,1,i-2))
1715         call matmat2(Dtilde(1,1,iti),Ug2der(1,1,i-2),DtUg2der(1,1,i-2))
1716         call matvec2(Ctilde(1,1,iti1),obrot_der(1,i-2),Ctobrder(1,i-2))
1717         call matvec2(Dtilde(1,1,iti),obrot2_der(1,i-2),Dtobr2der(1,i-2))
1718 c        print *,"itilde3 i iti iti1",i,iti,iti1
1719         do k=1,2
1720           muder(k,i-2)=Ub2der(k,i-2)
1721         enddo
1722         if (i.gt. nnt+1 .and. i.lt.nct+1) then
1723           if (itype(i-1).le.ntyp) then
1724             iti1 = itortyp(itype(i-1))
1725           else
1726             iti1=ntortyp+1
1727           endif
1728         else
1729           iti1=ntortyp+1
1730         endif
1731         do k=1,2
1732           mu(k,i-2)=Ub2(k,i-2)+b1(k,iti1)
1733         enddo
1734 C Vectors and matrices dependent on a single virtual-bond dihedral.
1735         call matvec2(DD(1,1,iti),b1tilde(1,iti1),auxvec(1))
1736         call matvec2(Ug2(1,1,i-2),auxvec(1),Ug2Db1t(1,i-2)) 
1737         call matvec2(Ug2der(1,1,i-2),auxvec(1),Ug2Db1tder(1,i-2)) 
1738         call matvec2(CC(1,1,iti1),Ub2(1,i-2),CUgb2(1,i-2))
1739         call matvec2(CC(1,1,iti1),Ub2der(1,i-2),CUgb2der(1,i-2))
1740         call matmat2(EUg(1,1,i-2),CC(1,1,iti1),EUgC(1,1,i-2))
1741         call matmat2(EUgder(1,1,i-2),CC(1,1,iti1),EUgCder(1,1,i-2))
1742         call matmat2(EUg(1,1,i-2),DD(1,1,iti1),EUgD(1,1,i-2))
1743         call matmat2(EUgder(1,1,i-2),DD(1,1,iti1),EUgDder(1,1,i-2))
1744 cd        write (iout,*) 'i',i,' mu ',(mu(k,i-2),k=1,2),
1745 cd     &  ' mu1',(b1(k,i-2),k=1,2),' mu2',(Ub2(k,i-2),k=1,2)
1746       enddo
1747 C Matrices dependent on two consecutive virtual-bond dihedrals.
1748 C The order of matrices is from left to right.
1749       do i=2,nres-1
1750         call matmat2(DtUg2(1,1,i-1),EUg(1,1,i),DtUg2EUg(1,1,i))
1751         call matmat2(DtUg2der(1,1,i-1),EUg(1,1,i),DtUg2EUgder(1,1,1,i))
1752         call matmat2(DtUg2(1,1,i-1),EUgder(1,1,i),DtUg2EUgder(1,1,2,i))
1753         call transpose2(DtUg2(1,1,i-1),auxmat(1,1))
1754         call matmat2(auxmat(1,1),EUg(1,1,i),Ug2DtEUg(1,1,i))
1755         call matmat2(auxmat(1,1),EUgder(1,1,i),Ug2DtEUgder(1,1,2,i))
1756         call transpose2(DtUg2der(1,1,i-1),auxmat(1,1))
1757         call matmat2(auxmat(1,1),EUg(1,1,i),Ug2DtEUgder(1,1,1,i))
1758       enddo
1759 cd      do i=1,nres
1760 cd        iti = itortyp(itype(i))
1761 cd        write (iout,*) i
1762 cd        do j=1,2
1763 cd        write (iout,'(2f10.5,5x,2f10.5,5x,2f10.5)') 
1764 cd     &  (EE(j,k,iti),k=1,2),(Ug(j,k,i),k=1,2),(EUg(j,k,i),k=1,2)
1765 cd        enddo
1766 cd      enddo
1767       return
1768       end
1769 C--------------------------------------------------------------------------
1770       subroutine eelec(ees,evdw1,eel_loc,eello_turn3,eello_turn4)
1771 C
1772 C This subroutine calculates the average interaction energy and its gradient
1773 C in the virtual-bond vectors between non-adjacent peptide groups, based on 
1774 C the potential described in Liwo et al., Protein Sci., 1993, 2, 1715. 
1775 C The potential depends both on the distance of peptide-group centers and on 
1776 C the orientation of the CA-CA virtual bonds.
1777
1778       implicit real*8 (a-h,o-z)
1779       include 'DIMENSIONS'
1780       include 'DIMENSIONS.ZSCOPT'
1781       include 'COMMON.CONTROL'
1782       include 'COMMON.IOUNITS'
1783       include 'COMMON.GEO'
1784       include 'COMMON.VAR'
1785       include 'COMMON.LOCAL'
1786       include 'COMMON.CHAIN'
1787       include 'COMMON.DERIV'
1788       include 'COMMON.INTERACT'
1789       include 'COMMON.CONTACTS'
1790       include 'COMMON.TORSION'
1791       include 'COMMON.VECTORS'
1792       include 'COMMON.FFIELD'
1793       dimension ggg(3),gggp(3),gggm(3),erij(3),dcosb(3),dcosg(3),
1794      &          erder(3,3),uryg(3,3),urzg(3,3),vryg(3,3),vrzg(3,3)
1795       double precision acipa(2,2),agg(3,4),aggi(3,4),aggi1(3,4),
1796      &    aggj(3,4),aggj1(3,4),a_temp(2,2),muij(4)
1797       common /locel/ a_temp,agg,aggi,aggi1,aggj,aggj1,j1
1798 c 4/26/02 - AL scaling factor for 1,4 repulsive VDW interactions
1799       double precision scal_el /0.5d0/
1800 C 12/13/98 
1801 C 13-go grudnia roku pamietnego... 
1802       double precision unmat(3,3) /1.0d0,0.0d0,0.0d0,
1803      &                   0.0d0,1.0d0,0.0d0,
1804      &                   0.0d0,0.0d0,1.0d0/
1805 cd      write(iout,*) 'In EELEC'
1806 cd      do i=1,nloctyp
1807 cd        write(iout,*) 'Type',i
1808 cd        write(iout,*) 'B1',B1(:,i)
1809 cd        write(iout,*) 'B2',B2(:,i)
1810 cd        write(iout,*) 'CC',CC(:,:,i)
1811 cd        write(iout,*) 'DD',DD(:,:,i)
1812 cd        write(iout,*) 'EE',EE(:,:,i)
1813 cd      enddo
1814 cd      call check_vecgrad
1815 cd      stop
1816       if (icheckgrad.eq.1) then
1817         do i=1,nres-1
1818           fac=1.0d0/dsqrt(scalar(dc(1,i),dc(1,i)))
1819           do k=1,3
1820             dc_norm(k,i)=dc(k,i)*fac
1821           enddo
1822 c          write (iout,*) 'i',i,' fac',fac
1823         enddo
1824       endif
1825       if (wel_loc.gt.0.0d0 .or. wcorr4.gt.0.0d0 .or. wcorr5.gt.0.0d0 
1826      &    .or. wcorr6.gt.0.0d0 .or. wturn3.gt.0.0d0 .or. 
1827      &    wturn4.gt.0.0d0 .or. wturn6.gt.0.0d0) then
1828 cd      if (wel_loc.gt.0.0d0) then
1829         if (icheckgrad.eq.1) then
1830         call vec_and_deriv_test
1831         else
1832         call vec_and_deriv
1833         endif
1834         call set_matrices
1835       endif
1836 cd      do i=1,nres-1
1837 cd        write (iout,*) 'i=',i
1838 cd        do k=1,3
1839 cd          write (iout,'(i5,2f10.5)') k,uy(k,i),uz(k,i)
1840 cd        enddo
1841 cd        do k=1,3
1842 cd          write (iout,'(f10.5,2x,3f10.5,2x,3f10.5)') 
1843 cd     &     uz(k,i),(uzgrad(k,l,1,i),l=1,3),(uzgrad(k,l,2,i),l=1,3)
1844 cd        enddo
1845 cd      enddo
1846       num_conti_hb=0
1847       ees=0.0D0
1848       evdw1=0.0D0
1849       eel_loc=0.0d0 
1850       eello_turn3=0.0d0
1851       eello_turn4=0.0d0
1852       ind=0
1853       do i=1,nres
1854         num_cont_hb(i)=0
1855       enddo
1856 cd      print '(a)','Enter EELEC'
1857 cd      write (iout,*) 'iatel_s=',iatel_s,' iatel_e=',iatel_e
1858       do i=1,nres
1859         gel_loc_loc(i)=0.0d0
1860         gcorr_loc(i)=0.0d0
1861       enddo
1862       do i=iatel_s,iatel_e
1863         if (itype(i).eq.ntyp1 .or. itype(i+1).eq.ntyp1) cycle
1864         if (itel(i).eq.0) goto 1215
1865         dxi=dc(1,i)
1866         dyi=dc(2,i)
1867         dzi=dc(3,i)
1868         dx_normi=dc_norm(1,i)
1869         dy_normi=dc_norm(2,i)
1870         dz_normi=dc_norm(3,i)
1871         xmedi=c(1,i)+0.5d0*dxi
1872         ymedi=c(2,i)+0.5d0*dyi
1873         zmedi=c(3,i)+0.5d0*dzi
1874         num_conti=0
1875 c        write (iout,*) 'i',i,' ielstart',ielstart(i),' ielend',ielend(i)
1876         do j=ielstart(i),ielend(i)
1877           if (itype(j).eq.ntyp1 .or. itype(j+1).eq.ntyp1) cycle
1878           if (itel(j).eq.0) goto 1216
1879           ind=ind+1
1880           iteli=itel(i)
1881           itelj=itel(j)
1882           if (j.eq.i+2 .and. itelj.eq.2) iteli=2
1883           aaa=app(iteli,itelj)
1884           bbb=bpp(iteli,itelj)
1885 C Diagnostics only!!!
1886 c         aaa=0.0D0
1887 c         bbb=0.0D0
1888 c         ael6i=0.0D0
1889 c         ael3i=0.0D0
1890 C End diagnostics
1891           ael6i=ael6(iteli,itelj)
1892           ael3i=ael3(iteli,itelj) 
1893           dxj=dc(1,j)
1894           dyj=dc(2,j)
1895           dzj=dc(3,j)
1896           dx_normj=dc_norm(1,j)
1897           dy_normj=dc_norm(2,j)
1898           dz_normj=dc_norm(3,j)
1899           xj=c(1,j)+0.5D0*dxj-xmedi
1900           yj=c(2,j)+0.5D0*dyj-ymedi
1901           zj=c(3,j)+0.5D0*dzj-zmedi
1902           rij=xj*xj+yj*yj+zj*zj
1903           rrmij=1.0D0/rij
1904           rij=dsqrt(rij)
1905           rmij=1.0D0/rij
1906           r3ij=rrmij*rmij
1907           r6ij=r3ij*r3ij  
1908           cosa=dx_normi*dx_normj+dy_normi*dy_normj+dz_normi*dz_normj
1909           cosb=(xj*dx_normi+yj*dy_normi+zj*dz_normi)*rmij
1910           cosg=(xj*dx_normj+yj*dy_normj+zj*dz_normj)*rmij
1911           fac=cosa-3.0D0*cosb*cosg
1912           ev1=aaa*r6ij*r6ij
1913 c 4/26/02 - AL scaling down 1,4 repulsive VDW interactions
1914           if (j.eq.i+2) ev1=scal_el*ev1
1915           ev2=bbb*r6ij
1916           fac3=ael6i*r6ij
1917           fac4=ael3i*r3ij
1918           evdwij=ev1+ev2
1919           el1=fac3*(4.0D0+fac*fac-3.0D0*(cosb*cosb+cosg*cosg))
1920           el2=fac4*fac       
1921           eesij=el1+el2
1922 c          write (iout,*) "i",i,iteli," j",j,itelj," eesij",eesij
1923 C 12/26/95 - for the evaluation of multi-body H-bonding interactions
1924           ees0ij=4.0D0+fac*fac-3.0D0*(cosb*cosb+cosg*cosg)
1925           ees=ees+eesij
1926           evdw1=evdw1+evdwij
1927 c             write (iout,'(a6,2i5,0pf7.3,2i5,2e11.3)') 
1928 c     &'evdw1',i,j,evdwij
1929 c     &,iteli,itelj,aaa,evdw1
1930
1931 c              write (iout,'(a6,2i5,0pf7.3)') 'ees',i,j,eesij
1932 c          write(iout,'(2(2i3,2x),7(1pd12.4)/2(3(1pd12.4),5x)/)')
1933 c     &      iteli,i,itelj,j,aaa,bbb,ael6i,ael3i,
1934 c     &      1.0D0/dsqrt(rrmij),evdwij,eesij,
1935 c     &      xmedi,ymedi,zmedi,xj,yj,zj
1936 C
1937 C Calculate contributions to the Cartesian gradient.
1938 C
1939 #ifdef SPLITELE
1940           facvdw=-6*rrmij*(ev1+evdwij) 
1941           facel=-3*rrmij*(el1+eesij)
1942           fac1=fac
1943           erij(1)=xj*rmij
1944           erij(2)=yj*rmij
1945           erij(3)=zj*rmij
1946           if (calc_grad) then
1947 *
1948 * Radial derivatives. First process both termini of the fragment (i,j)
1949
1950           ggg(1)=facel*xj
1951           ggg(2)=facel*yj
1952           ggg(3)=facel*zj
1953           do k=1,3
1954             ghalf=0.5D0*ggg(k)
1955             gelc(k,i)=gelc(k,i)+ghalf
1956             gelc(k,j)=gelc(k,j)+ghalf
1957           enddo
1958 *
1959 * Loop over residues i+1 thru j-1.
1960 *
1961           do k=i+1,j-1
1962             do l=1,3
1963               gelc(l,k)=gelc(l,k)+ggg(l)
1964             enddo
1965           enddo
1966           ggg(1)=facvdw*xj
1967           ggg(2)=facvdw*yj
1968           ggg(3)=facvdw*zj
1969           do k=1,3
1970             ghalf=0.5D0*ggg(k)
1971             gvdwpp(k,i)=gvdwpp(k,i)+ghalf
1972             gvdwpp(k,j)=gvdwpp(k,j)+ghalf
1973           enddo
1974 *
1975 * Loop over residues i+1 thru j-1.
1976 *
1977           do k=i+1,j-1
1978             do l=1,3
1979               gvdwpp(l,k)=gvdwpp(l,k)+ggg(l)
1980             enddo
1981           enddo
1982 #else
1983           facvdw=ev1+evdwij 
1984           facel=el1+eesij  
1985           fac1=fac
1986           fac=-3*rrmij*(facvdw+facvdw+facel)
1987           erij(1)=xj*rmij
1988           erij(2)=yj*rmij
1989           erij(3)=zj*rmij
1990           if (calc_grad) then
1991 *
1992 * Radial derivatives. First process both termini of the fragment (i,j)
1993
1994           ggg(1)=fac*xj
1995           ggg(2)=fac*yj
1996           ggg(3)=fac*zj
1997           do k=1,3
1998             ghalf=0.5D0*ggg(k)
1999             gelc(k,i)=gelc(k,i)+ghalf
2000             gelc(k,j)=gelc(k,j)+ghalf
2001           enddo
2002 *
2003 * Loop over residues i+1 thru j-1.
2004 *
2005           do k=i+1,j-1
2006             do l=1,3
2007               gelc(l,k)=gelc(l,k)+ggg(l)
2008             enddo
2009           enddo
2010 #endif
2011 *
2012 * Angular part
2013 *          
2014           ecosa=2.0D0*fac3*fac1+fac4
2015           fac4=-3.0D0*fac4
2016           fac3=-6.0D0*fac3
2017           ecosb=(fac3*(fac1*cosg+cosb)+cosg*fac4)
2018           ecosg=(fac3*(fac1*cosb+cosg)+cosb*fac4)
2019           do k=1,3
2020             dcosb(k)=rmij*(dc_norm(k,i)-erij(k)*cosb)
2021             dcosg(k)=rmij*(dc_norm(k,j)-erij(k)*cosg)
2022           enddo
2023 cd        print '(2i3,2(3(1pd14.5),3x))',i,j,(dcosb(k),k=1,3),
2024 cd   &          (dcosg(k),k=1,3)
2025           do k=1,3
2026             ggg(k)=ecosb*dcosb(k)+ecosg*dcosg(k) 
2027           enddo
2028           do k=1,3
2029             ghalf=0.5D0*ggg(k)
2030             gelc(k,i)=gelc(k,i)+ghalf
2031      &               +(ecosa*(dc_norm(k,j)-cosa*dc_norm(k,i))
2032      &               + ecosb*(erij(k)-cosb*dc_norm(k,i)))*vbld_inv(i+1)
2033             gelc(k,j)=gelc(k,j)+ghalf
2034      &               +(ecosa*(dc_norm(k,i)-cosa*dc_norm(k,j))
2035      &               + ecosg*(erij(k)-cosg*dc_norm(k,j)))*vbld_inv(j+1)
2036           enddo
2037           do k=i+1,j-1
2038             do l=1,3
2039               gelc(l,k)=gelc(l,k)+ggg(l)
2040             enddo
2041           enddo
2042           endif
2043
2044           IF (wel_loc.gt.0.0d0 .or. wcorr4.gt.0.0d0 .or. wcorr5.gt.0.0d0
2045      &        .or. wcorr6.gt.0.0d0 .or. wturn3.gt.0.0d0 
2046      &        .or. wturn4.gt.0.0d0 .or. wturn6.gt.0.0d0) THEN
2047 C
2048 C 9/25/99 Mixed third-order local-electrostatic terms. The local-interaction 
2049 C   energy of a peptide unit is assumed in the form of a second-order 
2050 C   Fourier series in the angles lambda1 and lambda2 (see Nishikawa et al.
2051 C   Macromolecules, 1974, 7, 797-806 for definition). This correlation terms
2052 C   are computed for EVERY pair of non-contiguous peptide groups.
2053 C
2054           if (j.lt.nres-1) then
2055             j1=j+1
2056             j2=j-1
2057           else
2058             j1=j-1
2059             j2=j-2
2060           endif
2061           kkk=0
2062           do k=1,2
2063             do l=1,2
2064               kkk=kkk+1
2065               muij(kkk)=mu(k,i)*mu(l,j)
2066             enddo
2067           enddo  
2068 cd         write (iout,*) 'EELEC: i',i,' j',j
2069 cd          write (iout,*) 'j',j,' j1',j1,' j2',j2
2070 cd          write(iout,*) 'muij',muij
2071           ury=scalar(uy(1,i),erij)
2072           urz=scalar(uz(1,i),erij)
2073           vry=scalar(uy(1,j),erij)
2074           vrz=scalar(uz(1,j),erij)
2075           a22=scalar(uy(1,i),uy(1,j))-3*ury*vry
2076           a23=scalar(uy(1,i),uz(1,j))-3*ury*vrz
2077           a32=scalar(uz(1,i),uy(1,j))-3*urz*vry
2078           a33=scalar(uz(1,i),uz(1,j))-3*urz*vrz
2079 C For diagnostics only
2080 cd          a22=1.0d0
2081 cd          a23=1.0d0
2082 cd          a32=1.0d0
2083 cd          a33=1.0d0
2084           fac=dsqrt(-ael6i)*r3ij
2085 cd          write (2,*) 'fac=',fac
2086 C For diagnostics only
2087 cd          fac=1.0d0
2088           a22=a22*fac
2089           a23=a23*fac
2090           a32=a32*fac
2091           a33=a33*fac
2092 cd          write (iout,'(4i5,4f10.5)')
2093 cd     &     i,itortyp(itype(i)),j,itortyp(itype(j)),a22,a23,a32,a33
2094 cd          write (iout,'(6f10.5)') (muij(k),k=1,4),fac,eel_loc_ij
2095 cd          write (iout,'(2(3f10.5,5x)/2(3f10.5,5x))') (uy(k,i),k=1,3),
2096 cd     &      (uz(k,i),k=1,3),(uy(k,j),k=1,3),(uz(k,j),k=1,3)
2097 cd          write (iout,'(4f10.5)') 
2098 cd     &      scalar(uy(1,i),uy(1,j)),scalar(uy(1,i),uz(1,j)),
2099 cd     &      scalar(uz(1,i),uy(1,j)),scalar(uz(1,i),uz(1,j))
2100 cd          write (iout,'(4f10.5)') ury,urz,vry,vrz
2101 cd           write (iout,'(2i3,9f10.5/)') i,j,
2102 cd     &      fac22,a22,fac23,a23,fac32,a32,fac33,a33,eel_loc_ij
2103           if (calc_grad) then
2104 C Derivatives of the elements of A in virtual-bond vectors
2105           call unormderiv(erij(1),unmat(1,1),rmij,erder(1,1))
2106 cd          do k=1,3
2107 cd            do l=1,3
2108 cd              erder(k,l)=0.0d0
2109 cd            enddo
2110 cd          enddo
2111           do k=1,3
2112             uryg(k,1)=scalar(erder(1,k),uy(1,i))
2113             uryg(k,2)=scalar(uygrad(1,k,1,i),erij(1))
2114             uryg(k,3)=scalar(uygrad(1,k,2,i),erij(1))
2115             urzg(k,1)=scalar(erder(1,k),uz(1,i))
2116             urzg(k,2)=scalar(uzgrad(1,k,1,i),erij(1))
2117             urzg(k,3)=scalar(uzgrad(1,k,2,i),erij(1))
2118             vryg(k,1)=scalar(erder(1,k),uy(1,j))
2119             vryg(k,2)=scalar(uygrad(1,k,1,j),erij(1))
2120             vryg(k,3)=scalar(uygrad(1,k,2,j),erij(1))
2121             vrzg(k,1)=scalar(erder(1,k),uz(1,j))
2122             vrzg(k,2)=scalar(uzgrad(1,k,1,j),erij(1))
2123             vrzg(k,3)=scalar(uzgrad(1,k,2,j),erij(1))
2124           enddo
2125 cd          do k=1,3
2126 cd            do l=1,3
2127 cd              uryg(k,l)=0.0d0
2128 cd              urzg(k,l)=0.0d0
2129 cd              vryg(k,l)=0.0d0
2130 cd              vrzg(k,l)=0.0d0
2131 cd            enddo
2132 cd          enddo
2133 C Compute radial contributions to the gradient
2134           facr=-3.0d0*rrmij
2135           a22der=a22*facr
2136           a23der=a23*facr
2137           a32der=a32*facr
2138           a33der=a33*facr
2139 cd          a22der=0.0d0
2140 cd          a23der=0.0d0
2141 cd          a32der=0.0d0
2142 cd          a33der=0.0d0
2143           agg(1,1)=a22der*xj
2144           agg(2,1)=a22der*yj
2145           agg(3,1)=a22der*zj
2146           agg(1,2)=a23der*xj
2147           agg(2,2)=a23der*yj
2148           agg(3,2)=a23der*zj
2149           agg(1,3)=a32der*xj
2150           agg(2,3)=a32der*yj
2151           agg(3,3)=a32der*zj
2152           agg(1,4)=a33der*xj
2153           agg(2,4)=a33der*yj
2154           agg(3,4)=a33der*zj
2155 C Add the contributions coming from er
2156           fac3=-3.0d0*fac
2157           do k=1,3
2158             agg(k,1)=agg(k,1)+fac3*(uryg(k,1)*vry+vryg(k,1)*ury)
2159             agg(k,2)=agg(k,2)+fac3*(uryg(k,1)*vrz+vrzg(k,1)*ury)
2160             agg(k,3)=agg(k,3)+fac3*(urzg(k,1)*vry+vryg(k,1)*urz)
2161             agg(k,4)=agg(k,4)+fac3*(urzg(k,1)*vrz+vrzg(k,1)*urz)
2162           enddo
2163           do k=1,3
2164 C Derivatives in DC(i) 
2165             ghalf1=0.5d0*agg(k,1)
2166             ghalf2=0.5d0*agg(k,2)
2167             ghalf3=0.5d0*agg(k,3)
2168             ghalf4=0.5d0*agg(k,4)
2169             aggi(k,1)=fac*(scalar(uygrad(1,k,1,i),uy(1,j))
2170      &      -3.0d0*uryg(k,2)*vry)+ghalf1
2171             aggi(k,2)=fac*(scalar(uygrad(1,k,1,i),uz(1,j))
2172      &      -3.0d0*uryg(k,2)*vrz)+ghalf2
2173             aggi(k,3)=fac*(scalar(uzgrad(1,k,1,i),uy(1,j))
2174      &      -3.0d0*urzg(k,2)*vry)+ghalf3
2175             aggi(k,4)=fac*(scalar(uzgrad(1,k,1,i),uz(1,j))
2176      &      -3.0d0*urzg(k,2)*vrz)+ghalf4
2177 C Derivatives in DC(i+1)
2178             aggi1(k,1)=fac*(scalar(uygrad(1,k,2,i),uy(1,j))
2179      &      -3.0d0*uryg(k,3)*vry)+agg(k,1)
2180             aggi1(k,2)=fac*(scalar(uygrad(1,k,2,i),uz(1,j))
2181      &      -3.0d0*uryg(k,3)*vrz)+agg(k,2)
2182             aggi1(k,3)=fac*(scalar(uzgrad(1,k,2,i),uy(1,j))
2183      &      -3.0d0*urzg(k,3)*vry)+agg(k,3)
2184             aggi1(k,4)=fac*(scalar(uzgrad(1,k,2,i),uz(1,j))
2185      &      -3.0d0*urzg(k,3)*vrz)+agg(k,4)
2186 C Derivatives in DC(j)
2187             aggj(k,1)=fac*(scalar(uygrad(1,k,1,j),uy(1,i))
2188      &      -3.0d0*vryg(k,2)*ury)+ghalf1
2189             aggj(k,2)=fac*(scalar(uzgrad(1,k,1,j),uy(1,i))
2190      &      -3.0d0*vrzg(k,2)*ury)+ghalf2
2191             aggj(k,3)=fac*(scalar(uygrad(1,k,1,j),uz(1,i))
2192      &      -3.0d0*vryg(k,2)*urz)+ghalf3
2193             aggj(k,4)=fac*(scalar(uzgrad(1,k,1,j),uz(1,i)) 
2194      &      -3.0d0*vrzg(k,2)*urz)+ghalf4
2195 C Derivatives in DC(j+1) or DC(nres-1)
2196             aggj1(k,1)=fac*(scalar(uygrad(1,k,2,j),uy(1,i))
2197      &      -3.0d0*vryg(k,3)*ury)
2198             aggj1(k,2)=fac*(scalar(uzgrad(1,k,2,j),uy(1,i))
2199      &      -3.0d0*vrzg(k,3)*ury)
2200             aggj1(k,3)=fac*(scalar(uygrad(1,k,2,j),uz(1,i))
2201      &      -3.0d0*vryg(k,3)*urz)
2202             aggj1(k,4)=fac*(scalar(uzgrad(1,k,2,j),uz(1,i)) 
2203      &      -3.0d0*vrzg(k,3)*urz)
2204 cd            aggi(k,1)=ghalf1
2205 cd            aggi(k,2)=ghalf2
2206 cd            aggi(k,3)=ghalf3
2207 cd            aggi(k,4)=ghalf4
2208 C Derivatives in DC(i+1)
2209 cd            aggi1(k,1)=agg(k,1)
2210 cd            aggi1(k,2)=agg(k,2)
2211 cd            aggi1(k,3)=agg(k,3)
2212 cd            aggi1(k,4)=agg(k,4)
2213 C Derivatives in DC(j)
2214 cd            aggj(k,1)=ghalf1
2215 cd            aggj(k,2)=ghalf2
2216 cd            aggj(k,3)=ghalf3
2217 cd            aggj(k,4)=ghalf4
2218 C Derivatives in DC(j+1)
2219 cd            aggj1(k,1)=0.0d0
2220 cd            aggj1(k,2)=0.0d0
2221 cd            aggj1(k,3)=0.0d0
2222 cd            aggj1(k,4)=0.0d0
2223             if (j.eq.nres-1 .and. i.lt.j-2) then
2224               do l=1,4
2225                 aggj1(k,l)=aggj1(k,l)+agg(k,l)
2226 cd                aggj1(k,l)=agg(k,l)
2227               enddo
2228             endif
2229           enddo
2230           endif
2231 c          goto 11111
2232 C Check the loc-el terms by numerical integration
2233           acipa(1,1)=a22
2234           acipa(1,2)=a23
2235           acipa(2,1)=a32
2236           acipa(2,2)=a33
2237           a22=-a22
2238           a23=-a23
2239           do l=1,2
2240             do k=1,3
2241               agg(k,l)=-agg(k,l)
2242               aggi(k,l)=-aggi(k,l)
2243               aggi1(k,l)=-aggi1(k,l)
2244               aggj(k,l)=-aggj(k,l)
2245               aggj1(k,l)=-aggj1(k,l)
2246             enddo
2247           enddo
2248           if (j.lt.nres-1) then
2249             a22=-a22
2250             a32=-a32
2251             do l=1,3,2
2252               do k=1,3
2253                 agg(k,l)=-agg(k,l)
2254                 aggi(k,l)=-aggi(k,l)
2255                 aggi1(k,l)=-aggi1(k,l)
2256                 aggj(k,l)=-aggj(k,l)
2257                 aggj1(k,l)=-aggj1(k,l)
2258               enddo
2259             enddo
2260           else
2261             a22=-a22
2262             a23=-a23
2263             a32=-a32
2264             a33=-a33
2265             do l=1,4
2266               do k=1,3
2267                 agg(k,l)=-agg(k,l)
2268                 aggi(k,l)=-aggi(k,l)
2269                 aggi1(k,l)=-aggi1(k,l)
2270                 aggj(k,l)=-aggj(k,l)
2271                 aggj1(k,l)=-aggj1(k,l)
2272               enddo
2273             enddo 
2274           endif    
2275           ENDIF ! WCORR
2276 11111     continue
2277           IF (wel_loc.gt.0.0d0) THEN
2278 C Contribution to the local-electrostatic energy coming from the i-j pair
2279           eel_loc_ij=a22*muij(1)+a23*muij(2)+a32*muij(3)
2280      &     +a33*muij(4)
2281 c          write (iout,*) 'i',i,' j',j,' eel_loc_ij',eel_loc_ij
2282 c          write (iout,'(a6,2i5,0pf7.3)')
2283 c     &            'eelloc',i,j,eel_loc_ij
2284 c          write (iout,*) a22,muij(1),a23,muij(2),a32,muij(3)
2285           eel_loc=eel_loc+eel_loc_ij
2286 C Partial derivatives in virtual-bond dihedral angles gamma
2287           if (calc_grad) then
2288           if (i.gt.1)
2289      &    gel_loc_loc(i-1)=gel_loc_loc(i-1)+ 
2290      &            a22*muder(1,i)*mu(1,j)+a23*muder(1,i)*mu(2,j)
2291      &           +a32*muder(2,i)*mu(1,j)+a33*muder(2,i)*mu(2,j)
2292           gel_loc_loc(j-1)=gel_loc_loc(j-1)+ 
2293      &            a22*mu(1,i)*muder(1,j)+a23*mu(1,i)*muder(2,j)
2294      &           +a32*mu(2,i)*muder(1,j)+a33*mu(2,i)*muder(2,j)
2295 cd          call checkint3(i,j,mu1,mu2,a22,a23,a32,a33,acipa,eel_loc_ij)
2296 cd          write(iout,*) 'agg  ',agg
2297 cd          write(iout,*) 'aggi ',aggi
2298 cd          write(iout,*) 'aggi1',aggi1
2299 cd          write(iout,*) 'aggj ',aggj
2300 cd          write(iout,*) 'aggj1',aggj1
2301
2302 C Derivatives of eello in DC(i+1) thru DC(j-1) or DC(nres-2)
2303           do l=1,3
2304             ggg(l)=agg(l,1)*muij(1)+
2305      &          agg(l,2)*muij(2)+agg(l,3)*muij(3)+agg(l,4)*muij(4)
2306           enddo
2307           do k=i+2,j2
2308             do l=1,3
2309               gel_loc(l,k)=gel_loc(l,k)+ggg(l)
2310             enddo
2311           enddo
2312 C Remaining derivatives of eello
2313           do l=1,3
2314             gel_loc(l,i)=gel_loc(l,i)+aggi(l,1)*muij(1)+
2315      &          aggi(l,2)*muij(2)+aggi(l,3)*muij(3)+aggi(l,4)*muij(4)
2316             gel_loc(l,i+1)=gel_loc(l,i+1)+aggi1(l,1)*muij(1)+
2317      &          aggi1(l,2)*muij(2)+aggi1(l,3)*muij(3)+aggi1(l,4)*muij(4)
2318             gel_loc(l,j)=gel_loc(l,j)+aggj(l,1)*muij(1)+
2319      &          aggj(l,2)*muij(2)+aggj(l,3)*muij(3)+aggj(l,4)*muij(4)
2320             gel_loc(l,j1)=gel_loc(l,j1)+aggj1(l,1)*muij(1)+
2321      &          aggj1(l,2)*muij(2)+aggj1(l,3)*muij(3)+aggj1(l,4)*muij(4)
2322           enddo
2323           endif
2324           ENDIF
2325           if (wturn3.gt.0.0d0 .or. wturn4.gt.0.0d0) then
2326 C Contributions from turns
2327             a_temp(1,1)=a22
2328             a_temp(1,2)=a23
2329             a_temp(2,1)=a32
2330             a_temp(2,2)=a33
2331             call eturn34(i,j,eello_turn3,eello_turn4)
2332           endif
2333 C Change 12/26/95 to calculate four-body contributions to H-bonding energy
2334           if (j.gt.i+1 .and. num_conti.le.maxconts) then
2335 C
2336 C Calculate the contact function. The ith column of the array JCONT will 
2337 C contain the numbers of atoms that make contacts with the atom I (of numbers
2338 C greater than I). The arrays FACONT and GACONT will contain the values of
2339 C the contact function and its derivative.
2340 c           r0ij=1.02D0*rpp(iteli,itelj)
2341 c           r0ij=1.11D0*rpp(iteli,itelj)
2342             r0ij=2.20D0*rpp(iteli,itelj)
2343 c           r0ij=1.55D0*rpp(iteli,itelj)
2344             call gcont(rij,r0ij,1.0D0,0.2d0*r0ij,fcont,fprimcont)
2345             if (fcont.gt.0.0D0) then
2346               num_conti=num_conti+1
2347               if (num_conti.gt.maxconts) then
2348                 write (iout,*) 'WARNING - max. # of contacts exceeded;',
2349      &                         ' will skip next contacts for this conf.'
2350               else
2351                 jcont_hb(num_conti,i)=j
2352                 IF (wcorr4.gt.0.0d0 .or. wcorr5.gt.0.0d0 .or. 
2353      &          wcorr6.gt.0.0d0 .or. wturn6.gt.0.0d0) THEN
2354 C 9/30/99 (AL) - store components necessary to evaluate higher-order loc-el
2355 C  terms.
2356                 d_cont(num_conti,i)=rij
2357 cd                write (2,'(3e15.5)') rij,r0ij+0.2d0*r0ij,rij
2358 C     --- Electrostatic-interaction matrix --- 
2359                 a_chuj(1,1,num_conti,i)=a22
2360                 a_chuj(1,2,num_conti,i)=a23
2361                 a_chuj(2,1,num_conti,i)=a32
2362                 a_chuj(2,2,num_conti,i)=a33
2363 C     --- Gradient of rij
2364                 do kkk=1,3
2365                   grij_hb_cont(kkk,num_conti,i)=erij(kkk)
2366                 enddo
2367 c             if (i.eq.1) then
2368 c                a_chuj(1,1,num_conti,i)=-0.61d0
2369 c                a_chuj(1,2,num_conti,i)= 0.4d0
2370 c                a_chuj(2,1,num_conti,i)= 0.65d0
2371 c                a_chuj(2,2,num_conti,i)= 0.50d0
2372 c             else if (i.eq.2) then
2373 c                a_chuj(1,1,num_conti,i)= 0.0d0
2374 c                a_chuj(1,2,num_conti,i)= 0.0d0
2375 c                a_chuj(2,1,num_conti,i)= 0.0d0
2376 c                a_chuj(2,2,num_conti,i)= 0.0d0
2377 c             endif
2378 C     --- and its gradients
2379 cd                write (iout,*) 'i',i,' j',j
2380 cd                do kkk=1,3
2381 cd                write (iout,*) 'iii 1 kkk',kkk
2382 cd                write (iout,*) agg(kkk,:)
2383 cd                enddo
2384 cd                do kkk=1,3
2385 cd                write (iout,*) 'iii 2 kkk',kkk
2386 cd                write (iout,*) aggi(kkk,:)
2387 cd                enddo
2388 cd                do kkk=1,3
2389 cd                write (iout,*) 'iii 3 kkk',kkk
2390 cd                write (iout,*) aggi1(kkk,:)
2391 cd                enddo
2392 cd                do kkk=1,3
2393 cd                write (iout,*) 'iii 4 kkk',kkk
2394 cd                write (iout,*) aggj(kkk,:)
2395 cd                enddo
2396 cd                do kkk=1,3
2397 cd                write (iout,*) 'iii 5 kkk',kkk
2398 cd                write (iout,*) aggj1(kkk,:)
2399 cd                enddo
2400                 kkll=0
2401                 do k=1,2
2402                   do l=1,2
2403                     kkll=kkll+1
2404                     do m=1,3
2405                       a_chuj_der(k,l,m,1,num_conti,i)=agg(m,kkll)
2406                       a_chuj_der(k,l,m,2,num_conti,i)=aggi(m,kkll)
2407                       a_chuj_der(k,l,m,3,num_conti,i)=aggi1(m,kkll)
2408                       a_chuj_der(k,l,m,4,num_conti,i)=aggj(m,kkll)
2409                       a_chuj_der(k,l,m,5,num_conti,i)=aggj1(m,kkll)
2410 c                      do mm=1,5
2411 c                      a_chuj_der(k,l,m,mm,num_conti,i)=0.0d0
2412 c                      enddo
2413                     enddo
2414                   enddo
2415                 enddo
2416                 ENDIF
2417                 IF (wcorr4.eq.0.0d0 .and. wcorr.gt.0.0d0) THEN
2418 C Calculate contact energies
2419                 cosa4=4.0D0*cosa
2420                 wij=cosa-3.0D0*cosb*cosg
2421                 cosbg1=cosb+cosg
2422                 cosbg2=cosb-cosg
2423 c               fac3=dsqrt(-ael6i)/r0ij**3     
2424                 fac3=dsqrt(-ael6i)*r3ij
2425                 ees0pij=dsqrt(4.0D0+cosa4+wij*wij-3.0D0*cosbg1*cosbg1)
2426                 ees0mij=dsqrt(4.0D0-cosa4+wij*wij-3.0D0*cosbg2*cosbg2)
2427 c               ees0mij=0.0D0
2428                 ees0p(num_conti,i)=0.5D0*fac3*(ees0pij+ees0mij)
2429                 ees0m(num_conti,i)=0.5D0*fac3*(ees0pij-ees0mij)
2430 C Diagnostics. Comment out or remove after debugging!
2431 c               ees0p(num_conti,i)=0.5D0*fac3*ees0pij
2432 c               ees0m(num_conti,i)=0.5D0*fac3*ees0mij
2433 c               ees0m(num_conti,i)=0.0D0
2434 C End diagnostics.
2435 c                write (iout,*) 'i=',i,' j=',j,' rij=',rij,' r0ij=',r0ij,
2436 c     & ' ees0ij=',ees0p(num_conti,i),ees0m(num_conti,i),' fcont=',fcont
2437                 facont_hb(num_conti,i)=fcont
2438                 if (calc_grad) then
2439 C Angular derivatives of the contact function
2440                 ees0pij1=fac3/ees0pij 
2441                 ees0mij1=fac3/ees0mij
2442                 fac3p=-3.0D0*fac3*rrmij
2443                 ees0pijp=0.5D0*fac3p*(ees0pij+ees0mij)
2444                 ees0mijp=0.5D0*fac3p*(ees0pij-ees0mij)
2445 c               ees0mij1=0.0D0
2446                 ecosa1=       ees0pij1*( 1.0D0+0.5D0*wij)
2447                 ecosb1=-1.5D0*ees0pij1*(wij*cosg+cosbg1)
2448                 ecosg1=-1.5D0*ees0pij1*(wij*cosb+cosbg1)
2449                 ecosa2=       ees0mij1*(-1.0D0+0.5D0*wij)
2450                 ecosb2=-1.5D0*ees0mij1*(wij*cosg+cosbg2) 
2451                 ecosg2=-1.5D0*ees0mij1*(wij*cosb-cosbg2)
2452                 ecosap=ecosa1+ecosa2
2453                 ecosbp=ecosb1+ecosb2
2454                 ecosgp=ecosg1+ecosg2
2455                 ecosam=ecosa1-ecosa2
2456                 ecosbm=ecosb1-ecosb2
2457                 ecosgm=ecosg1-ecosg2
2458 C Diagnostics
2459 c               ecosap=ecosa1
2460 c               ecosbp=ecosb1
2461 c               ecosgp=ecosg1
2462 c               ecosam=0.0D0
2463 c               ecosbm=0.0D0
2464 c               ecosgm=0.0D0
2465 C End diagnostics
2466                 fprimcont=fprimcont/rij
2467 cd              facont_hb(num_conti,i)=1.0D0
2468 C Following line is for diagnostics.
2469 cd              fprimcont=0.0D0
2470                 do k=1,3
2471                   dcosb(k)=rmij*(dc_norm(k,i)-erij(k)*cosb)
2472                   dcosg(k)=rmij*(dc_norm(k,j)-erij(k)*cosg)
2473                 enddo
2474                 do k=1,3
2475                   gggp(k)=ecosbp*dcosb(k)+ecosgp*dcosg(k)
2476                   gggm(k)=ecosbm*dcosb(k)+ecosgm*dcosg(k)
2477                 enddo
2478                 gggp(1)=gggp(1)+ees0pijp*xj
2479                 gggp(2)=gggp(2)+ees0pijp*yj
2480                 gggp(3)=gggp(3)+ees0pijp*zj
2481                 gggm(1)=gggm(1)+ees0mijp*xj
2482                 gggm(2)=gggm(2)+ees0mijp*yj
2483                 gggm(3)=gggm(3)+ees0mijp*zj
2484 C Derivatives due to the contact function
2485                 gacont_hbr(1,num_conti,i)=fprimcont*xj
2486                 gacont_hbr(2,num_conti,i)=fprimcont*yj
2487                 gacont_hbr(3,num_conti,i)=fprimcont*zj
2488                 do k=1,3
2489                   ghalfp=0.5D0*gggp(k)
2490                   ghalfm=0.5D0*gggm(k)
2491                   gacontp_hb1(k,num_conti,i)=ghalfp
2492      &              +(ecosap*(dc_norm(k,j)-cosa*dc_norm(k,i))
2493      &              + ecosbp*(erij(k)-cosb*dc_norm(k,i)))*vbld_inv(i+1)
2494                   gacontp_hb2(k,num_conti,i)=ghalfp
2495      &              +(ecosap*(dc_norm(k,i)-cosa*dc_norm(k,j))
2496      &              + ecosgp*(erij(k)-cosg*dc_norm(k,j)))*vbld_inv(j+1)
2497                   gacontp_hb3(k,num_conti,i)=gggp(k)
2498                   gacontm_hb1(k,num_conti,i)=ghalfm
2499      &              +(ecosam*(dc_norm(k,j)-cosa*dc_norm(k,i))
2500      &              + ecosbm*(erij(k)-cosb*dc_norm(k,i)))*vbld_inv(i+1)
2501                   gacontm_hb2(k,num_conti,i)=ghalfm
2502      &              +(ecosam*(dc_norm(k,i)-cosa*dc_norm(k,j))
2503      &              + ecosgm*(erij(k)-cosg*dc_norm(k,j)))*vbld_inv(j+1)
2504                   gacontm_hb3(k,num_conti,i)=gggm(k)
2505                 enddo
2506                 endif
2507 C Diagnostics. Comment out or remove after debugging!
2508 cdiag           do k=1,3
2509 cdiag             gacontp_hb1(k,num_conti,i)=0.0D0
2510 cdiag             gacontp_hb2(k,num_conti,i)=0.0D0
2511 cdiag             gacontp_hb3(k,num_conti,i)=0.0D0
2512 cdiag             gacontm_hb1(k,num_conti,i)=0.0D0
2513 cdiag             gacontm_hb2(k,num_conti,i)=0.0D0
2514 cdiag             gacontm_hb3(k,num_conti,i)=0.0D0
2515 cdiag           enddo
2516               ENDIF ! wcorr
2517               endif  ! num_conti.le.maxconts
2518             endif  ! fcont.gt.0
2519           endif    ! j.gt.i+1
2520  1216     continue
2521         enddo ! j
2522         num_cont_hb(i)=num_conti
2523  1215   continue
2524       enddo   ! i
2525 cd      do i=1,nres
2526 cd        write (iout,'(i3,3f10.5,5x,3f10.5)') 
2527 cd     &     i,(gel_loc(k,i),k=1,3),gel_loc_loc(i)
2528 cd      enddo
2529 c 12/7/99 Adam eello_turn3 will be considered as a separate energy term
2530 ccc      eel_loc=eel_loc+eello_turn3
2531       return
2532       end
2533 C-----------------------------------------------------------------------------
2534       subroutine eturn34(i,j,eello_turn3,eello_turn4)
2535 C Third- and fourth-order contributions from turns
2536       implicit real*8 (a-h,o-z)
2537       include 'DIMENSIONS'
2538       include 'DIMENSIONS.ZSCOPT'
2539       include 'COMMON.IOUNITS'
2540       include 'COMMON.GEO'
2541       include 'COMMON.VAR'
2542       include 'COMMON.LOCAL'
2543       include 'COMMON.CHAIN'
2544       include 'COMMON.DERIV'
2545       include 'COMMON.INTERACT'
2546       include 'COMMON.CONTACTS'
2547       include 'COMMON.TORSION'
2548       include 'COMMON.VECTORS'
2549       include 'COMMON.FFIELD'
2550       dimension ggg(3)
2551       double precision auxmat(2,2),auxmat1(2,2),auxmat2(2,2),pizda(2,2),
2552      &  e1t(2,2),e2t(2,2),e3t(2,2),e1tder(2,2),e2tder(2,2),e3tder(2,2),
2553      &  e1a(2,2),ae3(2,2),ae3e2(2,2),auxvec(2),auxvec1(2)
2554       double precision agg(3,4),aggi(3,4),aggi1(3,4),
2555      &    aggj(3,4),aggj1(3,4),a_temp(2,2)
2556       common /locel/ a_temp,agg,aggi,aggi1,aggj,aggj1,j1,j2
2557       if (j.eq.i+2) then
2558 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
2559 C
2560 C               Third-order contributions
2561 C        
2562 C                 (i+2)o----(i+3)
2563 C                      | |
2564 C                      | |
2565 C                 (i+1)o----i
2566 C
2567 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC   
2568 cd        call checkint_turn3(i,a_temp,eello_turn3_num)
2569         call matmat2(EUg(1,1,i+1),EUg(1,1,i+2),auxmat(1,1))
2570         call transpose2(auxmat(1,1),auxmat1(1,1))
2571         call matmat2(a_temp(1,1),auxmat1(1,1),pizda(1,1))
2572         eello_turn3=eello_turn3+0.5d0*(pizda(1,1)+pizda(2,2))
2573 cd        write (2,*) 'i,',i,' j',j,'eello_turn3',
2574 cd     &    0.5d0*(pizda(1,1)+pizda(2,2)),
2575 cd     &    ' eello_turn3_num',4*eello_turn3_num
2576         if (calc_grad) then
2577 C Derivatives in gamma(i)
2578         call matmat2(EUgder(1,1,i+1),EUg(1,1,i+2),auxmat2(1,1))
2579         call transpose2(auxmat2(1,1),pizda(1,1))
2580         call matmat2(a_temp(1,1),pizda(1,1),pizda(1,1))
2581         gel_loc_turn3(i)=gel_loc_turn3(i)+0.5d0*(pizda(1,1)+pizda(2,2))
2582 C Derivatives in gamma(i+1)
2583         call matmat2(EUg(1,1,i+1),EUgder(1,1,i+2),auxmat2(1,1))
2584         call transpose2(auxmat2(1,1),pizda(1,1))
2585         call matmat2(a_temp(1,1),pizda(1,1),pizda(1,1))
2586         gel_loc_turn3(i+1)=gel_loc_turn3(i+1)
2587      &    +0.5d0*(pizda(1,1)+pizda(2,2))
2588 C Cartesian derivatives
2589         do l=1,3
2590           a_temp(1,1)=aggi(l,1)
2591           a_temp(1,2)=aggi(l,2)
2592           a_temp(2,1)=aggi(l,3)
2593           a_temp(2,2)=aggi(l,4)
2594           call matmat2(a_temp(1,1),auxmat1(1,1),pizda(1,1))
2595           gcorr3_turn(l,i)=gcorr3_turn(l,i)
2596      &      +0.5d0*(pizda(1,1)+pizda(2,2))
2597           a_temp(1,1)=aggi1(l,1)
2598           a_temp(1,2)=aggi1(l,2)
2599           a_temp(2,1)=aggi1(l,3)
2600           a_temp(2,2)=aggi1(l,4)
2601           call matmat2(a_temp(1,1),auxmat1(1,1),pizda(1,1))
2602           gcorr3_turn(l,i+1)=gcorr3_turn(l,i+1)
2603      &      +0.5d0*(pizda(1,1)+pizda(2,2))
2604           a_temp(1,1)=aggj(l,1)
2605           a_temp(1,2)=aggj(l,2)
2606           a_temp(2,1)=aggj(l,3)
2607           a_temp(2,2)=aggj(l,4)
2608           call matmat2(a_temp(1,1),auxmat1(1,1),pizda(1,1))
2609           gcorr3_turn(l,j)=gcorr3_turn(l,j)
2610      &      +0.5d0*(pizda(1,1)+pizda(2,2))
2611           a_temp(1,1)=aggj1(l,1)
2612           a_temp(1,2)=aggj1(l,2)
2613           a_temp(2,1)=aggj1(l,3)
2614           a_temp(2,2)=aggj1(l,4)
2615           call matmat2(a_temp(1,1),auxmat1(1,1),pizda(1,1))
2616           gcorr3_turn(l,j1)=gcorr3_turn(l,j1)
2617      &      +0.5d0*(pizda(1,1)+pizda(2,2))
2618         enddo
2619         endif
2620       else if (j.eq.i+3 .and. itype(i+2).ne.ntyp1) then
2621 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
2622 C
2623 C               Fourth-order contributions
2624 C        
2625 C                 (i+3)o----(i+4)
2626 C                     /  |
2627 C               (i+2)o   |
2628 C                     \  |
2629 C                 (i+1)o----i
2630 C
2631 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC   
2632 cd        call checkint_turn4(i,a_temp,eello_turn4_num)
2633         iti1=itortyp(itype(i+1))
2634         iti2=itortyp(itype(i+2))
2635         iti3=itortyp(itype(i+3))
2636         call transpose2(EUg(1,1,i+1),e1t(1,1))
2637         call transpose2(Eug(1,1,i+2),e2t(1,1))
2638         call transpose2(Eug(1,1,i+3),e3t(1,1))
2639         call matmat2(e1t(1,1),a_temp(1,1),e1a(1,1))
2640         call matvec2(e1a(1,1),Ub2(1,i+3),auxvec(1))
2641         s1=scalar2(b1(1,iti2),auxvec(1))
2642         call matmat2(a_temp(1,1),e3t(1,1),ae3(1,1))
2643         call matvec2(ae3(1,1),Ub2(1,i+2),auxvec(1)) 
2644         s2=scalar2(b1(1,iti1),auxvec(1))
2645         call matmat2(ae3(1,1),e2t(1,1),ae3e2(1,1))
2646         call matmat2(ae3e2(1,1),e1t(1,1),pizda(1,1))
2647         s3=0.5d0*(pizda(1,1)+pizda(2,2))
2648         eello_turn4=eello_turn4-(s1+s2+s3)
2649 cd        write (2,*) 'i,',i,' j',j,'eello_turn4',-(s1+s2+s3),
2650 cd     &    ' eello_turn4_num',8*eello_turn4_num
2651 C Derivatives in gamma(i)
2652         if (calc_grad) then
2653         call transpose2(EUgder(1,1,i+1),e1tder(1,1))
2654         call matmat2(e1tder(1,1),a_temp(1,1),auxmat(1,1))
2655         call matvec2(auxmat(1,1),Ub2(1,i+3),auxvec(1))
2656         s1=scalar2(b1(1,iti2),auxvec(1))
2657         call matmat2(ae3e2(1,1),e1tder(1,1),pizda(1,1))
2658         s3=0.5d0*(pizda(1,1)+pizda(2,2))
2659         gel_loc_turn4(i)=gel_loc_turn4(i)-(s1+s3)
2660 C Derivatives in gamma(i+1)
2661         call transpose2(EUgder(1,1,i+2),e2tder(1,1))
2662         call matvec2(ae3(1,1),Ub2der(1,i+2),auxvec(1)) 
2663         s2=scalar2(b1(1,iti1),auxvec(1))
2664         call matmat2(ae3(1,1),e2tder(1,1),auxmat(1,1))
2665         call matmat2(auxmat(1,1),e1t(1,1),pizda(1,1))
2666         s3=0.5d0*(pizda(1,1)+pizda(2,2))
2667         gel_loc_turn4(i+1)=gel_loc_turn4(i+1)-(s2+s3)
2668 C Derivatives in gamma(i+2)
2669         call transpose2(EUgder(1,1,i+3),e3tder(1,1))
2670         call matvec2(e1a(1,1),Ub2der(1,i+3),auxvec(1))
2671         s1=scalar2(b1(1,iti2),auxvec(1))
2672         call matmat2(a_temp(1,1),e3tder(1,1),auxmat(1,1))
2673         call matvec2(auxmat(1,1),Ub2(1,i+2),auxvec(1)) 
2674         s2=scalar2(b1(1,iti1),auxvec(1))
2675         call matmat2(auxmat(1,1),e2t(1,1),auxmat(1,1))
2676         call matmat2(auxmat(1,1),e1t(1,1),pizda(1,1))
2677         s3=0.5d0*(pizda(1,1)+pizda(2,2))
2678         gel_loc_turn4(i+2)=gel_loc_turn4(i+2)-(s1+s2+s3)
2679 C Cartesian derivatives
2680 C Derivatives of this turn contributions in DC(i+2)
2681         if (j.lt.nres-1) then
2682           do l=1,3
2683             a_temp(1,1)=agg(l,1)
2684             a_temp(1,2)=agg(l,2)
2685             a_temp(2,1)=agg(l,3)
2686             a_temp(2,2)=agg(l,4)
2687             call matmat2(e1t(1,1),a_temp(1,1),e1a(1,1))
2688             call matvec2(e1a(1,1),Ub2(1,i+3),auxvec(1))
2689             s1=scalar2(b1(1,iti2),auxvec(1))
2690             call matmat2(a_temp(1,1),e3t(1,1),ae3(1,1))
2691             call matvec2(ae3(1,1),Ub2(1,i+2),auxvec(1)) 
2692             s2=scalar2(b1(1,iti1),auxvec(1))
2693             call matmat2(ae3(1,1),e2t(1,1),ae3e2(1,1))
2694             call matmat2(ae3e2(1,1),e1t(1,1),pizda(1,1))
2695             s3=0.5d0*(pizda(1,1)+pizda(2,2))
2696             ggg(l)=-(s1+s2+s3)
2697             gcorr4_turn(l,i+2)=gcorr4_turn(l,i+2)-(s1+s2+s3)
2698           enddo
2699         endif
2700 C Remaining derivatives of this turn contribution
2701         do l=1,3
2702           a_temp(1,1)=aggi(l,1)
2703           a_temp(1,2)=aggi(l,2)
2704           a_temp(2,1)=aggi(l,3)
2705           a_temp(2,2)=aggi(l,4)
2706           call matmat2(e1t(1,1),a_temp(1,1),e1a(1,1))
2707           call matvec2(e1a(1,1),Ub2(1,i+3),auxvec(1))
2708           s1=scalar2(b1(1,iti2),auxvec(1))
2709           call matmat2(a_temp(1,1),e3t(1,1),ae3(1,1))
2710           call matvec2(ae3(1,1),Ub2(1,i+2),auxvec(1)) 
2711           s2=scalar2(b1(1,iti1),auxvec(1))
2712           call matmat2(ae3(1,1),e2t(1,1),ae3e2(1,1))
2713           call matmat2(ae3e2(1,1),e1t(1,1),pizda(1,1))
2714           s3=0.5d0*(pizda(1,1)+pizda(2,2))
2715           gcorr4_turn(l,i)=gcorr4_turn(l,i)-(s1+s2+s3)
2716           a_temp(1,1)=aggi1(l,1)
2717           a_temp(1,2)=aggi1(l,2)
2718           a_temp(2,1)=aggi1(l,3)
2719           a_temp(2,2)=aggi1(l,4)
2720           call matmat2(e1t(1,1),a_temp(1,1),e1a(1,1))
2721           call matvec2(e1a(1,1),Ub2(1,i+3),auxvec(1))
2722           s1=scalar2(b1(1,iti2),auxvec(1))
2723           call matmat2(a_temp(1,1),e3t(1,1),ae3(1,1))
2724           call matvec2(ae3(1,1),Ub2(1,i+2),auxvec(1)) 
2725           s2=scalar2(b1(1,iti1),auxvec(1))
2726           call matmat2(ae3(1,1),e2t(1,1),ae3e2(1,1))
2727           call matmat2(ae3e2(1,1),e1t(1,1),pizda(1,1))
2728           s3=0.5d0*(pizda(1,1)+pizda(2,2))
2729           gcorr4_turn(l,i+1)=gcorr4_turn(l,i+1)-(s1+s2+s3)
2730           a_temp(1,1)=aggj(l,1)
2731           a_temp(1,2)=aggj(l,2)
2732           a_temp(2,1)=aggj(l,3)
2733           a_temp(2,2)=aggj(l,4)
2734           call matmat2(e1t(1,1),a_temp(1,1),e1a(1,1))
2735           call matvec2(e1a(1,1),Ub2(1,i+3),auxvec(1))
2736           s1=scalar2(b1(1,iti2),auxvec(1))
2737           call matmat2(a_temp(1,1),e3t(1,1),ae3(1,1))
2738           call matvec2(ae3(1,1),Ub2(1,i+2),auxvec(1)) 
2739           s2=scalar2(b1(1,iti1),auxvec(1))
2740           call matmat2(ae3(1,1),e2t(1,1),ae3e2(1,1))
2741           call matmat2(ae3e2(1,1),e1t(1,1),pizda(1,1))
2742           s3=0.5d0*(pizda(1,1)+pizda(2,2))
2743           gcorr4_turn(l,j)=gcorr4_turn(l,j)-(s1+s2+s3)
2744           a_temp(1,1)=aggj1(l,1)
2745           a_temp(1,2)=aggj1(l,2)
2746           a_temp(2,1)=aggj1(l,3)
2747           a_temp(2,2)=aggj1(l,4)
2748           call matmat2(e1t(1,1),a_temp(1,1),e1a(1,1))
2749           call matvec2(e1a(1,1),Ub2(1,i+3),auxvec(1))
2750           s1=scalar2(b1(1,iti2),auxvec(1))
2751           call matmat2(a_temp(1,1),e3t(1,1),ae3(1,1))
2752           call matvec2(ae3(1,1),Ub2(1,i+2),auxvec(1)) 
2753           s2=scalar2(b1(1,iti1),auxvec(1))
2754           call matmat2(ae3(1,1),e2t(1,1),ae3e2(1,1))
2755           call matmat2(ae3e2(1,1),e1t(1,1),pizda(1,1))
2756           s3=0.5d0*(pizda(1,1)+pizda(2,2))
2757           gcorr4_turn(l,j1)=gcorr4_turn(l,j1)-(s1+s2+s3)
2758         enddo
2759         endif
2760       endif          
2761       return
2762       end
2763 C-----------------------------------------------------------------------------
2764       subroutine vecpr(u,v,w)
2765       implicit real*8(a-h,o-z)
2766       dimension u(3),v(3),w(3)
2767       w(1)=u(2)*v(3)-u(3)*v(2)
2768       w(2)=-u(1)*v(3)+u(3)*v(1)
2769       w(3)=u(1)*v(2)-u(2)*v(1)
2770       return
2771       end
2772 C-----------------------------------------------------------------------------
2773       subroutine unormderiv(u,ugrad,unorm,ungrad)
2774 C This subroutine computes the derivatives of a normalized vector u, given
2775 C the derivatives computed without normalization conditions, ugrad. Returns
2776 C ungrad.
2777       implicit none
2778       double precision u(3),ugrad(3,3),unorm,ungrad(3,3)
2779       double precision vec(3)
2780       double precision scalar
2781       integer i,j
2782 c      write (2,*) 'ugrad',ugrad
2783 c      write (2,*) 'u',u
2784       do i=1,3
2785         vec(i)=scalar(ugrad(1,i),u(1))
2786       enddo
2787 c      write (2,*) 'vec',vec
2788       do i=1,3
2789         do j=1,3
2790           ungrad(j,i)=(ugrad(j,i)-u(j)*vec(i))*unorm
2791         enddo
2792       enddo
2793 c      write (2,*) 'ungrad',ungrad
2794       return
2795       end
2796 C-----------------------------------------------------------------------------
2797       subroutine escp(evdw2,evdw2_14)
2798 C
2799 C This subroutine calculates the excluded-volume interaction energy between
2800 C peptide-group centers and side chains and its gradient in virtual-bond and
2801 C side-chain vectors.
2802 C
2803       implicit real*8 (a-h,o-z)
2804       include 'DIMENSIONS'
2805       include 'DIMENSIONS.ZSCOPT'
2806       include 'COMMON.GEO'
2807       include 'COMMON.VAR'
2808       include 'COMMON.LOCAL'
2809       include 'COMMON.CHAIN'
2810       include 'COMMON.DERIV'
2811       include 'COMMON.INTERACT'
2812       include 'COMMON.FFIELD'
2813       include 'COMMON.IOUNITS'
2814       dimension ggg(3)
2815       evdw2=0.0D0
2816       evdw2_14=0.0d0
2817 cd    print '(a)','Enter ESCP'
2818 c      write (iout,*) 'iatscp_s=',iatscp_s,' iatscp_e=',iatscp_e,
2819 c     &  ' scal14',scal14
2820       do i=iatscp_s,iatscp_e
2821         if (itype(i).eq.ntyp1 .or. itype(i+1).eq.ntyp1) cycle
2822         iteli=itel(i)
2823 c        write (iout,*) "i",i," iteli",iteli," nscp_gr",nscp_gr(i),
2824 c     &   " iscp",(iscpstart(i,j),iscpend(i,j),j=1,nscp_gr(i))
2825         if (iteli.eq.0) goto 1225
2826         xi=0.5D0*(c(1,i)+c(1,i+1))
2827         yi=0.5D0*(c(2,i)+c(2,i+1))
2828         zi=0.5D0*(c(3,i)+c(3,i+1))
2829
2830         do iint=1,nscp_gr(i)
2831
2832         do j=iscpstart(i,iint),iscpend(i,iint)
2833           itypj=iabs(itype(j))
2834           if (itypj.eq.ntyp1) cycle
2835 C Uncomment following three lines for SC-p interactions
2836 c         xj=c(1,nres+j)-xi
2837 c         yj=c(2,nres+j)-yi
2838 c         zj=c(3,nres+j)-zi
2839 C Uncomment following three lines for Ca-p interactions
2840           xj=c(1,j)-xi
2841           yj=c(2,j)-yi
2842           zj=c(3,j)-zi
2843           rrij=1.0D0/(xj*xj+yj*yj+zj*zj)
2844           fac=rrij**expon2
2845           e1=fac*fac*aad(itypj,iteli)
2846           e2=fac*bad(itypj,iteli)
2847           if (iabs(j-i) .le. 2) then
2848             e1=scal14*e1
2849             e2=scal14*e2
2850             evdw2_14=evdw2_14+e1+e2
2851           endif
2852           evdwij=e1+e2
2853 c          write (iout,'(a6,2i5,0pf7.3,2i3,3e11.3)')
2854 c     &        'evdw2',i,j,evdwij,iteli,itypj,fac,aad(itypj,iteli),
2855 c     &       bad(itypj,iteli)
2856           evdw2=evdw2+evdwij
2857           if (calc_grad) then
2858 C
2859 C Calculate contributions to the gradient in the virtual-bond and SC vectors.
2860 C
2861           fac=-(evdwij+e1)*rrij
2862           ggg(1)=xj*fac
2863           ggg(2)=yj*fac
2864           ggg(3)=zj*fac
2865           if (j.lt.i) then
2866 cd          write (iout,*) 'j<i'
2867 C Uncomment following three lines for SC-p interactions
2868 c           do k=1,3
2869 c             gradx_scp(k,j)=gradx_scp(k,j)+ggg(k)
2870 c           enddo
2871           else
2872 cd          write (iout,*) 'j>i'
2873             do k=1,3
2874               ggg(k)=-ggg(k)
2875 C Uncomment following line for SC-p interactions
2876 c             gradx_scp(k,j)=gradx_scp(k,j)-ggg(k)
2877             enddo
2878           endif
2879           do k=1,3
2880             gvdwc_scp(k,i)=gvdwc_scp(k,i)-0.5D0*ggg(k)
2881           enddo
2882           kstart=min0(i+1,j)
2883           kend=max0(i-1,j-1)
2884 cd        write (iout,*) 'i=',i,' j=',j,' kstart=',kstart,' kend=',kend
2885 cd        write (iout,*) ggg(1),ggg(2),ggg(3)
2886           do k=kstart,kend
2887             do l=1,3
2888               gvdwc_scp(l,k)=gvdwc_scp(l,k)-ggg(l)
2889             enddo
2890           enddo
2891           endif
2892         enddo
2893         enddo ! iint
2894  1225   continue
2895       enddo ! i
2896       do i=1,nct
2897         do j=1,3
2898           gvdwc_scp(j,i)=expon*gvdwc_scp(j,i)
2899           gradx_scp(j,i)=expon*gradx_scp(j,i)
2900         enddo
2901       enddo
2902 C******************************************************************************
2903 C
2904 C                              N O T E !!!
2905 C
2906 C To save time the factor EXPON has been extracted from ALL components
2907 C of GVDWC and GRADX. Remember to multiply them by this factor before further 
2908 C use!
2909 C
2910 C******************************************************************************
2911       return
2912       end
2913 C--------------------------------------------------------------------------
2914       subroutine edis(ehpb)
2915
2916 C Evaluate bridge-strain energy and its gradient in virtual-bond and SC vectors.
2917 C
2918       implicit real*8 (a-h,o-z)
2919       include 'DIMENSIONS'
2920       include 'DIMENSIONS.ZSCOPT'
2921       include 'COMMON.SBRIDGE'
2922       include 'COMMON.CHAIN'
2923       include 'COMMON.DERIV'
2924       include 'COMMON.VAR'
2925       include 'COMMON.INTERACT'
2926       dimension ggg(3)
2927       ehpb=0.0D0
2928 cd    print *,'edis: nhpb=',nhpb,' fbr=',fbr
2929 cd    print *,'link_start=',link_start,' link_end=',link_end
2930       if (link_end.eq.0) return
2931       do i=link_start,link_end
2932 C If ihpb(i) and jhpb(i) > NRES, this is a SC-SC distance, otherwise a
2933 C CA-CA distance used in regularization of structure.
2934         ii=ihpb(i)
2935         jj=jhpb(i)
2936 C iii and jjj point to the residues for which the distance is assigned.
2937         if (ii.gt.nres) then
2938           iii=ii-nres
2939           jjj=jj-nres 
2940         else
2941           iii=ii
2942           jjj=jj
2943         endif
2944 C 24/11/03 AL: SS bridges handled separately because of introducing a specific
2945 C    distance and angle dependent SS bond potential.
2946         if (ii.gt.nres .and. iabs(itype(iii)).eq.1 .and. 
2947      & iabs(itype(jjj)).eq.1) then
2948           call ssbond_ene(iii,jjj,eij)
2949           ehpb=ehpb+2*eij
2950         else
2951 C Calculate the distance between the two points and its difference from the
2952 C target distance.
2953         dd=dist(ii,jj)
2954         rdis=dd-dhpb(i)
2955 C Get the force constant corresponding to this distance.
2956         waga=forcon(i)
2957 C Calculate the contribution to energy.
2958         ehpb=ehpb+waga*rdis*rdis
2959 C
2960 C Evaluate gradient.
2961 C
2962         fac=waga*rdis/dd
2963 cd      print *,'i=',i,' ii=',ii,' jj=',jj,' dhpb=',dhpb(i),' dd=',dd,
2964 cd   &   ' waga=',waga,' fac=',fac
2965         do j=1,3
2966           ggg(j)=fac*(c(j,jj)-c(j,ii))
2967         enddo
2968 cd      print '(i3,3(1pe14.5))',i,(ggg(j),j=1,3)
2969 C If this is a SC-SC distance, we need to calculate the contributions to the
2970 C Cartesian gradient in the SC vectors (ghpbx).
2971         if (iii.lt.ii) then
2972           do j=1,3
2973             ghpbx(j,iii)=ghpbx(j,iii)-ggg(j)
2974             ghpbx(j,jjj)=ghpbx(j,jjj)+ggg(j)
2975           enddo
2976         endif
2977         do j=iii,jjj-1
2978           do k=1,3
2979             ghpbc(k,j)=ghpbc(k,j)+ggg(k)
2980           enddo
2981         enddo
2982         endif
2983       enddo
2984       ehpb=0.5D0*ehpb
2985       return
2986       end
2987 C--------------------------------------------------------------------------
2988       subroutine ssbond_ene(i,j,eij)
2989
2990 C Calculate the distance and angle dependent SS-bond potential energy
2991 C using a free-energy function derived based on RHF/6-31G** ab initio
2992 C calculations of diethyl disulfide.
2993 C
2994 C A. Liwo and U. Kozlowska, 11/24/03
2995 C
2996       implicit real*8 (a-h,o-z)
2997       include 'DIMENSIONS'
2998       include 'DIMENSIONS.ZSCOPT'
2999       include 'COMMON.SBRIDGE'
3000       include 'COMMON.CHAIN'
3001       include 'COMMON.DERIV'
3002       include 'COMMON.LOCAL'
3003       include 'COMMON.INTERACT'
3004       include 'COMMON.VAR'
3005       include 'COMMON.IOUNITS'
3006       double precision erij(3),dcosom1(3),dcosom2(3),gg(3)
3007       itypi=iabs(itype(i))
3008       xi=c(1,nres+i)
3009       yi=c(2,nres+i)
3010       zi=c(3,nres+i)
3011       dxi=dc_norm(1,nres+i)
3012       dyi=dc_norm(2,nres+i)
3013       dzi=dc_norm(3,nres+i)
3014       dsci_inv=dsc_inv(itypi)
3015       itypj=iabs(itype(j))
3016       dscj_inv=dsc_inv(itypj)
3017       xj=c(1,nres+j)-xi
3018       yj=c(2,nres+j)-yi
3019       zj=c(3,nres+j)-zi
3020       dxj=dc_norm(1,nres+j)
3021       dyj=dc_norm(2,nres+j)
3022       dzj=dc_norm(3,nres+j)
3023       rrij=1.0D0/(xj*xj+yj*yj+zj*zj)
3024       rij=dsqrt(rrij)
3025       erij(1)=xj*rij
3026       erij(2)=yj*rij
3027       erij(3)=zj*rij
3028       om1=dxi*erij(1)+dyi*erij(2)+dzi*erij(3)
3029       om2=dxj*erij(1)+dyj*erij(2)+dzj*erij(3)
3030       om12=dxi*dxj+dyi*dyj+dzi*dzj
3031       do k=1,3
3032         dcosom1(k)=rij*(dc_norm(k,nres+i)-om1*erij(k))
3033         dcosom2(k)=rij*(dc_norm(k,nres+j)-om2*erij(k))
3034       enddo
3035       rij=1.0d0/rij
3036       deltad=rij-d0cm
3037       deltat1=1.0d0-om1
3038       deltat2=1.0d0+om2
3039       deltat12=om2-om1+2.0d0
3040       cosphi=om12-om1*om2
3041       eij=akcm*deltad*deltad+akth*(deltat1*deltat1+deltat2*deltat2)
3042      &  +akct*deltad*deltat12
3043      &  +v1ss*cosphi+v2ss*cosphi*cosphi+v3ss*cosphi*cosphi*cosphi
3044 c      write(iout,*) i,j,"rij",rij,"d0cm",d0cm," akcm",akcm," akth",akth,
3045 c     &  " akct",akct," deltad",deltad," deltat",deltat1,deltat2,
3046 c     &  " deltat12",deltat12," eij",eij 
3047       ed=2*akcm*deltad+akct*deltat12
3048       pom1=akct*deltad
3049       pom2=v1ss+2*v2ss*cosphi+3*v3ss*cosphi*cosphi
3050       eom1=-2*akth*deltat1-pom1-om2*pom2
3051       eom2= 2*akth*deltat2+pom1-om1*pom2
3052       eom12=pom2
3053       do k=1,3
3054         gg(k)=ed*erij(k)+eom1*dcosom1(k)+eom2*dcosom2(k)
3055       enddo
3056       do k=1,3
3057         ghpbx(k,i)=ghpbx(k,i)-gg(k)
3058      &            +(eom12*dc_norm(k,nres+j)+eom1*erij(k))*dsci_inv
3059         ghpbx(k,j)=ghpbx(k,j)+gg(k)
3060      &            +(eom12*dc_norm(k,nres+i)+eom2*erij(k))*dscj_inv
3061       enddo
3062 C
3063 C Calculate the components of the gradient in DC and X
3064 C
3065       do k=i,j-1
3066         do l=1,3
3067           ghpbc(l,k)=ghpbc(l,k)+gg(l)
3068         enddo
3069       enddo
3070       return
3071       end
3072 C--------------------------------------------------------------------------
3073       subroutine ebond(estr)
3074 c
3075 c Evaluate the energy of stretching of the CA-CA and CA-SC virtual bonds
3076 c
3077       implicit real*8 (a-h,o-z)
3078       include 'DIMENSIONS'
3079       include 'DIMENSIONS.ZSCOPT'
3080       include 'COMMON.LOCAL'
3081       include 'COMMON.GEO'
3082       include 'COMMON.INTERACT'
3083       include 'COMMON.DERIV'
3084       include 'COMMON.VAR'
3085       include 'COMMON.CHAIN'
3086       include 'COMMON.IOUNITS'
3087       include 'COMMON.NAMES'
3088       include 'COMMON.FFIELD'
3089       include 'COMMON.CONTROL'
3090       logical energy_dec /.false./
3091       double precision u(3),ud(3)
3092       estr=0.0d0
3093       estr1=0.0d0
3094 c      write (iout,*) "distchainmax",distchainmax
3095       do i=nnt+1,nct
3096         if (itype(i-1).eq.ntyp1 .or. itype(i).eq.ntyp1) then
3097           estr1=estr1+gnmr1(vbld(i),-1.0d0,distchainmax)
3098           do j=1,3
3099           gradb(j,i-1)=gnmr1prim(vbld(i),-1.0d0,distchainmax)
3100      &      *dc(j,i-1)/vbld(i)
3101           enddo
3102           if (energy_dec) write(iout,*)
3103      &       "estr1",i,vbld(i),distchainmax,
3104      &       gnmr1(vbld(i),-1.0d0,distchainmax)
3105         else
3106           diff = vbld(i)-vbldp0
3107 c          write (iout,*) i,vbld(i),vbldp0,diff,AKP*diff*diff
3108           estr=estr+diff*diff
3109           do j=1,3
3110             gradb(j,i-1)=AKP*diff*dc(j,i-1)/vbld(i)
3111           enddo
3112         endif
3113
3114       enddo
3115       estr=0.5d0*AKP*estr+estr1
3116 c
3117 c 09/18/07 AL: multimodal bond potential based on AM1 CA-SC PMF's included
3118 c
3119       do i=nnt,nct
3120         iti=iabs(itype(i))
3121         if (iti.ne.10 .and. iti.ne.ntyp1) then
3122           nbi=nbondterm(iti)
3123           if (nbi.eq.1) then
3124             diff=vbld(i+nres)-vbldsc0(1,iti)
3125 c            write (iout,*) i,iti,vbld(i+nres),vbldsc0(1,iti),diff,
3126 c     &      AKSC(1,iti),AKSC(1,iti)*diff*diff
3127             estr=estr+0.5d0*AKSC(1,iti)*diff*diff
3128             do j=1,3
3129               gradbx(j,i)=AKSC(1,iti)*diff*dc(j,i+nres)/vbld(i+nres)
3130             enddo
3131           else
3132             do j=1,nbi
3133               diff=vbld(i+nres)-vbldsc0(j,iti)
3134               ud(j)=aksc(j,iti)*diff
3135               u(j)=abond0(j,iti)+0.5d0*ud(j)*diff
3136             enddo
3137             uprod=u(1)
3138             do j=2,nbi
3139               uprod=uprod*u(j)
3140             enddo
3141             usum=0.0d0
3142             usumsqder=0.0d0
3143             do j=1,nbi
3144               uprod1=1.0d0
3145               uprod2=1.0d0
3146               do k=1,nbi
3147                 if (k.ne.j) then
3148                   uprod1=uprod1*u(k)
3149                   uprod2=uprod2*u(k)*u(k)
3150                 endif
3151               enddo
3152               usum=usum+uprod1
3153               usumsqder=usumsqder+ud(j)*uprod2
3154             enddo
3155 c            write (iout,*) i,iti,vbld(i+nres),(vbldsc0(j,iti),
3156 c     &      AKSC(j,iti),abond0(j,iti),u(j),j=1,nbi)
3157             estr=estr+uprod/usum
3158             do j=1,3
3159              gradbx(j,i)=usumsqder/(usum*usum)*dc(j,i+nres)/vbld(i+nres)
3160             enddo
3161           endif
3162         endif
3163       enddo
3164       return
3165       end
3166 #ifdef CRYST_THETA
3167 C--------------------------------------------------------------------------
3168       subroutine ebend(etheta)
3169 C
3170 C Evaluate the virtual-bond-angle energy given the virtual-bond dihedral
3171 C angles gamma and its derivatives in consecutive thetas and gammas.
3172 C
3173       implicit real*8 (a-h,o-z)
3174       include 'DIMENSIONS'
3175       include 'DIMENSIONS.ZSCOPT'
3176       include 'COMMON.LOCAL'
3177       include 'COMMON.GEO'
3178       include 'COMMON.INTERACT'
3179       include 'COMMON.DERIV'
3180       include 'COMMON.VAR'
3181       include 'COMMON.CHAIN'
3182       include 'COMMON.IOUNITS'
3183       include 'COMMON.NAMES'
3184       include 'COMMON.FFIELD'
3185       common /calcthet/ term1,term2,termm,diffak,ratak,
3186      & ak,aktc,termpre,termexp,sigc,sig0i,time11,time12,sigcsq,
3187      & delthe0,sig0inv,sigtc,sigsqtc,delthec,it
3188       double precision y(2),z(2)
3189       delta=0.02d0*pi
3190 c      time11=dexp(-2*time)
3191 c      time12=1.0d0
3192       etheta=0.0D0
3193 c      write (iout,*) "nres",nres
3194 c     write (*,'(a,i2)') 'EBEND ICG=',icg
3195 c      write (iout,*) ithet_start,ithet_end
3196       do i=ithet_start,ithet_end
3197         if (itype(i-1).eq.ntyp1) cycle
3198 C Zero the energy function and its derivative at 0 or pi.
3199         call splinthet(theta(i),0.5d0*delta,ss,ssd)
3200         it=itype(i-1)
3201         ichir1=isign(1,itype(i-2))
3202         ichir2=isign(1,itype(i))
3203          if (itype(i-2).eq.10) ichir1=isign(1,itype(i-1))
3204          if (itype(i).eq.10) ichir2=isign(1,itype(i-1))
3205          if (itype(i-1).eq.10) then
3206           itype1=isign(10,itype(i-2))
3207           ichir11=isign(1,itype(i-2))
3208           ichir12=isign(1,itype(i-2))
3209           itype2=isign(10,itype(i))
3210           ichir21=isign(1,itype(i))
3211           ichir22=isign(1,itype(i))
3212          endif
3213
3214         if (i.gt.3 .and. itype(i-2).ne.ntyp1) then
3215 #ifdef OSF
3216           phii=phi(i)
3217 c          icrc=0
3218 c          call proc_proc(phii,icrc)
3219           if (icrc.eq.1) phii=150.0
3220 #else
3221           phii=phi(i)
3222 #endif
3223           y(1)=dcos(phii)
3224           y(2)=dsin(phii)
3225         else
3226           y(1)=0.0D0
3227           y(2)=0.0D0
3228         endif
3229         if (i.lt.nres .and. itype(i).ne.ntyp1) then
3230 #ifdef OSF
3231           phii1=phi(i+1)
3232 c          icrc=0
3233 c          call proc_proc(phii1,icrc)
3234           if (icrc.eq.1) phii1=150.0
3235           phii1=pinorm(phii1)
3236           z(1)=cos(phii1)
3237 #else
3238           phii1=phi(i+1)
3239           z(1)=dcos(phii1)
3240 #endif
3241           z(2)=dsin(phii1)
3242         else
3243           z(1)=0.0D0
3244           z(2)=0.0D0
3245         endif
3246 C Calculate the "mean" value of theta from the part of the distribution
3247 C dependent on the adjacent virtual-bond-valence angles (gamma1 & gamma2).
3248 C In following comments this theta will be referred to as t_c.
3249         thet_pred_mean=0.0d0
3250         do k=1,2
3251             athetk=athet(k,it,ichir1,ichir2)
3252             bthetk=bthet(k,it,ichir1,ichir2)
3253           if (it.eq.10) then
3254              athetk=athet(k,itype1,ichir11,ichir12)
3255              bthetk=bthet(k,itype2,ichir21,ichir22)
3256           endif
3257           thet_pred_mean=thet_pred_mean+athetk*y(k)+bthetk*z(k)
3258         enddo
3259 c        write (iout,*) "thet_pred_mean",thet_pred_mean
3260         dthett=thet_pred_mean*ssd
3261         thet_pred_mean=thet_pred_mean*ss+a0thet(it)
3262 c        write (iout,*) "thet_pred_mean",thet_pred_mean
3263 C Derivatives of the "mean" values in gamma1 and gamma2.
3264         dthetg1=(-athet(1,it,ichir1,ichir2)*y(2)
3265      &+athet(2,it,ichir1,ichir2)*y(1))*ss
3266          dthetg2=(-bthet(1,it,ichir1,ichir2)*z(2)
3267      &          +bthet(2,it,ichir1,ichir2)*z(1))*ss
3268          if (it.eq.10) then
3269       dthetg1=(-athet(1,itype1,ichir11,ichir12)*y(2)
3270      &+athet(2,itype1,ichir11,ichir12)*y(1))*ss
3271         dthetg2=(-bthet(1,itype2,ichir21,ichir22)*z(2)
3272      &         +bthet(2,itype2,ichir21,ichir22)*z(1))*ss
3273          endif
3274         if (theta(i).gt.pi-delta) then
3275           call theteng(pi-delta,thet_pred_mean,theta0(it),f0,fprim0,
3276      &         E_tc0)
3277           call mixder(pi-delta,thet_pred_mean,theta0(it),fprim_tc0)
3278           call theteng(pi,thet_pred_mean,theta0(it),f1,fprim1,E_tc1)
3279           call spline1(theta(i),pi-delta,delta,f0,f1,fprim0,ethetai,
3280      &        E_theta)
3281           call spline2(theta(i),pi-delta,delta,E_tc0,E_tc1,fprim_tc0,
3282      &        E_tc)
3283         else if (theta(i).lt.delta) then
3284           call theteng(delta,thet_pred_mean,theta0(it),f0,fprim0,E_tc0)
3285           call theteng(0.0d0,thet_pred_mean,theta0(it),f1,fprim1,E_tc1)
3286           call spline1(theta(i),delta,-delta,f0,f1,fprim0,ethetai,
3287      &        E_theta)
3288           call mixder(delta,thet_pred_mean,theta0(it),fprim_tc0)
3289           call spline2(theta(i),delta,-delta,E_tc0,E_tc1,fprim_tc0,
3290      &        E_tc)
3291         else
3292           call theteng(theta(i),thet_pred_mean,theta0(it),ethetai,
3293      &        E_theta,E_tc)
3294         endif
3295         etheta=etheta+ethetai
3296 c        write (iout,'(2i3,3f8.3,f10.5)') i,it,rad2deg*theta(i),
3297 c     &    rad2deg*phii,rad2deg*phii1,ethetai
3298         if (i.gt.3) gloc(i-3,icg)=gloc(i-3,icg)+wang*E_tc*dthetg1
3299         if (i.lt.nres) gloc(i-2,icg)=gloc(i-2,icg)+wang*E_tc*dthetg2
3300         gloc(nphi+i-2,icg)=wang*(E_theta+E_tc*dthett)
3301 c 1215   continue
3302       enddo
3303 C Ufff.... We've done all this!!! 
3304       return
3305       end
3306 C---------------------------------------------------------------------------
3307       subroutine theteng(thetai,thet_pred_mean,theta0i,ethetai,E_theta,
3308      &     E_tc)
3309       implicit real*8 (a-h,o-z)
3310       include 'DIMENSIONS'
3311       include 'COMMON.LOCAL'
3312       include 'COMMON.IOUNITS'
3313       common /calcthet/ term1,term2,termm,diffak,ratak,
3314      & ak,aktc,termpre,termexp,sigc,sig0i,time11,time12,sigcsq,
3315      & delthe0,sig0inv,sigtc,sigsqtc,delthec,it
3316 C Calculate the contributions to both Gaussian lobes.
3317 C 6/6/97 - Deform the Gaussians using the factor of 1/(1+time)
3318 C The "polynomial part" of the "standard deviation" of this part of 
3319 C the distribution.
3320         sig=polthet(3,it)
3321         do j=2,0,-1
3322           sig=sig*thet_pred_mean+polthet(j,it)
3323         enddo
3324 C Derivative of the "interior part" of the "standard deviation of the" 
3325 C gamma-dependent Gaussian lobe in t_c.
3326         sigtc=3*polthet(3,it)
3327         do j=2,1,-1
3328           sigtc=sigtc*thet_pred_mean+j*polthet(j,it)
3329         enddo
3330         sigtc=sig*sigtc
3331 C Set the parameters of both Gaussian lobes of the distribution.
3332 C "Standard deviation" of the gamma-dependent Gaussian lobe (sigtc)
3333         fac=sig*sig+sigc0(it)
3334         sigcsq=fac+fac
3335         sigc=1.0D0/sigcsq
3336 C Following variable (sigsqtc) is -(1/2)d[sigma(t_c)**(-2))]/dt_c
3337         sigsqtc=-4.0D0*sigcsq*sigtc
3338 c       print *,i,sig,sigtc,sigsqtc
3339 C Following variable (sigtc) is d[sigma(t_c)]/dt_c
3340         sigtc=-sigtc/(fac*fac)
3341 C Following variable is sigma(t_c)**(-2)
3342         sigcsq=sigcsq*sigcsq
3343         sig0i=sig0(it)
3344         sig0inv=1.0D0/sig0i**2
3345         delthec=thetai-thet_pred_mean
3346         delthe0=thetai-theta0i
3347         term1=-0.5D0*sigcsq*delthec*delthec
3348         term2=-0.5D0*sig0inv*delthe0*delthe0
3349 C Following fuzzy logic is to avoid underflows in dexp and subsequent INFs and
3350 C NaNs in taking the logarithm. We extract the largest exponent which is added
3351 C to the energy (this being the log of the distribution) at the end of energy
3352 C term evaluation for this virtual-bond angle.
3353         if (term1.gt.term2) then
3354           termm=term1
3355           term2=dexp(term2-termm)
3356           term1=1.0d0
3357         else
3358           termm=term2
3359           term1=dexp(term1-termm)
3360           term2=1.0d0
3361         endif
3362 C The ratio between the gamma-independent and gamma-dependent lobes of
3363 C the distribution is a Gaussian function of thet_pred_mean too.
3364         diffak=gthet(2,it)-thet_pred_mean
3365         ratak=diffak/gthet(3,it)**2
3366         ak=dexp(gthet(1,it)-0.5D0*diffak*ratak)
3367 C Let's differentiate it in thet_pred_mean NOW.
3368         aktc=ak*ratak
3369 C Now put together the distribution terms to make complete distribution.
3370         termexp=term1+ak*term2
3371         termpre=sigc+ak*sig0i
3372 C Contribution of the bending energy from this theta is just the -log of
3373 C the sum of the contributions from the two lobes and the pre-exponential
3374 C factor. Simple enough, isn't it?
3375         ethetai=(-dlog(termexp)-termm+dlog(termpre))
3376 C NOW the derivatives!!!
3377 C 6/6/97 Take into account the deformation.
3378         E_theta=(delthec*sigcsq*term1
3379      &       +ak*delthe0*sig0inv*term2)/termexp
3380         E_tc=((sigtc+aktc*sig0i)/termpre
3381      &      -((delthec*sigcsq+delthec*delthec*sigsqtc)*term1+
3382      &       aktc*term2)/termexp)
3383       return
3384       end
3385 c-----------------------------------------------------------------------------
3386       subroutine mixder(thetai,thet_pred_mean,theta0i,E_tc_t)
3387       implicit real*8 (a-h,o-z)
3388       include 'DIMENSIONS'
3389       include 'COMMON.LOCAL'
3390       include 'COMMON.IOUNITS'
3391       common /calcthet/ term1,term2,termm,diffak,ratak,
3392      & ak,aktc,termpre,termexp,sigc,sig0i,time11,time12,sigcsq,
3393      & delthe0,sig0inv,sigtc,sigsqtc,delthec,it
3394       delthec=thetai-thet_pred_mean
3395       delthe0=thetai-theta0i
3396 C "Thank you" to MAPLE (probably spared one day of hand-differentiation).
3397       t3 = thetai-thet_pred_mean
3398       t6 = t3**2
3399       t9 = term1
3400       t12 = t3*sigcsq
3401       t14 = t12+t6*sigsqtc
3402       t16 = 1.0d0
3403       t21 = thetai-theta0i
3404       t23 = t21**2
3405       t26 = term2
3406       t27 = t21*t26
3407       t32 = termexp
3408       t40 = t32**2
3409       E_tc_t = -((sigcsq+2.D0*t3*sigsqtc)*t9-t14*sigcsq*t3*t16*t9
3410      & -aktc*sig0inv*t27)/t32+(t14*t9+aktc*t26)/t40
3411      & *(-t12*t9-ak*sig0inv*t27)
3412       return
3413       end
3414 #else
3415 C--------------------------------------------------------------------------
3416       subroutine ebend(etheta)
3417 C
3418 C Evaluate the virtual-bond-angle energy given the virtual-bond dihedral
3419 C angles gamma and its derivatives in consecutive thetas and gammas.
3420 C ab initio-derived potentials from 
3421 c Kozlowska et al., J. Phys.: Condens. Matter 19 (2007) 285203
3422 C
3423       implicit real*8 (a-h,o-z)
3424       include 'DIMENSIONS'
3425       include 'DIMENSIONS.ZSCOPT'
3426       include 'COMMON.LOCAL'
3427       include 'COMMON.GEO'
3428       include 'COMMON.INTERACT'
3429       include 'COMMON.DERIV'
3430       include 'COMMON.VAR'
3431       include 'COMMON.CHAIN'
3432       include 'COMMON.IOUNITS'
3433       include 'COMMON.NAMES'
3434       include 'COMMON.FFIELD'
3435       include 'COMMON.CONTROL'
3436       double precision coskt(mmaxtheterm),sinkt(mmaxtheterm),
3437      & cosph1(maxsingle),sinph1(maxsingle),cosph2(maxsingle),
3438      & sinph2(maxsingle),cosph1ph2(maxdouble,maxdouble),
3439      & sinph1ph2(maxdouble,maxdouble)
3440       logical lprn /.false./, lprn1 /.false./
3441       etheta=0.0D0
3442 c      write (iout,*) "ithetyp",(ithetyp(i),i=1,ntyp1)
3443       do i=ithet_start,ithet_end
3444 c        if (itype(i-1).eq.ntyp1) cycle
3445         if ((itype(i-1).eq.ntyp1).or.(itype(i-2).eq.ntyp1).or.
3446      &(itype(i).eq.ntyp1)) cycle
3447         if (iabs(itype(i+1)).eq.20) iblock=2
3448         if (iabs(itype(i+1)).ne.20) iblock=1
3449         dethetai=0.0d0
3450         dephii=0.0d0
3451         dephii1=0.0d0
3452         theti2=0.5d0*theta(i)
3453         ityp2=ithetyp((itype(i-1)))
3454         do k=1,nntheterm
3455           coskt(k)=dcos(k*theti2)
3456           sinkt(k)=dsin(k*theti2)
3457         enddo
3458         if (i.gt.3 .and. itype(i-3).ne.ntyp1) then
3459 #ifdef OSF
3460           phii=phi(i)
3461           if (phii.ne.phii) phii=150.0
3462 #else
3463           phii=phi(i)
3464 #endif
3465           ityp1=ithetyp((itype(i-2)))
3466           do k=1,nsingle
3467             cosph1(k)=dcos(k*phii)
3468             sinph1(k)=dsin(k*phii)
3469           enddo
3470         else
3471           phii=0.0d0
3472 c          ityp1=nthetyp+1
3473           do k=1,nsingle
3474             ityp1=ithetyp((itype(i-2)))
3475             cosph1(k)=0.0d0
3476             sinph1(k)=0.0d0
3477           enddo 
3478         endif
3479         if (i.lt.nres .and. itype(i+1).ne.ntyp1) then
3480 #ifdef OSF
3481           phii1=phi(i+1)
3482           if (phii1.ne.phii1) phii1=150.0
3483           phii1=pinorm(phii1)
3484 #else
3485           phii1=phi(i+1)
3486 #endif
3487           ityp3=ithetyp((itype(i)))
3488           do k=1,nsingle
3489             cosph2(k)=dcos(k*phii1)
3490             sinph2(k)=dsin(k*phii1)
3491           enddo
3492         else
3493           phii1=0.0d0
3494 c          ityp3=nthetyp+1
3495           ityp3=ithetyp((itype(i)))
3496           do k=1,nsingle
3497             cosph2(k)=0.0d0
3498             sinph2(k)=0.0d0
3499           enddo
3500         endif  
3501 c        write (iout,*) "i",i," ityp1",itype(i-2),ityp1,
3502 c     &   " ityp2",itype(i-1),ityp2," ityp3",itype(i),ityp3
3503 c        call flush(iout)
3504         ethetai=aa0thet(ityp1,ityp2,ityp3,iblock)
3505         do k=1,ndouble
3506           do l=1,k-1
3507             ccl=cosph1(l)*cosph2(k-l)
3508             ssl=sinph1(l)*sinph2(k-l)
3509             scl=sinph1(l)*cosph2(k-l)
3510             csl=cosph1(l)*sinph2(k-l)
3511             cosph1ph2(l,k)=ccl-ssl
3512             cosph1ph2(k,l)=ccl+ssl
3513             sinph1ph2(l,k)=scl+csl
3514             sinph1ph2(k,l)=scl-csl
3515           enddo
3516         enddo
3517         if (lprn) then
3518         write (iout,*) "i",i," ityp1",ityp1," ityp2",ityp2,
3519      &    " ityp3",ityp3," theti2",theti2," phii",phii," phii1",phii1
3520         write (iout,*) "coskt and sinkt"
3521         do k=1,nntheterm
3522           write (iout,*) k,coskt(k),sinkt(k)
3523         enddo
3524         endif
3525         do k=1,ntheterm
3526           ethetai=ethetai+aathet(k,ityp1,ityp2,ityp3,iblock)*sinkt(k)
3527           dethetai=dethetai+0.5d0*k*aathet(k,ityp1,ityp2,ityp3,iblock)
3528      &      *coskt(k)
3529           if (lprn)
3530      &    write (iout,*) "k",k,"
3531      &      aathet",aathet(k,ityp1,ityp2,ityp3,iblock),
3532      &     " ethetai",ethetai
3533         enddo
3534         if (lprn) then
3535         write (iout,*) "cosph and sinph"
3536         do k=1,nsingle
3537           write (iout,*) k,cosph1(k),sinph1(k),cosph2(k),sinph2(k)
3538         enddo
3539         write (iout,*) "cosph1ph2 and sinph2ph2"
3540         do k=2,ndouble
3541           do l=1,k-1
3542             write (iout,*) l,k,cosph1ph2(l,k),cosph1ph2(k,l),
3543      &         sinph1ph2(l,k),sinph1ph2(k,l) 
3544           enddo
3545         enddo
3546         write(iout,*) "ethetai",ethetai
3547         endif
3548         do m=1,ntheterm2
3549           do k=1,nsingle
3550             aux=bbthet(k,m,ityp1,ityp2,ityp3,iblock)*cosph1(k)
3551      &         +ccthet(k,m,ityp1,ityp2,ityp3,iblock)*sinph1(k)
3552      &         +ddthet(k,m,ityp1,ityp2,ityp3,iblock)*cosph2(k)
3553      &         +eethet(k,m,ityp1,ityp2,ityp3,iblock)*sinph2(k)
3554             ethetai=ethetai+sinkt(m)*aux
3555             dethetai=dethetai+0.5d0*m*aux*coskt(m)
3556             dephii=dephii+k*sinkt(m)*(
3557      &          ccthet(k,m,ityp1,ityp2,ityp3,iblock)*cosph1(k)-
3558      &          bbthet(k,m,ityp1,ityp2,ityp3,iblock)*sinph1(k))
3559             dephii1=dephii1+k*sinkt(m)*(
3560      &          eethet(k,m,ityp1,ityp2,ityp3,iblock)*cosph2(k)-
3561      &          ddthet(k,m,ityp1,ityp2,ityp3,iblock)*sinph2(k))
3562             if (lprn)
3563      &      write (iout,*) "m",m," k",k," bbthet",
3564      &         bbthet(k,m,ityp1,ityp2,ityp3,iblock)," ccthet",
3565      &         ccthet(k,m,ityp1,ityp2,ityp3,iblock)," ddthet",
3566      &         ddthet(k,m,ityp1,ityp2,ityp3,iblock)," eethet",
3567      &         eethet(k,m,ityp1,ityp2,ityp3,iblock)," ethetai",ethetai
3568           enddo
3569         enddo
3570         if (lprn)
3571      &  write(iout,*) "ethetai",ethetai
3572         do m=1,ntheterm3
3573           do k=2,ndouble
3574             do l=1,k-1
3575               aux=ffthet(l,k,m,ityp1,ityp2,ityp3,iblock)*cosph1ph2(l,k)+
3576      &            ffthet(k,l,m,ityp1,ityp2,ityp3,iblock)*cosph1ph2(k,l)+
3577      &            ggthet(l,k,m,ityp1,ityp2,ityp3,iblock)*sinph1ph2(l,k)+
3578      &            ggthet(k,l,m,ityp1,ityp2,ityp3,iblock)*sinph1ph2(k,l)
3579               ethetai=ethetai+sinkt(m)*aux
3580               dethetai=dethetai+0.5d0*m*coskt(m)*aux
3581               dephii=dephii+l*sinkt(m)*(
3582      &           -ffthet(l,k,m,ityp1,ityp2,ityp3,iblock)*sinph1ph2(l,k)-
3583      &            ffthet(k,l,m,ityp1,ityp2,ityp3,iblock)*sinph1ph2(k,l)+
3584      &            ggthet(l,k,m,ityp1,ityp2,ityp3,iblock)*cosph1ph2(l,k)+
3585      &            ggthet(k,l,m,ityp1,ityp2,ityp3,iblock)*cosph1ph2(k,l))
3586               dephii1=dephii1+(k-l)*sinkt(m)*(
3587      &           -ffthet(l,k,m,ityp1,ityp2,ityp3,iblock)*sinph1ph2(l,k)+
3588      &            ffthet(k,l,m,ityp1,ityp2,ityp3,iblock)*sinph1ph2(k,l)+
3589      &            ggthet(l,k,m,ityp1,ityp2,ityp3,iblock)*cosph1ph2(l,k)-
3590      &            ggthet(k,l,m,ityp1,ityp2,ityp3,iblock)*cosph1ph2(k,l))
3591               if (lprn) then
3592               write (iout,*) "m",m," k",k," l",l," ffthet",
3593      &            ffthet(l,k,m,ityp1,ityp2,ityp3,iblock),
3594      &            ffthet(k,l,m,ityp1,ityp2,ityp3,iblock)," ggthet",
3595      &            ggthet(l,k,m,ityp1,ityp2,ityp3,iblock),
3596      &            ggthet(k,l,m,ityp1,ityp2,ityp3,iblock),
3597      &            " ethetai",ethetai
3598               write (iout,*) cosph1ph2(l,k)*sinkt(m),
3599      &            cosph1ph2(k,l)*sinkt(m),
3600      &            sinph1ph2(l,k)*sinkt(m),sinph1ph2(k,l)*sinkt(m)
3601               endif
3602             enddo
3603           enddo
3604         enddo
3605 10      continue
3606         if (lprn1) write (iout,'(i2,3f8.1,9h ethetai ,f10.5)') 
3607      &   i,theta(i)*rad2deg,phii*rad2deg,
3608      &   phii1*rad2deg,ethetai
3609         etheta=etheta+ethetai
3610         if (i.gt.3) gloc(i-3,icg)=gloc(i-3,icg)+wang*dephii
3611         if (i.lt.nres) gloc(i-2,icg)=gloc(i-2,icg)+wang*dephii1
3612 c        gloc(nphi+i-2,icg)=wang*dethetai
3613         gloc(nphi+i-2,icg)=gloc(nphi+i-2,icg)+wang*dethetai
3614       enddo
3615       return
3616       end
3617 #endif
3618 #ifdef CRYST_SC
3619 c-----------------------------------------------------------------------------
3620       subroutine esc(escloc)
3621 C Calculate the local energy of a side chain and its derivatives in the
3622 C corresponding virtual-bond valence angles THETA and the spherical angles 
3623 C ALPHA and OMEGA.
3624       implicit real*8 (a-h,o-z)
3625       include 'DIMENSIONS'
3626       include 'DIMENSIONS.ZSCOPT'
3627       include 'COMMON.GEO'
3628       include 'COMMON.LOCAL'
3629       include 'COMMON.VAR'
3630       include 'COMMON.INTERACT'
3631       include 'COMMON.DERIV'
3632       include 'COMMON.CHAIN'
3633       include 'COMMON.IOUNITS'
3634       include 'COMMON.NAMES'
3635       include 'COMMON.FFIELD'
3636       double precision x(3),dersc(3),xemp(3),dersc0(3),dersc1(3),
3637      &     ddersc0(3),ddummy(3),xtemp(3),temp(3)
3638       common /sccalc/ time11,time12,time112,theti,it,nlobit
3639       delta=0.02d0*pi
3640       escloc=0.0D0
3641 c     write (iout,'(a)') 'ESC'
3642       do i=loc_start,loc_end
3643         it=itype(i)
3644         if (it.eq.ntyp1) cycle
3645         if (it.eq.10) goto 1
3646         nlobit=nlob(iabs(it))
3647 c       print *,'i=',i,' it=',it,' nlobit=',nlobit
3648 c       write (iout,*) 'i=',i,' ssa=',ssa,' ssad=',ssad
3649         theti=theta(i+1)-pipol
3650         x(1)=dtan(theti)
3651         x(2)=alph(i)
3652         x(3)=omeg(i)
3653 c        write (iout,*) "i",i," x",x(1),x(2),x(3)
3654
3655         if (x(2).gt.pi-delta) then
3656           xtemp(1)=x(1)
3657           xtemp(2)=pi-delta
3658           xtemp(3)=x(3)
3659           call enesc(xtemp,escloci0,dersc0,ddersc0,.true.)
3660           xtemp(2)=pi
3661           call enesc(xtemp,escloci1,dersc1,ddummy,.false.)
3662           call spline1(x(2),pi-delta,delta,escloci0,escloci1,dersc0(2),
3663      &        escloci,dersc(2))
3664           call spline2(x(2),pi-delta,delta,dersc0(1),dersc1(1),
3665      &        ddersc0(1),dersc(1))
3666           call spline2(x(2),pi-delta,delta,dersc0(3),dersc1(3),
3667      &        ddersc0(3),dersc(3))
3668           xtemp(2)=pi-delta
3669           call enesc_bound(xtemp,esclocbi0,dersc0,dersc12,.true.)
3670           xtemp(2)=pi
3671           call enesc_bound(xtemp,esclocbi1,dersc1,chuju,.false.)
3672           call spline1(x(2),pi-delta,delta,esclocbi0,esclocbi1,
3673      &            dersc0(2),esclocbi,dersc02)
3674           call spline2(x(2),pi-delta,delta,dersc0(1),dersc1(1),
3675      &            dersc12,dersc01)
3676           call splinthet(x(2),0.5d0*delta,ss,ssd)
3677           dersc0(1)=dersc01
3678           dersc0(2)=dersc02
3679           dersc0(3)=0.0d0
3680           do k=1,3
3681             dersc(k)=ss*dersc(k)+(1.0d0-ss)*dersc0(k)
3682           enddo
3683           dersc(2)=dersc(2)+ssd*(escloci-esclocbi)
3684 c         write (iout,*) 'i=',i,x(2)*rad2deg,escloci0,escloci,
3685 c    &             esclocbi,ss,ssd
3686           escloci=ss*escloci+(1.0d0-ss)*esclocbi
3687 c         escloci=esclocbi
3688 c         write (iout,*) escloci
3689         else if (x(2).lt.delta) then
3690           xtemp(1)=x(1)
3691           xtemp(2)=delta
3692           xtemp(3)=x(3)
3693           call enesc(xtemp,escloci0,dersc0,ddersc0,.true.)
3694           xtemp(2)=0.0d0
3695           call enesc(xtemp,escloci1,dersc1,ddummy,.false.)
3696           call spline1(x(2),delta,-delta,escloci0,escloci1,dersc0(2),
3697      &        escloci,dersc(2))
3698           call spline2(x(2),delta,-delta,dersc0(1),dersc1(1),
3699      &        ddersc0(1),dersc(1))
3700           call spline2(x(2),delta,-delta,dersc0(3),dersc1(3),
3701      &        ddersc0(3),dersc(3))
3702           xtemp(2)=delta
3703           call enesc_bound(xtemp,esclocbi0,dersc0,dersc12,.true.)
3704           xtemp(2)=0.0d0
3705           call enesc_bound(xtemp,esclocbi1,dersc1,chuju,.false.)
3706           call spline1(x(2),delta,-delta,esclocbi0,esclocbi1,
3707      &            dersc0(2),esclocbi,dersc02)
3708           call spline2(x(2),delta,-delta,dersc0(1),dersc1(1),
3709      &            dersc12,dersc01)
3710           dersc0(1)=dersc01
3711           dersc0(2)=dersc02
3712           dersc0(3)=0.0d0
3713           call splinthet(x(2),0.5d0*delta,ss,ssd)
3714           do k=1,3
3715             dersc(k)=ss*dersc(k)+(1.0d0-ss)*dersc0(k)
3716           enddo
3717           dersc(2)=dersc(2)+ssd*(escloci-esclocbi)
3718 c         write (iout,*) 'i=',i,x(2)*rad2deg,escloci0,escloci,
3719 c    &             esclocbi,ss,ssd
3720           escloci=ss*escloci+(1.0d0-ss)*esclocbi
3721 c         write (iout,*) escloci
3722         else
3723           call enesc(x,escloci,dersc,ddummy,.false.)
3724         endif
3725
3726         escloc=escloc+escloci
3727 c        write (iout,*) 'i=',i,' escloci=',escloci,' dersc=',dersc
3728
3729         gloc(nphi+i-1,icg)=gloc(nphi+i-1,icg)+
3730      &   wscloc*dersc(1)
3731         gloc(ialph(i,1),icg)=wscloc*dersc(2)
3732         gloc(ialph(i,1)+nside,icg)=wscloc*dersc(3)
3733     1   continue
3734       enddo
3735       return
3736       end
3737 C---------------------------------------------------------------------------
3738       subroutine enesc(x,escloci,dersc,ddersc,mixed)
3739       implicit real*8 (a-h,o-z)
3740       include 'DIMENSIONS'
3741       include 'COMMON.GEO'
3742       include 'COMMON.LOCAL'
3743       include 'COMMON.IOUNITS'
3744       common /sccalc/ time11,time12,time112,theti,it,nlobit
3745       double precision x(3),z(3),Ax(3,maxlob,-1:1),dersc(3),ddersc(3)
3746       double precision contr(maxlob,-1:1)
3747       logical mixed
3748 c       write (iout,*) 'it=',it,' nlobit=',nlobit
3749         escloc_i=0.0D0
3750         do j=1,3
3751           dersc(j)=0.0D0
3752           if (mixed) ddersc(j)=0.0d0
3753         enddo
3754         x3=x(3)
3755
3756 C Because of periodicity of the dependence of the SC energy in omega we have
3757 C to add up the contributions from x(3)-2*pi, x(3), and x(3+2*pi).
3758 C To avoid underflows, first compute & store the exponents.
3759
3760         do iii=-1,1
3761
3762           x(3)=x3+iii*dwapi
3763  
3764           do j=1,nlobit
3765             do k=1,3
3766               z(k)=x(k)-censc(k,j,it)
3767             enddo
3768             do k=1,3
3769               Axk=0.0D0
3770               do l=1,3
3771                 Axk=Axk+gaussc(l,k,j,it)*z(l)
3772               enddo
3773               Ax(k,j,iii)=Axk
3774             enddo 
3775             expfac=0.0D0 
3776             do k=1,3
3777               expfac=expfac+Ax(k,j,iii)*z(k)
3778             enddo
3779             contr(j,iii)=expfac
3780           enddo ! j
3781
3782         enddo ! iii
3783
3784         x(3)=x3
3785 C As in the case of ebend, we want to avoid underflows in exponentiation and
3786 C subsequent NaNs and INFs in energy calculation.
3787 C Find the largest exponent
3788         emin=contr(1,-1)
3789         do iii=-1,1
3790           do j=1,nlobit
3791             if (emin.gt.contr(j,iii)) emin=contr(j,iii)
3792           enddo 
3793         enddo
3794         emin=0.5D0*emin
3795 cd      print *,'it=',it,' emin=',emin
3796
3797 C Compute the contribution to SC energy and derivatives
3798         do iii=-1,1
3799
3800           do j=1,nlobit
3801             expfac=dexp(bsc(j,iabs(it))-0.5D0*contr(j,iii)+emin)
3802 cd          print *,'j=',j,' expfac=',expfac
3803             escloc_i=escloc_i+expfac
3804             do k=1,3
3805               dersc(k)=dersc(k)+Ax(k,j,iii)*expfac
3806             enddo
3807             if (mixed) then
3808               do k=1,3,2
3809                 ddersc(k)=ddersc(k)+(-Ax(2,j,iii)*Ax(k,j,iii)
3810      &            +gaussc(k,2,j,it))*expfac
3811               enddo
3812             endif
3813           enddo
3814
3815         enddo ! iii
3816
3817         dersc(1)=dersc(1)/cos(theti)**2
3818         ddersc(1)=ddersc(1)/cos(theti)**2
3819         ddersc(3)=ddersc(3)
3820
3821         escloci=-(dlog(escloc_i)-emin)
3822         do j=1,3
3823           dersc(j)=dersc(j)/escloc_i
3824         enddo
3825         if (mixed) then
3826           do j=1,3,2
3827             ddersc(j)=(ddersc(j)/escloc_i+dersc(2)*dersc(j))
3828           enddo
3829         endif
3830       return
3831       end
3832 C------------------------------------------------------------------------------
3833       subroutine enesc_bound(x,escloci,dersc,dersc12,mixed)
3834       implicit real*8 (a-h,o-z)
3835       include 'DIMENSIONS'
3836       include 'COMMON.GEO'
3837       include 'COMMON.LOCAL'
3838       include 'COMMON.IOUNITS'
3839       common /sccalc/ time11,time12,time112,theti,it,nlobit
3840       double precision x(3),z(3),Ax(3,maxlob),dersc(3)
3841       double precision contr(maxlob)
3842       logical mixed
3843
3844       escloc_i=0.0D0
3845
3846       do j=1,3
3847         dersc(j)=0.0D0
3848       enddo
3849
3850       do j=1,nlobit
3851         do k=1,2
3852           z(k)=x(k)-censc(k,j,it)
3853         enddo
3854         z(3)=dwapi
3855         do k=1,3
3856           Axk=0.0D0
3857           do l=1,3
3858             Axk=Axk+gaussc(l,k,j,it)*z(l)
3859           enddo
3860           Ax(k,j)=Axk
3861         enddo 
3862         expfac=0.0D0 
3863         do k=1,3
3864           expfac=expfac+Ax(k,j)*z(k)
3865         enddo
3866         contr(j)=expfac
3867       enddo ! j
3868
3869 C As in the case of ebend, we want to avoid underflows in exponentiation and
3870 C subsequent NaNs and INFs in energy calculation.
3871 C Find the largest exponent
3872       emin=contr(1)
3873       do j=1,nlobit
3874         if (emin.gt.contr(j)) emin=contr(j)
3875       enddo 
3876       emin=0.5D0*emin
3877  
3878 C Compute the contribution to SC energy and derivatives
3879
3880       dersc12=0.0d0
3881       do j=1,nlobit
3882         expfac=dexp(bsc(j,iabs(it))-0.5D0*contr(j)+emin)
3883         escloc_i=escloc_i+expfac
3884         do k=1,2
3885           dersc(k)=dersc(k)+Ax(k,j)*expfac
3886         enddo
3887         if (mixed) dersc12=dersc12+(-Ax(2,j)*Ax(1,j)
3888      &            +gaussc(1,2,j,it))*expfac
3889         dersc(3)=0.0d0
3890       enddo
3891
3892       dersc(1)=dersc(1)/cos(theti)**2
3893       dersc12=dersc12/cos(theti)**2
3894       escloci=-(dlog(escloc_i)-emin)
3895       do j=1,2
3896         dersc(j)=dersc(j)/escloc_i
3897       enddo
3898       if (mixed) dersc12=(dersc12/escloc_i+dersc(2)*dersc(1))
3899       return
3900       end
3901 #else
3902 c----------------------------------------------------------------------------------
3903       subroutine esc(escloc)
3904 C Calculate the local energy of a side chain and its derivatives in the
3905 C corresponding virtual-bond valence angles THETA and the spherical angles 
3906 C ALPHA and OMEGA derived from AM1 all-atom calculations.
3907 C added by Urszula Kozlowska. 07/11/2007
3908 C
3909       implicit real*8 (a-h,o-z)
3910       include 'DIMENSIONS'
3911       include 'DIMENSIONS.ZSCOPT'
3912       include 'COMMON.GEO'
3913       include 'COMMON.LOCAL'
3914       include 'COMMON.VAR'
3915       include 'COMMON.SCROT'
3916       include 'COMMON.INTERACT'
3917       include 'COMMON.DERIV'
3918       include 'COMMON.CHAIN'
3919       include 'COMMON.IOUNITS'
3920       include 'COMMON.NAMES'
3921       include 'COMMON.FFIELD'
3922       include 'COMMON.CONTROL'
3923       include 'COMMON.VECTORS'
3924       double precision x_prime(3),y_prime(3),z_prime(3)
3925      &    , sumene,dsc_i,dp2_i,x(65),
3926      &     xx,yy,zz,sumene1,sumene2,sumene3,sumene4,s1,s1_6,s2,s2_6,
3927      &    de_dxx,de_dyy,de_dzz,de_dt
3928       double precision s1_t,s1_6_t,s2_t,s2_6_t
3929       double precision 
3930      & dXX_Ci1(3),dYY_Ci1(3),dZZ_Ci1(3),dXX_Ci(3),
3931      & dYY_Ci(3),dZZ_Ci(3),dXX_XYZ(3),dYY_XYZ(3),dZZ_XYZ(3),
3932      & dt_dCi(3),dt_dCi1(3)
3933       common /sccalc/ time11,time12,time112,theti,it,nlobit
3934       delta=0.02d0*pi
3935       escloc=0.0D0
3936       do i=loc_start,loc_end
3937         if (itype(i).eq.ntyp1) cycle
3938         costtab(i+1) =dcos(theta(i+1))
3939         sinttab(i+1) =dsqrt(1-costtab(i+1)*costtab(i+1))
3940         cost2tab(i+1)=dsqrt(0.5d0*(1.0d0+costtab(i+1)))
3941         sint2tab(i+1)=dsqrt(0.5d0*(1.0d0-costtab(i+1)))
3942         cosfac2=0.5d0/(1.0d0+costtab(i+1))
3943         cosfac=dsqrt(cosfac2)
3944         sinfac2=0.5d0/(1.0d0-costtab(i+1))
3945         sinfac=dsqrt(sinfac2)
3946         it=iabs(itype(i))
3947         if (it.eq.10) goto 1
3948 c
3949 C  Compute the axes of tghe local cartesian coordinates system; store in
3950 c   x_prime, y_prime and z_prime 
3951 c
3952         do j=1,3
3953           x_prime(j) = 0.00
3954           y_prime(j) = 0.00
3955           z_prime(j) = 0.00
3956         enddo
3957 C        write(2,*) "dc_norm", dc_norm(1,i+nres),dc_norm(2,i+nres),
3958 C     &   dc_norm(3,i+nres)
3959         do j = 1,3
3960           x_prime(j) = (dc_norm(j,i) - dc_norm(j,i-1))*cosfac
3961           y_prime(j) = (dc_norm(j,i) + dc_norm(j,i-1))*sinfac
3962         enddo
3963         do j = 1,3
3964           z_prime(j) = -uz(j,i-1)*dsign(1.0d0,dfloat(itype(i)))
3965         enddo     
3966 c       write (2,*) "i",i
3967 c       write (2,*) "x_prime",(x_prime(j),j=1,3)
3968 c       write (2,*) "y_prime",(y_prime(j),j=1,3)
3969 c       write (2,*) "z_prime",(z_prime(j),j=1,3)
3970 c       write (2,*) "xx",scalar(x_prime(1),x_prime(1)),
3971 c      & " xy",scalar(x_prime(1),y_prime(1)),
3972 c      & " xz",scalar(x_prime(1),z_prime(1)),
3973 c      & " yy",scalar(y_prime(1),y_prime(1)),
3974 c      & " yz",scalar(y_prime(1),z_prime(1)),
3975 c      & " zz",scalar(z_prime(1),z_prime(1))
3976 c
3977 C Transform the unit vector of the ith side-chain centroid, dC_norm(*,i),
3978 C to local coordinate system. Store in xx, yy, zz.
3979 c
3980         xx=0.0d0
3981         yy=0.0d0
3982         zz=0.0d0
3983         do j = 1,3
3984           xx = xx + x_prime(j)*dc_norm(j,i+nres)
3985           yy = yy + y_prime(j)*dc_norm(j,i+nres)
3986           zz = zz + z_prime(j)*dc_norm(j,i+nres)
3987         enddo
3988
3989         xxtab(i)=xx
3990         yytab(i)=yy
3991         zztab(i)=zz
3992 C
3993 C Compute the energy of the ith side cbain
3994 C
3995 c        write (2,*) "xx",xx," yy",yy," zz",zz
3996         it=iabs(itype(i))
3997         do j = 1,65
3998           x(j) = sc_parmin(j,it) 
3999         enddo
4000 #ifdef CHECK_COORD
4001 Cc diagnostics - remove later
4002         xx1 = dcos(alph(2))
4003         yy1 = dsin(alph(2))*dcos(omeg(2))
4004         zz1 = -dsign(1.0d0,itype(i))*dsin(alph(2))*dsin(omeg(2))
4005         write(2,'(3f8.1,3f9.3,1x,3f9.3)') 
4006      &    alph(2)*rad2deg,omeg(2)*rad2deg,theta(3)*rad2deg,xx,yy,zz,
4007      &    xx1,yy1,zz1
4008 C,"  --- ", xx_w,yy_w,zz_w
4009 c end diagnostics
4010 #endif
4011         sumene1= x(1)+  x(2)*xx+  x(3)*yy+  x(4)*zz+  x(5)*xx**2
4012      &   + x(6)*yy**2+  x(7)*zz**2+  x(8)*xx*zz+  x(9)*xx*yy
4013      &   + x(10)*yy*zz
4014         sumene2=  x(11) + x(12)*xx + x(13)*yy + x(14)*zz + x(15)*xx**2
4015      & + x(16)*yy**2 + x(17)*zz**2 + x(18)*xx*zz + x(19)*xx*yy
4016      & + x(20)*yy*zz
4017         sumene3=  x(21) +x(22)*xx +x(23)*yy +x(24)*zz +x(25)*xx**2
4018      &  +x(26)*yy**2 +x(27)*zz**2 +x(28)*xx*zz +x(29)*xx*yy
4019      &  +x(30)*yy*zz +x(31)*xx**3 +x(32)*yy**3 +x(33)*zz**3
4020      &  +x(34)*(xx**2)*yy +x(35)*(xx**2)*zz +x(36)*(yy**2)*xx
4021      &  +x(37)*(yy**2)*zz +x(38)*(zz**2)*xx +x(39)*(zz**2)*yy
4022      &  +x(40)*xx*yy*zz
4023         sumene4= x(41) +x(42)*xx +x(43)*yy +x(44)*zz +x(45)*xx**2
4024      &  +x(46)*yy**2 +x(47)*zz**2 +x(48)*xx*zz +x(49)*xx*yy
4025      &  +x(50)*yy*zz +x(51)*xx**3 +x(52)*yy**3 +x(53)*zz**3
4026      &  +x(54)*(xx**2)*yy +x(55)*(xx**2)*zz +x(56)*(yy**2)*xx
4027      &  +x(57)*(yy**2)*zz +x(58)*(zz**2)*xx +x(59)*(zz**2)*yy
4028      &  +x(60)*xx*yy*zz
4029         dsc_i   = 0.743d0+x(61)
4030         dp2_i   = 1.9d0+x(62)
4031         dscp1=dsqrt(dsc_i**2+dp2_i**2-2*dsc_i*dp2_i
4032      &          *(xx*cost2tab(i+1)+yy*sint2tab(i+1)))
4033         dscp2=dsqrt(dsc_i**2+dp2_i**2-2*dsc_i*dp2_i
4034      &          *(xx*cost2tab(i+1)-yy*sint2tab(i+1)))
4035         s1=(1+x(63))/(0.1d0 + dscp1)
4036         s1_6=(1+x(64))/(0.1d0 + dscp1**6)
4037         s2=(1+x(65))/(0.1d0 + dscp2)
4038         s2_6=(1+x(65))/(0.1d0 + dscp2**6)
4039         sumene = ( sumene3*sint2tab(i+1) + sumene1)*(s1+s1_6)
4040      & + (sumene4*cost2tab(i+1) +sumene2)*(s2+s2_6)
4041 c        write(2,'(i2," sumene",7f9.3)') i,sumene1,sumene2,sumene3,
4042 c     &   sumene4,
4043 c     &   dscp1,dscp2,sumene
4044 c        sumene = enesc(x,xx,yy,zz,cost2tab(i+1),sint2tab(i+1))
4045         escloc = escloc + sumene
4046 c        write (2,*) "escloc",escloc
4047 c        write (2,*) "i",i," escloc",sumene,escloc,it,itype(i),
4048 c     &  zz,xx,yy
4049         if (.not. calc_grad) goto 1
4050 #ifdef DEBUG
4051 C
4052 C This section to check the numerical derivatives of the energy of ith side
4053 C chain in xx, yy, zz, and theta. Use the -DDEBUG compiler option or insert
4054 C #define DEBUG in the code to turn it on.
4055 C
4056         write (2,*) "sumene               =",sumene
4057         aincr=1.0d-7
4058         xxsave=xx
4059         xx=xx+aincr
4060         write (2,*) xx,yy,zz
4061         sumenep = enesc(x,xx,yy,zz,cost2tab(i+1),sint2tab(i+1))
4062         de_dxx_num=(sumenep-sumene)/aincr
4063         xx=xxsave
4064         write (2,*) "xx+ sumene from enesc=",sumenep
4065         yysave=yy
4066         yy=yy+aincr
4067         write (2,*) xx,yy,zz
4068         sumenep = enesc(x,xx,yy,zz,cost2tab(i+1),sint2tab(i+1))
4069         de_dyy_num=(sumenep-sumene)/aincr
4070         yy=yysave
4071         write (2,*) "yy+ sumene from enesc=",sumenep
4072         zzsave=zz
4073         zz=zz+aincr
4074         write (2,*) xx,yy,zz
4075         sumenep = enesc(x,xx,yy,zz,cost2tab(i+1),sint2tab(i+1))
4076         de_dzz_num=(sumenep-sumene)/aincr
4077         zz=zzsave
4078         write (2,*) "zz+ sumene from enesc=",sumenep
4079         costsave=cost2tab(i+1)
4080         sintsave=sint2tab(i+1)
4081         cost2tab(i+1)=dcos(0.5d0*(theta(i+1)+aincr))
4082         sint2tab(i+1)=dsin(0.5d0*(theta(i+1)+aincr))
4083         sumenep = enesc(x,xx,yy,zz,cost2tab(i+1),sint2tab(i+1))
4084         de_dt_num=(sumenep-sumene)/aincr
4085         write (2,*) " t+ sumene from enesc=",sumenep
4086         cost2tab(i+1)=costsave
4087         sint2tab(i+1)=sintsave
4088 C End of diagnostics section.
4089 #endif
4090 C        
4091 C Compute the gradient of esc
4092 C
4093         pom_s1=(1.0d0+x(63))/(0.1d0 + dscp1)**2
4094         pom_s16=6*(1.0d0+x(64))/(0.1d0 + dscp1**6)**2
4095         pom_s2=(1.0d0+x(65))/(0.1d0 + dscp2)**2
4096         pom_s26=6*(1.0d0+x(65))/(0.1d0 + dscp2**6)**2
4097         pom_dx=dsc_i*dp2_i*cost2tab(i+1)
4098         pom_dy=dsc_i*dp2_i*sint2tab(i+1)
4099         pom_dt1=-0.5d0*dsc_i*dp2_i*(xx*sint2tab(i+1)-yy*cost2tab(i+1))
4100         pom_dt2=-0.5d0*dsc_i*dp2_i*(xx*sint2tab(i+1)+yy*cost2tab(i+1))
4101         pom1=(sumene3*sint2tab(i+1)+sumene1)
4102      &     *(pom_s1/dscp1+pom_s16*dscp1**4)
4103         pom2=(sumene4*cost2tab(i+1)+sumene2)
4104      &     *(pom_s2/dscp2+pom_s26*dscp2**4)
4105         sumene1x=x(2)+2*x(5)*xx+x(8)*zz+ x(9)*yy
4106         sumene3x=x(22)+2*x(25)*xx+x(28)*zz+x(29)*yy+3*x(31)*xx**2
4107      &  +2*x(34)*xx*yy +2*x(35)*xx*zz +x(36)*(yy**2) +x(38)*(zz**2)
4108      &  +x(40)*yy*zz
4109         sumene2x=x(12)+2*x(15)*xx+x(18)*zz+ x(19)*yy
4110         sumene4x=x(42)+2*x(45)*xx +x(48)*zz +x(49)*yy +3*x(51)*xx**2
4111      &  +2*x(54)*xx*yy+2*x(55)*xx*zz+x(56)*(yy**2)+x(58)*(zz**2)
4112      &  +x(60)*yy*zz
4113         de_dxx =(sumene1x+sumene3x*sint2tab(i+1))*(s1+s1_6)
4114      &        +(sumene2x+sumene4x*cost2tab(i+1))*(s2+s2_6)
4115      &        +(pom1+pom2)*pom_dx
4116 #ifdef DEBUG
4117         write(2,*), "de_dxx = ", de_dxx,de_dxx_num
4118 #endif
4119 C
4120         sumene1y=x(3) + 2*x(6)*yy + x(9)*xx + x(10)*zz
4121         sumene3y=x(23) +2*x(26)*yy +x(29)*xx +x(30)*zz +3*x(32)*yy**2
4122      &  +x(34)*(xx**2) +2*x(36)*yy*xx +2*x(37)*yy*zz +x(39)*(zz**2)
4123      &  +x(40)*xx*zz
4124         sumene2y=x(13) + 2*x(16)*yy + x(19)*xx + x(20)*zz
4125         sumene4y=x(43)+2*x(46)*yy+x(49)*xx +x(50)*zz
4126      &  +3*x(52)*yy**2+x(54)*xx**2+2*x(56)*yy*xx +2*x(57)*yy*zz
4127      &  +x(59)*zz**2 +x(60)*xx*zz
4128         de_dyy =(sumene1y+sumene3y*sint2tab(i+1))*(s1+s1_6)
4129      &        +(sumene2y+sumene4y*cost2tab(i+1))*(s2+s2_6)
4130      &        +(pom1-pom2)*pom_dy
4131 #ifdef DEBUG
4132         write(2,*), "de_dyy = ", de_dyy,de_dyy_num
4133 #endif
4134 C
4135         de_dzz =(x(24) +2*x(27)*zz +x(28)*xx +x(30)*yy
4136      &  +3*x(33)*zz**2 +x(35)*xx**2 +x(37)*yy**2 +2*x(38)*zz*xx 
4137      &  +2*x(39)*zz*yy +x(40)*xx*yy)*sint2tab(i+1)*(s1+s1_6) 
4138      &  +(x(4) + 2*x(7)*zz+  x(8)*xx + x(10)*yy)*(s1+s1_6) 
4139      &  +(x(44)+2*x(47)*zz +x(48)*xx   +x(50)*yy  +3*x(53)*zz**2   
4140      &  +x(55)*xx**2 +x(57)*(yy**2)+2*x(58)*zz*xx +2*x(59)*zz*yy  
4141      &  +x(60)*xx*yy)*cost2tab(i+1)*(s2+s2_6)
4142      &  + ( x(14) + 2*x(17)*zz+  x(18)*xx + x(20)*yy)*(s2+s2_6)
4143 #ifdef DEBUG
4144         write(2,*), "de_dzz = ", de_dzz,de_dzz_num
4145 #endif
4146 C
4147         de_dt =  0.5d0*sumene3*cost2tab(i+1)*(s1+s1_6) 
4148      &  -0.5d0*sumene4*sint2tab(i+1)*(s2+s2_6)
4149      &  +pom1*pom_dt1+pom2*pom_dt2
4150 #ifdef DEBUG
4151         write(2,*), "de_dt = ", de_dt,de_dt_num
4152 #endif
4153
4154 C
4155        cossc=scalar(dc_norm(1,i),dc_norm(1,i+nres))
4156        cossc1=scalar(dc_norm(1,i-1),dc_norm(1,i+nres))
4157        cosfac2xx=cosfac2*xx
4158        sinfac2yy=sinfac2*yy
4159        do k = 1,3
4160          dt_dCi(k) = -(dc_norm(k,i-1)+costtab(i+1)*dc_norm(k,i))*
4161      &      vbld_inv(i+1)
4162          dt_dCi1(k)= -(dc_norm(k,i)+costtab(i+1)*dc_norm(k,i-1))*
4163      &      vbld_inv(i)
4164          pom=(dC_norm(k,i+nres)-cossc*dC_norm(k,i))*vbld_inv(i+1)
4165          pom1=(dC_norm(k,i+nres)-cossc1*dC_norm(k,i-1))*vbld_inv(i)
4166 c         write (iout,*) "i",i," k",k," pom",pom," pom1",pom1,
4167 c     &    " dt_dCi",dt_dCi(k)," dt_dCi1",dt_dCi1(k)
4168 c         write (iout,*) "dC_norm",(dC_norm(j,i),j=1,3),
4169 c     &   (dC_norm(j,i-1),j=1,3)," vbld_inv",vbld_inv(i+1),vbld_inv(i)
4170          dXX_Ci(k)=pom*cosfac-dt_dCi(k)*cosfac2xx
4171          dXX_Ci1(k)=-pom1*cosfac-dt_dCi1(k)*cosfac2xx
4172          dYY_Ci(k)=pom*sinfac+dt_dCi(k)*sinfac2yy
4173          dYY_Ci1(k)=pom1*sinfac+dt_dCi1(k)*sinfac2yy
4174          dZZ_Ci1(k)=0.0d0
4175          dZZ_Ci(k)=0.0d0
4176          do j=1,3
4177            dZZ_Ci(k)=dZZ_Ci(k)-uzgrad(j,k,2,i-1)
4178      & *dsign(1.0d0,dfloat(itype(i)))*dC_norm(j,i+nres)
4179            dZZ_Ci1(k)=dZZ_Ci1(k)-uzgrad(j,k,1,i-1)
4180      &  *dsign(1.0d0,dfloat(itype(i)))*dC_norm(j,i+nres)
4181          enddo
4182           
4183          dXX_XYZ(k)=vbld_inv(i+nres)*(x_prime(k)-xx*dC_norm(k,i+nres))
4184          dYY_XYZ(k)=vbld_inv(i+nres)*(y_prime(k)-yy*dC_norm(k,i+nres))
4185          dZZ_XYZ(k)=vbld_inv(i+nres)*(z_prime(k)-zz*dC_norm(k,i+nres))
4186 c
4187          dt_dCi(k) = -dt_dCi(k)/sinttab(i+1)
4188          dt_dCi1(k)= -dt_dCi1(k)/sinttab(i+1)
4189        enddo
4190
4191        do k=1,3
4192          dXX_Ctab(k,i)=dXX_Ci(k)
4193          dXX_C1tab(k,i)=dXX_Ci1(k)
4194          dYY_Ctab(k,i)=dYY_Ci(k)
4195          dYY_C1tab(k,i)=dYY_Ci1(k)
4196          dZZ_Ctab(k,i)=dZZ_Ci(k)
4197          dZZ_C1tab(k,i)=dZZ_Ci1(k)
4198          dXX_XYZtab(k,i)=dXX_XYZ(k)
4199          dYY_XYZtab(k,i)=dYY_XYZ(k)
4200          dZZ_XYZtab(k,i)=dZZ_XYZ(k)
4201        enddo
4202
4203        do k = 1,3
4204 c         write (iout,*) "k",k," dxx_ci1",dxx_ci1(k)," dyy_ci1",
4205 c     &    dyy_ci1(k)," dzz_ci1",dzz_ci1(k)
4206 c         write (iout,*) "k",k," dxx_ci",dxx_ci(k)," dyy_ci",
4207 c     &    dyy_ci(k)," dzz_ci",dzz_ci(k)
4208 c         write (iout,*) "k",k," dt_dci",dt_dci(k)," dt_dci",
4209 c     &    dt_dci(k)
4210 c         write (iout,*) "k",k," dxx_XYZ",dxx_XYZ(k)," dyy_XYZ",
4211 c     &    dyy_XYZ(k)," dzz_XYZ",dzz_XYZ(k) 
4212          gscloc(k,i-1)=gscloc(k,i-1)+de_dxx*dxx_ci1(k)
4213      &    +de_dyy*dyy_ci1(k)+de_dzz*dzz_ci1(k)+de_dt*dt_dCi1(k)
4214          gscloc(k,i)=gscloc(k,i)+de_dxx*dxx_Ci(k)
4215      &    +de_dyy*dyy_Ci(k)+de_dzz*dzz_Ci(k)+de_dt*dt_dCi(k)
4216          gsclocx(k,i)=                 de_dxx*dxx_XYZ(k)
4217      &    +de_dyy*dyy_XYZ(k)+de_dzz*dzz_XYZ(k)
4218        enddo
4219 c       write(iout,*) "ENERGY GRAD = ", (gscloc(k,i-1),k=1,3),
4220 c     &  (gscloc(k,i),k=1,3),(gsclocx(k,i),k=1,3)  
4221
4222 C to check gradient call subroutine check_grad
4223
4224     1 continue
4225       enddo
4226       return
4227       end
4228 #endif
4229 c------------------------------------------------------------------------------
4230       subroutine gcont(rij,r0ij,eps0ij,delta,fcont,fprimcont)
4231 C
4232 C This procedure calculates two-body contact function g(rij) and its derivative:
4233 C
4234 C           eps0ij                                     !       x < -1
4235 C g(rij) =  esp0ij*(-0.9375*x+0.625*x**3-0.1875*x**5)  ! -1 =< x =< 1
4236 C            0                                         !       x > 1
4237 C
4238 C where x=(rij-r0ij)/delta
4239 C
4240 C rij - interbody distance, r0ij - contact distance, eps0ij - contact energy
4241 C
4242       implicit none
4243       double precision rij,r0ij,eps0ij,fcont,fprimcont
4244       double precision x,x2,x4,delta
4245 c     delta=0.02D0*r0ij
4246 c      delta=0.2D0*r0ij
4247       x=(rij-r0ij)/delta
4248       if (x.lt.-1.0D0) then
4249         fcont=eps0ij
4250         fprimcont=0.0D0
4251       else if (x.le.1.0D0) then  
4252         x2=x*x
4253         x4=x2*x2
4254         fcont=eps0ij*(x*(-0.9375D0+0.6250D0*x2-0.1875D0*x4)+0.5D0)
4255         fprimcont=eps0ij * (-0.9375D0+1.8750D0*x2-0.9375D0*x4)/delta
4256       else
4257         fcont=0.0D0
4258         fprimcont=0.0D0
4259       endif
4260       return
4261       end
4262 c------------------------------------------------------------------------------
4263       subroutine splinthet(theti,delta,ss,ssder)
4264       implicit real*8 (a-h,o-z)
4265       include 'DIMENSIONS'
4266       include 'DIMENSIONS.ZSCOPT'
4267       include 'COMMON.VAR'
4268       include 'COMMON.GEO'
4269       thetup=pi-delta
4270       thetlow=delta
4271       if (theti.gt.pipol) then
4272         call gcont(theti,thetup,1.0d0,delta,ss,ssder)
4273       else
4274         call gcont(-theti,-thetlow,1.0d0,delta,ss,ssder)
4275         ssder=-ssder
4276       endif
4277       return
4278       end
4279 c------------------------------------------------------------------------------
4280       subroutine spline1(x,x0,delta,f0,f1,fprim0,f,fprim)
4281       implicit none
4282       double precision x,x0,delta,f0,f1,fprim0,f,fprim
4283       double precision ksi,ksi2,ksi3,a1,a2,a3
4284       a1=fprim0*delta/(f1-f0)
4285       a2=3.0d0-2.0d0*a1
4286       a3=a1-2.0d0
4287       ksi=(x-x0)/delta
4288       ksi2=ksi*ksi
4289       ksi3=ksi2*ksi  
4290       f=f0+(f1-f0)*ksi*(a1+ksi*(a2+a3*ksi))
4291       fprim=(f1-f0)/delta*(a1+ksi*(2*a2+3*ksi*a3))
4292       return
4293       end
4294 c------------------------------------------------------------------------------
4295       subroutine spline2(x,x0,delta,f0x,f1x,fprim0x,fx)
4296       implicit none
4297       double precision x,x0,delta,f0x,f1x,fprim0x,fx
4298       double precision ksi,ksi2,ksi3,a1,a2,a3
4299       ksi=(x-x0)/delta  
4300       ksi2=ksi*ksi
4301       ksi3=ksi2*ksi
4302       a1=fprim0x*delta
4303       a2=3*(f1x-f0x)-2*fprim0x*delta
4304       a3=fprim0x*delta-2*(f1x-f0x)
4305       fx=f0x+a1*ksi+a2*ksi2+a3*ksi3
4306       return
4307       end
4308 C-----------------------------------------------------------------------------
4309 #ifdef CRYST_TOR
4310 C-----------------------------------------------------------------------------
4311       subroutine etor(etors,edihcnstr,fact)
4312       implicit real*8 (a-h,o-z)
4313       include 'DIMENSIONS'
4314       include 'DIMENSIONS.ZSCOPT'
4315       include 'COMMON.VAR'
4316       include 'COMMON.GEO'
4317       include 'COMMON.LOCAL'
4318       include 'COMMON.TORSION'
4319       include 'COMMON.INTERACT'
4320       include 'COMMON.DERIV'
4321       include 'COMMON.CHAIN'
4322       include 'COMMON.NAMES'
4323       include 'COMMON.IOUNITS'
4324       include 'COMMON.FFIELD'
4325       include 'COMMON.TORCNSTR'
4326       logical lprn
4327 C Set lprn=.true. for debugging
4328       lprn=.false.
4329 c      lprn=.true.
4330       etors=0.0D0
4331       do i=iphi_start,iphi_end
4332         if (itype(i-2).eq.ntyp1 .or. itype(i-1).eq.ntyp1
4333      &      .or. itype(i).eq.ntyp1) cycle
4334         itori=itortyp(itype(i-2))
4335         itori1=itortyp(itype(i-1))
4336         phii=phi(i)
4337         gloci=0.0D0
4338 C Proline-Proline pair is a special case...
4339         if (itori.eq.3 .and. itori1.eq.3) then
4340           if (phii.gt.-dwapi3) then
4341             cosphi=dcos(3*phii)
4342             fac=1.0D0/(1.0D0-cosphi)
4343             etorsi=v1(1,3,3)*fac
4344             etorsi=etorsi+etorsi
4345             etors=etors+etorsi-v1(1,3,3)
4346             gloci=gloci-3*fac*etorsi*dsin(3*phii)
4347           endif
4348           do j=1,3
4349             v1ij=v1(j+1,itori,itori1)
4350             v2ij=v2(j+1,itori,itori1)
4351             cosphi=dcos(j*phii)
4352             sinphi=dsin(j*phii)
4353             etors=etors+v1ij*cosphi+v2ij*sinphi+dabs(v1ij)+dabs(v2ij)
4354             gloci=gloci+j*(v2ij*cosphi-v1ij*sinphi)
4355           enddo
4356         else 
4357           do j=1,nterm_old
4358             v1ij=v1(j,itori,itori1)
4359             v2ij=v2(j,itori,itori1)
4360             cosphi=dcos(j*phii)
4361             sinphi=dsin(j*phii)
4362             etors=etors+v1ij*cosphi+v2ij*sinphi+dabs(v1ij)+dabs(v2ij)
4363             gloci=gloci+j*(v2ij*cosphi-v1ij*sinphi)
4364           enddo
4365         endif
4366         if (lprn)
4367      &  write (iout,'(2(a3,2x,i3,2x),2i3,6f8.3/26x,6f8.3/)')
4368      &  restyp(itype(i-2)),i-2,restyp(itype(i-1)),i-1,itori,itori1,
4369      &  (v1(j,itori,itori1),j=1,6),(v2(j,itori,itori1),j=1,6)
4370         gloc(i-3,icg)=gloc(i-3,icg)+wtor*fact*gloci
4371 c       write (iout,*) 'i=',i,' gloc=',gloc(i-3,icg)
4372       enddo
4373 ! 6/20/98 - dihedral angle constraints
4374       edihcnstr=0.0d0
4375       do i=1,ndih_constr
4376         itori=idih_constr(i)
4377         phii=phi(itori)
4378         difi=phii-phi0(i)
4379         if (difi.gt.drange(i)) then
4380           difi=difi-drange(i)
4381           edihcnstr=edihcnstr+0.25d0*ftors*difi**4
4382           gloc(itori-3,icg)=gloc(itori-3,icg)+ftors*difi**3
4383         else if (difi.lt.-drange(i)) then
4384           difi=difi+drange(i)
4385           edihcnstr=edihcnstr+0.25d0*ftors*difi**4
4386           gloc(itori-3,icg)=gloc(itori-3,icg)+ftors*difi**3
4387         endif
4388 !        write (iout,'(2i5,2f8.3,2e14.5)') i,itori,rad2deg*phii,
4389 !     &    rad2deg*difi,0.25d0*ftors*difi**4,gloc(itori-3,icg)
4390       enddo
4391 !      write (iout,*) 'edihcnstr',edihcnstr
4392       return
4393       end
4394 c------------------------------------------------------------------------------
4395 #else
4396       subroutine etor(etors,edihcnstr,fact)
4397       implicit real*8 (a-h,o-z)
4398       include 'DIMENSIONS'
4399       include 'DIMENSIONS.ZSCOPT'
4400       include 'COMMON.VAR'
4401       include 'COMMON.GEO'
4402       include 'COMMON.LOCAL'
4403       include 'COMMON.TORSION'
4404       include 'COMMON.INTERACT'
4405       include 'COMMON.DERIV'
4406       include 'COMMON.CHAIN'
4407       include 'COMMON.NAMES'
4408       include 'COMMON.IOUNITS'
4409       include 'COMMON.FFIELD'
4410       include 'COMMON.TORCNSTR'
4411       logical lprn
4412 C Set lprn=.true. for debugging
4413       lprn=.false.
4414 c      lprn=.true.
4415       etors=0.0D0
4416       do i=iphi_start,iphi_end
4417         if (itype(i-2).eq.ntyp1 .or. itype(i-1).eq.ntyp1
4418      &       .or. itype(i).eq.ntyp1) cycle
4419         if (itel(i-2).eq.0 .or. itel(i-1).eq.0) goto 1215
4420          if (iabs(itype(i)).eq.20) then
4421          iblock=2
4422          else
4423          iblock=1
4424          endif
4425         itori=itortyp(itype(i-2))
4426         itori1=itortyp(itype(i-1))
4427         phii=phi(i)
4428         gloci=0.0D0
4429 C Regular cosine and sine terms
4430         do j=1,nterm(itori,itori1,iblock)
4431           v1ij=v1(j,itori,itori1,iblock)
4432           v2ij=v2(j,itori,itori1,iblock)
4433           cosphi=dcos(j*phii)
4434           sinphi=dsin(j*phii)
4435           etors=etors+v1ij*cosphi+v2ij*sinphi
4436           gloci=gloci+j*(v2ij*cosphi-v1ij*sinphi)
4437         enddo
4438 C Lorentz terms
4439 C                         v1
4440 C  E = SUM ----------------------------------- - v1
4441 C          [v2 cos(phi/2)+v3 sin(phi/2)]^2 + 1
4442 C
4443         cosphi=dcos(0.5d0*phii)
4444         sinphi=dsin(0.5d0*phii)
4445         do j=1,nlor(itori,itori1,iblock)
4446           vl1ij=vlor1(j,itori,itori1)
4447           vl2ij=vlor2(j,itori,itori1)
4448           vl3ij=vlor3(j,itori,itori1)
4449           pom=vl2ij*cosphi+vl3ij*sinphi
4450           pom1=1.0d0/(pom*pom+1.0d0)
4451           etors=etors+vl1ij*pom1
4452 c          if (energy_dec) etors_ii=etors_ii+
4453 c     &                vl1ij*pom1
4454           pom=-pom*pom1*pom1
4455           gloci=gloci+vl1ij*(vl3ij*cosphi-vl2ij*sinphi)*pom
4456         enddo
4457 C Subtract the constant term
4458         etors=etors-v0(itori,itori1,iblock)
4459         if (lprn)
4460      &  write (iout,'(2(a3,2x,i3,2x),2i3,6f8.3/26x,6f8.3/)')
4461      &  restyp(itype(i-2)),i-2,restyp(itype(i-1)),i-1,itori,itori1,
4462      &  (v1(j,itori,itori1,1),j=1,6),(v2(j,itori,itori1,1),j=1,6)
4463         gloc(i-3,icg)=gloc(i-3,icg)+wtor*fact*gloci
4464 c       write (iout,*) 'i=',i,' gloc=',gloc(i-3,icg)
4465  1215   continue
4466       enddo
4467 ! 6/20/98 - dihedral angle constraints
4468       edihcnstr=0.0d0
4469       do i=1,ndih_constr
4470         itori=idih_constr(i)
4471         phii=phi(itori)
4472         difi=pinorm(phii-phi0(i))
4473         edihi=0.0d0
4474         if (difi.gt.drange(i)) then
4475           difi=difi-drange(i)
4476           edihcnstr=edihcnstr+0.25d0*ftors*difi**4
4477           gloc(itori-3,icg)=gloc(itori-3,icg)+ftors*difi**3
4478           edihi=0.25d0*ftors*difi**4
4479         else if (difi.lt.-drange(i)) then
4480           difi=difi+drange(i)
4481           edihcnstr=edihcnstr+0.25d0*ftors*difi**4
4482           gloc(itori-3,icg)=gloc(itori-3,icg)+ftors*difi**3
4483           edihi=0.25d0*ftors*difi**4
4484         else
4485           difi=0.0d0
4486         endif
4487 c        write (iout,'(2i5,4f10.5,e15.5)') i,itori,phii,phi0(i),difi,
4488 c     &    drange(i),edihi
4489 !        write (iout,'(2i5,2f8.3,2e14.5)') i,itori,rad2deg*phii,
4490 !     &    rad2deg*difi,0.25d0*ftors*difi**4,gloc(itori-3,icg)
4491       enddo
4492 !      write (iout,*) 'edihcnstr',edihcnstr
4493       return
4494       end
4495 c----------------------------------------------------------------------------
4496       subroutine etor_d(etors_d,fact2)
4497 C 6/23/01 Compute double torsional energy
4498       implicit real*8 (a-h,o-z)
4499       include 'DIMENSIONS'
4500       include 'DIMENSIONS.ZSCOPT'
4501       include 'COMMON.VAR'
4502       include 'COMMON.GEO'
4503       include 'COMMON.LOCAL'
4504       include 'COMMON.TORSION'
4505       include 'COMMON.INTERACT'
4506       include 'COMMON.DERIV'
4507       include 'COMMON.CHAIN'
4508       include 'COMMON.NAMES'
4509       include 'COMMON.IOUNITS'
4510       include 'COMMON.FFIELD'
4511       include 'COMMON.TORCNSTR'
4512       logical lprn
4513 C Set lprn=.true. for debugging
4514       lprn=.false.
4515 c     lprn=.true.
4516       etors_d=0.0D0
4517       do i=iphi_start,iphi_end-1
4518         if (itype(i-2).eq.ntyp1.or. itype(i-1).eq.ntyp1
4519      &      .or. itype(i).eq.ntyp1 .or. itype(i+1).eq.ntyp1) cycle
4520         if (itel(i-2).eq.0 .or. itel(i-1).eq.0 .or. itel(i).eq.0) 
4521      &     goto 1215
4522         itori=itortyp(itype(i-2))
4523         itori1=itortyp(itype(i-1))
4524         itori2=itortyp(itype(i))
4525         phii=phi(i)
4526         phii1=phi(i+1)
4527         gloci1=0.0D0
4528         gloci2=0.0D0
4529         iblock=1
4530         if (iabs(itype(i+1)).eq.20) iblock=2
4531 C Regular cosine and sine terms
4532         do j=1,ntermd_1(itori,itori1,itori2,iblock)
4533           v1cij=v1c(1,j,itori,itori1,itori2,iblock)
4534           v1sij=v1s(1,j,itori,itori1,itori2,iblock)
4535           v2cij=v1c(2,j,itori,itori1,itori2,iblock)
4536           v2sij=v1s(2,j,itori,itori1,itori2,iblock)
4537           cosphi1=dcos(j*phii)
4538           sinphi1=dsin(j*phii)
4539           cosphi2=dcos(j*phii1)
4540           sinphi2=dsin(j*phii1)
4541           etors_d=etors_d+v1cij*cosphi1+v1sij*sinphi1+
4542      &     v2cij*cosphi2+v2sij*sinphi2
4543           gloci1=gloci1+j*(v1sij*cosphi1-v1cij*sinphi1)
4544           gloci2=gloci2+j*(v2sij*cosphi2-v2cij*sinphi2)
4545         enddo
4546         do k=2,ntermd_2(itori,itori1,itori2,iblock)
4547           do l=1,k-1
4548             v1cdij = v2c(k,l,itori,itori1,itori2,iblock)
4549             v2cdij = v2c(l,k,itori,itori1,itori2,iblock)
4550             v1sdij = v2s(k,l,itori,itori1,itori2,iblock)
4551             v2sdij = v2s(l,k,itori,itori1,itori2,iblock)
4552             cosphi1p2=dcos(l*phii+(k-l)*phii1)
4553             cosphi1m2=dcos(l*phii-(k-l)*phii1)
4554             sinphi1p2=dsin(l*phii+(k-l)*phii1)
4555             sinphi1m2=dsin(l*phii-(k-l)*phii1)
4556             etors_d=etors_d+v1cdij*cosphi1p2+v2cdij*cosphi1m2+
4557      &        v1sdij*sinphi1p2+v2sdij*sinphi1m2
4558             gloci1=gloci1+l*(v1sdij*cosphi1p2+v2sdij*cosphi1m2
4559      &        -v1cdij*sinphi1p2-v2cdij*sinphi1m2)
4560             gloci2=gloci2+(k-l)*(v1sdij*cosphi1p2-v2sdij*cosphi1m2
4561      &        -v1cdij*sinphi1p2+v2cdij*sinphi1m2)
4562           enddo
4563         enddo
4564         gloc(i-3,icg)=gloc(i-3,icg)+wtor_d*fact2*gloci1
4565         gloc(i-2,icg)=gloc(i-2,icg)+wtor_d*fact2*gloci2
4566  1215   continue
4567       enddo
4568       return
4569       end
4570 #endif
4571 c------------------------------------------------------------------------------
4572       subroutine eback_sc_corr(esccor)
4573 c 7/21/2007 Correlations between the backbone-local and side-chain-local
4574 c        conformational states; temporarily implemented as differences
4575 c        between UNRES torsional potentials (dependent on three types of
4576 c        residues) and the torsional potentials dependent on all 20 types
4577 c        of residues computed from AM1 energy surfaces of terminally-blocked
4578 c        amino-acid residues.
4579       implicit real*8 (a-h,o-z)
4580       include 'DIMENSIONS'
4581       include 'DIMENSIONS.ZSCOPT'
4582       include 'COMMON.VAR'
4583       include 'COMMON.GEO'
4584       include 'COMMON.LOCAL'
4585       include 'COMMON.TORSION'
4586       include 'COMMON.SCCOR'
4587       include 'COMMON.INTERACT'
4588       include 'COMMON.DERIV'
4589       include 'COMMON.CHAIN'
4590       include 'COMMON.NAMES'
4591       include 'COMMON.IOUNITS'
4592       include 'COMMON.FFIELD'
4593       include 'COMMON.CONTROL'
4594       logical lprn
4595 C Set lprn=.true. for debugging
4596       lprn=.false.
4597 c      lprn=.true.
4598 c      write (iout,*) "EBACK_SC_COR",iphi_start,iphi_end,nterm_sccor
4599       esccor=0.0D0
4600       do i=itau_start,itau_end
4601         if ((itype(i-2).eq.ntyp1).or.(itype(i-1).eq.ntyp1)) cycle
4602         esccor_ii=0.0D0
4603         isccori=isccortyp(itype(i-2))
4604         isccori1=isccortyp(itype(i-1))
4605         phii=phi(i)
4606         do intertyp=1,3 !intertyp
4607 cc Added 09 May 2012 (Adasko)
4608 cc  Intertyp means interaction type of backbone mainchain correlation: 
4609 c   1 = SC...Ca...Ca...Ca
4610 c   2 = Ca...Ca...Ca...SC
4611 c   3 = SC...Ca...Ca...SCi
4612         gloci=0.0D0
4613         if (((intertyp.eq.3).and.((itype(i-2).eq.10).or.
4614      &      (itype(i-1).eq.10).or.(itype(i-2).eq.ntyp1).or.
4615      &      (itype(i-1).eq.ntyp1)))
4616      &    .or. ((intertyp.eq.1).and.((itype(i-2).eq.10)
4617      &     .or.(itype(i-2).eq.ntyp1).or.(itype(i-1).eq.ntyp1)
4618      &     .or.(itype(i).eq.ntyp1)))
4619      &    .or.((intertyp.eq.2).and.((itype(i-1).eq.10).or.
4620      &      (itype(i-1).eq.ntyp1).or.(itype(i-2).eq.ntyp1).or.
4621      &      (itype(i-3).eq.ntyp1)))) cycle
4622         if ((intertyp.eq.2).and.(i.eq.4).and.(itype(1).eq.ntyp1)) cycle
4623         if ((intertyp.eq.1).and.(i.eq.nres).and.(itype(nres).eq.ntyp1))
4624      & cycle
4625        do j=1,nterm_sccor(isccori,isccori1)
4626           v1ij=v1sccor(j,intertyp,isccori,isccori1)
4627           v2ij=v2sccor(j,intertyp,isccori,isccori1)
4628           cosphi=dcos(j*tauangle(intertyp,i))
4629           sinphi=dsin(j*tauangle(intertyp,i))
4630            esccor=esccor+v1ij*cosphi+v2ij*sinphi
4631            gloci=gloci+j*(v2ij*cosphi-v1ij*sinphi)
4632          enddo
4633 c      write (iout,*) "EBACK_SC_COR",i,v1ij*cosphi+v2ij*sinphi,intertyp,
4634 c     & nterm_sccor(isccori,isccori1),isccori,isccori1
4635 c        gloc_sc(intertyp,i-3,icg)=gloc_sc(intertyp,i-3,icg)+wsccor*gloci
4636         if (lprn)
4637      &  write (iout,'(2(a3,2x,i3,2x),2i3,6f8.3/26x,6f8.3/)')
4638      &  restyp(itype(i-2)),i-2,restyp(itype(i-1)),i-1,itori,itori1,
4639      &  (v1sccor(j,1,itori,itori1),j=1,6)
4640      &  ,(v2sccor(j,1,itori,itori1),j=1,6)
4641 c        gsccor_loc(i-3)=gloci
4642        enddo !intertyp
4643       enddo
4644       return
4645       end
4646 c------------------------------------------------------------------------------
4647       subroutine multibody(ecorr)
4648 C This subroutine calculates multi-body contributions to energy following
4649 C the idea of Skolnick et al. If side chains I and J make a contact and
4650 C at the same time side chains I+1 and J+1 make a contact, an extra 
4651 C contribution equal to sqrt(eps(i,j)*eps(i+1,j+1)) is added.
4652       implicit real*8 (a-h,o-z)
4653       include 'DIMENSIONS'
4654       include 'COMMON.IOUNITS'
4655       include 'COMMON.DERIV'
4656       include 'COMMON.INTERACT'
4657       include 'COMMON.CONTACTS'
4658       double precision gx(3),gx1(3)
4659       logical lprn
4660
4661 C Set lprn=.true. for debugging
4662       lprn=.false.
4663
4664       if (lprn) then
4665         write (iout,'(a)') 'Contact function values:'
4666         do i=nnt,nct-2
4667           write (iout,'(i2,20(1x,i2,f10.5))') 
4668      &        i,(jcont(j,i),facont(j,i),j=1,num_cont(i))
4669         enddo
4670       endif
4671       ecorr=0.0D0
4672       do i=nnt,nct
4673         do j=1,3
4674           gradcorr(j,i)=0.0D0
4675           gradxorr(j,i)=0.0D0
4676         enddo
4677       enddo
4678       do i=nnt,nct-2
4679
4680         DO ISHIFT = 3,4
4681
4682         i1=i+ishift
4683         num_conti=num_cont(i)
4684         num_conti1=num_cont(i1)
4685         do jj=1,num_conti
4686           j=jcont(jj,i)
4687           do kk=1,num_conti1
4688             j1=jcont(kk,i1)
4689             if (j1.eq.j+ishift .or. j1.eq.j-ishift) then
4690 cd          write(iout,*)'i=',i,' j=',j,' i1=',i1,' j1=',j1,
4691 cd   &                   ' ishift=',ishift
4692 C Contacts I--J and I+ISHIFT--J+-ISHIFT1 occur simultaneously. 
4693 C The system gains extra energy.
4694               ecorr=ecorr+esccorr(i,j,i1,j1,jj,kk)
4695             endif   ! j1==j+-ishift
4696           enddo     ! kk  
4697         enddo       ! jj
4698
4699         ENDDO ! ISHIFT
4700
4701       enddo         ! i
4702       return
4703       end
4704 c------------------------------------------------------------------------------
4705       double precision function esccorr(i,j,k,l,jj,kk)
4706       implicit real*8 (a-h,o-z)
4707       include 'DIMENSIONS'
4708       include 'COMMON.IOUNITS'
4709       include 'COMMON.DERIV'
4710       include 'COMMON.INTERACT'
4711       include 'COMMON.CONTACTS'
4712       double precision gx(3),gx1(3)
4713       logical lprn
4714       lprn=.false.
4715       eij=facont(jj,i)
4716       ekl=facont(kk,k)
4717 cd    write (iout,'(4i5,3f10.5)') i,j,k,l,eij,ekl,-eij*ekl
4718 C Calculate the multi-body contribution to energy.
4719 C Calculate multi-body contributions to the gradient.
4720 cd    write (iout,'(2(2i3,3f10.5))')i,j,(gacont(m,jj,i),m=1,3),
4721 cd   & k,l,(gacont(m,kk,k),m=1,3)
4722       do m=1,3
4723         gx(m) =ekl*gacont(m,jj,i)
4724         gx1(m)=eij*gacont(m,kk,k)
4725         gradxorr(m,i)=gradxorr(m,i)-gx(m)
4726         gradxorr(m,j)=gradxorr(m,j)+gx(m)
4727         gradxorr(m,k)=gradxorr(m,k)-gx1(m)
4728         gradxorr(m,l)=gradxorr(m,l)+gx1(m)
4729       enddo
4730       do m=i,j-1
4731         do ll=1,3
4732           gradcorr(ll,m)=gradcorr(ll,m)+gx(ll)
4733         enddo
4734       enddo
4735       do m=k,l-1
4736         do ll=1,3
4737           gradcorr(ll,m)=gradcorr(ll,m)+gx1(ll)
4738         enddo
4739       enddo 
4740       esccorr=-eij*ekl
4741       return
4742       end
4743 c------------------------------------------------------------------------------
4744 #ifdef MPL
4745       subroutine pack_buffer(dimen1,dimen2,atom,indx,buffer)
4746       implicit real*8 (a-h,o-z)
4747       include 'DIMENSIONS' 
4748       integer dimen1,dimen2,atom,indx
4749       double precision buffer(dimen1,dimen2)
4750       double precision zapas 
4751       common /contacts_hb/ zapas(3,ntyp,maxres,7),
4752      &   facont_hb(ntyp,maxres),ees0p(ntyp,maxres),ees0m(ntyp,maxres),
4753      &         num_cont_hb(maxres),jcont_hb(ntyp,maxres)
4754       num_kont=num_cont_hb(atom)
4755       do i=1,num_kont
4756         do k=1,7
4757           do j=1,3
4758             buffer(i,indx+(k-1)*3+j)=zapas(j,i,atom,k)
4759           enddo ! j
4760         enddo ! k
4761         buffer(i,indx+22)=facont_hb(i,atom)
4762         buffer(i,indx+23)=ees0p(i,atom)
4763         buffer(i,indx+24)=ees0m(i,atom)
4764         buffer(i,indx+25)=dfloat(jcont_hb(i,atom))
4765       enddo ! i
4766       buffer(1,indx+26)=dfloat(num_kont)
4767       return
4768       end
4769 c------------------------------------------------------------------------------
4770       subroutine unpack_buffer(dimen1,dimen2,atom,indx,buffer)
4771       implicit real*8 (a-h,o-z)
4772       include 'DIMENSIONS' 
4773       integer dimen1,dimen2,atom,indx
4774       double precision buffer(dimen1,dimen2)
4775       double precision zapas 
4776       common /contacts_hb/ zapas(3,ntyp,maxres,7),
4777      &         facont_hb(ntyp,maxres),ees0p(ntyp,maxres),
4778      &         ees0m(ntyp,maxres),
4779      &         num_cont_hb(maxres),jcont_hb(ntyp,maxres)
4780       num_kont=buffer(1,indx+26)
4781       num_kont_old=num_cont_hb(atom)
4782       num_cont_hb(atom)=num_kont+num_kont_old
4783       do i=1,num_kont
4784         ii=i+num_kont_old
4785         do k=1,7    
4786           do j=1,3
4787             zapas(j,ii,atom,k)=buffer(i,indx+(k-1)*3+j)
4788           enddo ! j 
4789         enddo ! k 
4790         facont_hb(ii,atom)=buffer(i,indx+22)
4791         ees0p(ii,atom)=buffer(i,indx+23)
4792         ees0m(ii,atom)=buffer(i,indx+24)
4793         jcont_hb(ii,atom)=buffer(i,indx+25)
4794       enddo ! i
4795       return
4796       end
4797 c------------------------------------------------------------------------------
4798 #endif
4799       subroutine multibody_hb(ecorr,ecorr5,ecorr6,n_corr,n_corr1)
4800 C This subroutine calculates multi-body contributions to hydrogen-bonding 
4801       implicit real*8 (a-h,o-z)
4802       include 'DIMENSIONS'
4803       include 'DIMENSIONS.ZSCOPT'
4804       include 'COMMON.IOUNITS'
4805 #ifdef MPL
4806       include 'COMMON.INFO'
4807 #endif
4808       include 'COMMON.FFIELD'
4809       include 'COMMON.DERIV'
4810       include 'COMMON.INTERACT'
4811       include 'COMMON.CONTACTS'
4812 #ifdef MPL
4813       parameter (max_cont=maxconts)
4814       parameter (max_dim=2*(8*3+2))
4815       parameter (msglen1=max_cont*max_dim*4)
4816       parameter (msglen2=2*msglen1)
4817       integer source,CorrelType,CorrelID,Error
4818       double precision buffer(max_cont,max_dim)
4819 #endif
4820       double precision gx(3),gx1(3)
4821       logical lprn,ldone
4822
4823 C Set lprn=.true. for debugging
4824       lprn=.false.
4825 #ifdef MPL
4826       n_corr=0
4827       n_corr1=0
4828       if (fgProcs.le.1) goto 30
4829       if (lprn) then
4830         write (iout,'(a)') 'Contact function values:'
4831         do i=nnt,nct-2
4832           write (iout,'(2i3,50(1x,i2,f5.2))') 
4833      &    i,num_cont_hb(i),(jcont_hb(j,i),facont_hb(j,i),
4834      &    j=1,num_cont_hb(i))
4835         enddo
4836       endif
4837 C Caution! Following code assumes that electrostatic interactions concerning
4838 C a given atom are split among at most two processors!
4839       CorrelType=477
4840       CorrelID=MyID+1
4841       ldone=.false.
4842       do i=1,max_cont
4843         do j=1,max_dim
4844           buffer(i,j)=0.0D0
4845         enddo
4846       enddo
4847       mm=mod(MyRank,2)
4848 cd    write (iout,*) 'MyRank',MyRank,' mm',mm
4849       if (mm) 20,20,10 
4850    10 continue
4851 cd    write (iout,*) 'Sending: MyRank',MyRank,' mm',mm,' ldone',ldone
4852       if (MyRank.gt.0) then
4853 C Send correlation contributions to the preceding processor
4854         msglen=msglen1
4855         nn=num_cont_hb(iatel_s)
4856         call pack_buffer(max_cont,max_dim,iatel_s,0,buffer)
4857 cd      write (iout,*) 'The BUFFER array:'
4858 cd      do i=1,nn
4859 cd        write (iout,'(i2,9(3f8.3,2x))') i,(buffer(i,j),j=1,26)
4860 cd      enddo
4861         if (ielstart(iatel_s).gt.iatel_s+ispp) then
4862           msglen=msglen2
4863             call pack_buffer(max_cont,max_dim,iatel_s+1,26,buffer)
4864 C Clear the contacts of the atom passed to the neighboring processor
4865         nn=num_cont_hb(iatel_s+1)
4866 cd      do i=1,nn
4867 cd        write (iout,'(i2,9(3f8.3,2x))') i,(buffer(i,j+26),j=1,26)
4868 cd      enddo
4869             num_cont_hb(iatel_s)=0
4870         endif 
4871 cd      write (iout,*) 'Processor ',MyID,MyRank,
4872 cd   & ' is sending correlation contribution to processor',MyID-1,
4873 cd   & ' msglen=',msglen
4874 cd      write (*,*) 'Processor ',MyID,MyRank,
4875 cd   & ' is sending correlation contribution to processor',MyID-1,
4876 cd   & ' msglen=',msglen,' CorrelType=',CorrelType
4877         call mp_bsend(buffer,msglen,MyID-1,CorrelType,CorrelID)
4878 cd      write (iout,*) 'Processor ',MyID,
4879 cd   & ' has sent correlation contribution to processor',MyID-1,
4880 cd   & ' msglen=',msglen,' CorrelID=',CorrelID
4881 cd      write (*,*) 'Processor ',MyID,
4882 cd   & ' has sent correlation contribution to processor',MyID-1,
4883 cd   & ' msglen=',msglen,' CorrelID=',CorrelID
4884         msglen=msglen1
4885       endif ! (MyRank.gt.0)
4886       if (ldone) goto 30
4887       ldone=.true.
4888    20 continue
4889 cd    write (iout,*) 'Receiving: MyRank',MyRank,' mm',mm,' ldone',ldone
4890       if (MyRank.lt.fgProcs-1) then
4891 C Receive correlation contributions from the next processor
4892         msglen=msglen1
4893         if (ielend(iatel_e).lt.nct-1) msglen=msglen2
4894 cd      write (iout,*) 'Processor',MyID,
4895 cd   & ' is receiving correlation contribution from processor',MyID+1,
4896 cd   & ' msglen=',msglen,' CorrelType=',CorrelType
4897 cd      write (*,*) 'Processor',MyID,
4898 cd   & ' is receiving correlation contribution from processor',MyID+1,
4899 cd   & ' msglen=',msglen,' CorrelType=',CorrelType
4900         nbytes=-1
4901         do while (nbytes.le.0)
4902           call mp_probe(MyID+1,CorrelType,nbytes)
4903         enddo
4904 cd      print *,'Processor',MyID,' msglen',msglen,' nbytes',nbytes
4905         call mp_brecv(buffer,msglen,MyID+1,CorrelType,nbytes)
4906 cd      write (iout,*) 'Processor',MyID,
4907 cd   & ' has received correlation contribution from processor',MyID+1,
4908 cd   & ' msglen=',msglen,' nbytes=',nbytes
4909 cd      write (iout,*) 'The received BUFFER array:'
4910 cd      do i=1,max_cont
4911 cd        write (iout,'(i2,9(3f8.3,2x))') i,(buffer(i,j),j=1,52)
4912 cd      enddo
4913         if (msglen.eq.msglen1) then
4914           call unpack_buffer(max_cont,max_dim,iatel_e+1,0,buffer)
4915         else if (msglen.eq.msglen2)  then
4916           call unpack_buffer(max_cont,max_dim,iatel_e,0,buffer) 
4917           call unpack_buffer(max_cont,max_dim,iatel_e+1,26,buffer) 
4918         else
4919           write (iout,*) 
4920      & 'ERROR!!!! message length changed while processing correlations.'
4921           write (*,*) 
4922      & 'ERROR!!!! message length changed while processing correlations.'
4923           call mp_stopall(Error)
4924         endif ! msglen.eq.msglen1
4925       endif ! MyRank.lt.fgProcs-1
4926       if (ldone) goto 30
4927       ldone=.true.
4928       goto 10
4929    30 continue
4930 #endif
4931       if (lprn) then
4932         write (iout,'(a)') 'Contact function values:'
4933         do i=nnt,nct-2
4934           write (iout,'(2i3,50(1x,i2,f5.2))') 
4935      &    i,num_cont_hb(i),(jcont_hb(j,i),facont_hb(j,i),
4936      &    j=1,num_cont_hb(i))
4937         enddo
4938       endif
4939       ecorr=0.0D0
4940 C Remove the loop below after debugging !!!
4941       do i=nnt,nct
4942         do j=1,3
4943           gradcorr(j,i)=0.0D0
4944           gradxorr(j,i)=0.0D0
4945         enddo
4946       enddo
4947 C Calculate the local-electrostatic correlation terms
4948       do i=iatel_s,iatel_e+1
4949         i1=i+1
4950         num_conti=num_cont_hb(i)
4951         num_conti1=num_cont_hb(i+1)
4952         do jj=1,num_conti
4953           j=jcont_hb(jj,i)
4954           do kk=1,num_conti1
4955             j1=jcont_hb(kk,i1)
4956 c            write (iout,*) 'i=',i,' j=',j,' i1=',i1,' j1=',j1,
4957 c     &         ' jj=',jj,' kk=',kk
4958             if (j1.eq.j+1 .or. j1.eq.j-1) then
4959 C Contacts I-J and (I+1)-(J+1) or (I+1)-(J-1) occur simultaneously. 
4960 C The system gains extra energy.
4961               ecorr=ecorr+ehbcorr(i,j,i+1,j1,jj,kk,0.72D0,0.32D0)
4962               n_corr=n_corr+1
4963             else if (j1.eq.j) then
4964 C Contacts I-J and I-(J+1) occur simultaneously. 
4965 C The system loses extra energy.
4966 c             ecorr=ecorr+ehbcorr(i,j,i+1,j,jj,kk,0.60D0,-0.40D0) 
4967             endif
4968           enddo ! kk
4969           do kk=1,num_conti
4970             j1=jcont_hb(kk,i)
4971 c           write (iout,*) 'i=',i,' j=',j,' i1=',i1,' j1=',j1,
4972 c    &         ' jj=',jj,' kk=',kk
4973             if (j1.eq.j+1) then
4974 C Contacts I-J and (I+1)-J occur simultaneously. 
4975 C The system loses extra energy.
4976 c             ecorr=ecorr+ehbcorr(i,j,i,j+1,jj,kk,0.60D0,-0.40D0)
4977             endif ! j1==j+1
4978           enddo ! kk
4979         enddo ! jj
4980       enddo ! i
4981       return
4982       end
4983 c------------------------------------------------------------------------------
4984       subroutine multibody_eello(ecorr,ecorr5,ecorr6,eturn6,n_corr,
4985      &  n_corr1)
4986 C This subroutine calculates multi-body contributions to hydrogen-bonding 
4987       implicit real*8 (a-h,o-z)
4988       include 'DIMENSIONS'
4989       include 'DIMENSIONS.ZSCOPT'
4990       include 'COMMON.IOUNITS'
4991 #ifdef MPL
4992       include 'COMMON.INFO'
4993 #endif
4994       include 'COMMON.FFIELD'
4995       include 'COMMON.DERIV'
4996       include 'COMMON.INTERACT'
4997       include 'COMMON.CONTACTS'
4998 #ifdef MPL
4999       parameter (max_cont=maxconts)
5000       parameter (max_dim=2*(8*3+2))
5001       parameter (msglen1=max_cont*max_dim*4)
5002       parameter (msglen2=2*msglen1)
5003       integer source,CorrelType,CorrelID,Error
5004       double precision buffer(max_cont,max_dim)
5005 #endif
5006       double precision gx(3),gx1(3)
5007       logical lprn,ldone
5008
5009 C Set lprn=.true. for debugging
5010       lprn=.false.
5011       eturn6=0.0d0
5012 #ifdef MPL
5013       n_corr=0
5014       n_corr1=0
5015       if (fgProcs.le.1) goto 30
5016       if (lprn) then
5017         write (iout,'(a)') 'Contact function values:'
5018         do i=nnt,nct-2
5019           write (iout,'(2i3,50(1x,i2,f5.2))') 
5020      &    i,num_cont_hb(i),(jcont_hb(j,i),facont_hb(j,i),
5021      &    j=1,num_cont_hb(i))
5022         enddo
5023       endif
5024 C Caution! Following code assumes that electrostatic interactions concerning
5025 C a given atom are split among at most two processors!
5026       CorrelType=477
5027       CorrelID=MyID+1
5028       ldone=.false.
5029       do i=1,max_cont
5030         do j=1,max_dim
5031           buffer(i,j)=0.0D0
5032         enddo
5033       enddo
5034       mm=mod(MyRank,2)
5035 cd    write (iout,*) 'MyRank',MyRank,' mm',mm
5036       if (mm) 20,20,10 
5037    10 continue
5038 cd    write (iout,*) 'Sending: MyRank',MyRank,' mm',mm,' ldone',ldone
5039       if (MyRank.gt.0) then
5040 C Send correlation contributions to the preceding processor
5041         msglen=msglen1
5042         nn=num_cont_hb(iatel_s)
5043         call pack_buffer(max_cont,max_dim,iatel_s,0,buffer)
5044 cd      write (iout,*) 'The BUFFER array:'
5045 cd      do i=1,nn
5046 cd        write (iout,'(i2,9(3f8.3,2x))') i,(buffer(i,j),j=1,26)
5047 cd      enddo
5048         if (ielstart(iatel_s).gt.iatel_s+ispp) then
5049           msglen=msglen2
5050             call pack_buffer(max_cont,max_dim,iatel_s+1,26,buffer)
5051 C Clear the contacts of the atom passed to the neighboring processor
5052         nn=num_cont_hb(iatel_s+1)
5053 cd      do i=1,nn
5054 cd        write (iout,'(i2,9(3f8.3,2x))') i,(buffer(i,j+26),j=1,26)
5055 cd      enddo
5056             num_cont_hb(iatel_s)=0
5057         endif 
5058 cd      write (iout,*) 'Processor ',MyID,MyRank,
5059 cd   & ' is sending correlation contribution to processor',MyID-1,
5060 cd   & ' msglen=',msglen
5061 cd      write (*,*) 'Processor ',MyID,MyRank,
5062 cd   & ' is sending correlation contribution to processor',MyID-1,
5063 cd   & ' msglen=',msglen,' CorrelType=',CorrelType
5064         call mp_bsend(buffer,msglen,MyID-1,CorrelType,CorrelID)
5065 cd      write (iout,*) 'Processor ',MyID,
5066 cd   & ' has sent correlation contribution to processor',MyID-1,
5067 cd   & ' msglen=',msglen,' CorrelID=',CorrelID
5068 cd      write (*,*) 'Processor ',MyID,
5069 cd   & ' has sent correlation contribution to processor',MyID-1,
5070 cd   & ' msglen=',msglen,' CorrelID=',CorrelID
5071         msglen=msglen1
5072       endif ! (MyRank.gt.0)
5073       if (ldone) goto 30
5074       ldone=.true.
5075    20 continue
5076 cd    write (iout,*) 'Receiving: MyRank',MyRank,' mm',mm,' ldone',ldone
5077       if (MyRank.lt.fgProcs-1) then
5078 C Receive correlation contributions from the next processor
5079         msglen=msglen1
5080         if (ielend(iatel_e).lt.nct-1) msglen=msglen2
5081 cd      write (iout,*) 'Processor',MyID,
5082 cd   & ' is receiving correlation contribution from processor',MyID+1,
5083 cd   & ' msglen=',msglen,' CorrelType=',CorrelType
5084 cd      write (*,*) 'Processor',MyID,
5085 cd   & ' is receiving correlation contribution from processor',MyID+1,
5086 cd   & ' msglen=',msglen,' CorrelType=',CorrelType
5087         nbytes=-1
5088         do while (nbytes.le.0)
5089           call mp_probe(MyID+1,CorrelType,nbytes)
5090         enddo
5091 cd      print *,'Processor',MyID,' msglen',msglen,' nbytes',nbytes
5092         call mp_brecv(buffer,msglen,MyID+1,CorrelType,nbytes)
5093 cd      write (iout,*) 'Processor',MyID,
5094 cd   & ' has received correlation contribution from processor',MyID+1,
5095 cd   & ' msglen=',msglen,' nbytes=',nbytes
5096 cd      write (iout,*) 'The received BUFFER array:'
5097 cd      do i=1,max_cont
5098 cd        write (iout,'(i2,9(3f8.3,2x))') i,(buffer(i,j),j=1,52)
5099 cd      enddo
5100         if (msglen.eq.msglen1) then
5101           call unpack_buffer(max_cont,max_dim,iatel_e+1,0,buffer)
5102         else if (msglen.eq.msglen2)  then
5103           call unpack_buffer(max_cont,max_dim,iatel_e,0,buffer) 
5104           call unpack_buffer(max_cont,max_dim,iatel_e+1,26,buffer) 
5105         else
5106           write (iout,*) 
5107      & 'ERROR!!!! message length changed while processing correlations.'
5108           write (*,*) 
5109      & 'ERROR!!!! message length changed while processing correlations.'
5110           call mp_stopall(Error)
5111         endif ! msglen.eq.msglen1
5112       endif ! MyRank.lt.fgProcs-1
5113       if (ldone) goto 30
5114       ldone=.true.
5115       goto 10
5116    30 continue
5117 #endif
5118       if (lprn) then
5119         write (iout,'(a)') 'Contact function values:'
5120         do i=nnt,nct-2
5121           write (iout,'(2i3,50(1x,i2,f5.2))') 
5122      &    i,num_cont_hb(i),(jcont_hb(j,i),facont_hb(j,i),
5123      &    j=1,num_cont_hb(i))
5124         enddo
5125       endif
5126       ecorr=0.0D0
5127       ecorr5=0.0d0
5128       ecorr6=0.0d0
5129 C Remove the loop below after debugging !!!
5130       do i=nnt,nct
5131         do j=1,3
5132           gradcorr(j,i)=0.0D0
5133           gradxorr(j,i)=0.0D0
5134         enddo
5135       enddo
5136 C Calculate the dipole-dipole interaction energies
5137       if (wcorr6.gt.0.0d0 .or. wturn6.gt.0.0d0) then
5138       do i=iatel_s,iatel_e+1
5139         num_conti=num_cont_hb(i)
5140         do jj=1,num_conti
5141           j=jcont_hb(jj,i)
5142           call dipole(i,j,jj)
5143         enddo
5144       enddo
5145       endif
5146 C Calculate the local-electrostatic correlation terms
5147       do i=iatel_s,iatel_e+1
5148         i1=i+1
5149         num_conti=num_cont_hb(i)
5150         num_conti1=num_cont_hb(i+1)
5151         do jj=1,num_conti
5152           j=jcont_hb(jj,i)
5153           do kk=1,num_conti1
5154             j1=jcont_hb(kk,i1)
5155 c            write (*,*) 'i=',i,' j=',j,' i1=',i1,' j1=',j1,
5156 c     &         ' jj=',jj,' kk=',kk
5157             if (j1.eq.j+1 .or. j1.eq.j-1) then
5158 C Contacts I-J and (I+1)-(J+1) or (I+1)-(J-1) occur simultaneously. 
5159 C The system gains extra energy.
5160               n_corr=n_corr+1
5161               sqd1=dsqrt(d_cont(jj,i))
5162               sqd2=dsqrt(d_cont(kk,i1))
5163               sred_geom = sqd1*sqd2
5164               IF (sred_geom.lt.cutoff_corr) THEN
5165                 call gcont(sred_geom,r0_corr,1.0D0,delt_corr,
5166      &            ekont,fprimcont)
5167 c               write (*,*) 'i=',i,' j=',j,' i1=',i1,' j1=',j1,
5168 c     &         ' jj=',jj,' kk=',kk
5169                 fac_prim1=0.5d0*sqd2/sqd1*fprimcont
5170                 fac_prim2=0.5d0*sqd1/sqd2*fprimcont
5171                 do l=1,3
5172                   g_contij(l,1)=fac_prim1*grij_hb_cont(l,jj,i)
5173                   g_contij(l,2)=fac_prim2*grij_hb_cont(l,kk,i1)
5174                 enddo
5175                 n_corr1=n_corr1+1
5176 cd               write (iout,*) 'sred_geom=',sred_geom,
5177 cd     &          ' ekont=',ekont,' fprim=',fprimcont
5178                 call calc_eello(i,j,i+1,j1,jj,kk)
5179                 if (wcorr4.gt.0.0d0) 
5180      &            ecorr=ecorr+eello4(i,j,i+1,j1,jj,kk)
5181                 if (wcorr5.gt.0.0d0)
5182      &            ecorr5=ecorr5+eello5(i,j,i+1,j1,jj,kk)
5183 c                print *,"wcorr5",ecorr5
5184 cd                write(2,*)'wcorr6',wcorr6,' wturn6',wturn6
5185 cd                write(2,*)'ijkl',i,j,i+1,j1 
5186                 if (wcorr6.gt.0.0d0 .and. (j.ne.i+4 .or. j1.ne.i+3
5187      &               .or. wturn6.eq.0.0d0))then
5188 cd                  write (iout,*) '******ecorr6: i,j,i+1,j1',i,j,i+1,j1
5189                   ecorr6=ecorr6+eello6(i,j,i+1,j1,jj,kk)
5190 cd                write (iout,*) 'ecorr',ecorr,' ecorr5=',ecorr5,
5191 cd     &            'ecorr6=',ecorr6
5192 cd                write (iout,'(4e15.5)') sred_geom,
5193 cd     &          dabs(eello4(i,j,i+1,j1,jj,kk)),
5194 cd     &          dabs(eello5(i,j,i+1,j1,jj,kk)),
5195 cd     &          dabs(eello6(i,j,i+1,j1,jj,kk))
5196                 else if (wturn6.gt.0.0d0
5197      &            .and. (j.eq.i+4 .and. j1.eq.i+3)) then
5198 cd                  write (iout,*) '******eturn6: i,j,i+1,j1',i,j,i+1,j1
5199                   eturn6=eturn6+eello_turn6(i,jj,kk)
5200 cd                  write (2,*) 'multibody_eello:eturn6',eturn6
5201                 endif
5202               ENDIF
5203 1111          continue
5204             else if (j1.eq.j) then
5205 C Contacts I-J and I-(J+1) occur simultaneously. 
5206 C The system loses extra energy.
5207 c             ecorr=ecorr+ehbcorr(i,j,i+1,j,jj,kk,0.60D0,-0.40D0) 
5208             endif
5209           enddo ! kk
5210           do kk=1,num_conti
5211             j1=jcont_hb(kk,i)
5212 c           write (iout,*) 'i=',i,' j=',j,' i1=',i1,' j1=',j1,
5213 c    &         ' jj=',jj,' kk=',kk
5214             if (j1.eq.j+1) then
5215 C Contacts I-J and (I+1)-J occur simultaneously. 
5216 C The system loses extra energy.
5217 c             ecorr=ecorr+ehbcorr(i,j,i,j+1,jj,kk,0.60D0,-0.40D0)
5218             endif ! j1==j+1
5219           enddo ! kk
5220         enddo ! jj
5221       enddo ! i
5222       return
5223       end
5224 c------------------------------------------------------------------------------
5225       double precision function ehbcorr(i,j,k,l,jj,kk,coeffp,coeffm)
5226       implicit real*8 (a-h,o-z)
5227       include 'DIMENSIONS'
5228       include 'COMMON.IOUNITS'
5229       include 'COMMON.DERIV'
5230       include 'COMMON.INTERACT'
5231       include 'COMMON.CONTACTS'
5232       double precision gx(3),gx1(3)
5233       logical lprn
5234       lprn=.false.
5235       eij=facont_hb(jj,i)
5236       ekl=facont_hb(kk,k)
5237       ees0pij=ees0p(jj,i)
5238       ees0pkl=ees0p(kk,k)
5239       ees0mij=ees0m(jj,i)
5240       ees0mkl=ees0m(kk,k)
5241       ekont=eij*ekl
5242       ees=-(coeffp*ees0pij*ees0pkl+coeffm*ees0mij*ees0mkl)
5243 cd    ees=-(coeffp*ees0pkl+coeffm*ees0mkl)
5244 C Following 4 lines for diagnostics.
5245 cd    ees0pkl=0.0D0
5246 cd    ees0pij=1.0D0
5247 cd    ees0mkl=0.0D0
5248 cd    ees0mij=1.0D0
5249 c     write (iout,*)'Contacts have occurred for peptide groups',i,j,
5250 c    &   ' and',k,l
5251 c     write (iout,*)'Contacts have occurred for peptide groups',
5252 c    &  i,j,' fcont:',eij,' eij',' eesij',ees0pij,ees0mij,' and ',k,l
5253 c    & ,' fcont ',ekl,' eeskl',ees0pkl,ees0mkl,' ees=',ees
5254 C Calculate the multi-body contribution to energy.
5255       ecorr=ecorr+ekont*ees
5256       if (calc_grad) then
5257 C Calculate multi-body contributions to the gradient.
5258       do ll=1,3
5259         ghalf=0.5D0*ees*ekl*gacont_hbr(ll,jj,i)
5260         gradcorr(ll,i)=gradcorr(ll,i)+ghalf
5261      &  -ekont*(coeffp*ees0pkl*gacontp_hb1(ll,jj,i)+
5262      &  coeffm*ees0mkl*gacontm_hb1(ll,jj,i))
5263         gradcorr(ll,j)=gradcorr(ll,j)+ghalf
5264      &  -ekont*(coeffp*ees0pkl*gacontp_hb2(ll,jj,i)+
5265      &  coeffm*ees0mkl*gacontm_hb2(ll,jj,i))
5266         ghalf=0.5D0*ees*eij*gacont_hbr(ll,kk,k)
5267         gradcorr(ll,k)=gradcorr(ll,k)+ghalf
5268      &  -ekont*(coeffp*ees0pij*gacontp_hb1(ll,kk,k)+
5269      &  coeffm*ees0mij*gacontm_hb1(ll,kk,k))
5270         gradcorr(ll,l)=gradcorr(ll,l)+ghalf
5271      &  -ekont*(coeffp*ees0pij*gacontp_hb2(ll,kk,k)+
5272      &  coeffm*ees0mij*gacontm_hb2(ll,kk,k))
5273       enddo
5274       do m=i+1,j-1
5275         do ll=1,3
5276           gradcorr(ll,m)=gradcorr(ll,m)+
5277      &     ees*ekl*gacont_hbr(ll,jj,i)-
5278      &     ekont*(coeffp*ees0pkl*gacontp_hb3(ll,jj,i)+
5279      &     coeffm*ees0mkl*gacontm_hb3(ll,jj,i))
5280         enddo
5281       enddo
5282       do m=k+1,l-1
5283         do ll=1,3
5284           gradcorr(ll,m)=gradcorr(ll,m)+
5285      &     ees*eij*gacont_hbr(ll,kk,k)-
5286      &     ekont*(coeffp*ees0pij*gacontp_hb3(ll,kk,k)+
5287      &     coeffm*ees0mij*gacontm_hb3(ll,kk,k))
5288         enddo
5289       enddo 
5290       endif
5291       ehbcorr=ekont*ees
5292       return
5293       end
5294 C---------------------------------------------------------------------------
5295       subroutine dipole(i,j,jj)
5296       implicit real*8 (a-h,o-z)
5297       include 'DIMENSIONS'
5298       include 'DIMENSIONS.ZSCOPT'
5299       include 'COMMON.IOUNITS'
5300       include 'COMMON.CHAIN'
5301       include 'COMMON.FFIELD'
5302       include 'COMMON.DERIV'
5303       include 'COMMON.INTERACT'
5304       include 'COMMON.CONTACTS'
5305       include 'COMMON.TORSION'
5306       include 'COMMON.VAR'
5307       include 'COMMON.GEO'
5308       dimension dipi(2,2),dipj(2,2),dipderi(2),dipderj(2),auxvec(2),
5309      &  auxmat(2,2)
5310       iti1 = itortyp(itype(i+1))
5311       if (j.lt.nres-1) then
5312         if (itype(j).le.ntyp) then
5313           itj1 = itortyp(itype(j+1))
5314         else
5315           itj=ntortyp+1 
5316         endif
5317       else
5318         itj1=ntortyp+1
5319       endif
5320       do iii=1,2
5321         dipi(iii,1)=Ub2(iii,i)
5322         dipderi(iii)=Ub2der(iii,i)
5323         dipi(iii,2)=b1(iii,iti1)
5324         dipj(iii,1)=Ub2(iii,j)
5325         dipderj(iii)=Ub2der(iii,j)
5326         dipj(iii,2)=b1(iii,itj1)
5327       enddo
5328       kkk=0
5329       do iii=1,2
5330         call matvec2(a_chuj(1,1,jj,i),dipj(1,iii),auxvec(1)) 
5331         do jjj=1,2
5332           kkk=kkk+1
5333           dip(kkk,jj,i)=scalar2(dipi(1,jjj),auxvec(1))
5334         enddo
5335       enddo
5336       if (.not.calc_grad) return
5337       do kkk=1,5
5338         do lll=1,3
5339           mmm=0
5340           do iii=1,2
5341             call matvec2(a_chuj_der(1,1,lll,kkk,jj,i),dipj(1,iii),
5342      &        auxvec(1))
5343             do jjj=1,2
5344               mmm=mmm+1
5345               dipderx(lll,kkk,mmm,jj,i)=scalar2(dipi(1,jjj),auxvec(1))
5346             enddo
5347           enddo
5348         enddo
5349       enddo
5350       call transpose2(a_chuj(1,1,jj,i),auxmat(1,1))
5351       call matvec2(auxmat(1,1),dipderi(1),auxvec(1))
5352       do iii=1,2
5353         dipderg(iii,jj,i)=scalar2(auxvec(1),dipj(1,iii))
5354       enddo
5355       call matvec2(a_chuj(1,1,jj,i),dipderj(1),auxvec(1))
5356       do iii=1,2
5357         dipderg(iii+2,jj,i)=scalar2(auxvec(1),dipi(1,iii))
5358       enddo
5359       return
5360       end
5361 C---------------------------------------------------------------------------
5362       subroutine calc_eello(i,j,k,l,jj,kk)
5363
5364 C This subroutine computes matrices and vectors needed to calculate 
5365 C the fourth-, fifth-, and sixth-order local-electrostatic terms.
5366 C
5367       implicit real*8 (a-h,o-z)
5368       include 'DIMENSIONS'
5369       include 'DIMENSIONS.ZSCOPT'
5370       include 'COMMON.IOUNITS'
5371       include 'COMMON.CHAIN'
5372       include 'COMMON.DERIV'
5373       include 'COMMON.INTERACT'
5374       include 'COMMON.CONTACTS'
5375       include 'COMMON.TORSION'
5376       include 'COMMON.VAR'
5377       include 'COMMON.GEO'
5378       include 'COMMON.FFIELD'
5379       double precision aa1(2,2),aa2(2,2),aa1t(2,2),aa2t(2,2),
5380      &  aa1tder(2,2,3,5),aa2tder(2,2,3,5),auxmat(2,2)
5381       logical lprn
5382       common /kutas/ lprn
5383 cd      write (iout,*) 'calc_eello: i=',i,' j=',j,' k=',k,' l=',l,
5384 cd     & ' jj=',jj,' kk=',kk
5385 cd      if (i.ne.2 .or. j.ne.4 .or. k.ne.3 .or. l.ne.5) return
5386       do iii=1,2
5387         do jjj=1,2
5388           aa1(iii,jjj)=a_chuj(iii,jjj,jj,i)
5389           aa2(iii,jjj)=a_chuj(iii,jjj,kk,k)
5390         enddo
5391       enddo
5392       call transpose2(aa1(1,1),aa1t(1,1))
5393       call transpose2(aa2(1,1),aa2t(1,1))
5394       do kkk=1,5
5395         do lll=1,3
5396           call transpose2(a_chuj_der(1,1,lll,kkk,jj,i),
5397      &      aa1tder(1,1,lll,kkk))
5398           call transpose2(a_chuj_der(1,1,lll,kkk,kk,k),
5399      &      aa2tder(1,1,lll,kkk))
5400         enddo
5401       enddo 
5402       if (l.eq.j+1) then
5403 C parallel orientation of the two CA-CA-CA frames.
5404         if (i.gt.1 .and. itype(i).le.ntyp) then
5405           iti=itortyp(itype(i))
5406         else
5407           iti=ntortyp+1
5408         endif
5409         itk1=itortyp(itype(k+1))
5410         itj=itortyp(itype(j))
5411         if (l.lt.nres-1 .and. itype(l+1).le.ntyp) then
5412           itl1=itortyp(itype(l+1))
5413         else
5414           itl1=ntortyp+1
5415         endif
5416 C A1 kernel(j+1) A2T
5417 cd        do iii=1,2
5418 cd          write (iout,'(3f10.5,5x,3f10.5)') 
5419 cd     &     (EUg(iii,jjj,k),jjj=1,2),(EUg(iii,jjj,l),jjj=1,2)
5420 cd        enddo
5421         call kernel(aa1(1,1),aa2t(1,1),a_chuj_der(1,1,1,1,jj,i),
5422      &   aa2tder(1,1,1,1),1,.false.,EUg(1,1,l),EUgder(1,1,l),
5423      &   AEA(1,1,1),AEAderg(1,1,1),AEAderx(1,1,1,1,1,1))
5424 C Following matrices are needed only for 6-th order cumulants
5425         IF (wcorr6.gt.0.0d0) THEN
5426         call kernel(aa1(1,1),aa2t(1,1),a_chuj_der(1,1,1,1,jj,i),
5427      &   aa2tder(1,1,1,1),1,.false.,EUgC(1,1,l),EUgCder(1,1,l),
5428      &   AECA(1,1,1),AECAderg(1,1,1),AECAderx(1,1,1,1,1,1))
5429         call kernel(aa1(1,1),aa2t(1,1),a_chuj_der(1,1,1,1,jj,i),
5430      &   aa2tder(1,1,1,1),2,.false.,Ug2DtEUg(1,1,l),
5431      &   Ug2DtEUgder(1,1,1,l),ADtEA(1,1,1),ADtEAderg(1,1,1,1),
5432      &   ADtEAderx(1,1,1,1,1,1))
5433         lprn=.false.
5434         call kernel(aa1(1,1),aa2t(1,1),a_chuj_der(1,1,1,1,jj,i),
5435      &   aa2tder(1,1,1,1),2,.false.,DtUg2EUg(1,1,l),
5436      &   DtUg2EUgder(1,1,1,l),ADtEA1(1,1,1),ADtEA1derg(1,1,1,1),
5437      &   ADtEA1derx(1,1,1,1,1,1))
5438         ENDIF
5439 C End 6-th order cumulants
5440 cd        lprn=.false.
5441 cd        if (lprn) then
5442 cd        write (2,*) 'In calc_eello6'
5443 cd        do iii=1,2
5444 cd          write (2,*) 'iii=',iii
5445 cd          do kkk=1,5
5446 cd            write (2,*) 'kkk=',kkk
5447 cd            do jjj=1,2
5448 cd              write (2,'(3(2f10.5),5x)') 
5449 cd     &        ((ADtEA1derx(jjj,mmm,lll,kkk,iii,1),mmm=1,2),lll=1,3)
5450 cd            enddo
5451 cd          enddo
5452 cd        enddo
5453 cd        endif
5454         call transpose2(EUgder(1,1,k),auxmat(1,1))
5455         call matmat2(auxmat(1,1),AEA(1,1,1),EAEAderg(1,1,1,1))
5456         call transpose2(EUg(1,1,k),auxmat(1,1))
5457         call matmat2(auxmat(1,1),AEA(1,1,1),EAEA(1,1,1))
5458         call matmat2(auxmat(1,1),AEAderg(1,1,1),EAEAderg(1,1,2,1))
5459         do iii=1,2
5460           do kkk=1,5
5461             do lll=1,3
5462               call matmat2(auxmat(1,1),AEAderx(1,1,lll,kkk,iii,1),
5463      &          EAEAderx(1,1,lll,kkk,iii,1))
5464             enddo
5465           enddo
5466         enddo
5467 C A1T kernel(i+1) A2
5468         call kernel(aa1t(1,1),aa2(1,1),aa1tder(1,1,1,1),
5469      &   a_chuj_der(1,1,1,1,kk,k),1,.false.,EUg(1,1,k),EUgder(1,1,k),
5470      &   AEA(1,1,2),AEAderg(1,1,2),AEAderx(1,1,1,1,1,2))
5471 C Following matrices are needed only for 6-th order cumulants
5472         IF (wcorr6.gt.0.0d0) THEN
5473         call kernel(aa1t(1,1),aa2(1,1),aa1tder(1,1,1,1),
5474      &   a_chuj_der(1,1,1,1,kk,k),1,.false.,EUgC(1,1,k),EUgCder(1,1,k),
5475      &   AECA(1,1,2),AECAderg(1,1,2),AECAderx(1,1,1,1,1,2))
5476         call kernel(aa1t(1,1),aa2(1,1),aa1tder(1,1,1,1),
5477      &   a_chuj_der(1,1,1,1,kk,k),2,.false.,Ug2DtEUg(1,1,k),
5478      &   Ug2DtEUgder(1,1,1,k),ADtEA(1,1,2),ADtEAderg(1,1,1,2),
5479      &   ADtEAderx(1,1,1,1,1,2))
5480         call kernel(aa1t(1,1),aa2(1,1),aa1tder(1,1,1,1),
5481      &   a_chuj_der(1,1,1,1,kk,k),2,.false.,DtUg2EUg(1,1,k),
5482      &   DtUg2EUgder(1,1,1,k),ADtEA1(1,1,2),ADtEA1derg(1,1,1,2),
5483      &   ADtEA1derx(1,1,1,1,1,2))
5484         ENDIF
5485 C End 6-th order cumulants
5486         call transpose2(EUgder(1,1,l),auxmat(1,1))
5487         call matmat2(auxmat(1,1),AEA(1,1,2),EAEAderg(1,1,1,2))
5488         call transpose2(EUg(1,1,l),auxmat(1,1))
5489         call matmat2(auxmat(1,1),AEA(1,1,2),EAEA(1,1,2))
5490         call matmat2(auxmat(1,1),AEAderg(1,1,2),EAEAderg(1,1,2,2))
5491         do iii=1,2
5492           do kkk=1,5
5493             do lll=1,3
5494               call matmat2(auxmat(1,1),AEAderx(1,1,lll,kkk,iii,2),
5495      &          EAEAderx(1,1,lll,kkk,iii,2))
5496             enddo
5497           enddo
5498         enddo
5499 C AEAb1 and AEAb2
5500 C Calculate the vectors and their derivatives in virtual-bond dihedral angles.
5501 C They are needed only when the fifth- or the sixth-order cumulants are
5502 C indluded.
5503         IF (wcorr5.gt.0.0d0 .or. wcorr6.gt.0.0d0) THEN
5504         call transpose2(AEA(1,1,1),auxmat(1,1))
5505         call matvec2(auxmat(1,1),b1(1,iti),AEAb1(1,1,1))
5506         call matvec2(auxmat(1,1),Ub2(1,i),AEAb2(1,1,1))
5507         call matvec2(auxmat(1,1),Ub2der(1,i),AEAb2derg(1,2,1,1))
5508         call transpose2(AEAderg(1,1,1),auxmat(1,1))
5509         call matvec2(auxmat(1,1),b1(1,iti),AEAb1derg(1,1,1))
5510         call matvec2(auxmat(1,1),Ub2(1,i),AEAb2derg(1,1,1,1))
5511         call matvec2(AEA(1,1,1),b1(1,itk1),AEAb1(1,2,1))
5512         call matvec2(AEAderg(1,1,1),b1(1,itk1),AEAb1derg(1,2,1))
5513         call matvec2(AEA(1,1,1),Ub2(1,k+1),AEAb2(1,2,1))
5514         call matvec2(AEAderg(1,1,1),Ub2(1,k+1),AEAb2derg(1,1,2,1))
5515         call matvec2(AEA(1,1,1),Ub2der(1,k+1),AEAb2derg(1,2,2,1))
5516         call transpose2(AEA(1,1,2),auxmat(1,1))
5517         call matvec2(auxmat(1,1),b1(1,itj),AEAb1(1,1,2))
5518         call matvec2(auxmat(1,1),Ub2(1,j),AEAb2(1,1,2))
5519         call matvec2(auxmat(1,1),Ub2der(1,j),AEAb2derg(1,2,1,2))
5520         call transpose2(AEAderg(1,1,2),auxmat(1,1))
5521         call matvec2(auxmat(1,1),b1(1,itj),AEAb1derg(1,1,2))
5522         call matvec2(auxmat(1,1),Ub2(1,j),AEAb2derg(1,1,1,2))
5523         call matvec2(AEA(1,1,2),b1(1,itl1),AEAb1(1,2,2))
5524         call matvec2(AEAderg(1,1,2),b1(1,itl1),AEAb1derg(1,2,2))
5525         call matvec2(AEA(1,1,2),Ub2(1,l+1),AEAb2(1,2,2))
5526         call matvec2(AEAderg(1,1,2),Ub2(1,l+1),AEAb2derg(1,1,2,2))
5527         call matvec2(AEA(1,1,2),Ub2der(1,l+1),AEAb2derg(1,2,2,2))
5528 C Calculate the Cartesian derivatives of the vectors.
5529         do iii=1,2
5530           do kkk=1,5
5531             do lll=1,3
5532               call transpose2(AEAderx(1,1,lll,kkk,iii,1),auxmat(1,1))
5533               call matvec2(auxmat(1,1),b1(1,iti),
5534      &          AEAb1derx(1,lll,kkk,iii,1,1))
5535               call matvec2(auxmat(1,1),Ub2(1,i),
5536      &          AEAb2derx(1,lll,kkk,iii,1,1))
5537               call matvec2(AEAderx(1,1,lll,kkk,iii,1),b1(1,itk1),
5538      &          AEAb1derx(1,lll,kkk,iii,2,1))
5539               call matvec2(AEAderx(1,1,lll,kkk,iii,1),Ub2(1,k+1),
5540      &          AEAb2derx(1,lll,kkk,iii,2,1))
5541               call transpose2(AEAderx(1,1,lll,kkk,iii,2),auxmat(1,1))
5542               call matvec2(auxmat(1,1),b1(1,itj),
5543      &          AEAb1derx(1,lll,kkk,iii,1,2))
5544               call matvec2(auxmat(1,1),Ub2(1,j),
5545      &          AEAb2derx(1,lll,kkk,iii,1,2))
5546               call matvec2(AEAderx(1,1,lll,kkk,iii,2),b1(1,itl1),
5547      &          AEAb1derx(1,lll,kkk,iii,2,2))
5548               call matvec2(AEAderx(1,1,lll,kkk,iii,2),Ub2(1,l+1),
5549      &          AEAb2derx(1,lll,kkk,iii,2,2))
5550             enddo
5551           enddo
5552         enddo
5553         ENDIF
5554 C End vectors
5555       else
5556 C Antiparallel orientation of the two CA-CA-CA frames.
5557         if (i.gt.1 .and. itype(i).le.ntyp) then
5558           iti=itortyp(itype(i))
5559         else
5560           iti=ntortyp+1
5561         endif
5562         itk1=itortyp(itype(k+1))
5563         itl=itortyp(itype(l))
5564         itj=itortyp(itype(j))
5565         if (j.lt.nres-1 .and. itype(j+1).le.ntyp) then
5566           itj1=itortyp(itype(j+1))
5567         else 
5568           itj1=ntortyp+1
5569         endif
5570 C A2 kernel(j-1)T A1T
5571         call kernel(aa1(1,1),aa2t(1,1),a_chuj_der(1,1,1,1,jj,i),
5572      &   aa2tder(1,1,1,1),1,.true.,EUg(1,1,j),EUgder(1,1,j),
5573      &   AEA(1,1,1),AEAderg(1,1,1),AEAderx(1,1,1,1,1,1))
5574 C Following matrices are needed only for 6-th order cumulants
5575         IF (wcorr6.gt.0.0d0 .or. (wturn6.gt.0.0d0 .and.
5576      &     j.eq.i+4 .and. l.eq.i+3)) THEN
5577         call kernel(aa1(1,1),aa2t(1,1),a_chuj_der(1,1,1,1,jj,i),
5578      &   aa2tder(1,1,1,1),1,.true.,EUgC(1,1,j),EUgCder(1,1,j),
5579      &   AECA(1,1,1),AECAderg(1,1,1),AECAderx(1,1,1,1,1,1))
5580         call kernel(aa2(1,1),aa2t(1,1),a_chuj_der(1,1,1,1,jj,i),
5581      &   aa2tder(1,1,1,1),2,.true.,Ug2DtEUg(1,1,j),
5582      &   Ug2DtEUgder(1,1,1,j),ADtEA(1,1,1),ADtEAderg(1,1,1,1),
5583      &   ADtEAderx(1,1,1,1,1,1))
5584         call kernel(aa1(1,1),aa2t(1,1),a_chuj_der(1,1,1,1,jj,i),
5585      &   aa2tder(1,1,1,1),2,.true.,DtUg2EUg(1,1,j),
5586      &   DtUg2EUgder(1,1,1,j),ADtEA1(1,1,1),ADtEA1derg(1,1,1,1),
5587      &   ADtEA1derx(1,1,1,1,1,1))
5588         ENDIF
5589 C End 6-th order cumulants
5590         call transpose2(EUgder(1,1,k),auxmat(1,1))
5591         call matmat2(auxmat(1,1),AEA(1,1,1),EAEAderg(1,1,1,1))
5592         call transpose2(EUg(1,1,k),auxmat(1,1))
5593         call matmat2(auxmat(1,1),AEA(1,1,1),EAEA(1,1,1))
5594         call matmat2(auxmat(1,1),AEAderg(1,1,1),EAEAderg(1,1,2,1))
5595         do iii=1,2
5596           do kkk=1,5
5597             do lll=1,3
5598               call matmat2(auxmat(1,1),AEAderx(1,1,lll,kkk,iii,1),
5599      &          EAEAderx(1,1,lll,kkk,iii,1))
5600             enddo
5601           enddo
5602         enddo
5603 C A2T kernel(i+1)T A1
5604         call kernel(aa2t(1,1),aa1(1,1),aa2tder(1,1,1,1),
5605      &   a_chuj_der(1,1,1,1,jj,i),1,.true.,EUg(1,1,k),EUgder(1,1,k),
5606      &   AEA(1,1,2),AEAderg(1,1,2),AEAderx(1,1,1,1,1,2))
5607 C Following matrices are needed only for 6-th order cumulants
5608         IF (wcorr6.gt.0.0d0 .or. (wturn6.gt.0.0d0 .and.
5609      &     j.eq.i+4 .and. l.eq.i+3)) THEN
5610         call kernel(aa2t(1,1),aa1(1,1),aa2tder(1,1,1,1),
5611      &   a_chuj_der(1,1,1,1,jj,i),1,.true.,EUgC(1,1,k),EUgCder(1,1,k),
5612      &   AECA(1,1,2),AECAderg(1,1,2),AECAderx(1,1,1,1,1,2))
5613         call kernel(aa2t(1,1),aa1(1,1),aa2tder(1,1,1,1),
5614      &   a_chuj_der(1,1,1,1,jj,i),2,.true.,Ug2DtEUg(1,1,k),
5615      &   Ug2DtEUgder(1,1,1,k),ADtEA(1,1,2),ADtEAderg(1,1,1,2),
5616      &   ADtEAderx(1,1,1,1,1,2))
5617         call kernel(aa2t(1,1),aa1(1,1),aa2tder(1,1,1,1),
5618      &   a_chuj_der(1,1,1,1,jj,i),2,.true.,DtUg2EUg(1,1,k),
5619      &   DtUg2EUgder(1,1,1,k),ADtEA1(1,1,2),ADtEA1derg(1,1,1,2),
5620      &   ADtEA1derx(1,1,1,1,1,2))
5621         ENDIF
5622 C End 6-th order cumulants
5623         call transpose2(EUgder(1,1,j),auxmat(1,1))
5624         call matmat2(auxmat(1,1),AEA(1,1,1),EAEAderg(1,1,2,2))
5625         call transpose2(EUg(1,1,j),auxmat(1,1))
5626         call matmat2(auxmat(1,1),AEA(1,1,2),EAEA(1,1,2))
5627         call matmat2(auxmat(1,1),AEAderg(1,1,2),EAEAderg(1,1,2,2))
5628         do iii=1,2
5629           do kkk=1,5
5630             do lll=1,3
5631               call matmat2(auxmat(1,1),AEAderx(1,1,lll,kkk,iii,2),
5632      &          EAEAderx(1,1,lll,kkk,iii,2))
5633             enddo
5634           enddo
5635         enddo
5636 C AEAb1 and AEAb2
5637 C Calculate the vectors and their derivatives in virtual-bond dihedral angles.
5638 C They are needed only when the fifth- or the sixth-order cumulants are
5639 C indluded.
5640         IF (wcorr5.gt.0.0d0 .or. wcorr6.gt.0.0d0 .or.
5641      &    (wturn6.gt.0.0d0 .and. j.eq.i+4 .and. l.eq.i+3)) THEN
5642         call transpose2(AEA(1,1,1),auxmat(1,1))
5643         call matvec2(auxmat(1,1),b1(1,iti),AEAb1(1,1,1))
5644         call matvec2(auxmat(1,1),Ub2(1,i),AEAb2(1,1,1))
5645         call matvec2(auxmat(1,1),Ub2der(1,i),AEAb2derg(1,2,1,1))
5646         call transpose2(AEAderg(1,1,1),auxmat(1,1))
5647         call matvec2(auxmat(1,1),b1(1,iti),AEAb1derg(1,1,1))
5648         call matvec2(auxmat(1,1),Ub2(1,i),AEAb2derg(1,1,1,1))
5649         call matvec2(AEA(1,1,1),b1(1,itk1),AEAb1(1,2,1))
5650         call matvec2(AEAderg(1,1,1),b1(1,itk1),AEAb1derg(1,2,1))
5651         call matvec2(AEA(1,1,1),Ub2(1,k+1),AEAb2(1,2,1))
5652         call matvec2(AEAderg(1,1,1),Ub2(1,k+1),AEAb2derg(1,1,2,1))
5653         call matvec2(AEA(1,1,1),Ub2der(1,k+1),AEAb2derg(1,2,2,1))
5654         call transpose2(AEA(1,1,2),auxmat(1,1))
5655         call matvec2(auxmat(1,1),b1(1,itj1),AEAb1(1,1,2))
5656         call matvec2(auxmat(1,1),Ub2(1,l),AEAb2(1,1,2))
5657         call matvec2(auxmat(1,1),Ub2der(1,l),AEAb2derg(1,2,1,2))
5658         call transpose2(AEAderg(1,1,2),auxmat(1,1))
5659         call matvec2(auxmat(1,1),b1(1,itl),AEAb1(1,1,2))
5660         call matvec2(auxmat(1,1),Ub2(1,l),AEAb2derg(1,1,1,2))
5661         call matvec2(AEA(1,1,2),b1(1,itj1),AEAb1(1,2,2))
5662         call matvec2(AEAderg(1,1,2),b1(1,itj1),AEAb1derg(1,2,2))
5663         call matvec2(AEA(1,1,2),Ub2(1,j),AEAb2(1,2,2))
5664         call matvec2(AEAderg(1,1,2),Ub2(1,j),AEAb2derg(1,1,2,2))
5665         call matvec2(AEA(1,1,2),Ub2der(1,j),AEAb2derg(1,2,2,2))
5666 C Calculate the Cartesian derivatives of the vectors.
5667         do iii=1,2
5668           do kkk=1,5
5669             do lll=1,3
5670               call transpose2(AEAderx(1,1,lll,kkk,iii,1),auxmat(1,1))
5671               call matvec2(auxmat(1,1),b1(1,iti),
5672      &          AEAb1derx(1,lll,kkk,iii,1,1))
5673               call matvec2(auxmat(1,1),Ub2(1,i),
5674      &          AEAb2derx(1,lll,kkk,iii,1,1))
5675               call matvec2(AEAderx(1,1,lll,kkk,iii,1),b1(1,itk1),
5676      &          AEAb1derx(1,lll,kkk,iii,2,1))
5677               call matvec2(AEAderx(1,1,lll,kkk,iii,1),Ub2(1,k+1),
5678      &          AEAb2derx(1,lll,kkk,iii,2,1))
5679               call transpose2(AEAderx(1,1,lll,kkk,iii,2),auxmat(1,1))
5680               call matvec2(auxmat(1,1),b1(1,itl),
5681      &          AEAb1derx(1,lll,kkk,iii,1,2))
5682               call matvec2(auxmat(1,1),Ub2(1,l),
5683      &          AEAb2derx(1,lll,kkk,iii,1,2))
5684               call matvec2(AEAderx(1,1,lll,kkk,iii,2),b1(1,itj1),
5685      &          AEAb1derx(1,lll,kkk,iii,2,2))
5686               call matvec2(AEAderx(1,1,lll,kkk,iii,2),Ub2(1,j),
5687      &          AEAb2derx(1,lll,kkk,iii,2,2))
5688             enddo
5689           enddo
5690         enddo
5691         ENDIF
5692 C End vectors
5693       endif
5694       return
5695       end
5696 C---------------------------------------------------------------------------
5697       subroutine kernel(aa1,aa2t,aa1derx,aa2tderx,nderg,transp,
5698      &  KK,KKderg,AKA,AKAderg,AKAderx)
5699       implicit none
5700       integer nderg
5701       logical transp
5702       double precision aa1(2,2),aa2t(2,2),aa1derx(2,2,3,5),
5703      &  aa2tderx(2,2,3,5),KK(2,2),KKderg(2,2,nderg),AKA(2,2),
5704      &  AKAderg(2,2,nderg),AKAderx(2,2,3,5,2)
5705       integer iii,kkk,lll
5706       integer jjj,mmm
5707       logical lprn
5708       common /kutas/ lprn
5709       call prodmat3(aa1(1,1),aa2t(1,1),KK(1,1),transp,AKA(1,1))
5710       do iii=1,nderg 
5711         call prodmat3(aa1(1,1),aa2t(1,1),KKderg(1,1,iii),transp,
5712      &    AKAderg(1,1,iii))
5713       enddo
5714 cd      if (lprn) write (2,*) 'In kernel'
5715       do kkk=1,5
5716 cd        if (lprn) write (2,*) 'kkk=',kkk
5717         do lll=1,3
5718           call prodmat3(aa1derx(1,1,lll,kkk),aa2t(1,1),
5719      &      KK(1,1),transp,AKAderx(1,1,lll,kkk,1))
5720 cd          if (lprn) then
5721 cd            write (2,*) 'lll=',lll
5722 cd            write (2,*) 'iii=1'
5723 cd            do jjj=1,2
5724 cd              write (2,'(3(2f10.5),5x)') 
5725 cd     &        (AKAderx(jjj,mmm,lll,kkk,1),mmm=1,2)
5726 cd            enddo
5727 cd          endif
5728           call prodmat3(aa1(1,1),aa2tderx(1,1,lll,kkk),
5729      &      KK(1,1),transp,AKAderx(1,1,lll,kkk,2))
5730 cd          if (lprn) then
5731 cd            write (2,*) 'lll=',lll
5732 cd            write (2,*) 'iii=2'
5733 cd            do jjj=1,2
5734 cd              write (2,'(3(2f10.5),5x)') 
5735 cd     &        (AKAderx(jjj,mmm,lll,kkk,2),mmm=1,2)
5736 cd            enddo
5737 cd          endif
5738         enddo
5739       enddo
5740       return
5741       end
5742 C---------------------------------------------------------------------------
5743       double precision function eello4(i,j,k,l,jj,kk)
5744       implicit real*8 (a-h,o-z)
5745       include 'DIMENSIONS'
5746       include 'DIMENSIONS.ZSCOPT'
5747       include 'COMMON.IOUNITS'
5748       include 'COMMON.CHAIN'
5749       include 'COMMON.DERIV'
5750       include 'COMMON.INTERACT'
5751       include 'COMMON.CONTACTS'
5752       include 'COMMON.TORSION'
5753       include 'COMMON.VAR'
5754       include 'COMMON.GEO'
5755       double precision pizda(2,2),ggg1(3),ggg2(3)
5756 cd      if (i.ne.1 .or. j.ne.5 .or. k.ne.2 .or.l.ne.4) then
5757 cd        eello4=0.0d0
5758 cd        return
5759 cd      endif
5760 cd      print *,'eello4:',i,j,k,l,jj,kk
5761 cd      write (2,*) 'i',i,' j',j,' k',k,' l',l
5762 cd      call checkint4(i,j,k,l,jj,kk,eel4_num)
5763 cold      eij=facont_hb(jj,i)
5764 cold      ekl=facont_hb(kk,k)
5765 cold      ekont=eij*ekl
5766       eel4=-EAEA(1,1,1)-EAEA(2,2,1)
5767       if (calc_grad) then
5768 cd      eel41=-EAEA(1,1,2)-EAEA(2,2,2)
5769       gcorr_loc(k-1)=gcorr_loc(k-1)
5770      &   -ekont*(EAEAderg(1,1,1,1)+EAEAderg(2,2,1,1))
5771       if (l.eq.j+1) then
5772         gcorr_loc(l-1)=gcorr_loc(l-1)
5773      &     -ekont*(EAEAderg(1,1,2,1)+EAEAderg(2,2,2,1))
5774       else
5775         gcorr_loc(j-1)=gcorr_loc(j-1)
5776      &     -ekont*(EAEAderg(1,1,2,1)+EAEAderg(2,2,2,1))
5777       endif
5778       do iii=1,2
5779         do kkk=1,5
5780           do lll=1,3
5781             derx(lll,kkk,iii)=-EAEAderx(1,1,lll,kkk,iii,1)
5782      &                        -EAEAderx(2,2,lll,kkk,iii,1)
5783 cd            derx(lll,kkk,iii)=0.0d0
5784           enddo
5785         enddo
5786       enddo
5787 cd      gcorr_loc(l-1)=0.0d0
5788 cd      gcorr_loc(j-1)=0.0d0
5789 cd      gcorr_loc(k-1)=0.0d0
5790 cd      eel4=1.0d0
5791 cd      write (iout,*)'Contacts have occurred for peptide groups',
5792 cd     &  i,j,' fcont:',eij,' eij',' and ',k,l,
5793 cd     &  ' fcont ',ekl,' eel4=',eel4,' eel4_num',16*eel4_num
5794       if (j.lt.nres-1) then
5795         j1=j+1
5796         j2=j-1
5797       else
5798         j1=j-1
5799         j2=j-2
5800       endif
5801       if (l.lt.nres-1) then
5802         l1=l+1
5803         l2=l-1
5804       else
5805         l1=l-1
5806         l2=l-2
5807       endif
5808       do ll=1,3
5809 cold        ghalf=0.5d0*eel4*ekl*gacont_hbr(ll,jj,i)
5810         ggg1(ll)=eel4*g_contij(ll,1)
5811         ggg2(ll)=eel4*g_contij(ll,2)
5812         ghalf=0.5d0*ggg1(ll)
5813 cd        ghalf=0.0d0
5814         gradcorr(ll,i)=gradcorr(ll,i)+ghalf+ekont*derx(ll,2,1)
5815         gradcorr(ll,i+1)=gradcorr(ll,i+1)+ekont*derx(ll,3,1)
5816         gradcorr(ll,j)=gradcorr(ll,j)+ghalf+ekont*derx(ll,4,1)
5817         gradcorr(ll,j1)=gradcorr(ll,j1)+ekont*derx(ll,5,1)
5818 cold        ghalf=0.5d0*eel4*eij*gacont_hbr(ll,kk,k)
5819         ghalf=0.5d0*ggg2(ll)
5820 cd        ghalf=0.0d0
5821         gradcorr(ll,k)=gradcorr(ll,k)+ghalf+ekont*derx(ll,2,2)
5822         gradcorr(ll,k+1)=gradcorr(ll,k+1)+ekont*derx(ll,3,2)
5823         gradcorr(ll,l)=gradcorr(ll,l)+ghalf+ekont*derx(ll,4,2)
5824         gradcorr(ll,l1)=gradcorr(ll,l1)+ekont*derx(ll,5,2)
5825       enddo
5826 cd      goto 1112
5827       do m=i+1,j-1
5828         do ll=1,3
5829 cold          gradcorr(ll,m)=gradcorr(ll,m)+eel4*ekl*gacont_hbr(ll,jj,i)
5830           gradcorr(ll,m)=gradcorr(ll,m)+ggg1(ll)
5831         enddo
5832       enddo
5833       do m=k+1,l-1
5834         do ll=1,3
5835 cold          gradcorr(ll,m)=gradcorr(ll,m)+eel4*eij*gacont_hbr(ll,kk,k)
5836           gradcorr(ll,m)=gradcorr(ll,m)+ggg2(ll)
5837         enddo
5838       enddo
5839 1112  continue
5840       do m=i+2,j2
5841         do ll=1,3
5842           gradcorr(ll,m)=gradcorr(ll,m)+ekont*derx(ll,1,1)
5843         enddo
5844       enddo
5845       do m=k+2,l2
5846         do ll=1,3
5847           gradcorr(ll,m)=gradcorr(ll,m)+ekont*derx(ll,1,2)
5848         enddo
5849       enddo 
5850 cd      do iii=1,nres-3
5851 cd        write (2,*) iii,gcorr_loc(iii)
5852 cd      enddo
5853       endif
5854       eello4=ekont*eel4
5855 cd      write (2,*) 'ekont',ekont
5856 cd      write (iout,*) 'eello4',ekont*eel4
5857       return
5858       end
5859 C---------------------------------------------------------------------------
5860       double precision function eello5(i,j,k,l,jj,kk)
5861       implicit real*8 (a-h,o-z)
5862       include 'DIMENSIONS'
5863       include 'DIMENSIONS.ZSCOPT'
5864       include 'COMMON.IOUNITS'
5865       include 'COMMON.CHAIN'
5866       include 'COMMON.DERIV'
5867       include 'COMMON.INTERACT'
5868       include 'COMMON.CONTACTS'
5869       include 'COMMON.TORSION'
5870       include 'COMMON.VAR'
5871       include 'COMMON.GEO'
5872       double precision pizda(2,2),auxmat(2,2),auxmat1(2,2),vv(2)
5873       double precision ggg1(3),ggg2(3)
5874 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
5875 C                                                                              C
5876 C                            Parallel chains                                   C
5877 C                                                                              C
5878 C          o             o                   o             o                   C
5879 C         /l\           / \             \   / \           / \   /              C
5880 C        /   \         /   \             \ /   \         /   \ /               C
5881 C       j| o |l1       | o |              o| o |         | o |o                C
5882 C     \  |/k\|         |/ \|  /            |/ \|         |/ \|                 C
5883 C      \i/   \         /   \ /             /   \         /   \                 C
5884 C       o    k1             o                                                  C
5885 C         (I)          (II)                (III)          (IV)                 C
5886 C                                                                              C
5887 C      eello5_1        eello5_2            eello5_3       eello5_4             C
5888 C                                                                              C
5889 C                            Antiparallel chains                               C
5890 C                                                                              C
5891 C          o             o                   o             o                   C
5892 C         /j\           / \             \   / \           / \   /              C
5893 C        /   \         /   \             \ /   \         /   \ /               C
5894 C      j1| o |l        | o |              o| o |         | o |o                C
5895 C     \  |/k\|         |/ \|  /            |/ \|         |/ \|                 C
5896 C      \i/   \         /   \ /             /   \         /   \                 C
5897 C       o     k1            o                                                  C
5898 C         (I)          (II)                (III)          (IV)                 C
5899 C                                                                              C
5900 C      eello5_1        eello5_2            eello5_3       eello5_4             C
5901 C                                                                              C
5902 C o denotes a local interaction, vertical lines an electrostatic interaction.  C
5903 C                                                                              C
5904 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
5905 cd      if (i.ne.2 .or. j.ne.6 .or. k.ne.3 .or. l.ne.5) then
5906 cd        eello5=0.0d0
5907 cd        return
5908 cd      endif
5909 cd      write (iout,*)
5910 cd     &   'EELLO5: Contacts have occurred for peptide groups',i,j,
5911 cd     &   ' and',k,l
5912       itk=itortyp(itype(k))
5913       itl=itortyp(itype(l))
5914       itj=itortyp(itype(j))
5915       eello5_1=0.0d0
5916       eello5_2=0.0d0
5917       eello5_3=0.0d0
5918       eello5_4=0.0d0
5919 cd      call checkint5(i,j,k,l,jj,kk,eel5_1_num,eel5_2_num,
5920 cd     &   eel5_3_num,eel5_4_num)
5921       do iii=1,2
5922         do kkk=1,5
5923           do lll=1,3
5924             derx(lll,kkk,iii)=0.0d0
5925           enddo
5926         enddo
5927       enddo
5928 cd      eij=facont_hb(jj,i)
5929 cd      ekl=facont_hb(kk,k)
5930 cd      ekont=eij*ekl
5931 cd      write (iout,*)'Contacts have occurred for peptide groups',
5932 cd     &  i,j,' fcont:',eij,' eij',' and ',k,l
5933 cd      goto 1111
5934 C Contribution from the graph I.
5935 cd      write (2,*) 'AEA  ',AEA(1,1,1),AEA(2,1,1),AEA(1,2,1),AEA(2,2,1)
5936 cd      write (2,*) 'AEAb2',AEAb2(1,1,1),AEAb2(2,1,1)
5937       call transpose2(EUg(1,1,k),auxmat(1,1))
5938       call matmat2(AEA(1,1,1),auxmat(1,1),pizda(1,1))
5939       vv(1)=pizda(1,1)-pizda(2,2)
5940       vv(2)=pizda(1,2)+pizda(2,1)
5941       eello5_1=scalar2(AEAb2(1,1,1),Ub2(1,k))
5942      & +0.5d0*scalar2(vv(1),Dtobr2(1,i))
5943       if (calc_grad) then
5944 C Explicit gradient in virtual-dihedral angles.
5945       if (i.gt.1) g_corr5_loc(i-1)=g_corr5_loc(i-1)
5946      & +ekont*(scalar2(AEAb2derg(1,2,1,1),Ub2(1,k))
5947      & +0.5d0*scalar2(vv(1),Dtobr2der(1,i)))
5948       call transpose2(EUgder(1,1,k),auxmat1(1,1))
5949       call matmat2(AEA(1,1,1),auxmat1(1,1),pizda(1,1))
5950       vv(1)=pizda(1,1)-pizda(2,2)
5951       vv(2)=pizda(1,2)+pizda(2,1)
5952       g_corr5_loc(k-1)=g_corr5_loc(k-1)
5953      & +ekont*(scalar2(AEAb2(1,1,1),Ub2der(1,k))
5954      & +0.5d0*scalar2(vv(1),Dtobr2(1,i)))
5955       call matmat2(AEAderg(1,1,1),auxmat(1,1),pizda(1,1))
5956       vv(1)=pizda(1,1)-pizda(2,2)
5957       vv(2)=pizda(1,2)+pizda(2,1)
5958       if (l.eq.j+1) then
5959         if (l.lt.nres-1) g_corr5_loc(l-1)=g_corr5_loc(l-1)
5960      &   +ekont*(scalar2(AEAb2derg(1,1,1,1),Ub2(1,k))
5961      &   +0.5d0*scalar2(vv(1),Dtobr2(1,i)))
5962       else
5963         if (j.lt.nres-1) g_corr5_loc(j-1)=g_corr5_loc(j-1)
5964      &   +ekont*(scalar2(AEAb2derg(1,1,1,1),Ub2(1,k))
5965      &   +0.5d0*scalar2(vv(1),Dtobr2(1,i)))
5966       endif 
5967 C Cartesian gradient
5968       do iii=1,2
5969         do kkk=1,5
5970           do lll=1,3
5971             call matmat2(AEAderx(1,1,lll,kkk,iii,1),auxmat(1,1),
5972      &        pizda(1,1))
5973             vv(1)=pizda(1,1)-pizda(2,2)
5974             vv(2)=pizda(1,2)+pizda(2,1)
5975             derx(lll,kkk,iii)=derx(lll,kkk,iii)
5976      &       +scalar2(AEAb2derx(1,lll,kkk,iii,1,1),Ub2(1,k))
5977      &       +0.5d0*scalar2(vv(1),Dtobr2(1,i))
5978           enddo
5979         enddo
5980       enddo
5981 c      goto 1112
5982       endif
5983 c1111  continue
5984 C Contribution from graph II 
5985       call transpose2(EE(1,1,itk),auxmat(1,1))
5986       call matmat2(auxmat(1,1),AEA(1,1,1),pizda(1,1))
5987       vv(1)=pizda(1,1)+pizda(2,2)
5988       vv(2)=pizda(2,1)-pizda(1,2)
5989       eello5_2=scalar2(AEAb1(1,2,1),b1(1,itk))
5990      & -0.5d0*scalar2(vv(1),Ctobr(1,k))
5991       if (calc_grad) then
5992 C Explicit gradient in virtual-dihedral angles.
5993       g_corr5_loc(k-1)=g_corr5_loc(k-1)
5994      & -0.5d0*ekont*scalar2(vv(1),Ctobrder(1,k))
5995       call matmat2(auxmat(1,1),AEAderg(1,1,1),pizda(1,1))
5996       vv(1)=pizda(1,1)+pizda(2,2)
5997       vv(2)=pizda(2,1)-pizda(1,2)
5998       if (l.eq.j+1) then
5999         g_corr5_loc(l-1)=g_corr5_loc(l-1)
6000      &   +ekont*(scalar2(AEAb1derg(1,2,1),b1(1,itk))
6001      &   -0.5d0*scalar2(vv(1),Ctobr(1,k)))
6002       else
6003         g_corr5_loc(j-1)=g_corr5_loc(j-1)
6004      &   +ekont*(scalar2(AEAb1derg(1,2,1),b1(1,itk))
6005      &   -0.5d0*scalar2(vv(1),Ctobr(1,k)))
6006       endif
6007 C Cartesian gradient
6008       do iii=1,2
6009         do kkk=1,5
6010           do lll=1,3
6011             call matmat2(auxmat(1,1),AEAderx(1,1,lll,kkk,iii,1),
6012      &        pizda(1,1))
6013             vv(1)=pizda(1,1)+pizda(2,2)
6014             vv(2)=pizda(2,1)-pizda(1,2)
6015             derx(lll,kkk,iii)=derx(lll,kkk,iii)
6016      &       +scalar2(AEAb1derx(1,lll,kkk,iii,2,1),b1(1,itk))
6017      &       -0.5d0*scalar2(vv(1),Ctobr(1,k))
6018           enddo
6019         enddo
6020       enddo
6021 cd      goto 1112
6022       endif
6023 cd1111  continue
6024       if (l.eq.j+1) then
6025 cd        goto 1110
6026 C Parallel orientation
6027 C Contribution from graph III
6028         call transpose2(EUg(1,1,l),auxmat(1,1))
6029         call matmat2(AEA(1,1,2),auxmat(1,1),pizda(1,1))
6030         vv(1)=pizda(1,1)-pizda(2,2)
6031         vv(2)=pizda(1,2)+pizda(2,1)
6032         eello5_3=scalar2(AEAb2(1,1,2),Ub2(1,l))
6033      &   +0.5d0*scalar2(vv(1),Dtobr2(1,j))
6034         if (calc_grad) then
6035 C Explicit gradient in virtual-dihedral angles.
6036         g_corr5_loc(j-1)=g_corr5_loc(j-1)
6037      &   +ekont*(scalar2(AEAb2derg(1,2,1,2),Ub2(1,l))
6038      &   +0.5d0*scalar2(vv(1),Dtobr2der(1,j)))
6039         call matmat2(AEAderg(1,1,2),auxmat(1,1),pizda(1,1))
6040         vv(1)=pizda(1,1)-pizda(2,2)
6041         vv(2)=pizda(1,2)+pizda(2,1)
6042         g_corr5_loc(k-1)=g_corr5_loc(k-1)
6043      &   +ekont*(scalar2(AEAb2derg(1,1,1,2),Ub2(1,l))
6044      &   +0.5d0*scalar2(vv(1),Dtobr2(1,j)))
6045         call transpose2(EUgder(1,1,l),auxmat1(1,1))
6046         call matmat2(AEA(1,1,2),auxmat1(1,1),pizda(1,1))
6047         vv(1)=pizda(1,1)-pizda(2,2)
6048         vv(2)=pizda(1,2)+pizda(2,1)
6049         g_corr5_loc(l-1)=g_corr5_loc(l-1)
6050      &   +ekont*(scalar2(AEAb2(1,1,2),Ub2der(1,l))
6051      &   +0.5d0*scalar2(vv(1),Dtobr2(1,j)))
6052 C Cartesian gradient
6053         do iii=1,2
6054           do kkk=1,5
6055             do lll=1,3
6056               call matmat2(AEAderx(1,1,lll,kkk,iii,2),auxmat(1,1),
6057      &          pizda(1,1))
6058               vv(1)=pizda(1,1)-pizda(2,2)
6059               vv(2)=pizda(1,2)+pizda(2,1)
6060               derx(lll,kkk,iii)=derx(lll,kkk,iii)
6061      &         +scalar2(AEAb2derx(1,lll,kkk,iii,1,2),Ub2(1,l))
6062      &         +0.5d0*scalar2(vv(1),Dtobr2(1,j))
6063             enddo
6064           enddo
6065         enddo
6066 cd        goto 1112
6067         endif
6068 C Contribution from graph IV
6069 cd1110    continue
6070         call transpose2(EE(1,1,itl),auxmat(1,1))
6071         call matmat2(auxmat(1,1),AEA(1,1,2),pizda(1,1))
6072         vv(1)=pizda(1,1)+pizda(2,2)
6073         vv(2)=pizda(2,1)-pizda(1,2)
6074         eello5_4=scalar2(AEAb1(1,2,2),b1(1,itl))
6075      &   -0.5d0*scalar2(vv(1),Ctobr(1,l))
6076         if (calc_grad) then
6077 C Explicit gradient in virtual-dihedral angles.
6078         g_corr5_loc(l-1)=g_corr5_loc(l-1)
6079      &   -0.5d0*ekont*scalar2(vv(1),Ctobrder(1,l))
6080         call matmat2(auxmat(1,1),AEAderg(1,1,2),pizda(1,1))
6081         vv(1)=pizda(1,1)+pizda(2,2)
6082         vv(2)=pizda(2,1)-pizda(1,2)
6083         g_corr5_loc(k-1)=g_corr5_loc(k-1)
6084      &   +ekont*(scalar2(AEAb1derg(1,2,2),b1(1,itl))
6085      &   -0.5d0*scalar2(vv(1),Ctobr(1,l)))
6086 C Cartesian gradient
6087         do iii=1,2
6088           do kkk=1,5
6089             do lll=1,3
6090               call matmat2(auxmat(1,1),AEAderx(1,1,lll,kkk,iii,2),
6091      &          pizda(1,1))
6092               vv(1)=pizda(1,1)+pizda(2,2)
6093               vv(2)=pizda(2,1)-pizda(1,2)
6094               derx(lll,kkk,iii)=derx(lll,kkk,iii)
6095      &         +scalar2(AEAb1derx(1,lll,kkk,iii,2,2),b1(1,itl))
6096      &         -0.5d0*scalar2(vv(1),Ctobr(1,l))
6097             enddo
6098           enddo
6099         enddo
6100         endif
6101       else
6102 C Antiparallel orientation
6103 C Contribution from graph III
6104 c        goto 1110
6105         call transpose2(EUg(1,1,j),auxmat(1,1))
6106         call matmat2(AEA(1,1,2),auxmat(1,1),pizda(1,1))
6107         vv(1)=pizda(1,1)-pizda(2,2)
6108         vv(2)=pizda(1,2)+pizda(2,1)
6109         eello5_3=scalar2(AEAb2(1,1,2),Ub2(1,j))
6110      &   +0.5d0*scalar2(vv(1),Dtobr2(1,l))
6111         if (calc_grad) then
6112 C Explicit gradient in virtual-dihedral angles.
6113         g_corr5_loc(l-1)=g_corr5_loc(l-1)
6114      &   +ekont*(scalar2(AEAb2derg(1,2,1,2),Ub2(1,j))
6115      &   +0.5d0*scalar2(vv(1),Dtobr2der(1,l)))
6116         call matmat2(AEAderg(1,1,2),auxmat(1,1),pizda(1,1))
6117         vv(1)=pizda(1,1)-pizda(2,2)
6118         vv(2)=pizda(1,2)+pizda(2,1)
6119         g_corr5_loc(k-1)=g_corr5_loc(k-1)
6120      &   +ekont*(scalar2(AEAb2derg(1,1,1,2),Ub2(1,j))
6121      &   +0.5d0*scalar2(vv(1),Dtobr2(1,l)))
6122         call transpose2(EUgder(1,1,j),auxmat1(1,1))
6123         call matmat2(AEA(1,1,2),auxmat1(1,1),pizda(1,1))
6124         vv(1)=pizda(1,1)-pizda(2,2)
6125         vv(2)=pizda(1,2)+pizda(2,1)
6126         g_corr5_loc(j-1)=g_corr5_loc(j-1)
6127      &   +ekont*(scalar2(AEAb2(1,1,2),Ub2der(1,j))
6128      &   +0.5d0*scalar2(vv(1),Dtobr2(1,l)))
6129 C Cartesian gradient
6130         do iii=1,2
6131           do kkk=1,5
6132             do lll=1,3
6133               call matmat2(AEAderx(1,1,lll,kkk,iii,2),auxmat(1,1),
6134      &          pizda(1,1))
6135               vv(1)=pizda(1,1)-pizda(2,2)
6136               vv(2)=pizda(1,2)+pizda(2,1)
6137               derx(lll,kkk,3-iii)=derx(lll,kkk,3-iii)
6138      &         +scalar2(AEAb2derx(1,lll,kkk,iii,1,2),Ub2(1,j))
6139      &         +0.5d0*scalar2(vv(1),Dtobr2(1,l))
6140             enddo
6141           enddo
6142         enddo
6143 cd        goto 1112
6144         endif
6145 C Contribution from graph IV
6146 1110    continue
6147         call transpose2(EE(1,1,itj),auxmat(1,1))
6148         call matmat2(auxmat(1,1),AEA(1,1,2),pizda(1,1))
6149         vv(1)=pizda(1,1)+pizda(2,2)
6150         vv(2)=pizda(2,1)-pizda(1,2)
6151         eello5_4=scalar2(AEAb1(1,2,2),b1(1,itj))
6152      &   -0.5d0*scalar2(vv(1),Ctobr(1,j))
6153         if (calc_grad) then
6154 C Explicit gradient in virtual-dihedral angles.
6155         g_corr5_loc(j-1)=g_corr5_loc(j-1)
6156      &   -0.5d0*ekont*scalar2(vv(1),Ctobrder(1,j))
6157         call matmat2(auxmat(1,1),AEAderg(1,1,2),pizda(1,1))
6158         vv(1)=pizda(1,1)+pizda(2,2)
6159         vv(2)=pizda(2,1)-pizda(1,2)
6160         g_corr5_loc(k-1)=g_corr5_loc(k-1)
6161      &   +ekont*(scalar2(AEAb1derg(1,2,2),b1(1,itj))
6162      &   -0.5d0*scalar2(vv(1),Ctobr(1,j)))
6163 C Cartesian gradient
6164         do iii=1,2
6165           do kkk=1,5
6166             do lll=1,3
6167               call matmat2(auxmat(1,1),AEAderx(1,1,lll,kkk,iii,2),
6168      &          pizda(1,1))
6169               vv(1)=pizda(1,1)+pizda(2,2)
6170               vv(2)=pizda(2,1)-pizda(1,2)
6171               derx(lll,kkk,3-iii)=derx(lll,kkk,3-iii)
6172      &         +scalar2(AEAb1derx(1,lll,kkk,iii,2,2),b1(1,itj))
6173      &         -0.5d0*scalar2(vv(1),Ctobr(1,j))
6174             enddo
6175           enddo
6176         enddo
6177       endif
6178       endif
6179 1112  continue
6180       eel5=eello5_1+eello5_2+eello5_3+eello5_4
6181 cd      if (i.eq.2 .and. j.eq.8 .and. k.eq.3 .and. l.eq.7) then
6182 cd        write (2,*) 'ijkl',i,j,k,l
6183 cd        write (2,*) 'eello5_1',eello5_1,' eello5_2',eello5_2,
6184 cd     &     ' eello5_3',eello5_3,' eello5_4',eello5_4
6185 cd      endif
6186 cd      write(iout,*) 'eello5_1',eello5_1,' eel5_1_num',16*eel5_1_num
6187 cd      write(iout,*) 'eello5_2',eello5_2,' eel5_2_num',16*eel5_2_num
6188 cd      write(iout,*) 'eello5_3',eello5_3,' eel5_3_num',16*eel5_3_num
6189 cd      write(iout,*) 'eello5_4',eello5_4,' eel5_4_num',16*eel5_4_num
6190       if (calc_grad) then
6191       if (j.lt.nres-1) then
6192         j1=j+1
6193         j2=j-1
6194       else
6195         j1=j-1
6196         j2=j-2
6197       endif
6198       if (l.lt.nres-1) then
6199         l1=l+1
6200         l2=l-1
6201       else
6202         l1=l-1
6203         l2=l-2
6204       endif
6205 cd      eij=1.0d0
6206 cd      ekl=1.0d0
6207 cd      ekont=1.0d0
6208 cd      write (2,*) 'eij',eij,' ekl',ekl,' ekont',ekont
6209       do ll=1,3
6210         ggg1(ll)=eel5*g_contij(ll,1)
6211         ggg2(ll)=eel5*g_contij(ll,2)
6212 cold        ghalf=0.5d0*eel5*ekl*gacont_hbr(ll,jj,i)
6213         ghalf=0.5d0*ggg1(ll)
6214 cd        ghalf=0.0d0
6215         gradcorr5(ll,i)=gradcorr5(ll,i)+ghalf+ekont*derx(ll,2,1)
6216         gradcorr5(ll,i+1)=gradcorr5(ll,i+1)+ekont*derx(ll,3,1)
6217         gradcorr5(ll,j)=gradcorr5(ll,j)+ghalf+ekont*derx(ll,4,1)
6218         gradcorr5(ll,j1)=gradcorr5(ll,j1)+ekont*derx(ll,5,1)
6219 cold        ghalf=0.5d0*eel5*eij*gacont_hbr(ll,kk,k)
6220         ghalf=0.5d0*ggg2(ll)
6221 cd        ghalf=0.0d0
6222         gradcorr5(ll,k)=gradcorr5(ll,k)+ghalf+ekont*derx(ll,2,2)
6223         gradcorr5(ll,k+1)=gradcorr5(ll,k+1)+ekont*derx(ll,3,2)
6224         gradcorr5(ll,l)=gradcorr5(ll,l)+ghalf+ekont*derx(ll,4,2)
6225         gradcorr5(ll,l1)=gradcorr5(ll,l1)+ekont*derx(ll,5,2)
6226       enddo
6227 cd      goto 1112
6228       do m=i+1,j-1
6229         do ll=1,3
6230 cold          gradcorr5(ll,m)=gradcorr5(ll,m)+eel5*ekl*gacont_hbr(ll,jj,i)
6231           gradcorr5(ll,m)=gradcorr5(ll,m)+ggg1(ll)
6232         enddo
6233       enddo
6234       do m=k+1,l-1
6235         do ll=1,3
6236 cold          gradcorr5(ll,m)=gradcorr5(ll,m)+eel5*eij*gacont_hbr(ll,kk,k)
6237           gradcorr5(ll,m)=gradcorr5(ll,m)+ggg2(ll)
6238         enddo
6239       enddo
6240 c1112  continue
6241       do m=i+2,j2
6242         do ll=1,3
6243           gradcorr5(ll,m)=gradcorr5(ll,m)+ekont*derx(ll,1,1)
6244         enddo
6245       enddo
6246       do m=k+2,l2
6247         do ll=1,3
6248           gradcorr5(ll,m)=gradcorr5(ll,m)+ekont*derx(ll,1,2)
6249         enddo
6250       enddo 
6251 cd      do iii=1,nres-3
6252 cd        write (2,*) iii,g_corr5_loc(iii)
6253 cd      enddo
6254       endif
6255       eello5=ekont*eel5
6256 cd      write (2,*) 'ekont',ekont
6257 cd      write (iout,*) 'eello5',ekont*eel5
6258       return
6259       end
6260 c--------------------------------------------------------------------------
6261       double precision function eello6(i,j,k,l,jj,kk)
6262       implicit real*8 (a-h,o-z)
6263       include 'DIMENSIONS'
6264       include 'DIMENSIONS.ZSCOPT'
6265       include 'COMMON.IOUNITS'
6266       include 'COMMON.CHAIN'
6267       include 'COMMON.DERIV'
6268       include 'COMMON.INTERACT'
6269       include 'COMMON.CONTACTS'
6270       include 'COMMON.TORSION'
6271       include 'COMMON.VAR'
6272       include 'COMMON.GEO'
6273       include 'COMMON.FFIELD'
6274       double precision ggg1(3),ggg2(3)
6275 cd      if (i.ne.1 .or. j.ne.3 .or. k.ne.2 .or. l.ne.4) then
6276 cd        eello6=0.0d0
6277 cd        return
6278 cd      endif
6279 cd      write (iout,*)
6280 cd     &   'EELLO6: Contacts have occurred for peptide groups',i,j,
6281 cd     &   ' and',k,l
6282       eello6_1=0.0d0
6283       eello6_2=0.0d0
6284       eello6_3=0.0d0
6285       eello6_4=0.0d0
6286       eello6_5=0.0d0
6287       eello6_6=0.0d0
6288 cd      call checkint6(i,j,k,l,jj,kk,eel6_1_num,eel6_2_num,
6289 cd     &   eel6_3_num,eel6_4_num,eel6_5_num,eel6_6_num)
6290       do iii=1,2
6291         do kkk=1,5
6292           do lll=1,3
6293             derx(lll,kkk,iii)=0.0d0
6294           enddo
6295         enddo
6296       enddo
6297 cd      eij=facont_hb(jj,i)
6298 cd      ekl=facont_hb(kk,k)
6299 cd      ekont=eij*ekl
6300 cd      eij=1.0d0
6301 cd      ekl=1.0d0
6302 cd      ekont=1.0d0
6303       if (l.eq.j+1) then
6304         eello6_1=eello6_graph1(i,j,k,l,1,.false.)
6305         eello6_2=eello6_graph1(j,i,l,k,2,.false.)
6306         eello6_3=eello6_graph2(i,j,k,l,jj,kk,.false.)
6307         eello6_4=eello6_graph4(i,j,k,l,jj,kk,1,.false.)
6308         eello6_5=eello6_graph4(j,i,l,k,jj,kk,2,.false.)
6309         eello6_6=eello6_graph3(i,j,k,l,jj,kk,.false.)
6310       else
6311         eello6_1=eello6_graph1(i,j,k,l,1,.false.)
6312         eello6_2=eello6_graph1(l,k,j,i,2,.true.)
6313         eello6_3=eello6_graph2(i,l,k,j,jj,kk,.true.)
6314         eello6_4=eello6_graph4(i,j,k,l,jj,kk,1,.false.)
6315         if (wturn6.eq.0.0d0 .or. j.ne.i+4) then
6316           eello6_5=eello6_graph4(l,k,j,i,kk,jj,2,.true.)
6317         else
6318           eello6_5=0.0d0
6319         endif
6320         eello6_6=eello6_graph3(i,l,k,j,jj,kk,.true.)
6321       endif
6322 C If turn contributions are considered, they will be handled separately.
6323       eel6=eello6_1+eello6_2+eello6_3+eello6_4+eello6_5+eello6_6
6324 cd      write(iout,*) 'eello6_1',eello6_1,' eel6_1_num',16*eel6_1_num
6325 cd      write(iout,*) 'eello6_2',eello6_2,' eel6_2_num',16*eel6_2_num
6326 cd      write(iout,*) 'eello6_3',eello6_3,' eel6_3_num',16*eel6_3_num
6327 cd      write(iout,*) 'eello6_4',eello6_4,' eel6_4_num',16*eel6_4_num
6328 cd      write(iout,*) 'eello6_5',eello6_5,' eel6_5_num',16*eel6_5_num
6329 cd      write(iout,*) 'eello6_6',eello6_6,' eel6_6_num',16*eel6_6_num
6330 cd      goto 1112
6331       if (calc_grad) then
6332       if (j.lt.nres-1) then
6333         j1=j+1
6334         j2=j-1
6335       else
6336         j1=j-1
6337         j2=j-2
6338       endif
6339       if (l.lt.nres-1) then
6340         l1=l+1
6341         l2=l-1
6342       else
6343         l1=l-1
6344         l2=l-2
6345       endif
6346       do ll=1,3
6347         ggg1(ll)=eel6*g_contij(ll,1)
6348         ggg2(ll)=eel6*g_contij(ll,2)
6349 cold        ghalf=0.5d0*eel6*ekl*gacont_hbr(ll,jj,i)
6350         ghalf=0.5d0*ggg1(ll)
6351 cd        ghalf=0.0d0
6352         gradcorr6(ll,i)=gradcorr6(ll,i)+ghalf+ekont*derx(ll,2,1)
6353         gradcorr6(ll,i+1)=gradcorr6(ll,i+1)+ekont*derx(ll,3,1)
6354         gradcorr6(ll,j)=gradcorr6(ll,j)+ghalf+ekont*derx(ll,4,1)
6355         gradcorr6(ll,j1)=gradcorr6(ll,j1)+ekont*derx(ll,5,1)
6356         ghalf=0.5d0*ggg2(ll)
6357 cold        ghalf=0.5d0*eel6*eij*gacont_hbr(ll,kk,k)
6358 cd        ghalf=0.0d0
6359         gradcorr6(ll,k)=gradcorr6(ll,k)+ghalf+ekont*derx(ll,2,2)
6360         gradcorr6(ll,k+1)=gradcorr6(ll,k+1)+ekont*derx(ll,3,2)
6361         gradcorr6(ll,l)=gradcorr6(ll,l)+ghalf+ekont*derx(ll,4,2)
6362         gradcorr6(ll,l1)=gradcorr6(ll,l1)+ekont*derx(ll,5,2)
6363       enddo
6364 cd      goto 1112
6365       do m=i+1,j-1
6366         do ll=1,3
6367 cold          gradcorr6(ll,m)=gradcorr6(ll,m)+eel6*ekl*gacont_hbr(ll,jj,i)
6368           gradcorr6(ll,m)=gradcorr6(ll,m)+ggg1(ll)
6369         enddo
6370       enddo
6371       do m=k+1,l-1
6372         do ll=1,3
6373 cold          gradcorr6(ll,m)=gradcorr6(ll,m)+eel6*eij*gacont_hbr(ll,kk,k)
6374           gradcorr6(ll,m)=gradcorr6(ll,m)+ggg2(ll)
6375         enddo
6376       enddo
6377 1112  continue
6378       do m=i+2,j2
6379         do ll=1,3
6380           gradcorr6(ll,m)=gradcorr6(ll,m)+ekont*derx(ll,1,1)
6381         enddo
6382       enddo
6383       do m=k+2,l2
6384         do ll=1,3
6385           gradcorr6(ll,m)=gradcorr6(ll,m)+ekont*derx(ll,1,2)
6386         enddo
6387       enddo 
6388 cd      do iii=1,nres-3
6389 cd        write (2,*) iii,g_corr6_loc(iii)
6390 cd      enddo
6391       endif
6392       eello6=ekont*eel6
6393 cd      write (2,*) 'ekont',ekont
6394 cd      write (iout,*) 'eello6',ekont*eel6
6395       return
6396       end
6397 c--------------------------------------------------------------------------
6398       double precision function eello6_graph1(i,j,k,l,imat,swap)
6399       implicit real*8 (a-h,o-z)
6400       include 'DIMENSIONS'
6401       include 'DIMENSIONS.ZSCOPT'
6402       include 'COMMON.IOUNITS'
6403       include 'COMMON.CHAIN'
6404       include 'COMMON.DERIV'
6405       include 'COMMON.INTERACT'
6406       include 'COMMON.CONTACTS'
6407       include 'COMMON.TORSION'
6408       include 'COMMON.VAR'
6409       include 'COMMON.GEO'
6410       double precision vv(2),vv1(2),pizda(2,2),auxmat(2,2),pizda1(2,2)
6411       logical swap
6412       logical lprn
6413       common /kutas/ lprn
6414 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
6415 C                                                                              C 
6416 C      Parallel       Antiparallel                                             C
6417 C                                                                              C
6418 C          o             o                                                     C
6419 C         /l\           /j\                                                    C
6420 C        /   \         /   \                                                   C
6421 C       /| o |         | o |\                                                  C
6422 C     \ j|/k\|  /   \  |/k\|l /                                                C
6423 C      \ /   \ /     \ /   \ /                                                 C
6424 C       o     o       o     o                                                  C
6425 C       i             i                                                        C
6426 C                                                                              C
6427 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
6428       itk=itortyp(itype(k))
6429       s1= scalar2(AEAb1(1,2,imat),CUgb2(1,i))
6430       s2=-scalar2(AEAb2(1,1,imat),Ug2Db1t(1,k))
6431       s3= scalar2(AEAb2(1,1,imat),CUgb2(1,k))
6432       call transpose2(EUgC(1,1,k),auxmat(1,1))
6433       call matmat2(AEA(1,1,imat),auxmat(1,1),pizda1(1,1))
6434       vv1(1)=pizda1(1,1)-pizda1(2,2)
6435       vv1(2)=pizda1(1,2)+pizda1(2,1)
6436       s4=0.5d0*scalar2(vv1(1),Dtobr2(1,i))
6437       vv(1)=AEAb1(1,2,imat)*b1(1,itk)-AEAb1(2,2,imat)*b1(2,itk)
6438       vv(2)=AEAb1(1,2,imat)*b1(2,itk)+AEAb1(2,2,imat)*b1(1,itk)
6439       s5=scalar2(vv(1),Dtobr2(1,i))
6440 cd      write (2,*) 's1',s1,' s2',s2,' s3',s3,' s4', s4,' s5',s5
6441       eello6_graph1=-0.5d0*(s1+s2+s3+s4+s5)
6442       if (.not. calc_grad) return
6443       if (i.gt.1) g_corr6_loc(i-1)=g_corr6_loc(i-1)
6444      & -0.5d0*ekont*(scalar2(AEAb1(1,2,imat),CUgb2der(1,i))
6445      & -scalar2(AEAb2derg(1,2,1,imat),Ug2Db1t(1,k))
6446      & +scalar2(AEAb2derg(1,2,1,imat),CUgb2(1,k))
6447      & +0.5d0*scalar2(vv1(1),Dtobr2der(1,i))
6448      & +scalar2(vv(1),Dtobr2der(1,i)))
6449       call matmat2(AEAderg(1,1,imat),auxmat(1,1),pizda1(1,1))
6450       vv1(1)=pizda1(1,1)-pizda1(2,2)
6451       vv1(2)=pizda1(1,2)+pizda1(2,1)
6452       vv(1)=AEAb1derg(1,2,imat)*b1(1,itk)-AEAb1derg(2,2,imat)*b1(2,itk)
6453       vv(2)=AEAb1derg(1,2,imat)*b1(2,itk)+AEAb1derg(2,2,imat)*b1(1,itk)
6454       if (l.eq.j+1) then
6455         g_corr6_loc(l-1)=g_corr6_loc(l-1)
6456      & +ekont*(-0.5d0*(scalar2(AEAb1derg(1,2,imat),CUgb2(1,i))
6457      & -scalar2(AEAb2derg(1,1,1,imat),Ug2Db1t(1,k))
6458      & +scalar2(AEAb2derg(1,1,1,imat),CUgb2(1,k))
6459      & +0.5d0*scalar2(vv1(1),Dtobr2(1,i))+scalar2(vv(1),Dtobr2(1,i))))
6460       else
6461         g_corr6_loc(j-1)=g_corr6_loc(j-1)
6462      & +ekont*(-0.5d0*(scalar2(AEAb1derg(1,2,imat),CUgb2(1,i))
6463      & -scalar2(AEAb2derg(1,1,1,imat),Ug2Db1t(1,k))
6464      & +scalar2(AEAb2derg(1,1,1,imat),CUgb2(1,k))
6465      & +0.5d0*scalar2(vv1(1),Dtobr2(1,i))+scalar2(vv(1),Dtobr2(1,i))))
6466       endif
6467       call transpose2(EUgCder(1,1,k),auxmat(1,1))
6468       call matmat2(AEA(1,1,imat),auxmat(1,1),pizda1(1,1))
6469       vv1(1)=pizda1(1,1)-pizda1(2,2)
6470       vv1(2)=pizda1(1,2)+pizda1(2,1)
6471       if (k.gt.1) g_corr6_loc(k-1)=g_corr6_loc(k-1)
6472      & +ekont*(-0.5d0*(-scalar2(AEAb2(1,1,imat),Ug2Db1tder(1,k))
6473      & +scalar2(AEAb2(1,1,imat),CUgb2der(1,k))
6474      & +0.5d0*scalar2(vv1(1),Dtobr2(1,i))))
6475       do iii=1,2
6476         if (swap) then
6477           ind=3-iii
6478         else
6479           ind=iii
6480         endif
6481         do kkk=1,5
6482           do lll=1,3
6483             s1= scalar2(AEAb1derx(1,lll,kkk,iii,2,imat),CUgb2(1,i))
6484             s2=-scalar2(AEAb2derx(1,lll,kkk,iii,1,imat),Ug2Db1t(1,k))
6485             s3= scalar2(AEAb2derx(1,lll,kkk,iii,1,imat),CUgb2(1,k))
6486             call transpose2(EUgC(1,1,k),auxmat(1,1))
6487             call matmat2(AEAderx(1,1,lll,kkk,iii,imat),auxmat(1,1),
6488      &        pizda1(1,1))
6489             vv1(1)=pizda1(1,1)-pizda1(2,2)
6490             vv1(2)=pizda1(1,2)+pizda1(2,1)
6491             s4=0.5d0*scalar2(vv1(1),Dtobr2(1,i))
6492             vv(1)=AEAb1derx(1,lll,kkk,iii,2,imat)*b1(1,itk)
6493      &       -AEAb1derx(2,lll,kkk,iii,2,imat)*b1(2,itk)
6494             vv(2)=AEAb1derx(1,lll,kkk,iii,2,imat)*b1(2,itk)
6495      &       +AEAb1derx(2,lll,kkk,iii,2,imat)*b1(1,itk)
6496             s5=scalar2(vv(1),Dtobr2(1,i))
6497             derx(lll,kkk,ind)=derx(lll,kkk,ind)-0.5d0*(s1+s2+s3+s4+s5)
6498           enddo
6499         enddo
6500       enddo
6501       return
6502       end
6503 c----------------------------------------------------------------------------
6504       double precision function eello6_graph2(i,j,k,l,jj,kk,swap)
6505       implicit real*8 (a-h,o-z)
6506       include 'DIMENSIONS'
6507       include 'DIMENSIONS.ZSCOPT'
6508       include 'COMMON.IOUNITS'
6509       include 'COMMON.CHAIN'
6510       include 'COMMON.DERIV'
6511       include 'COMMON.INTERACT'
6512       include 'COMMON.CONTACTS'
6513       include 'COMMON.TORSION'
6514       include 'COMMON.VAR'
6515       include 'COMMON.GEO'
6516       logical swap
6517       double precision vv(2),pizda(2,2),auxmat(2,2),auxvec(2),
6518      & auxvec1(2),auxvec2(2),auxmat1(2,2)
6519       logical lprn
6520       common /kutas/ lprn
6521 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
6522 C                                                                              C
6523 C      Parallel       Antiparallel                                             C
6524 C                                                                              C
6525 C          o             o                                                     C
6526 C     \   /l\           /j\   /                                                C
6527 C      \ /   \         /   \ /                                                 C
6528 C       o| o |         | o |o                                                  C
6529 C     \ j|/k\|      \  |/k\|l                                                  C
6530 C      \ /   \       \ /   \                                                   C
6531 C       o             o                                                        C
6532 C       i             i                                                        C
6533 C                                                                              C
6534 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
6535 cd      write (2,*) 'eello6_graph2: i,',i,' j',j,' k',k,' l',l
6536 C AL 7/4/01 s1 would occur in the sixth-order moment, 
6537 C           but not in a cluster cumulant
6538 #ifdef MOMENT
6539       s1=dip(1,jj,i)*dip(1,kk,k)
6540 #endif
6541       call matvec2(ADtEA1(1,1,1),Ub2(1,k),auxvec(1))
6542       s2=-0.5d0*scalar2(Ub2(1,i),auxvec(1))
6543       call matvec2(ADtEA(1,1,2),Ub2(1,l),auxvec1(1))
6544       s3=-0.5d0*scalar2(Ub2(1,j),auxvec1(1))
6545       call transpose2(EUg(1,1,k),auxmat(1,1))
6546       call matmat2(ADtEA1(1,1,1),auxmat(1,1),pizda(1,1))
6547       vv(1)=pizda(1,1)-pizda(2,2)
6548       vv(2)=pizda(1,2)+pizda(2,1)
6549       s4=-0.25d0*scalar2(vv(1),Dtobr2(1,i))
6550 cd      write (2,*) 'eello6_graph2:','s1',s1,' s2',s2,' s3',s3,' s4',s4
6551 #ifdef MOMENT
6552       eello6_graph2=-(s1+s2+s3+s4)
6553 #else
6554       eello6_graph2=-(s2+s3+s4)
6555 #endif
6556 c      eello6_graph2=-s3
6557       if (.not. calc_grad) return
6558 C Derivatives in gamma(i-1)
6559       if (i.gt.1) then
6560 #ifdef MOMENT
6561         s1=dipderg(1,jj,i)*dip(1,kk,k)
6562 #endif
6563         s2=-0.5d0*scalar2(Ub2der(1,i),auxvec(1))
6564         call matvec2(ADtEAderg(1,1,1,2),Ub2(1,l),auxvec2(1))
6565         s3=-0.5d0*scalar2(Ub2(1,j),auxvec2(1))
6566         s4=-0.25d0*scalar2(vv(1),Dtobr2der(1,i))
6567 #ifdef MOMENT
6568         g_corr6_loc(i-1)=g_corr6_loc(i-1)-ekont*(s1+s2+s3+s4)
6569 #else
6570         g_corr6_loc(i-1)=g_corr6_loc(i-1)-ekont*(s2+s3+s4)
6571 #endif
6572 c        g_corr6_loc(i-1)=g_corr6_loc(i-1)-s3
6573       endif
6574 C Derivatives in gamma(k-1)
6575 #ifdef MOMENT
6576       s1=dip(1,jj,i)*dipderg(1,kk,k)
6577 #endif
6578       call matvec2(ADtEA1(1,1,1),Ub2der(1,k),auxvec2(1))
6579       s2=-0.5d0*scalar2(Ub2(1,i),auxvec2(1))
6580       call matvec2(ADtEAderg(1,1,2,2),Ub2(1,l),auxvec2(1))
6581       s3=-0.5d0*scalar2(Ub2(1,j),auxvec2(1))
6582       call transpose2(EUgder(1,1,k),auxmat1(1,1))
6583       call matmat2(ADtEA1(1,1,1),auxmat1(1,1),pizda(1,1))
6584       vv(1)=pizda(1,1)-pizda(2,2)
6585       vv(2)=pizda(1,2)+pizda(2,1)
6586       s4=-0.25d0*scalar2(vv(1),Dtobr2(1,i))
6587 #ifdef MOMENT
6588       g_corr6_loc(k-1)=g_corr6_loc(k-1)-ekont*(s1+s2+s3+s4)
6589 #else
6590       g_corr6_loc(k-1)=g_corr6_loc(k-1)-ekont*(s2+s3+s4)
6591 #endif
6592 c      g_corr6_loc(k-1)=g_corr6_loc(k-1)-s3
6593 C Derivatives in gamma(j-1) or gamma(l-1)
6594       if (j.gt.1) then
6595 #ifdef MOMENT
6596         s1=dipderg(3,jj,i)*dip(1,kk,k) 
6597 #endif
6598         call matvec2(ADtEA1derg(1,1,1,1),Ub2(1,k),auxvec2(1))
6599         s2=-0.5d0*scalar2(Ub2(1,i),auxvec2(1))
6600         s3=-0.5d0*scalar2(Ub2der(1,j),auxvec1(1))
6601         call matmat2(ADtEA1derg(1,1,1,1),auxmat(1,1),pizda(1,1))
6602         vv(1)=pizda(1,1)-pizda(2,2)
6603         vv(2)=pizda(1,2)+pizda(2,1)
6604         s4=-0.25d0*scalar2(vv(1),Dtobr2(1,i))
6605 #ifdef MOMENT
6606         if (swap) then
6607           g_corr6_loc(l-1)=g_corr6_loc(l-1)-ekont*s1
6608         else
6609           g_corr6_loc(j-1)=g_corr6_loc(j-1)-ekont*s1
6610         endif
6611 #endif
6612         g_corr6_loc(j-1)=g_corr6_loc(j-1)-ekont*(s2+s3+s4)
6613 c        g_corr6_loc(j-1)=g_corr6_loc(j-1)-s3
6614       endif
6615 C Derivatives in gamma(l-1) or gamma(j-1)
6616       if (l.gt.1) then 
6617 #ifdef MOMENT
6618         s1=dip(1,jj,i)*dipderg(3,kk,k)
6619 #endif
6620         call matvec2(ADtEA1derg(1,1,2,1),Ub2(1,k),auxvec2(1))
6621         s2=-0.5d0*scalar2(Ub2(1,i),auxvec2(1))
6622         call matvec2(ADtEA(1,1,2),Ub2der(1,l),auxvec2(1))
6623         s3=-0.5d0*scalar2(Ub2(1,j),auxvec2(1))
6624         call matmat2(ADtEA1derg(1,1,2,1),auxmat(1,1),pizda(1,1))
6625         vv(1)=pizda(1,1)-pizda(2,2)
6626         vv(2)=pizda(1,2)+pizda(2,1)
6627         s4=-0.25d0*scalar2(vv(1),Dtobr2(1,i))
6628 #ifdef MOMENT
6629         if (swap) then
6630           g_corr6_loc(j-1)=g_corr6_loc(j-1)-ekont*s1
6631         else
6632           g_corr6_loc(l-1)=g_corr6_loc(l-1)-ekont*s1
6633         endif
6634 #endif
6635         g_corr6_loc(l-1)=g_corr6_loc(l-1)-ekont*(s2+s3+s4)
6636 c        g_corr6_loc(l-1)=g_corr6_loc(l-1)-s3
6637       endif
6638 C Cartesian derivatives.
6639       if (lprn) then
6640         write (2,*) 'In eello6_graph2'
6641         do iii=1,2
6642           write (2,*) 'iii=',iii
6643           do kkk=1,5
6644             write (2,*) 'kkk=',kkk
6645             do jjj=1,2
6646               write (2,'(3(2f10.5),5x)') 
6647      &        ((ADtEA1derx(jjj,mmm,lll,kkk,iii,1),mmm=1,2),lll=1,3)
6648             enddo
6649           enddo
6650         enddo
6651       endif
6652       do iii=1,2
6653         do kkk=1,5
6654           do lll=1,3
6655 #ifdef MOMENT
6656             if (iii.eq.1) then
6657               s1=dipderx(lll,kkk,1,jj,i)*dip(1,kk,k)
6658             else
6659               s1=dip(1,jj,i)*dipderx(lll,kkk,1,kk,k)
6660             endif
6661 #endif
6662             call matvec2(ADtEA1derx(1,1,lll,kkk,iii,1),Ub2(1,k),
6663      &        auxvec(1))
6664             s2=-0.5d0*scalar2(Ub2(1,i),auxvec(1))
6665             call matvec2(ADtEAderx(1,1,lll,kkk,iii,2),Ub2(1,l),
6666      &        auxvec(1))
6667             s3=-0.5d0*scalar2(Ub2(1,j),auxvec(1))
6668             call transpose2(EUg(1,1,k),auxmat(1,1))
6669             call matmat2(ADtEA1derx(1,1,lll,kkk,iii,1),auxmat(1,1),
6670      &        pizda(1,1))
6671             vv(1)=pizda(1,1)-pizda(2,2)
6672             vv(2)=pizda(1,2)+pizda(2,1)
6673             s4=-0.25d0*scalar2(vv(1),Dtobr2(1,i))
6674 cd            write (2,*) 's1',s1,' s2',s2,' s3',s3,' s4',s4
6675 #ifdef MOMENT
6676             derx(lll,kkk,iii)=derx(lll,kkk,iii)-(s1+s2+s4)
6677 #else
6678             derx(lll,kkk,iii)=derx(lll,kkk,iii)-(s2+s4)
6679 #endif
6680             if (swap) then
6681               derx(lll,kkk,3-iii)=derx(lll,kkk,3-iii)-s3
6682             else
6683               derx(lll,kkk,iii)=derx(lll,kkk,iii)-s3
6684             endif
6685           enddo
6686         enddo
6687       enddo
6688       return
6689       end
6690 c----------------------------------------------------------------------------
6691       double precision function eello6_graph3(i,j,k,l,jj,kk,swap)
6692       implicit real*8 (a-h,o-z)
6693       include 'DIMENSIONS'
6694       include 'DIMENSIONS.ZSCOPT'
6695       include 'COMMON.IOUNITS'
6696       include 'COMMON.CHAIN'
6697       include 'COMMON.DERIV'
6698       include 'COMMON.INTERACT'
6699       include 'COMMON.CONTACTS'
6700       include 'COMMON.TORSION'
6701       include 'COMMON.VAR'
6702       include 'COMMON.GEO'
6703       double precision vv(2),pizda(2,2),auxmat(2,2),auxvec(2)
6704       logical swap
6705 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
6706 C                                                                              C 
6707 C      Parallel       Antiparallel                                             C
6708 C                                                                              C
6709 C          o             o                                                     C
6710 C         /l\   /   \   /j\                                                    C
6711 C        /   \ /     \ /   \                                                   C
6712 C       /| o |o       o| o |\                                                  C
6713 C       j|/k\|  /      |/k\|l /                                                C
6714 C        /   \ /       /   \ /                                                 C
6715 C       /     o       /     o                                                  C
6716 C       i             i                                                        C
6717 C                                                                              C
6718 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
6719 C
6720 C 4/7/01 AL Component s1 was removed, because it pertains to the respective 
6721 C           energy moment and not to the cluster cumulant.
6722       iti=itortyp(itype(i))
6723       if (j.lt.nres-1 .and. itype(j+1).le.ntyp) then
6724         itj1=itortyp(itype(j+1))
6725       else
6726         itj1=ntortyp+1
6727       endif
6728       itk=itortyp(itype(k))
6729       itk1=itortyp(itype(k+1))
6730       if (l.lt.nres-1 .and. itype(l+1).le.ntyp) then
6731         itl1=itortyp(itype(l+1))
6732       else
6733         itl1=ntortyp+1
6734       endif
6735 #ifdef MOMENT
6736       s1=dip(4,jj,i)*dip(4,kk,k)
6737 #endif
6738       call matvec2(AECA(1,1,1),b1(1,itk1),auxvec(1))
6739       s2=0.5d0*scalar2(b1(1,itk),auxvec(1))
6740       call matvec2(AECA(1,1,2),b1(1,itl1),auxvec(1))
6741       s3=0.5d0*scalar2(b1(1,itj1),auxvec(1))
6742       call transpose2(EE(1,1,itk),auxmat(1,1))
6743       call matmat2(auxmat(1,1),AECA(1,1,1),pizda(1,1))
6744       vv(1)=pizda(1,1)+pizda(2,2)
6745       vv(2)=pizda(2,1)-pizda(1,2)
6746       s4=-0.25d0*scalar2(vv(1),Ctobr(1,k))
6747 cd      write (2,*) 'eello6_graph3:','s1',s1,' s2',s2,' s3',s3,' s4',s4
6748 #ifdef MOMENT
6749       eello6_graph3=-(s1+s2+s3+s4)
6750 #else
6751       eello6_graph3=-(s2+s3+s4)
6752 #endif
6753 c      eello6_graph3=-s4
6754       if (.not. calc_grad) return
6755 C Derivatives in gamma(k-1)
6756       call matvec2(AECAderg(1,1,2),b1(1,itl1),auxvec(1))
6757       s3=0.5d0*scalar2(b1(1,itj1),auxvec(1))
6758       s4=-0.25d0*scalar2(vv(1),Ctobrder(1,k))
6759       g_corr6_loc(k-1)=g_corr6_loc(k-1)-ekont*(s3+s4)
6760 C Derivatives in gamma(l-1)
6761       call matvec2(AECAderg(1,1,1),b1(1,itk1),auxvec(1))
6762       s2=0.5d0*scalar2(b1(1,itk),auxvec(1))
6763       call matmat2(auxmat(1,1),AECAderg(1,1,1),pizda(1,1))
6764       vv(1)=pizda(1,1)+pizda(2,2)
6765       vv(2)=pizda(2,1)-pizda(1,2)
6766       s4=-0.25d0*scalar2(vv(1),Ctobr(1,k))
6767       g_corr6_loc(l-1)=g_corr6_loc(l-1)-ekont*(s2+s4) 
6768 C Cartesian derivatives.
6769       do iii=1,2
6770         do kkk=1,5
6771           do lll=1,3
6772 #ifdef MOMENT
6773             if (iii.eq.1) then
6774               s1=dipderx(lll,kkk,4,jj,i)*dip(4,kk,k)
6775             else
6776               s1=dip(4,jj,i)*dipderx(lll,kkk,4,kk,k)
6777             endif
6778 #endif
6779             call matvec2(AECAderx(1,1,lll,kkk,iii,1),b1(1,itk1),
6780      &        auxvec(1))
6781             s2=0.5d0*scalar2(b1(1,itk),auxvec(1))
6782             call matvec2(AECAderx(1,1,lll,kkk,iii,2),b1(1,itl1),
6783      &        auxvec(1))
6784             s3=0.5d0*scalar2(b1(1,itj1),auxvec(1))
6785             call matmat2(auxmat(1,1),AECAderx(1,1,lll,kkk,iii,1),
6786      &        pizda(1,1))
6787             vv(1)=pizda(1,1)+pizda(2,2)
6788             vv(2)=pizda(2,1)-pizda(1,2)
6789             s4=-0.25d0*scalar2(vv(1),Ctobr(1,k))
6790 #ifdef MOMENT
6791             derx(lll,kkk,iii)=derx(lll,kkk,iii)-(s1+s2+s4)
6792 #else
6793             derx(lll,kkk,iii)=derx(lll,kkk,iii)-(s2+s4)
6794 #endif
6795             if (swap) then
6796               derx(lll,kkk,3-iii)=derx(lll,kkk,3-iii)-s3
6797             else
6798               derx(lll,kkk,iii)=derx(lll,kkk,iii)-s3
6799             endif
6800 c            derx(lll,kkk,iii)=derx(lll,kkk,iii)-s4
6801           enddo
6802         enddo
6803       enddo
6804       return
6805       end
6806 c----------------------------------------------------------------------------
6807       double precision function eello6_graph4(i,j,k,l,jj,kk,imat,swap)
6808       implicit real*8 (a-h,o-z)
6809       include 'DIMENSIONS'
6810       include 'DIMENSIONS.ZSCOPT'
6811       include 'COMMON.IOUNITS'
6812       include 'COMMON.CHAIN'
6813       include 'COMMON.DERIV'
6814       include 'COMMON.INTERACT'
6815       include 'COMMON.CONTACTS'
6816       include 'COMMON.TORSION'
6817       include 'COMMON.VAR'
6818       include 'COMMON.GEO'
6819       include 'COMMON.FFIELD'
6820       double precision vv(2),pizda(2,2),auxmat(2,2),auxvec(2),
6821      & auxvec1(2),auxmat1(2,2)
6822       logical swap
6823 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
6824 C                                                                              C 
6825 C      Parallel       Antiparallel                                             C
6826 C                                                                              C
6827 C          o             o                                                     C
6828 C         /l\   /   \   /j\                                                    C
6829 C        /   \ /     \ /   \                                                   C
6830 C       /| o |o       o| o |\                                                  C
6831 C     \ j|/k\|      \  |/k\|l                                                  C
6832 C      \ /   \       \ /   \                                                   C
6833 C       o     \       o     \                                                  C
6834 C       i             i                                                        C
6835 C                                                                              C
6836 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
6837 C
6838 C 4/7/01 AL Component s1 was removed, because it pertains to the respective 
6839 C           energy moment and not to the cluster cumulant.
6840 cd      write (2,*) 'eello_graph4: wturn6',wturn6
6841       iti=itortyp(itype(i))
6842       itj=itortyp(itype(j))
6843       if (j.lt.nres-1 .and. itype(j+1).le.ntyp) then
6844         itj1=itortyp(itype(j+1))
6845       else
6846         itj1=ntortyp+1
6847       endif
6848       itk=itortyp(itype(k))
6849       if (k.lt.nres-1 .and. itype(k+1).le.ntyp) then
6850         itk1=itortyp(itype(k+1))
6851       else
6852         itk1=ntortyp+1
6853       endif
6854       itl=itortyp(itype(l))
6855       if (l.lt.nres-1) then
6856         itl1=itortyp(itype(l+1))
6857       else
6858         itl1=ntortyp+1
6859       endif
6860 cd      write (2,*) 'eello6_graph4:','i',i,' j',j,' k',k,' l',l
6861 cd      write (2,*) 'iti',iti,' itj',itj,' itj1',itj1,' itk',itk,
6862 cd     & ' itl',itl,' itl1',itl1
6863 #ifdef MOMENT
6864       if (imat.eq.1) then
6865         s1=dip(3,jj,i)*dip(3,kk,k)
6866       else
6867         s1=dip(2,jj,j)*dip(2,kk,l)
6868       endif
6869 #endif
6870       call matvec2(AECA(1,1,imat),Ub2(1,k),auxvec(1))
6871       s2=0.5d0*scalar2(Ub2(1,i),auxvec(1))
6872       if (j.eq.l+1) then
6873         call matvec2(ADtEA1(1,1,3-imat),b1(1,itj1),auxvec1(1))
6874         s3=-0.5d0*scalar2(b1(1,itj),auxvec1(1))
6875       else
6876         call matvec2(ADtEA1(1,1,3-imat),b1(1,itl1),auxvec1(1))
6877         s3=-0.5d0*scalar2(b1(1,itl),auxvec1(1))
6878       endif
6879       call transpose2(EUg(1,1,k),auxmat(1,1))
6880       call matmat2(AECA(1,1,imat),auxmat(1,1),pizda(1,1))
6881       vv(1)=pizda(1,1)-pizda(2,2)
6882       vv(2)=pizda(2,1)+pizda(1,2)
6883       s4=0.25d0*scalar2(vv(1),Dtobr2(1,i))
6884 cd      write (2,*) 'eello6_graph4:','s1',s1,' s2',s2,' s3',s3,' s4',s4
6885 #ifdef MOMENT
6886       eello6_graph4=-(s1+s2+s3+s4)
6887 #else
6888       eello6_graph4=-(s2+s3+s4)
6889 #endif
6890       if (.not. calc_grad) return
6891 C Derivatives in gamma(i-1)
6892       if (i.gt.1) then
6893 #ifdef MOMENT
6894         if (imat.eq.1) then
6895           s1=dipderg(2,jj,i)*dip(3,kk,k)
6896         else
6897           s1=dipderg(4,jj,j)*dip(2,kk,l)
6898         endif
6899 #endif
6900         s2=0.5d0*scalar2(Ub2der(1,i),auxvec(1))
6901         if (j.eq.l+1) then
6902           call matvec2(ADtEA1derg(1,1,1,3-imat),b1(1,itj1),auxvec1(1))
6903           s3=-0.5d0*scalar2(b1(1,itj),auxvec1(1))
6904         else
6905           call matvec2(ADtEA1derg(1,1,1,3-imat),b1(1,itl1),auxvec1(1))
6906           s3=-0.5d0*scalar2(b1(1,itl),auxvec1(1))
6907         endif
6908         s4=0.25d0*scalar2(vv(1),Dtobr2der(1,i))
6909         if (wturn6.gt.0.0d0 .and. k.eq.l+4 .and. i.eq.j+2) then
6910 cd          write (2,*) 'turn6 derivatives'
6911 #ifdef MOMENT
6912           gel_loc_turn6(i-1)=gel_loc_turn6(i-1)-ekont*(s1+s2+s3+s4)
6913 #else
6914           gel_loc_turn6(i-1)=gel_loc_turn6(i-1)-ekont*(s2+s3+s4)
6915 #endif
6916         else
6917 #ifdef MOMENT
6918           g_corr6_loc(i-1)=g_corr6_loc(i-1)-ekont*(s1+s2+s3+s4)
6919 #else
6920           g_corr6_loc(i-1)=g_corr6_loc(i-1)-ekont*(s2+s3+s4)
6921 #endif
6922         endif
6923       endif
6924 C Derivatives in gamma(k-1)
6925 #ifdef MOMENT
6926       if (imat.eq.1) then
6927         s1=dip(3,jj,i)*dipderg(2,kk,k)
6928       else
6929         s1=dip(2,jj,j)*dipderg(4,kk,l)
6930       endif
6931 #endif
6932       call matvec2(AECA(1,1,imat),Ub2der(1,k),auxvec1(1))
6933       s2=0.5d0*scalar2(Ub2(1,i),auxvec1(1))
6934       if (j.eq.l+1) then
6935         call matvec2(ADtEA1derg(1,1,2,3-imat),b1(1,itj1),auxvec1(1))
6936         s3=-0.5d0*scalar2(b1(1,itj),auxvec1(1))
6937       else
6938         call matvec2(ADtEA1derg(1,1,2,3-imat),b1(1,itl1),auxvec1(1))
6939         s3=-0.5d0*scalar2(b1(1,itl),auxvec1(1))
6940       endif
6941       call transpose2(EUgder(1,1,k),auxmat1(1,1))
6942       call matmat2(AECA(1,1,imat),auxmat1(1,1),pizda(1,1))
6943       vv(1)=pizda(1,1)-pizda(2,2)
6944       vv(2)=pizda(2,1)+pizda(1,2)
6945       s4=0.25d0*scalar2(vv(1),Dtobr2(1,i))
6946       if (wturn6.gt.0.0d0 .and. k.eq.l+4 .and. i.eq.j+2) then
6947 #ifdef MOMENT
6948         gel_loc_turn6(k-1)=gel_loc_turn6(k-1)-ekont*(s1+s2+s3+s4)
6949 #else
6950         gel_loc_turn6(k-1)=gel_loc_turn6(k-1)-ekont*(s2+s3+s4)
6951 #endif
6952       else
6953 #ifdef MOMENT
6954         g_corr6_loc(k-1)=g_corr6_loc(k-1)-ekont*(s1+s2+s3+s4)
6955 #else
6956         g_corr6_loc(k-1)=g_corr6_loc(k-1)-ekont*(s2+s3+s4)
6957 #endif
6958       endif
6959 C Derivatives in gamma(j-1) or gamma(l-1)
6960       if (l.eq.j+1 .and. l.gt.1) then
6961         call matvec2(AECAderg(1,1,imat),Ub2(1,k),auxvec(1))
6962         s2=0.5d0*scalar2(Ub2(1,i),auxvec(1))
6963         call matmat2(AECAderg(1,1,imat),auxmat(1,1),pizda(1,1))
6964         vv(1)=pizda(1,1)-pizda(2,2)
6965         vv(2)=pizda(2,1)+pizda(1,2)
6966         s4=0.25d0*scalar2(vv(1),Dtobr2(1,i))
6967         g_corr6_loc(l-1)=g_corr6_loc(l-1)-ekont*(s2+s4)
6968       else if (j.gt.1) then
6969         call matvec2(AECAderg(1,1,imat),Ub2(1,k),auxvec(1))
6970         s2=0.5d0*scalar2(Ub2(1,i),auxvec(1))
6971         call matmat2(AECAderg(1,1,imat),auxmat(1,1),pizda(1,1))
6972         vv(1)=pizda(1,1)-pizda(2,2)
6973         vv(2)=pizda(2,1)+pizda(1,2)
6974         s4=0.25d0*scalar2(vv(1),Dtobr2(1,i))
6975         if (wturn6.gt.0.0d0 .and. k.eq.l+4 .and. i.eq.j+2) then
6976           gel_loc_turn6(j-1)=gel_loc_turn6(j-1)-ekont*(s2+s4)
6977         else
6978           g_corr6_loc(j-1)=g_corr6_loc(j-1)-ekont*(s2+s4)
6979         endif
6980       endif
6981 C Cartesian derivatives.
6982       do iii=1,2
6983         do kkk=1,5
6984           do lll=1,3
6985 #ifdef MOMENT
6986             if (iii.eq.1) then
6987               if (imat.eq.1) then
6988                 s1=dipderx(lll,kkk,3,jj,i)*dip(3,kk,k)
6989               else
6990                 s1=dipderx(lll,kkk,2,jj,j)*dip(2,kk,l)
6991               endif
6992             else
6993               if (imat.eq.1) then
6994                 s1=dip(3,jj,i)*dipderx(lll,kkk,3,kk,k)
6995               else
6996                 s1=dip(2,jj,j)*dipderx(lll,kkk,2,kk,l)
6997               endif
6998             endif
6999 #endif
7000             call matvec2(AECAderx(1,1,lll,kkk,iii,imat),Ub2(1,k),
7001      &        auxvec(1))
7002             s2=0.5d0*scalar2(Ub2(1,i),auxvec(1))
7003             if (j.eq.l+1) then
7004               call matvec2(ADtEA1derx(1,1,lll,kkk,iii,3-imat),
7005      &          b1(1,itj1),auxvec(1))
7006               s3=-0.5d0*scalar2(b1(1,itj),auxvec(1))
7007             else
7008               call matvec2(ADtEA1derx(1,1,lll,kkk,iii,3-imat),
7009      &          b1(1,itl1),auxvec(1))
7010               s3=-0.5d0*scalar2(b1(1,itl),auxvec(1))
7011             endif
7012             call matmat2(AECAderx(1,1,lll,kkk,iii,imat),auxmat(1,1),
7013      &        pizda(1,1))
7014             vv(1)=pizda(1,1)-pizda(2,2)
7015             vv(2)=pizda(2,1)+pizda(1,2)
7016             s4=0.25d0*scalar2(vv(1),Dtobr2(1,i))
7017             if (swap) then
7018               if (wturn6.gt.0.0d0 .and. k.eq.l+4 .and. i.eq.j+2) then
7019 #ifdef MOMENT
7020                 derx_turn(lll,kkk,3-iii)=derx_turn(lll,kkk,3-iii)
7021      &             -(s1+s2+s4)
7022 #else
7023                 derx_turn(lll,kkk,3-iii)=derx_turn(lll,kkk,3-iii)
7024      &             -(s2+s4)
7025 #endif
7026                 derx_turn(lll,kkk,iii)=derx_turn(lll,kkk,iii)-s3
7027               else
7028 #ifdef MOMENT
7029                 derx(lll,kkk,3-iii)=derx(lll,kkk,3-iii)-(s1+s2+s4)
7030 #else
7031                 derx(lll,kkk,3-iii)=derx(lll,kkk,3-iii)-(s2+s4)
7032 #endif
7033                 derx(lll,kkk,iii)=derx(lll,kkk,iii)-s3
7034               endif
7035             else
7036 #ifdef MOMENT
7037               derx(lll,kkk,iii)=derx(lll,kkk,iii)-(s1+s2+s4)
7038 #else
7039               derx(lll,kkk,iii)=derx(lll,kkk,iii)-(s2+s4)
7040 #endif
7041               if (l.eq.j+1) then
7042                 derx(lll,kkk,iii)=derx(lll,kkk,iii)-s3
7043               else 
7044                 derx(lll,kkk,3-iii)=derx(lll,kkk,3-iii)-s3
7045               endif
7046             endif 
7047           enddo
7048         enddo
7049       enddo
7050       return
7051       end
7052 c----------------------------------------------------------------------------
7053       double precision function eello_turn6(i,jj,kk)
7054       implicit real*8 (a-h,o-z)
7055       include 'DIMENSIONS'
7056       include 'DIMENSIONS.ZSCOPT'
7057       include 'COMMON.IOUNITS'
7058       include 'COMMON.CHAIN'
7059       include 'COMMON.DERIV'
7060       include 'COMMON.INTERACT'
7061       include 'COMMON.CONTACTS'
7062       include 'COMMON.TORSION'
7063       include 'COMMON.VAR'
7064       include 'COMMON.GEO'
7065       double precision vtemp1(2),vtemp2(2),vtemp3(2),vtemp4(2),
7066      &  atemp(2,2),auxmat(2,2),achuj_temp(2,2),gtemp(2,2),gvec(2),
7067      &  ggg1(3),ggg2(3)
7068       double precision vtemp1d(2),vtemp2d(2),vtemp3d(2),vtemp4d(2),
7069      &  atempd(2,2),auxmatd(2,2),achuj_tempd(2,2),gtempd(2,2),gvecd(2)
7070 C 4/7/01 AL Components s1, s8, and s13 were removed, because they pertain to
7071 C           the respective energy moment and not to the cluster cumulant.
7072       eello_turn6=0.0d0
7073       j=i+4
7074       k=i+1
7075       l=i+3
7076       iti=itortyp(itype(i))
7077       itk=itortyp(itype(k))
7078       itk1=itortyp(itype(k+1))
7079       itl=itortyp(itype(l))
7080       itj=itortyp(itype(j))
7081 cd      write (2,*) 'itk',itk,' itk1',itk1,' itl',itl,' itj',itj
7082 cd      write (2,*) 'i',i,' k',k,' j',j,' l',l
7083 cd      if (i.ne.1 .or. j.ne.3 .or. k.ne.2 .or. l.ne.4) then
7084 cd        eello6=0.0d0
7085 cd        return
7086 cd      endif
7087 cd      write (iout,*)
7088 cd     &   'EELLO6: Contacts have occurred for peptide groups',i,j,
7089 cd     &   ' and',k,l
7090 cd      call checkint_turn6(i,jj,kk,eel_turn6_num)
7091       do iii=1,2
7092         do kkk=1,5
7093           do lll=1,3
7094             derx_turn(lll,kkk,iii)=0.0d0
7095           enddo
7096         enddo
7097       enddo
7098 cd      eij=1.0d0
7099 cd      ekl=1.0d0
7100 cd      ekont=1.0d0
7101       eello6_5=eello6_graph4(l,k,j,i,kk,jj,2,.true.)
7102 cd      eello6_5=0.0d0
7103 cd      write (2,*) 'eello6_5',eello6_5
7104 #ifdef MOMENT
7105       call transpose2(AEA(1,1,1),auxmat(1,1))
7106       call matmat2(EUg(1,1,i+1),auxmat(1,1),auxmat(1,1))
7107       ss1=scalar2(Ub2(1,i+2),b1(1,itl))
7108       s1 = (auxmat(1,1)+auxmat(2,2))*ss1
7109 #else
7110       s1 = 0.0d0
7111 #endif
7112       call matvec2(EUg(1,1,i+2),b1(1,itl),vtemp1(1))
7113       call matvec2(AEA(1,1,1),vtemp1(1),vtemp1(1))
7114       s2 = scalar2(b1(1,itk),vtemp1(1))
7115 #ifdef MOMENT
7116       call transpose2(AEA(1,1,2),atemp(1,1))
7117       call matmat2(atemp(1,1),EUg(1,1,i+4),atemp(1,1))
7118       call matvec2(Ug2(1,1,i+2),dd(1,1,itk1),vtemp2(1))
7119       s8 = -(atemp(1,1)+atemp(2,2))*scalar2(cc(1,1,itl),vtemp2(1))
7120 #else
7121       s8=0.0d0
7122 #endif
7123       call matmat2(EUg(1,1,i+3),AEA(1,1,2),auxmat(1,1))
7124       call matvec2(auxmat(1,1),Ub2(1,i+4),vtemp3(1))
7125       s12 = scalar2(Ub2(1,i+2),vtemp3(1))
7126 #ifdef MOMENT
7127       call transpose2(a_chuj(1,1,kk,i+1),achuj_temp(1,1))
7128       call matmat2(achuj_temp(1,1),EUg(1,1,i+2),gtemp(1,1))
7129       call matmat2(gtemp(1,1),EUg(1,1,i+3),gtemp(1,1)) 
7130       call matvec2(a_chuj(1,1,jj,i),Ub2(1,i+4),vtemp4(1)) 
7131       ss13 = scalar2(b1(1,itk),vtemp4(1))
7132       s13 = (gtemp(1,1)+gtemp(2,2))*ss13
7133 #else
7134       s13=0.0d0
7135 #endif
7136 c      write (2,*) 's1,s2,s8,s12,s13',s1,s2,s8,s12,s13
7137 c      s1=0.0d0
7138 c      s2=0.0d0
7139 c      s8=0.0d0
7140 c      s12=0.0d0
7141 c      s13=0.0d0
7142       eel_turn6 = eello6_5 - 0.5d0*(s1+s2+s12+s8+s13)
7143       if (calc_grad) then
7144 C Derivatives in gamma(i+2)
7145 #ifdef MOMENT
7146       call transpose2(AEA(1,1,1),auxmatd(1,1))
7147       call matmat2(EUgder(1,1,i+1),auxmatd(1,1),auxmatd(1,1))
7148       s1d = (auxmatd(1,1)+auxmatd(2,2))*ss1
7149       call transpose2(AEAderg(1,1,2),atempd(1,1))
7150       call matmat2(atempd(1,1),EUg(1,1,i+4),atempd(1,1))
7151       s8d = -(atempd(1,1)+atempd(2,2))*scalar2(cc(1,1,itl),vtemp2(1))
7152 #else
7153       s8d=0.0d0
7154 #endif
7155       call matmat2(EUg(1,1,i+3),AEAderg(1,1,2),auxmatd(1,1))
7156       call matvec2(auxmatd(1,1),Ub2(1,i+4),vtemp3d(1))
7157       s12d = scalar2(Ub2(1,i+2),vtemp3d(1))
7158 c      s1d=0.0d0
7159 c      s2d=0.0d0
7160 c      s8d=0.0d0
7161 c      s12d=0.0d0
7162 c      s13d=0.0d0
7163       gel_loc_turn6(i)=gel_loc_turn6(i)-0.5d0*ekont*(s1d+s8d+s12d)
7164 C Derivatives in gamma(i+3)
7165 #ifdef MOMENT
7166       call transpose2(AEA(1,1,1),auxmatd(1,1))
7167       call matmat2(EUg(1,1,i+1),auxmatd(1,1),auxmatd(1,1))
7168       ss1d=scalar2(Ub2der(1,i+2),b1(1,itl))
7169       s1d = (auxmatd(1,1)+auxmatd(2,2))*ss1d
7170 #else
7171       s1d=0.0d0
7172 #endif
7173       call matvec2(EUgder(1,1,i+2),b1(1,itl),vtemp1d(1))
7174       call matvec2(AEA(1,1,1),vtemp1d(1),vtemp1d(1))
7175       s2d = scalar2(b1(1,itk),vtemp1d(1))
7176 #ifdef MOMENT
7177       call matvec2(Ug2der(1,1,i+2),dd(1,1,itk1),vtemp2d(1))
7178       s8d = -(atemp(1,1)+atemp(2,2))*scalar2(cc(1,1,itl),vtemp2d(1))
7179 #endif
7180       s12d = scalar2(Ub2der(1,i+2),vtemp3(1))
7181 #ifdef MOMENT
7182       call matmat2(achuj_temp(1,1),EUgder(1,1,i+2),gtempd(1,1))
7183       call matmat2(gtempd(1,1),EUg(1,1,i+3),gtempd(1,1)) 
7184       s13d = (gtempd(1,1)+gtempd(2,2))*ss13
7185 #else
7186       s13d=0.0d0
7187 #endif
7188 c      s1d=0.0d0
7189 c      s2d=0.0d0
7190 c      s8d=0.0d0
7191 c      s12d=0.0d0
7192 c      s13d=0.0d0
7193 #ifdef MOMENT
7194       gel_loc_turn6(i+1)=gel_loc_turn6(i+1)
7195      &               -0.5d0*ekont*(s1d+s2d+s8d+s12d+s13d)
7196 #else
7197       gel_loc_turn6(i+1)=gel_loc_turn6(i+1)
7198      &               -0.5d0*ekont*(s2d+s12d)
7199 #endif
7200 C Derivatives in gamma(i+4)
7201       call matmat2(EUgder(1,1,i+3),AEA(1,1,2),auxmatd(1,1))
7202       call matvec2(auxmatd(1,1),Ub2(1,i+4),vtemp3d(1))
7203       s12d = scalar2(Ub2(1,i+2),vtemp3d(1))
7204 #ifdef MOMENT
7205       call matmat2(achuj_temp(1,1),EUg(1,1,i+2),gtempd(1,1))
7206       call matmat2(gtempd(1,1),EUgder(1,1,i+3),gtempd(1,1)) 
7207       s13d = (gtempd(1,1)+gtempd(2,2))*ss13
7208 #else
7209       s13d = 0.0d0
7210 #endif
7211 c      s1d=0.0d0
7212 c      s2d=0.0d0
7213 c      s8d=0.0d0
7214 C      s12d=0.0d0
7215 c      s13d=0.0d0
7216 #ifdef MOMENT
7217       gel_loc_turn6(i+2)=gel_loc_turn6(i+2)-0.5d0*ekont*(s12d+s13d)
7218 #else
7219       gel_loc_turn6(i+2)=gel_loc_turn6(i+2)-0.5d0*ekont*(s12d)
7220 #endif
7221 C Derivatives in gamma(i+5)
7222 #ifdef MOMENT
7223       call transpose2(AEAderg(1,1,1),auxmatd(1,1))
7224       call matmat2(EUg(1,1,i+1),auxmatd(1,1),auxmatd(1,1))
7225       s1d = (auxmatd(1,1)+auxmatd(2,2))*ss1
7226 #else
7227       s1d = 0.0d0
7228 #endif
7229       call matvec2(EUg(1,1,i+2),b1(1,itl),vtemp1d(1))
7230       call matvec2(AEAderg(1,1,1),vtemp1d(1),vtemp1d(1))
7231       s2d = scalar2(b1(1,itk),vtemp1d(1))
7232 #ifdef MOMENT
7233       call transpose2(AEA(1,1,2),atempd(1,1))
7234       call matmat2(atempd(1,1),EUgder(1,1,i+4),atempd(1,1))
7235       s8d = -(atempd(1,1)+atempd(2,2))*scalar2(cc(1,1,itl),vtemp2(1))
7236 #else
7237       s8d = 0.0d0
7238 #endif
7239       call matvec2(auxmat(1,1),Ub2der(1,i+4),vtemp3d(1))
7240       s12d = scalar2(Ub2(1,i+2),vtemp3d(1))
7241 #ifdef MOMENT
7242       call matvec2(a_chuj(1,1,jj,i),Ub2der(1,i+4),vtemp4d(1)) 
7243       ss13d = scalar2(b1(1,itk),vtemp4d(1))
7244       s13d = (gtemp(1,1)+gtemp(2,2))*ss13d
7245 #else
7246       s13d = 0.0d0
7247 #endif
7248 c      s1d=0.0d0
7249 c      s2d=0.0d0
7250 c      s8d=0.0d0
7251 c      s12d=0.0d0
7252 c      s13d=0.0d0
7253 #ifdef MOMENT
7254       gel_loc_turn6(i+3)=gel_loc_turn6(i+3)
7255      &               -0.5d0*ekont*(s1d+s2d+s8d+s12d+s13d)
7256 #else
7257       gel_loc_turn6(i+3)=gel_loc_turn6(i+3)
7258      &               -0.5d0*ekont*(s2d+s12d)
7259 #endif
7260 C Cartesian derivatives
7261       do iii=1,2
7262         do kkk=1,5
7263           do lll=1,3
7264 #ifdef MOMENT
7265             call transpose2(AEAderx(1,1,lll,kkk,iii,1),auxmatd(1,1))
7266             call matmat2(EUg(1,1,i+1),auxmatd(1,1),auxmatd(1,1))
7267             s1d = (auxmatd(1,1)+auxmatd(2,2))*ss1
7268 #else
7269             s1d = 0.0d0
7270 #endif
7271             call matvec2(EUg(1,1,i+2),b1(1,itl),vtemp1(1))
7272             call matvec2(AEAderx(1,1,lll,kkk,iii,1),vtemp1(1),
7273      &          vtemp1d(1))
7274             s2d = scalar2(b1(1,itk),vtemp1d(1))
7275 #ifdef MOMENT
7276             call transpose2(AEAderx(1,1,lll,kkk,iii,2),atempd(1,1))
7277             call matmat2(atempd(1,1),EUg(1,1,i+4),atempd(1,1))
7278             s8d = -(atempd(1,1)+atempd(2,2))*
7279      &           scalar2(cc(1,1,itl),vtemp2(1))
7280 #else
7281             s8d = 0.0d0
7282 #endif
7283             call matmat2(EUg(1,1,i+3),AEAderx(1,1,lll,kkk,iii,2),
7284      &           auxmatd(1,1))
7285             call matvec2(auxmatd(1,1),Ub2(1,i+4),vtemp3d(1))
7286             s12d = scalar2(Ub2(1,i+2),vtemp3d(1))
7287 c      s1d=0.0d0
7288 c      s2d=0.0d0
7289 c      s8d=0.0d0
7290 c      s12d=0.0d0
7291 c      s13d=0.0d0
7292 #ifdef MOMENT
7293             derx_turn(lll,kkk,iii) = derx_turn(lll,kkk,iii) 
7294      &        - 0.5d0*(s1d+s2d)
7295 #else
7296             derx_turn(lll,kkk,iii) = derx_turn(lll,kkk,iii) 
7297      &        - 0.5d0*s2d
7298 #endif
7299 #ifdef MOMENT
7300             derx_turn(lll,kkk,3-iii) = derx_turn(lll,kkk,3-iii) 
7301      &        - 0.5d0*(s8d+s12d)
7302 #else
7303             derx_turn(lll,kkk,3-iii) = derx_turn(lll,kkk,3-iii) 
7304      &        - 0.5d0*s12d
7305 #endif
7306           enddo
7307         enddo
7308       enddo
7309 #ifdef MOMENT
7310       do kkk=1,5
7311         do lll=1,3
7312           call transpose2(a_chuj_der(1,1,lll,kkk,kk,i+1),
7313      &      achuj_tempd(1,1))
7314           call matmat2(achuj_tempd(1,1),EUg(1,1,i+2),gtempd(1,1))
7315           call matmat2(gtempd(1,1),EUg(1,1,i+3),gtempd(1,1)) 
7316           s13d=(gtempd(1,1)+gtempd(2,2))*ss13
7317           derx_turn(lll,kkk,2) = derx_turn(lll,kkk,2)-0.5d0*s13d
7318           call matvec2(a_chuj_der(1,1,lll,kkk,jj,i),Ub2(1,i+4),
7319      &      vtemp4d(1)) 
7320           ss13d = scalar2(b1(1,itk),vtemp4d(1))
7321           s13d = (gtemp(1,1)+gtemp(2,2))*ss13d
7322           derx_turn(lll,kkk,1) = derx_turn(lll,kkk,1)-0.5d0*s13d
7323         enddo
7324       enddo
7325 #endif
7326 cd      write(iout,*) 'eel6_turn6',eel_turn6,' eel_turn6_num',
7327 cd     &  16*eel_turn6_num
7328 cd      goto 1112
7329       if (j.lt.nres-1) then
7330         j1=j+1
7331         j2=j-1
7332       else
7333         j1=j-1
7334         j2=j-2
7335       endif
7336       if (l.lt.nres-1) then
7337         l1=l+1
7338         l2=l-1
7339       else
7340         l1=l-1
7341         l2=l-2
7342       endif
7343       do ll=1,3
7344         ggg1(ll)=eel_turn6*g_contij(ll,1)
7345         ggg2(ll)=eel_turn6*g_contij(ll,2)
7346         ghalf=0.5d0*ggg1(ll)
7347 cd        ghalf=0.0d0
7348         gcorr6_turn(ll,i)=gcorr6_turn(ll,i)+ghalf
7349      &    +ekont*derx_turn(ll,2,1)
7350         gcorr6_turn(ll,i+1)=gcorr6_turn(ll,i+1)+ekont*derx_turn(ll,3,1)
7351         gcorr6_turn(ll,j)=gcorr6_turn(ll,j)+ghalf
7352      &    +ekont*derx_turn(ll,4,1)
7353         gcorr6_turn(ll,j1)=gcorr6_turn(ll,j1)+ekont*derx_turn(ll,5,1)
7354         ghalf=0.5d0*ggg2(ll)
7355 cd        ghalf=0.0d0
7356         gcorr6_turn(ll,k)=gcorr6_turn(ll,k)+ghalf
7357      &    +ekont*derx_turn(ll,2,2)
7358         gcorr6_turn(ll,k+1)=gcorr6_turn(ll,k+1)+ekont*derx_turn(ll,3,2)
7359         gcorr6_turn(ll,l)=gcorr6_turn(ll,l)+ghalf
7360      &    +ekont*derx_turn(ll,4,2)
7361         gcorr6_turn(ll,l1)=gcorr6_turn(ll,l1)+ekont*derx_turn(ll,5,2)
7362       enddo
7363 cd      goto 1112
7364       do m=i+1,j-1
7365         do ll=1,3
7366           gcorr6_turn(ll,m)=gcorr6_turn(ll,m)+ggg1(ll)
7367         enddo
7368       enddo
7369       do m=k+1,l-1
7370         do ll=1,3
7371           gcorr6_turn(ll,m)=gcorr6_turn(ll,m)+ggg2(ll)
7372         enddo
7373       enddo
7374 1112  continue
7375       do m=i+2,j2
7376         do ll=1,3
7377           gcorr6_turn(ll,m)=gcorr6_turn(ll,m)+ekont*derx_turn(ll,1,1)
7378         enddo
7379       enddo
7380       do m=k+2,l2
7381         do ll=1,3
7382           gcorr6_turn(ll,m)=gcorr6_turn(ll,m)+ekont*derx_turn(ll,1,2)
7383         enddo
7384       enddo 
7385 cd      do iii=1,nres-3
7386 cd        write (2,*) iii,g_corr6_loc(iii)
7387 cd      enddo
7388       endif
7389       eello_turn6=ekont*eel_turn6
7390 cd      write (2,*) 'ekont',ekont
7391 cd      write (2,*) 'eel_turn6',ekont*eel_turn6
7392       return
7393       end
7394 crc-------------------------------------------------
7395       SUBROUTINE MATVEC2(A1,V1,V2)
7396       implicit real*8 (a-h,o-z)
7397       include 'DIMENSIONS'
7398       DIMENSION A1(2,2),V1(2),V2(2)
7399 c      DO 1 I=1,2
7400 c        VI=0.0
7401 c        DO 3 K=1,2
7402 c    3     VI=VI+A1(I,K)*V1(K)
7403 c        Vaux(I)=VI
7404 c    1 CONTINUE
7405
7406       vaux1=a1(1,1)*v1(1)+a1(1,2)*v1(2)
7407       vaux2=a1(2,1)*v1(1)+a1(2,2)*v1(2)
7408
7409       v2(1)=vaux1
7410       v2(2)=vaux2
7411       END
7412 C---------------------------------------
7413       SUBROUTINE MATMAT2(A1,A2,A3)
7414       implicit real*8 (a-h,o-z)
7415       include 'DIMENSIONS'
7416       DIMENSION A1(2,2),A2(2,2),A3(2,2)
7417 c      DIMENSION AI3(2,2)
7418 c        DO  J=1,2
7419 c          A3IJ=0.0
7420 c          DO K=1,2
7421 c           A3IJ=A3IJ+A1(I,K)*A2(K,J)
7422 c          enddo
7423 c          A3(I,J)=A3IJ
7424 c       enddo
7425 c      enddo
7426
7427       ai3_11=a1(1,1)*a2(1,1)+a1(1,2)*a2(2,1)
7428       ai3_12=a1(1,1)*a2(1,2)+a1(1,2)*a2(2,2)
7429       ai3_21=a1(2,1)*a2(1,1)+a1(2,2)*a2(2,1)
7430       ai3_22=a1(2,1)*a2(1,2)+a1(2,2)*a2(2,2)
7431
7432       A3(1,1)=AI3_11
7433       A3(2,1)=AI3_21
7434       A3(1,2)=AI3_12
7435       A3(2,2)=AI3_22
7436       END
7437
7438 c-------------------------------------------------------------------------
7439       double precision function scalar2(u,v)
7440       implicit none
7441       double precision u(2),v(2)
7442       double precision sc
7443       integer i
7444       scalar2=u(1)*v(1)+u(2)*v(2)
7445       return
7446       end
7447
7448 C-----------------------------------------------------------------------------
7449
7450       subroutine transpose2(a,at)
7451       implicit none
7452       double precision a(2,2),at(2,2)
7453       at(1,1)=a(1,1)
7454       at(1,2)=a(2,1)
7455       at(2,1)=a(1,2)
7456       at(2,2)=a(2,2)
7457       return
7458       end
7459 c--------------------------------------------------------------------------
7460       subroutine transpose(n,a,at)
7461       implicit none
7462       integer n,i,j
7463       double precision a(n,n),at(n,n)
7464       do i=1,n
7465         do j=1,n
7466           at(j,i)=a(i,j)
7467         enddo
7468       enddo
7469       return
7470       end
7471 C---------------------------------------------------------------------------
7472       subroutine prodmat3(a1,a2,kk,transp,prod)
7473       implicit none
7474       integer i,j
7475       double precision a1(2,2),a2(2,2),a2t(2,2),kk(2,2),prod(2,2)
7476       logical transp
7477 crc      double precision auxmat(2,2),prod_(2,2)
7478
7479       if (transp) then
7480 crc        call transpose2(kk(1,1),auxmat(1,1))
7481 crc        call matmat2(a1(1,1),auxmat(1,1),auxmat(1,1))
7482 crc        call matmat2(auxmat(1,1),a2(1,1),prod_(1,1)) 
7483         
7484            prod(1,1)=(a1(1,1)*kk(1,1)+a1(1,2)*kk(1,2))*a2(1,1)
7485      & +(a1(1,1)*kk(2,1)+a1(1,2)*kk(2,2))*a2(2,1)
7486            prod(1,2)=(a1(1,1)*kk(1,1)+a1(1,2)*kk(1,2))*a2(1,2)
7487      & +(a1(1,1)*kk(2,1)+a1(1,2)*kk(2,2))*a2(2,2)
7488            prod(2,1)=(a1(2,1)*kk(1,1)+a1(2,2)*kk(1,2))*a2(1,1)
7489      & +(a1(2,1)*kk(2,1)+a1(2,2)*kk(2,2))*a2(2,1)
7490            prod(2,2)=(a1(2,1)*kk(1,1)+a1(2,2)*kk(1,2))*a2(1,2)
7491      & +(a1(2,1)*kk(2,1)+a1(2,2)*kk(2,2))*a2(2,2)
7492
7493       else
7494 crc        call matmat2(a1(1,1),kk(1,1),auxmat(1,1))
7495 crc        call matmat2(auxmat(1,1),a2(1,1),prod_(1,1))
7496
7497            prod(1,1)=(a1(1,1)*kk(1,1)+a1(1,2)*kk(2,1))*a2(1,1)
7498      &  +(a1(1,1)*kk(1,2)+a1(1,2)*kk(2,2))*a2(2,1)
7499            prod(1,2)=(a1(1,1)*kk(1,1)+a1(1,2)*kk(2,1))*a2(1,2)
7500      &  +(a1(1,1)*kk(1,2)+a1(1,2)*kk(2,2))*a2(2,2)
7501            prod(2,1)=(a1(2,1)*kk(1,1)+a1(2,2)*kk(2,1))*a2(1,1)
7502      &  +(a1(2,1)*kk(1,2)+a1(2,2)*kk(2,2))*a2(2,1)
7503            prod(2,2)=(a1(2,1)*kk(1,1)+a1(2,2)*kk(2,1))*a2(1,2)
7504      &  +(a1(2,1)*kk(1,2)+a1(2,2)*kk(2,2))*a2(2,2)
7505
7506       endif
7507 c      call transpose2(a2(1,1),a2t(1,1))
7508
7509 crc      print *,transp
7510 crc      print *,((prod_(i,j),i=1,2),j=1,2)
7511 crc      print *,((prod(i,j),i=1,2),j=1,2)
7512
7513       return
7514       end
7515 C-----------------------------------------------------------------------------
7516       double precision function scalar(u,v)
7517       implicit none
7518       double precision u(3),v(3)
7519       double precision sc
7520       integer i
7521       sc=0.0d0
7522       do i=1,3
7523         sc=sc+u(i)*v(i)
7524       enddo
7525       scalar=sc
7526       return
7527       end
7528