bugfix in shield FGPROC>1
[unres.git] / source / unres / src_MD-M / energy_split-sep.F
1       subroutine etotal_long(energia)
2       implicit real*8 (a-h,o-z)
3       include 'DIMENSIONS'
4 c
5 c Compute the long-range slow-varying contributions to the energy
6 c
7 #ifndef ISNAN
8       external proc_proc
9 #ifdef WINPGI
10 cMS$ATTRIBUTES C ::  proc_proc
11 #endif
12 #endif
13 #ifdef MPI
14       include "mpif.h"
15       double precision weights_(n_ene)
16 #endif
17       include 'COMMON.SETUP'
18       include 'COMMON.IOUNITS'
19       double precision energia(0:n_ene)
20       include 'COMMON.FFIELD'
21       include 'COMMON.DERIV'
22       include 'COMMON.INTERACT'
23       include 'COMMON.SBRIDGE'
24       include 'COMMON.CHAIN'
25       include 'COMMON.VAR'
26       include 'COMMON.LOCAL'
27       include 'COMMON.MD'
28       include 'COMMON.CONTROL'
29       include 'COMMON.SHIELD'
30 c      write(iout,'(a,i2)')'Calling etotal_long ipot=',ipot
31       if (modecalc.eq.12.or.modecalc.eq.14) then
32 #ifdef MPI
33 c        if (fg_rank.eq.0) call int_from_cart1(.false.)
34 #else
35         call int_from_cart1(.false.)
36 #endif
37       endif
38 #ifdef MPI      
39 c      write(iout,*) "ETOTAL_LONG Processor",fg_rank,
40 c     & " absolute rank",myrank," nfgtasks",nfgtasks
41       call flush(iout)
42       if (nfgtasks.gt.1) then
43         time00=MPI_Wtime()
44 C FG slaves call the following matching MPI_Bcast in ERGASTULUM
45         if (fg_rank.eq.0) then
46           call MPI_Bcast(3,1,MPI_INTEGER,king,FG_COMM,IERROR)
47 c          write (iout,*) "Processor",myrank," BROADCAST iorder"
48 c          call flush(iout)
49 C FG master sets up the WEIGHTS_ array which will be broadcast to the 
50 C FG slaves as WEIGHTS array.
51           weights_(1)=wsc
52           weights_(2)=wscp
53           weights_(3)=welec
54           weights_(4)=wcorr
55           weights_(5)=wcorr5
56           weights_(6)=wcorr6
57           weights_(7)=wel_loc
58           weights_(8)=wturn3
59           weights_(9)=wturn4
60           weights_(10)=wturn6
61           weights_(11)=wang
62           weights_(12)=wscloc
63           weights_(13)=wtor
64           weights_(14)=wtor_d
65           weights_(15)=wstrain
66           weights_(16)=wvdwpp
67           weights_(17)=wbond
68           weights_(18)=scal14
69           weights_(21)=wsccor
70 C FG Master broadcasts the WEIGHTS_ array
71           call MPI_Bcast(weights_(1),n_ene,
72      &        MPI_DOUBLE_PRECISION,king,FG_COMM,IERROR)
73         else
74 C FG slaves receive the WEIGHTS array
75           call MPI_Bcast(weights(1),n_ene,
76      &        MPI_DOUBLE_PRECISION,king,FG_COMM,IERROR)
77           wsc=weights(1)
78           wscp=weights(2)
79           welec=weights(3)
80           wcorr=weights(4)
81           wcorr5=weights(5)
82           wcorr6=weights(6)
83           wel_loc=weights(7)
84           wturn3=weights(8)
85           wturn4=weights(9)
86           wturn6=weights(10)
87           wang=weights(11)
88           wscloc=weights(12)
89           wtor=weights(13)
90           wtor_d=weights(14)
91           wstrain=weights(15)
92           wvdwpp=weights(16)
93           wbond=weights(17)
94           scal14=weights(18)
95           wsccor=weights(21)
96         endif
97         call MPI_Bcast(dc(1,1),6*nres,MPI_DOUBLE_PRECISION,
98      &    king,FG_COMM,IERR)
99          time_Bcast=time_Bcast+MPI_Wtime()-time00
100          time_Bcastw=time_Bcastw+MPI_Wtime()-time00
101 c        call chainbuild_cart
102 c        call int_from_cart1(.false.)
103       endif
104 c      write (iout,*) 'Processor',myrank,
105 c     &  ' calling etotal_short ipot=',ipot
106 c      call flush(iout)
107 c      print *,'Processor',myrank,' nnt=',nnt,' nct=',nct
108 #endif     
109 cd    print *,'nnt=',nnt,' nct=',nct
110 C
111 C Compute the side-chain and electrostatic interaction energy
112 C
113       goto (101,102,103,104,105,106) ipot
114 C Lennard-Jones potential.
115   101 call elj_long(evdw)
116 cd    print '(a)','Exit ELJ'
117       goto 107
118 C Lennard-Jones-Kihara potential (shifted).
119   102 call eljk_long(evdw)
120       goto 107
121 C Berne-Pechukas potential (dilated LJ, angular dependence).
122   103 call ebp_long(evdw)
123       goto 107
124 C Gay-Berne potential (shifted LJ, angular dependence).
125   104 call egb_long(evdw)
126       goto 107
127 C Gay-Berne-Vorobjev potential (shifted LJ, angular dependence).
128   105 call egbv_long(evdw)
129       goto 107
130 C Soft-sphere potential
131   106 call e_softsphere(evdw)
132 C
133 C Calculate electrostatic (H-bonding) energy of the main chain.
134 C
135   107 continue
136       call vec_and_deriv
137       if (shield_mode.eq.1) then
138        call set_shield_fac
139       else if  (shield_mode.eq.2) then
140        call set_shield_fac2
141       if (nfgtasks.gt.1) then
142 C#define DEBUG
143 #ifdef DEBUG
144        write(iout,*) "befor reduce fac_shield reduce"
145        do i=1,nres
146         write(2,*) "fac",itype(i),fac_shield(i),grad_shield(1,i)
147         write(2,*) "list", shield_list(1,i),ishield_list(i),
148      &  grad_shield_side(1,1,i),grad_shield_loc(1,1,i)
149        enddo
150 #endif
151        call MPI_Allgatherv(fac_shield(ivec_start),ivec_count(fg_rank1),
152      &  MPI_DOUBLE_PRECISION,fac_shield(1),ivec_count(0),ivec_displ(0),
153      &  MPI_DOUBLE_PRECISION,FG_COMM,IERR)
154        call MPI_Allgatherv(shield_list(1,ivec_start),
155      &  ivec_count(fg_rank1),
156      &  MPI_I50,shield_list(1,1),ivec_count(0),
157      &  ivec_displ(0),
158      &  MPI_I50,FG_COMM,IERR)
159        call MPI_Allgatherv(ishield_list(ivec_start),
160      &  ivec_count(fg_rank1),
161      &  MPI_INTEGER,ishield_list(1),ivec_count(0),
162      &  ivec_displ(0),
163      &  MPI_INTEGER,FG_COMM,IERR)
164        call MPI_Allgatherv(grad_shield(1,ivec_start),
165      &  ivec_count(fg_rank1),
166      &  MPI_UYZ,grad_shield(1,1),ivec_count(0),
167      &  ivec_displ(0),
168      &  MPI_UYZ,FG_COMM,IERR)
169        call MPI_Allgatherv(grad_shield_side(1,1,ivec_start),
170      &  ivec_count(fg_rank1),
171      &  MPI_SHI,grad_shield_side(1,1,1),ivec_count(0),
172      &  ivec_displ(0),
173      &  MPI_SHI,FG_COMM,IERR)
174        call MPI_Allgatherv(grad_shield_loc(1,1,ivec_start),
175      &  ivec_count(fg_rank1),
176      &  MPI_SHI,grad_shield_loc(1,1,1),ivec_count(0),
177      &  ivec_displ(0),
178      &  MPI_SHI,FG_COMM,IERR)
179 #ifdef DEBUG
180        write(iout,*) "after reduce fac_shield reduce"
181        do i=1,nres
182         write(2,*) "fac",itype(i),fac_shield(i),grad_shield(1,i)
183         write(2,*) "list", shield_list(1,i),ishield_list(i),
184      &  grad_shield_side(1,1,i),grad_shield_loc(1,1,i)
185        enddo
186 #endif
187 C#undef DEBUG
188       endif
189
190       endif
191       if (ipot.lt.6) then
192 #ifdef SPLITELE
193          if (welec.gt.0d0.or.wvdwpp.gt.0d0.or.wel_loc.gt.0d0.or.
194      &       wturn3.gt.0d0.or.wturn4.gt.0d0 .or. wcorr.gt.0.0d0
195      &       .or. wcorr4.gt.0.0d0 .or. wcorr5.gt.0.d0
196      &       .or. wcorr6.gt.0.0d0 .or. wturn6.gt.0.0d0 ) then
197 #else
198          if (welec.gt.0d0.or.wel_loc.gt.0d0.or.
199      &       wturn3.gt.0d0.or.wturn4.gt.0d0 .or. wcorr.gt.0.0d0
200      &       .or. wcorr4.gt.0.0d0 .or. wcorr5.gt.0.d0
201      &       .or. wcorr6.gt.0.0d0 .or. wturn6.gt.0.0d0 ) then
202 #endif
203            call eelec_scale(ees,evdw1,eel_loc,eello_turn3,eello_turn4)
204          else
205             ees=0
206             evdw1=0
207             eel_loc=0
208             eello_turn3=0
209             eello_turn4=0
210          endif
211       else
212 c        write (iout,*) "Soft-spheer ELEC potential"
213         call eelec_soft_sphere(ees,evdw1,eel_loc,eello_turn3,
214      &   eello_turn4)
215       endif
216 C
217 C Calculate excluded-volume interaction energy between peptide groups
218 C and side chains.
219 C
220       if (ipot.lt.6) then
221        if(wscp.gt.0d0) then
222         call escp_long(evdw2,evdw2_14)
223        else
224         evdw2=0
225         evdw2_14=0
226        endif
227       else
228         call escp_soft_sphere(evdw2,evdw2_14)
229       endif
230
231 C 12/1/95 Multi-body terms
232 C
233       n_corr=0
234       n_corr1=0
235       if ((wcorr4.gt.0.0d0 .or. wcorr5.gt.0.0d0 .or. wcorr6.gt.0.0d0 
236      &    .or. wturn6.gt.0.0d0) .and. ipot.lt.6) then
237          call multibody_eello(ecorr,ecorr5,ecorr6,eturn6,n_corr,n_corr1)
238 c         write (2,*) 'n_corr=',n_corr,' n_corr1=',n_corr1,
239 c     &" ecorr",ecorr," ecorr5",ecorr5," ecorr6",ecorr6," eturn6",eturn6
240       else
241          ecorr=0.0d0
242          ecorr5=0.0d0
243          ecorr6=0.0d0
244          eturn6=0.0d0
245       endif
246       if ((wcorr4.eq.0.0d0 .and. wcorr.gt.0.0d0) .and. ipot.lt.6) then
247          call multibody_hb(ecorr,ecorr5,ecorr6,n_corr,n_corr1)
248       endif
249 C
250 C      if (wliptran.gt.0) then
251 C        call Eliptransfer(eliptran)
252 C      else
253 C       eliptran=0.0d0
254 C      endif 
255 C If performing constraint dynamics, call the constraint energy
256 C  after the equilibration time
257       if(usampl.and.totT.gt.eq_time) then
258          call EconstrQ   
259          call Econstr_back
260       else
261          Uconst=0.0d0
262          Uconst_back=0.0d0
263       endif
264
265 C Sum the energies
266 C
267       do i=1,n_ene
268         energia(i)=0.0d0
269       enddo
270       energia(1)=evdw
271 #ifdef SCP14
272       energia(2)=evdw2-evdw2_14
273       energia(18)=evdw2_14
274 #else
275       energia(2)=evdw2
276       energia(18)=0.0d0
277 #endif
278 #ifdef SPLITELE
279       energia(3)=ees
280       energia(16)=evdw1
281 #else
282       energia(3)=ees+evdw1
283       energia(16)=0.0d0
284 #endif
285       energia(4)=ecorr
286       energia(5)=ecorr5
287       energia(6)=ecorr6
288       energia(7)=eel_loc
289       energia(8)=eello_turn3
290       energia(9)=eello_turn4
291       energia(10)=eturn6
292       energia(20)=Uconst+Uconst_back
293       call sum_energy(energia,.true.)
294 C      call enerprint
295 C      write (iout,*) "Exit ETOTAL_LONG"
296 C      call flush(iout)
297       return
298       end
299 c------------------------------------------------------------------------------
300       subroutine etotal_short(energia)
301       implicit real*8 (a-h,o-z)
302       include 'DIMENSIONS'
303 c
304 c Compute the short-range fast-varying contributions to the energy
305 c
306 #ifndef ISNAN
307       external proc_proc
308 #ifdef WINPGI
309 cMS$ATTRIBUTES C ::  proc_proc
310 #endif
311 #endif
312 #ifdef MPI
313       include "mpif.h"
314       double precision weights_(n_ene)
315 #endif
316       include 'COMMON.SETUP'
317       include 'COMMON.IOUNITS'
318       double precision energia(0:n_ene)
319       include 'COMMON.FFIELD'
320       include 'COMMON.DERIV'
321       include 'COMMON.INTERACT'
322       include 'COMMON.SBRIDGE'
323       include 'COMMON.CHAIN'
324       include 'COMMON.VAR'
325       include 'COMMON.LOCAL'
326       include 'COMMON.CONTROL'
327       include 'COMMON.SHIELD'
328 c      write(iout,'(a,i2)')'Calling etotal_short ipot=',ipot
329 c      call flush(iout)
330       if (modecalc.eq.12.or.modecalc.eq.14) then
331 #ifdef MPI
332         if (fg_rank.eq.0) call int_from_cart1(.false.)
333 #else
334         call int_from_cart1(.false.)
335 #endif
336       endif
337 #ifdef MPI      
338 c      write(iout,*) "ETOTAL_SHORT Processor",fg_rank,
339 c     & " absolute rank",myrank," nfgtasks",nfgtasks
340 c      call flush(iout)
341       if (nfgtasks.gt.1) then
342         time00=MPI_Wtime()
343 C FG slaves call the following matching MPI_Bcast in ERGASTULUM
344         if (fg_rank.eq.0) then
345           call MPI_Bcast(2,1,MPI_INTEGER,king,FG_COMM,IERROR)
346 c          write (iout,*) "Processor",myrank," BROADCAST iorder"
347 c          call flush(iout)
348 C FG master sets up the WEIGHTS_ array which will be broadcast to the 
349 C FG slaves as WEIGHTS array.
350           weights_(1)=wsc
351           weights_(2)=wscp
352           weights_(3)=welec
353           weights_(4)=wcorr
354           weights_(5)=wcorr5
355           weights_(6)=wcorr6
356           weights_(7)=wel_loc
357           weights_(8)=wturn3
358           weights_(9)=wturn4
359           weights_(10)=wturn6
360           weights_(11)=wang
361           weights_(12)=wscloc
362           weights_(13)=wtor
363           weights_(14)=wtor_d
364           weights_(15)=wstrain
365           weights_(16)=wvdwpp
366           weights_(17)=wbond
367           weights_(18)=scal14
368           weights_(21)=wsccor
369 C FG Master broadcasts the WEIGHTS_ array
370           call MPI_Bcast(weights_(1),n_ene,
371      &        MPI_DOUBLE_PRECISION,king,FG_COMM,IERROR)
372         else
373 C FG slaves receive the WEIGHTS array
374           call MPI_Bcast(weights(1),n_ene,
375      &        MPI_DOUBLE_PRECISION,king,FG_COMM,IERROR)
376           wsc=weights(1)
377           wscp=weights(2)
378           welec=weights(3)
379           wcorr=weights(4)
380           wcorr5=weights(5)
381           wcorr6=weights(6)
382           wel_loc=weights(7)
383           wturn3=weights(8)
384           wturn4=weights(9)
385           wturn6=weights(10)
386           wang=weights(11)
387           wscloc=weights(12)
388           wtor=weights(13)
389           wtor_d=weights(14)
390           wstrain=weights(15)
391           wvdwpp=weights(16)
392           wbond=weights(17)
393           scal14=weights(18)
394           wsccor=weights(21)
395         endif
396 c        write (iout,*),"Processor",myrank," BROADCAST weights"
397         call MPI_Bcast(c(1,1),maxres6,MPI_DOUBLE_PRECISION,
398      &    king,FG_COMM,IERR)
399 c        write (iout,*) "Processor",myrank," BROADCAST c"
400         call MPI_Bcast(dc(1,1),maxres6,MPI_DOUBLE_PRECISION,
401      &    king,FG_COMM,IERR)
402 c        write (iout,*) "Processor",myrank," BROADCAST dc"
403         call MPI_Bcast(dc_norm(1,1),maxres6,MPI_DOUBLE_PRECISION,
404      &    king,FG_COMM,IERR)
405 c        write (iout,*) "Processor",myrank," BROADCAST dc_norm"
406         call MPI_Bcast(theta(1),nres,MPI_DOUBLE_PRECISION,
407      &    king,FG_COMM,IERR)
408 c        write (iout,*) "Processor",myrank," BROADCAST theta"
409         call MPI_Bcast(phi(1),nres,MPI_DOUBLE_PRECISION,
410      &    king,FG_COMM,IERR)
411 c        write (iout,*) "Processor",myrank," BROADCAST phi"
412         call MPI_Bcast(alph(1),nres,MPI_DOUBLE_PRECISION,
413      &    king,FG_COMM,IERR)
414 c        write (iout,*) "Processor",myrank," BROADCAST alph"
415         call MPI_Bcast(omeg(1),nres,MPI_DOUBLE_PRECISION,
416      &    king,FG_COMM,IERR)
417 c        write (iout,*) "Processor",myrank," BROADCAST omeg"
418         call MPI_Bcast(vbld(1),2*nres,MPI_DOUBLE_PRECISION,
419      &    king,FG_COMM,IERR)
420 c        write (iout,*) "Processor",myrank," BROADCAST vbld"
421         call MPI_Bcast(vbld_inv(1),2*nres,MPI_DOUBLE_PRECISION,
422      &    king,FG_COMM,IERR)
423          time_Bcast=time_Bcast+MPI_Wtime()-time00
424 c        write (iout,*) "Processor",myrank," BROADCAST vbld_inv"
425       endif
426 c      write (iout,*) 'Processor',myrank,
427 c     &  ' calling etotal_short ipot=',ipot
428 c      call flush(iout)
429 c      print *,'Processor',myrank,' nnt=',nnt,' nct=',nct
430 #endif     
431 c      call int_from_cart1(.false.)
432 C
433 C Compute the side-chain and electrostatic interaction energy
434 C
435       goto (101,102,103,104,105,106) ipot
436 C Lennard-Jones potential.
437   101 call elj_short(evdw)
438 cd    print '(a)','Exit ELJ'
439       goto 107
440 C Lennard-Jones-Kihara potential (shifted).
441   102 call eljk_short(evdw)
442       goto 107
443 C Berne-Pechukas potential (dilated LJ, angular dependence).
444   103 call ebp_short(evdw)
445       goto 107
446 C Gay-Berne potential (shifted LJ, angular dependence).
447   104 call egb_short(evdw)
448       goto 107
449 C Gay-Berne-Vorobjev potential (shifted LJ, angular dependence).
450   105 call egbv_short(evdw)
451       goto 107
452 C Soft-sphere potential - already dealt with in the long-range part
453   106 evdw=0.0d0
454 c  106 call e_softsphere_short(evdw)
455 C
456 C Calculate electrostatic (H-bonding) energy of the main chain.
457 C
458   107 continue
459 C     call vec_and_deriv
460 C     if (shield_mode.eq.1) then
461 C       call set_shield_fac
462 C      else if  (shield_mode.eq.2) then
463 C       call set_shield_fac2
464 C      if (nfgtasks.gt.1) then
465 C#define DEBUG
466 C#ifdef DEBUG
467 C       write(iout,*) "befor reduce fac_shield reduce"
468 C       do i=1,nres
469 C        write(2,*) "fac",itype(i),fac_shield(i),grad_shield(1,i)
470 C        write(2,*) "list", shield_list(1,i),ishield_list(i),
471 C     &  grad_shield_side(1,1,i),grad_shield_loc(1,1,i)
472 C       enddo
473 C#endif
474 C       call MPI_Allgatherv(fac_shield(ivec_start),ivec_count(fg_rank1),
475 C     &  MPI_DOUBLE_PRECISION,fac_shield(1),ivec_count(0),ivec_displ(0),
476 C     &  MPI_DOUBLE_PRECISION,FG_COMM,IERR)
477 C       call MPI_Allgatherv(shield_list(1,ivec_start),
478 C     &  ivec_count(fg_rank1),
479 C     &  MPI_I50,shield_list(1,1),ivec_count(0),
480 C     &  ivec_displ(0),
481 C     &  MPI_I50,FG_COMM,IERR)
482 C       call MPI_Allgatherv(ishield_list(ivec_start),
483 C     &  ivec_count(fg_rank1),
484 C     &  MPI_INTEGER,ishield_list(1),ivec_count(0),
485 C     &  ivec_displ(0),
486 C     &  MPI_INTEGER,FG_COMM,IERR)
487 C       call MPI_Allgatherv(grad_shield(1,ivec_start),
488 C     &  ivec_count(fg_rank1),
489 C     &  MPI_UYZ,grad_shield(1,1),ivec_count(0),
490 C     &  ivec_displ(0),
491 C     &  MPI_UYZ,FG_COMM,IERR)
492 C       call MPI_Allgatherv(grad_shield_side(1,1,ivec_start),
493 C     &  ivec_count(fg_rank1),
494 C     &  MPI_SHI,grad_shield_side(1,1,1),ivec_count(0),
495 C     &  ivec_displ(0),
496 C     &  MPI_SHI,FG_COMM,IERR)
497 C       call MPI_Allgatherv(grad_shield_loc(1,1,ivec_start),
498 C     &  ivec_count(fg_rank1),
499 C     &  MPI_SHI,grad_shield_loc(1,1,1),ivec_count(0),
500 C     &  ivec_displ(0),
501 C     &  MPI_SHI,FG_COMM,IERR)
502 C#ifdef DEBUG
503 C       write(iout,*) "after reduce fac_shield reduce"
504 C       do i=1,nres
505 C        write(2,*) "fac",itype(i),fac_shield(i),grad_shield(1,i)
506 C        write(2,*) "list", shield_list(1,i),ishield_list(i),
507 C     &  grad_shield_side(1,1,i),grad_shield_loc(1,1,i)
508 C       enddo
509 C#endif
510 C#undef DEBUG
511 C      endif
512 C
513 C      endif
514 C      if (ipot.lt.6) then
515 C#ifdef SPLITELE
516 C         if (welec.gt.0d0.or.wvdwpp.gt.0d0.or.wel_loc.gt.0d0.or.
517 C     &       wturn3.gt.0d0.or.wturn4.gt.0d0 .or. wcorr.gt.0.0d0
518 C     &       .or. wcorr4.gt.0.0d0 .or. wcorr5.gt.0.d0
519 C     &       .or. wcorr6.gt.0.0d0 .or. wturn6.gt.0.0d0 ) then
520 C#else
521 C         if (welec.gt.0d0.or.wel_loc.gt.0d0.or.
522 C     &       wturn3.gt.0d0.or.wturn4.gt.0d0 .or. wcorr.gt.0.0d0
523 C     &       .or. wcorr4.gt.0.0d0 .or. wcorr5.gt.0.d0
524 C     &       .or. wcorr6.gt.0.0d0 .or. wturn6.gt.0.0d0 ) then
525 C#endif
526 C           call eelec_scale(ees,evdw1,eel_loc,eello_turn3,eello_turn4)
527 C         else
528 C            ees=0
529 C            evdw1=0
530 C            eel_loc=0
531 C            eello_turn3=0
532 C            eello_turn4=0
533 C         endif
534 C      else
535 c        write (iout,*) "Soft-spheer ELEC potential"
536 C        call eelec_soft_sphere(ees,evdw1,eel_loc,eello_turn3,
537 C     &   eello_turn4)
538 C      endif
539
540 c
541 c Calculate the short-range part of Evdwpp
542 c
543       call evdwpp_short(evdw1)
544 c
545 c Calculate the short-range part of ESCp
546 c
547       if (ipot.lt.6) then
548         call escp_short(evdw2,evdw2_14)
549       endif
550 c
551 c Calculate the bond-stretching energy
552 c
553       call ebond(estr)
554
555 C Calculate the disulfide-bridge and other energy and the contributions
556 C from other distance constraints.
557       call edis(ehpb)
558 C
559 C Calculate the virtual-bond-angle energy.
560 C
561       call ebend(ebe,ethetcnstr)
562 C
563 C Calculate the SC local energy.
564 C
565       call vec_and_deriv
566       call esc(escloc)
567 C
568 C Calculate the virtual-bond torsional energy.
569 C
570       call etor(etors,edihcnstr)
571 C
572 C 6/23/01 Calculate double-torsional energy
573 C
574       call etor_d(etors_d)
575 C
576 C 21/5/07 Calculate local sicdechain correlation energy
577 C
578       if (wsccor.gt.0.0d0) then
579         call eback_sc_corr(esccor)
580       else
581         esccor=0.0d0
582       endif
583       if (wliptran.gt.0) then
584         call Eliptransfer(eliptran)
585       else
586        eliptran=0.0d0
587       endif
588 C       print *,eliptran,wliptran
589 C
590 C Put energy components into an array
591 C
592       do i=1,n_ene
593         energia(i)=0.0d0
594       enddo
595       energia(1)=evdw
596 #ifdef SCP14
597       energia(2)=evdw2-evdw2_14
598       energia(18)=evdw2_14
599 #else
600       energia(2)=evdw2
601       energia(18)=0.0d0
602 #endif
603 #ifdef SPLITELE
604       energia(3)=ees
605       energia(16)=evdw1
606 #else
607       energia(3)=ees+evdw1
608       energia(16)=0.0d0
609 #endif
610       energia(4)=ecorr
611       energia(5)=ecorr5
612       energia(6)=ecorr6
613       energia(7)=eel_loc
614       energia(8)=eello_turn3
615       energia(9)=eello_turn4
616
617       energia(11)=ebe
618       energia(12)=escloc
619       energia(13)=etors
620       energia(14)=etors_d
621       energia(15)=ehpb
622       energia(17)=estr
623       energia(19)=edihcnstr
624       energia(21)=esccor
625       energia(22)=eliptran
626 C      write (iout,*) "ETOTAL_SHORT before SUM_ENERGY"
627       call flush(iout)
628       call sum_energy(energia,.true.)
629 C      write (iout,*) "Exit ETOTAL_SHORT"
630 C      call flush(iout)
631       return
632       end