respa cleaning (my bug) and shielding MPI and previous bug cleaning
[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
250 C If performing constraint dynamics, call the constraint energy
251 C  after the equilibration time
252       if(usampl.and.totT.gt.eq_time) then
253          call EconstrQ   
254          call Econstr_back
255       else
256          Uconst=0.0d0
257          Uconst_back=0.0d0
258       endif
259
260 C Sum the energies
261 C
262       do i=1,n_ene
263         energia(i)=0.0d0
264       enddo
265       energia(1)=evdw
266 #ifdef SCP14
267       energia(2)=evdw2-evdw2_14
268       energia(18)=evdw2_14
269 #else
270       energia(2)=evdw2
271       energia(18)=0.0d0
272 #endif
273 #ifdef SPLITELE
274       energia(3)=ees
275       energia(16)=evdw1
276 #else
277       energia(3)=ees+evdw1
278       energia(16)=0.0d0
279 #endif
280       energia(4)=ecorr
281       energia(5)=ecorr5
282       energia(6)=ecorr6
283       energia(7)=eel_loc
284       energia(8)=eello_turn3
285       energia(9)=eello_turn4
286       energia(10)=eturn6
287       energia(20)=Uconst+Uconst_back
288       call sum_energy(energia,.true.)
289 C      call enerprint
290       write (iout,*) "Exit ETOTAL_LONG"
291       call flush(iout)
292       return
293       end
294 c------------------------------------------------------------------------------
295       subroutine etotal_short(energia)
296       implicit real*8 (a-h,o-z)
297       include 'DIMENSIONS'
298 c
299 c Compute the short-range fast-varying contributions to the energy
300 c
301 #ifndef ISNAN
302       external proc_proc
303 #ifdef WINPGI
304 cMS$ATTRIBUTES C ::  proc_proc
305 #endif
306 #endif
307 #ifdef MPI
308       include "mpif.h"
309       double precision weights_(n_ene)
310 #endif
311       include 'COMMON.SETUP'
312       include 'COMMON.IOUNITS'
313       double precision energia(0:n_ene)
314       include 'COMMON.FFIELD'
315       include 'COMMON.DERIV'
316       include 'COMMON.INTERACT'
317       include 'COMMON.SBRIDGE'
318       include 'COMMON.CHAIN'
319       include 'COMMON.VAR'
320       include 'COMMON.LOCAL'
321       include 'COMMON.CONTROL'
322       include 'COMMON.SHIELD'
323 c      write(iout,'(a,i2)')'Calling etotal_short ipot=',ipot
324 c      call flush(iout)
325       if (modecalc.eq.12.or.modecalc.eq.14) then
326 #ifdef MPI
327         if (fg_rank.eq.0) call int_from_cart1(.false.)
328 #else
329         call int_from_cart1(.false.)
330 #endif
331       endif
332 #ifdef MPI      
333 c      write(iout,*) "ETOTAL_SHORT Processor",fg_rank,
334 c     & " absolute rank",myrank," nfgtasks",nfgtasks
335 c      call flush(iout)
336       if (nfgtasks.gt.1) then
337         time00=MPI_Wtime()
338 C FG slaves call the following matching MPI_Bcast in ERGASTULUM
339         if (fg_rank.eq.0) then
340           call MPI_Bcast(2,1,MPI_INTEGER,king,FG_COMM,IERROR)
341 c          write (iout,*) "Processor",myrank," BROADCAST iorder"
342 c          call flush(iout)
343 C FG master sets up the WEIGHTS_ array which will be broadcast to the 
344 C FG slaves as WEIGHTS array.
345           weights_(1)=wsc
346           weights_(2)=wscp
347           weights_(3)=welec
348           weights_(4)=wcorr
349           weights_(5)=wcorr5
350           weights_(6)=wcorr6
351           weights_(7)=wel_loc
352           weights_(8)=wturn3
353           weights_(9)=wturn4
354           weights_(10)=wturn6
355           weights_(11)=wang
356           weights_(12)=wscloc
357           weights_(13)=wtor
358           weights_(14)=wtor_d
359           weights_(15)=wstrain
360           weights_(16)=wvdwpp
361           weights_(17)=wbond
362           weights_(18)=scal14
363           weights_(21)=wsccor
364 C FG Master broadcasts the WEIGHTS_ array
365           call MPI_Bcast(weights_(1),n_ene,
366      &        MPI_DOUBLE_PRECISION,king,FG_COMM,IERROR)
367         else
368 C FG slaves receive the WEIGHTS array
369           call MPI_Bcast(weights(1),n_ene,
370      &        MPI_DOUBLE_PRECISION,king,FG_COMM,IERROR)
371           wsc=weights(1)
372           wscp=weights(2)
373           welec=weights(3)
374           wcorr=weights(4)
375           wcorr5=weights(5)
376           wcorr6=weights(6)
377           wel_loc=weights(7)
378           wturn3=weights(8)
379           wturn4=weights(9)
380           wturn6=weights(10)
381           wang=weights(11)
382           wscloc=weights(12)
383           wtor=weights(13)
384           wtor_d=weights(14)
385           wstrain=weights(15)
386           wvdwpp=weights(16)
387           wbond=weights(17)
388           scal14=weights(18)
389           wsccor=weights(21)
390         endif
391 c        write (iout,*),"Processor",myrank," BROADCAST weights"
392         call MPI_Bcast(c(1,1),maxres6,MPI_DOUBLE_PRECISION,
393      &    king,FG_COMM,IERR)
394 c        write (iout,*) "Processor",myrank," BROADCAST c"
395         call MPI_Bcast(dc(1,1),maxres6,MPI_DOUBLE_PRECISION,
396      &    king,FG_COMM,IERR)
397 c        write (iout,*) "Processor",myrank," BROADCAST dc"
398         call MPI_Bcast(dc_norm(1,1),maxres6,MPI_DOUBLE_PRECISION,
399      &    king,FG_COMM,IERR)
400 c        write (iout,*) "Processor",myrank," BROADCAST dc_norm"
401         call MPI_Bcast(theta(1),nres,MPI_DOUBLE_PRECISION,
402      &    king,FG_COMM,IERR)
403 c        write (iout,*) "Processor",myrank," BROADCAST theta"
404         call MPI_Bcast(phi(1),nres,MPI_DOUBLE_PRECISION,
405      &    king,FG_COMM,IERR)
406 c        write (iout,*) "Processor",myrank," BROADCAST phi"
407         call MPI_Bcast(alph(1),nres,MPI_DOUBLE_PRECISION,
408      &    king,FG_COMM,IERR)
409 c        write (iout,*) "Processor",myrank," BROADCAST alph"
410         call MPI_Bcast(omeg(1),nres,MPI_DOUBLE_PRECISION,
411      &    king,FG_COMM,IERR)
412 c        write (iout,*) "Processor",myrank," BROADCAST omeg"
413         call MPI_Bcast(vbld(1),2*nres,MPI_DOUBLE_PRECISION,
414      &    king,FG_COMM,IERR)
415 c        write (iout,*) "Processor",myrank," BROADCAST vbld"
416         call MPI_Bcast(vbld_inv(1),2*nres,MPI_DOUBLE_PRECISION,
417      &    king,FG_COMM,IERR)
418          time_Bcast=time_Bcast+MPI_Wtime()-time00
419 c        write (iout,*) "Processor",myrank," BROADCAST vbld_inv"
420       endif
421 c      write (iout,*) 'Processor',myrank,
422 c     &  ' calling etotal_short ipot=',ipot
423 c      call flush(iout)
424 c      print *,'Processor',myrank,' nnt=',nnt,' nct=',nct
425 #endif     
426 c      call int_from_cart1(.false.)
427 C
428 C Compute the side-chain and electrostatic interaction energy
429 C
430       goto (101,102,103,104,105,106) ipot
431 C Lennard-Jones potential.
432   101 call elj_short(evdw)
433 cd    print '(a)','Exit ELJ'
434       goto 107
435 C Lennard-Jones-Kihara potential (shifted).
436   102 call eljk_short(evdw)
437       goto 107
438 C Berne-Pechukas potential (dilated LJ, angular dependence).
439   103 call ebp_short(evdw)
440       goto 107
441 C Gay-Berne potential (shifted LJ, angular dependence).
442   104 call egb_short(evdw)
443       goto 107
444 C Gay-Berne-Vorobjev potential (shifted LJ, angular dependence).
445   105 call egbv_short(evdw)
446       goto 107
447 C Soft-sphere potential - already dealt with in the long-range part
448   106 evdw=0.0d0
449 c  106 call e_softsphere_short(evdw)
450 C
451 C Calculate electrostatic (H-bonding) energy of the main chain.
452 C
453   107 continue
454 C     call vec_and_deriv
455 C     if (shield_mode.eq.1) then
456 C       call set_shield_fac
457 C      else if  (shield_mode.eq.2) then
458 C       call set_shield_fac2
459 C      if (nfgtasks.gt.1) then
460 C#define DEBUG
461 C#ifdef DEBUG
462 C       write(iout,*) "befor reduce fac_shield reduce"
463 C       do i=1,nres
464 C        write(2,*) "fac",itype(i),fac_shield(i),grad_shield(1,i)
465 C        write(2,*) "list", shield_list(1,i),ishield_list(i),
466 C     &  grad_shield_side(1,1,i),grad_shield_loc(1,1,i)
467 C       enddo
468 C#endif
469 C       call MPI_Allgatherv(fac_shield(ivec_start),ivec_count(fg_rank1),
470 C     &  MPI_DOUBLE_PRECISION,fac_shield(1),ivec_count(0),ivec_displ(0),
471 C     &  MPI_DOUBLE_PRECISION,FG_COMM,IERR)
472 C       call MPI_Allgatherv(shield_list(1,ivec_start),
473 C     &  ivec_count(fg_rank1),
474 C     &  MPI_I50,shield_list(1,1),ivec_count(0),
475 C     &  ivec_displ(0),
476 C     &  MPI_I50,FG_COMM,IERR)
477 C       call MPI_Allgatherv(ishield_list(ivec_start),
478 C     &  ivec_count(fg_rank1),
479 C     &  MPI_INTEGER,ishield_list(1),ivec_count(0),
480 C     &  ivec_displ(0),
481 C     &  MPI_INTEGER,FG_COMM,IERR)
482 C       call MPI_Allgatherv(grad_shield(1,ivec_start),
483 C     &  ivec_count(fg_rank1),
484 C     &  MPI_UYZ,grad_shield(1,1),ivec_count(0),
485 C     &  ivec_displ(0),
486 C     &  MPI_UYZ,FG_COMM,IERR)
487 C       call MPI_Allgatherv(grad_shield_side(1,1,ivec_start),
488 C     &  ivec_count(fg_rank1),
489 C     &  MPI_SHI,grad_shield_side(1,1,1),ivec_count(0),
490 C     &  ivec_displ(0),
491 C     &  MPI_SHI,FG_COMM,IERR)
492 C       call MPI_Allgatherv(grad_shield_loc(1,1,ivec_start),
493 C     &  ivec_count(fg_rank1),
494 C     &  MPI_SHI,grad_shield_loc(1,1,1),ivec_count(0),
495 C     &  ivec_displ(0),
496 C     &  MPI_SHI,FG_COMM,IERR)
497 C#ifdef DEBUG
498 C       write(iout,*) "after reduce fac_shield reduce"
499 C       do i=1,nres
500 C        write(2,*) "fac",itype(i),fac_shield(i),grad_shield(1,i)
501 C        write(2,*) "list", shield_list(1,i),ishield_list(i),
502 C     &  grad_shield_side(1,1,i),grad_shield_loc(1,1,i)
503 C       enddo
504 C#endif
505 C#undef DEBUG
506 C      endif
507 C
508 C      endif
509 C      if (ipot.lt.6) then
510 C#ifdef SPLITELE
511 C         if (welec.gt.0d0.or.wvdwpp.gt.0d0.or.wel_loc.gt.0d0.or.
512 C     &       wturn3.gt.0d0.or.wturn4.gt.0d0 .or. wcorr.gt.0.0d0
513 C     &       .or. wcorr4.gt.0.0d0 .or. wcorr5.gt.0.d0
514 C     &       .or. wcorr6.gt.0.0d0 .or. wturn6.gt.0.0d0 ) then
515 C#else
516 C         if (welec.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#endif
521 C           call eelec_scale(ees,evdw1,eel_loc,eello_turn3,eello_turn4)
522 C         else
523 C            ees=0
524 C            evdw1=0
525 C            eel_loc=0
526 C            eello_turn3=0
527 C            eello_turn4=0
528 C         endif
529 C      else
530 c        write (iout,*) "Soft-spheer ELEC potential"
531 C        call eelec_soft_sphere(ees,evdw1,eel_loc,eello_turn3,
532 C     &   eello_turn4)
533 C      endif
534
535 c
536 c Calculate the short-range part of Evdwpp
537 c
538       call evdwpp_short(evdw1)
539 c
540 c Calculate the short-range part of ESCp
541 c
542       if (ipot.lt.6) then
543         call escp_short(evdw2,evdw2_14)
544       endif
545 c
546 c Calculate the bond-stretching energy
547 c
548       call ebond(estr)
549
550 C Calculate the disulfide-bridge and other energy and the contributions
551 C from other distance constraints.
552       call edis(ehpb)
553 C
554 C Calculate the virtual-bond-angle energy.
555 C
556       call ebend(ebe,ethetcnstr)
557 C
558 C Calculate the SC local energy.
559 C
560       call vec_and_deriv
561       call esc(escloc)
562 C
563 C Calculate the virtual-bond torsional energy.
564 C
565       call etor(etors,edihcnstr)
566 C
567 C 6/23/01 Calculate double-torsional energy
568 C
569       call etor_d(etors_d)
570 C
571 C 21/5/07 Calculate local sicdechain correlation energy
572 C
573       if (wsccor.gt.0.0d0) then
574         call eback_sc_corr(esccor)
575       else
576         esccor=0.0d0
577       endif
578 C
579 C Put energy components into an array
580 C
581       do i=1,n_ene
582         energia(i)=0.0d0
583       enddo
584       energia(1)=evdw
585 #ifdef SCP14
586       energia(2)=evdw2-evdw2_14
587       energia(18)=evdw2_14
588 #else
589       energia(2)=evdw2
590       energia(18)=0.0d0
591 #endif
592 #ifdef SPLITELE
593       energia(3)=ees
594       energia(16)=evdw1
595 #else
596       energia(3)=ees+evdw1
597       energia(16)=0.0d0
598 #endif
599       energia(4)=ecorr
600       energia(5)=ecorr5
601       energia(6)=ecorr6
602       energia(7)=eel_loc
603       energia(8)=eello_turn3
604       energia(9)=eello_turn4
605
606       energia(11)=ebe
607       energia(12)=escloc
608       energia(13)=etors
609       energia(14)=etors_d
610       energia(15)=ehpb
611       energia(17)=estr
612       energia(19)=edihcnstr
613       energia(21)=esccor
614       write (iout,*) "ETOTAL_SHORT before SUM_ENERGY"
615       call flush(iout)
616       call sum_energy(energia,.true.)
617       write (iout,*) "Exit ETOTAL_SHORT"
618       call flush(iout)
619       return
620       end