update new files
[unres.git] / source / maxlik / src-Fmatch_safe / sum_gradient.F
1       subroutine sum_gradient
2       implicit none
3       include 'DIMENSIONS'
4       include "DIMENSIONS.ZSCOPT"
5 #ifndef ISNAN
6       external proc_proc
7 #ifdef WINPGI
8 cMS$ATTRIBUTES C ::  proc_proc
9 #endif
10 #endif
11 #ifdef MPI
12       include 'mpif.h'
13       include "COMMON.SETUP"
14       integer ierr,ierror
15 #endif
16       double precision gradbufc(3,-1:maxres),gradbufx(3,-1:maxres),
17      & glocbuf(4*maxres),gradbufc_sum(3,-1:maxres)
18      & ,gloc_scbuf(3,-1:maxres)
19       include 'COMMON.IOUNITS'
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.CONTROL'
27       include 'COMMON.TIME1'
28       include 'COMMON.SCCOR'
29       integer i,j,k
30       double precision time00,time_barrier,time_barrier_g,time_reduce
31 #ifdef TIMING
32       time01=MPI_Wtime()
33 #endif
34 #ifdef DEBUG
35       write (iout,*) "sum_gradient gvdwc, gvdwx"
36       do i=1,nres
37         write (iout,'(i3,3f10.5,5x,3f10.5,5x,f10.5)') 
38      &   i,(gvdwx(j,i),j=1,3),(gvdwc(j,i),j=1,3)
39       enddo
40       call flush(iout)
41 #endif
42 #ifdef MPI
43 C FG slaves call the following matching MPI_Bcast in ERGASTULUM
44         if (nfgtasks.gt.1 .and. fg_rank.eq.0) 
45      &    call MPI_Bcast(1,1,MPI_INTEGER,king,FG_COMM,IERROR)
46 #endif
47 C
48 C 9/29/08 AL Transform parts of gradients in site coordinates to the gradient
49 C            in virtual-bond-vector coordinates
50 C
51 #ifdef DEBUG
52 c      write (iout,*) "gel_loc gel_loc_long and gel_loc_loc"
53 c      do i=1,nres-1
54 c        write (iout,'(i5,3f10.5,2x,3f10.5,2x,f10.5)') 
55 c     &   i,(gel_loc(j,i),j=1,3),(gel_loc_long(j,i),j=1,3),gel_loc_loc(i)
56 c      enddo
57 c      write (iout,*) "gel_loc_tur3 gel_loc_turn4"
58 c      do i=1,nres-1
59 c        write (iout,'(i5,3f10.5,2x,f10.5)') 
60 c     &  i,(gcorr4_turn(j,i),j=1,3),gel_loc_turn4(i)
61 c      enddo
62       write (iout,*) "gradcorr5 gradcorr5_long gradcorr5_loc"
63       do i=1,nres
64         write (iout,'(i3,3f10.5,5x,3f10.5,5x,f10.5)') 
65      &   i,(gradcorr5(j,i),j=1,3),(gradcorr5_long(j,i),j=1,3),
66      &   g_corr5_loc(i)
67       enddo
68       call flush(iout)
69 #endif
70 #ifdef SPLITELE
71       do i=0,nct
72         do j=1,3
73           gradbufc(j,i)=wsc*gvdwc(j,i)+
74      &                wscp*(gvdwc_scp(j,i)+gvdwc_scpp(j,i))+
75      &                welec*gelc_long(j,i)+wvdwpp*gvdwpp(j,i)+
76      &                wel_loc*gel_loc_long(j,i)+
77      &                wcorr*gradcorr_long(j,i)+
78      &                wcorr5*gradcorr5_long(j,i)+
79      &                wcorr6*gradcorr6_long(j,i)+
80      &                wturn6*gcorr6_turn_long(j,i)+
81      &                wstrain*ghpbc(j,i)
82      &                +wliptran*gliptranc(j,i)
83 c     &                +gradafm(j,i)
84      &                 +welec*gshieldc(j,i)
85      &                 +wcorr*gshieldc_ec(j,i)
86      &                 +wturn3*gshieldc_t3(j,i)
87      &                 +wturn4*gshieldc_t4(j,i)
88      &                 +wel_loc*gshieldc_ll(j,i)
89 c     &                +wtube*gg_tube(j,i)
90
91
92
93         enddo
94       enddo 
95 #else
96       do i=0,nct
97         do j=1,3
98           gradbufc(j,i)=wsc*gvdwc(j,i)+
99      &                wscp*(gvdwc_scp(j,i)+gvdwc_scpp(j,i))+
100      &                welec*gelc_long(j,i)+
101      &                wbond*gradb(j,i)+
102      &                wel_loc*gel_loc_long(j,i)+
103      &                wcorr*gradcorr_long(j,i)+
104      &                wcorr5*gradcorr5_long(j,i)+
105      &                wcorr6*gradcorr6_long(j,i)+
106      &                wturn6*gcorr6_turn_long(j,i)+
107      &                wstrain*ghpbc(j,i)
108      &                +wliptran*gliptranc(j,i)
109 c     &                +gradafm(j,i)
110      &                 +welec*gshieldc(j,i)
111      &                 +wcorr*gshieldc_ec(j,i)
112      &                 +wturn4*gshieldc_t4(j,i)
113      &                 +wel_loc*gshieldc_ll(j,i)
114 c     &                +wtube*gg_tube(j,i)
115
116
117
118         enddo
119       enddo 
120 #endif
121 #ifdef MPI
122       if (nfgtasks.gt.1) then
123       time00=MPI_Wtime()
124 #ifdef DEBUG
125       write (iout,*) "gradbufc before allreduce"
126       do i=1,nres
127         write (iout,'(i3,3f10.5)') i,(gradbufc(j,i),j=1,3)
128       enddo
129       call flush(iout)
130 #endif
131       do i=0,nres
132         do j=1,3
133           gradbufc_sum(j,i)=gradbufc(j,i)
134         enddo
135       enddo
136 c      call MPI_AllReduce(gradbufc(1,1),gradbufc_sum(1,1),3*nres,
137 c     &    MPI_DOUBLE_PRECISION,MPI_SUM,FG_COMM,IERR)
138 c      time_reduce=time_reduce+MPI_Wtime()-time00
139 #ifdef DEBUG
140 c      write (iout,*) "gradbufc_sum after allreduce"
141 c      do i=1,nres
142 c        write (iout,'(i3,3f10.5)') i,(gradbufc_sum(j,i),j=1,3)
143 c      enddo
144 c      call flush(iout)
145 #endif
146 #ifdef TIMING
147 c      time_allreduce=time_allreduce+MPI_Wtime()-time00
148 #endif
149       do i=nnt,nres
150         do k=1,3
151           gradbufc(k,i)=0.0d0
152         enddo
153       enddo
154 #ifdef DEBUG
155       write (iout,*) "igrad_start",igrad_start," igrad_end",igrad_end
156       write (iout,*) (i," jgrad_start",jgrad_start(i),
157      &                  " jgrad_end  ",jgrad_end(i),
158      &                  i=igrad_start,igrad_end)
159 #endif
160 c
161 c Obsolete and inefficient code; we can make the effort O(n) and, therefore,
162 c do not parallelize this part.
163 c
164 c      do i=igrad_start,igrad_end
165 c        do j=jgrad_start(i),jgrad_end(i)
166 c          do k=1,3
167 c            gradbufc(k,i)=gradbufc(k,i)+gradbufc_sum(k,j)
168 c          enddo
169 c        enddo
170 c      enddo
171       do j=1,3
172         gradbufc(j,nres-1)=gradbufc_sum(j,nres)
173       enddo
174       do i=nres-2,-1,-1
175         do j=1,3
176           gradbufc(j,i)=gradbufc(j,i+1)+gradbufc_sum(j,i+1)
177         enddo
178       enddo
179 #ifdef DEBUG
180       write (iout,*) "gradbufc after summing"
181       do i=1,nres
182         write (iout,'(i3,3f10.5)') i,(gradbufc(j,i),j=1,3)
183       enddo
184       call flush(iout)
185 #endif
186       else
187 #endif
188 #ifdef DEBUG
189       write (iout,*) "gradbufc"
190       do i=1,nres
191         write (iout,'(i3,3f10.5)') i,(gradbufc(j,i),j=1,3)
192       enddo
193       call flush(iout)
194 #endif
195       do i=-1,nres
196         do j=1,3
197           gradbufc_sum(j,i)=gradbufc(j,i)
198           gradbufc(j,i)=0.0d0
199         enddo
200       enddo
201       do j=1,3
202         gradbufc(j,nres-1)=gradbufc_sum(j,nres)
203       enddo
204       do i=nres-2,-1,-1
205         do j=1,3
206           gradbufc(j,i)=gradbufc(j,i+1)+gradbufc_sum(j,i+1)
207         enddo
208       enddo
209 c      do i=nnt,nres-1
210 c        do k=1,3
211 c          gradbufc(k,i)=0.0d0
212 c        enddo
213 c        do j=i+1,nres
214 c          do k=1,3
215 c            gradbufc(k,i)=gradbufc(k,i)+gradbufc(k,j)
216 c          enddo
217 c        enddo
218 c      enddo
219 #ifdef DEBUG
220       write (iout,*) "gradbufc after summing"
221       do i=1,nres
222         write (iout,'(i3,3f10.5)') i,(gradbufc(j,i),j=1,3)
223       enddo
224       call flush(iout)
225 #endif
226 #ifdef MPI
227       endif
228 #endif
229       do k=1,3
230         gradbufc(k,nres)=0.0d0
231       enddo
232       do i=-1,nct
233         do j=1,3
234 #ifdef SPLITELE
235 C          print *,gradbufc(1,13)
236 C          print *,welec*gelc(1,13)
237 C          print *,wel_loc*gel_loc(1,13)
238 C          print *,0.5d0*(wscp*gvdwc_scpp(1,13))
239 C          print *,welec*gelc_long(1,13)+wvdwpp*gvdwpp(1,13)
240 C          print *,wel_loc*gel_loc_long(1,13)
241 C          print *,gradafm(1,13),"AFM"
242           gradc(j,i,icg)=gradbufc(j,i)+welec*gelc(j,i)+
243      &                wel_loc*gel_loc(j,i)+
244      &                0.5d0*(wscp*gvdwc_scpp(j,i)+
245      &                welec*gelc_long(j,i)+wvdwpp*gvdwpp(j,i)+
246      &                wel_loc*gel_loc_long(j,i)+
247      &                wcorr*gradcorr_long(j,i)+
248      &                wcorr5*gradcorr5_long(j,i)+
249      &                wcorr6*gradcorr6_long(j,i)+
250      &                wturn6*gcorr6_turn_long(j,i))+
251      &                wbond*gradb(j,i)+
252      &                wcorr*gradcorr(j,i)+
253      &                wturn3*gcorr3_turn(j,i)+
254      &                wturn4*gcorr4_turn(j,i)+
255      &                wcorr5*gradcorr5(j,i)+
256      &                wcorr6*gradcorr6(j,i)+
257      &                wturn6*gcorr6_turn(j,i)+
258      &                wsccor*gsccorc(j,i)
259      &               +wscloc*gscloc(j,i)
260      &               +wliptran*gliptranc(j,i)
261 c     &                +gradafm(j,i)
262      &                 +welec*gshieldc(j,i)
263      &                 +welec*gshieldc_loc(j,i)
264      &                 +wcorr*gshieldc_ec(j,i)
265      &                 +wcorr*gshieldc_loc_ec(j,i)
266      &                 +wturn3*gshieldc_t3(j,i)
267      &                 +wturn3*gshieldc_loc_t3(j,i)
268      &                 +wturn4*gshieldc_t4(j,i)
269      &                 +wturn4*gshieldc_loc_t4(j,i)
270      &                 +wel_loc*gshieldc_ll(j,i)
271      &                 +wel_loc*gshieldc_loc_ll(j,i)
272 c     &                +wtube*gg_tube(j,i)
273
274 #else
275           gradc(j,i,icg)=gradbufc(j,i)+welec*gelc(j,i)+
276      &                wel_loc*gel_loc(j,i)+
277      &                0.5d0*(wscp*gvdwc_scpp(j,i)+
278      &                welec*gelc_long(j,i)+
279      &                wel_loc*gel_loc_long(j,i)+
280      &                wcorr*gcorr_long(j,i)+
281      &                wcorr5*gradcorr5_long(j,i)+
282      &                wcorr6*gradcorr6_long(j,i)+
283      &                wturn6*gcorr6_turn_long(j,i))+
284      &                wbond*gradb(j,i)+
285      &                wcorr*gradcorr(j,i)+
286      &                wturn3*gcorr3_turn(j,i)+
287      &                wturn4*gcorr4_turn(j,i)+
288      &                wcorr5*gradcorr5(j,i)+
289      &                wcorr6*gradcorr6(j,i)+
290      &                wturn6*gcorr6_turn(j,i)+
291      &                wsccor*gsccorc(j,i)
292      &               +wscloc*gscloc(j,i)
293      &               +wliptran*gliptranc(j,i)
294 c     &                +gradafm(j,i)
295      &                 +welec*gshieldc(j,i)
296      &                 +welec*gshieldc_loc(j,i)
297      &                 +wcorr*gshieldc_ec(j,i)
298      &                 +wcorr*gshieldc_loc_ec(j,i)
299      &                 +wturn3*gshieldc_t3(j,i)
300      &                 +wturn3*gshieldc_loc_t3(j,i)
301      &                 +wturn4*gshieldc_t4(j,i)
302      &                 +wturn4*gshieldc_loc_t4(j,i)
303      &                 +wel_loc*gshieldc_ll(j,i)
304      &                 +wel_loc*gshieldc_loc_ll(j,i)
305 c     &                +wtube*gg_tube(j,i)
306
307
308 #endif
309           gradx(j,i,icg)=wsc*gvdwx(j,i)+wscp*gradx_scp(j,i)+
310      &                  wbond*gradbx(j,i)+
311      &                  wstrain*ghpbx(j,i)+wcorr*gradxorr(j,i)+
312      &                  wsccor*gsccorx(j,i)
313      &                 +wscloc*gsclocx(j,i)
314      &                 +wliptran*gliptranx(j,i)
315      &                 +welec*gshieldx(j,i)
316      &                 +wcorr*gshieldx_ec(j,i)
317      &                 +wturn3*gshieldx_t3(j,i)
318      &                 +wturn4*gshieldx_t4(j,i)
319      &                 +wel_loc*gshieldx_ll(j,i)
320 c     &                 +wtube*gg_tube_sc(j,i)
321
322
323
324         enddo
325       enddo 
326 #ifdef DEBUG
327       write (iout,*) "gloc before adding corr"
328       do i=1,4*nres
329         write (iout,*) i,gloc(i,icg)
330       enddo
331 #endif
332       do i=1,nres-3
333         gloc(i,icg)=gloc(i,icg)+wcorr*gcorr_loc(i)
334      &   +wcorr5*g_corr5_loc(i)
335      &   +wcorr6*g_corr6_loc(i)
336      &   +wturn4*gel_loc_turn4(i)
337      &   +wturn3*gel_loc_turn3(i)
338      &   +wturn6*gel_loc_turn6(i)
339      &   +wel_loc*gel_loc_loc(i)
340       enddo
341 #ifdef DEBUG
342       write (iout,*) "gloc after adding corr"
343       do i=1,4*nres
344         write (iout,*) i,gloc(i,icg)
345       enddo
346 #endif
347 #ifdef MPI
348       if (nfgtasks.gt.1) then
349         do j=1,3
350           do i=1,nres
351             gradbufc(j,i)=gradc(j,i,icg)
352             gradbufx(j,i)=gradx(j,i,icg)
353           enddo
354         enddo
355         do i=1,4*nres
356           glocbuf(i)=gloc(i,icg)
357         enddo
358 c#define DEBUG
359 #ifdef DEBUG
360       write (iout,*) "gloc_sc before reduce"
361       do i=1,nres
362        do j=1,1
363         write (iout,*) i,j,gloc_sc(j,i,icg)
364        enddo
365       enddo
366 #endif
367 c#undef DEBUG
368         do i=1,nres
369          do j=1,3
370           gloc_scbuf(j,i)=gloc_sc(j,i,icg)
371          enddo
372         enddo
373         time00=MPI_Wtime()
374 c        write (iout,*) "GRAD: MPI_BARRIER"
375 c        call flush(iout)
376         call MPI_Barrier(FG_COMM,IERR)
377         time_barrier_g=time_barrier_g+MPI_Wtime()-time00
378         time00=MPI_Wtime()
379         call MPI_Reduce(gradbufc(1,1),gradc(1,1,icg),3*nres,
380      &    MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR)
381         call MPI_Reduce(gradbufx(1,1),gradx(1,1,icg),3*nres,
382      &    MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR)
383         call MPI_Reduce(glocbuf(1),gloc(1,icg),4*nres,
384      &    MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR)
385         time_reduce=time_reduce+MPI_Wtime()-time00
386         call MPI_Reduce(gloc_scbuf(1,1),gloc_sc(1,1,icg),3*nres,
387      &    MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR)
388         time_reduce=time_reduce+MPI_Wtime()-time00
389 c#define DEBUG
390 #ifdef DEBUG
391       write (iout,*) "gloc_sc after reduce"
392       do i=1,nres
393        do j=1,1
394         write (iout,*) i,j,gloc_sc(j,i,icg)
395        enddo
396       enddo
397 #endif
398 c#undef DEBUG
399 #ifdef DEBUG
400       write (iout,*) "gloc after reduce"
401       do i=1,4*nres
402         write (iout,*) i,gloc(i,icg)
403       enddo
404 #endif
405       endif
406 #endif
407 #ifdef DEBUG
408       write (iout,*) "gradc gradx gloc"
409       do i=1,nres
410         write (iout,'(i5,3f10.5,5x,3f10.5,5x,f10.5)') 
411      &   i,(gradc(j,i,icg),j=1,3),(gradx(j,i,icg),j=1,3),gloc(i,icg)
412       enddo 
413 #endif
414 #ifdef TIMING
415       time_sumgradient=time_sumgradient+MPI_Wtime()-time01
416 #endif
417       return
418       end
419 c---------------------------------------------------------------------------
420       subroutine sum_gradient_compon
421       implicit none
422       include 'DIMENSIONS'
423       include "DIMENSIONS.ZSCOPT"
424 #ifndef ISNAN
425       external proc_proc
426 #ifdef WINPGI
427 cMS$ATTRIBUTES C ::  proc_proc
428 #endif
429 #endif
430 #ifdef MPI
431       include 'mpif.h'
432       include "COMMON.SETUP"
433       integer ierr,ierror
434 #endif
435       double precision gradbufc(3,-1:maxres),gradbufc_sum(3,-1:maxres),
436      &  glocbuf(max_ene,4*maxres)
437       include 'COMMON.IOUNITS'
438       include 'COMMON.FFIELD'
439       include 'COMMON.DERIV'
440       include 'COMMON.INTERACT'
441       include 'COMMON.SBRIDGE'
442       include 'COMMON.CHAIN'
443       include 'COMMON.VAR'
444       include 'COMMON.CONTROL'
445       include 'COMMON.TIME1'
446       include 'COMMON.SCCOR'
447       include "COMMON.NAMES"
448       include "COMMON.ENERGIES"
449       integer i,j,k,ii,iene
450       double precision time_reduce,time00
451 #ifdef TIMING
452       time01=MPI_Wtime()
453 #endif
454 C Sum up the components of the Cartesian gradient.
455 C
456       gcompon=0.0d0
457       gcompon_loc=0.0d0
458       gcomponx=0.0d0
459       do i=1,nct
460         do j=1,3
461 c Part in coordinates
462           gcompon(1,j,i)=gvdwc(j,i)
463           gcompon(2,j,i)=gvdwc_scp(j,i)+gvdwc_scpp(j,i)
464           gcompon(3,j,i)=gelc_long(j,i)
465           gcompon(16,j,i)=gvdwpp(j,i)
466           gcompon(7,j,i)=gel_loc_long(j,i)
467           gcompon(4,j,i)=gradcorr_long(j,i)
468           gcompon(5,j,i)=gradcorr5_long(j,i)
469           gcompon(6,j,i)=gradcorr6_long(j,i)
470           gcompon(10,j,i)=gcorr6_turn_long(j,i)
471           gcompon(15,j,i)=ghpbc(j,i)
472           gcompon(22,j,i)=gliptranc(j,i)
473 c Part in vectors
474           gcompon_loc(2,j,i)=0.5d0*gvdwc_scpp(j,i)
475           gcompon_loc(3,j,i)=0.5d0*gelc_long(j,i)+gelc(j,i)
476           gcompon_loc(16,j,i)=0.5d0*gvdwpp(j,i)
477           gcompon_loc(4,j,i)=0.5d0*gradcorr_long(j,i)+gradcorr(j,i)
478           gcompon_loc(5,j,i)=0.5d0*gradcorr5_long(j,i)+gradcorr5(j,i)
479           gcompon_loc(6,j,i)=0.5d0*gradcorr6_long(j,i)+gradcorr6(j,i)
480           gcompon_loc(7,j,i)=0.5d0*gel_loc_long(j,i)+gel_loc(j,i)
481           gcompon_loc(8,j,i)=gcorr3_turn(j,i)
482           gcompon_loc(9,j,i)=gcorr4_turn(j,i)
483           gcompon_loc(10,j,i)=0.5d0*gcorr6_turn_long(j,i)
484      &        +gcorr6_turn(j,i)
485           gcompon_loc(12,j,i)=gscloc(j,i)
486           gcompon_loc(17,j,i)=gradb(j,i)
487           gcompon_loc(19,j,i)=gsccorc(j,i)
488           gcompon_loc(22,j,i)=gliptranc(j,i)
489 c sidechain components
490           gcomponx(1,j,i)=gvdwx(j,i)
491           gcomponx(2,j,i)=gradx_scp(j,i)
492           gcomponx(12,j,i)=gsclocx(j,i)
493           gcomponx(17,j,i)=gradbx(j,i)
494           gcomponx(19,j,i)=gsccorx(j,i)
495           gcomponx(22,j,i)=gliptranx(j,i)
496           if (shield_mode.gt.0) then
497 c Part in coordinates
498             gcompon(3,j,i)=gcompon(3,j,i)+gshieldc(j,i)
499             gcompon(4,j,i)=gcompon(4,j,i)+gshieldc_ec(j,i)
500             gcompon(7,j,i)=gcompon(7,j,i)+gshieldc_ll(j,i)
501             gcompon(8,j,i)=gcompon(8,j,i)+gshieldc_t3(j,i)
502             gcompon(9,j,i)=gcompon(9,j,i)+gshieldc_t4(j,i)
503 c part in vectors
504             gcompon_loc(3,j,i)=gshieldc(j,i)+gshieldc_loc(j,i) 
505             gcompon_loc(4,j,i)=gcompon_loc(4,j,i)+gshieldc_ec(j,i)
506      &         +gshieldc_loc_ec(j,i)
507             gcompon_loc(7,j,i)=gcompon_loc(7,j,i)+gshieldc_ll(j,i)
508      &         +gshieldc_loc_ll(j,i)
509             gcompon_loc(8,j,i)=gcompon_loc(8,j,i)+gshieldc_t3(j,i)
510      &         +gshieldc_loc_t3(j,i)
511             gcompon_loc(9,j,i)=gcompon_loc(9,j,i)+gshieldc_t4(j,i)
512      &         +gshieldc_loc_t4(j,i)
513 c Part in sidechains
514             gcomponx(3,j,i)=gshieldx(j,i)
515             gcomponx(4,j,i)=gshieldx_ec(j,i)
516             gcomponx(7,j,i)=gshieldx_ec(j,i)
517             gcomponx(8,j,i)=gshieldx_t3(j,i)
518             gcomponx(9,j,i)=gshieldx_t4(j,i)
519           endif
520 #ifndef SPLITELE
521           gcompon(3,j,i)=gcompon(3,j,i)+gcompon(16,j,i)
522           gcompon_loc(3,j,i)=gcompon_loc(3,j,i)+gcompon_loc(16,j,i) 
523 #endif
524         enddo
525       enddo
526
527       do i=1,nres-3
528 c        gloc(i,icg)=gloc(i,icg)+wcorr*gcorr_loc(i)
529 c     &   +wcorr5*g_corr5_loc(i)
530 c     &   +wcorr6*g_corr6_loc(i)
531 c     &   +wturn4*gel_loc_turn4(i)
532 c     &   +wturn3*gel_loc_turn3(i)
533 c     &   +wturn6*gel_loc_turn6(i)
534 c     &   +wel_loc*gel_loc_loc(i)
535         gloc_compon(4,i)=gcorr_loc(i)
536         gloc_compon(5,i)=g_corr5_loc(i)
537         gloc_compon(6,i)=g_corr6_loc(i)
538         gloc_compon(7,i)=gel_loc_loc(i)
539         gloc_compon(8,i)=gel_loc_turn3(i)
540         gloc_compon(9,i)=gel_loc_turn4(i) 
541         gloc_compon(10,i)=gel_loc_turn6(i)
542 c     &   +wsccor*gsccor_loc(i)
543 c BYLA ROZNICA Z CLUSTER< OSTATNIA LINIA DODANA
544       enddo
545 #ifdef DEBUG
546       write (iout,*) "Gradient components right after assignment"
547       do iene=1,n_ene
548
549         write (iout,'(a,i3,1x,a)') "Component",iene,ename(iene)
550
551         do i=1,nres
552           write (iout,'(a4,i5,3e15.5,5x,e15.5,5x,3e15.5,5x,e15.5)') 
553      &      restyp(itype(i)),i,(gcompon(iene,j,i),j=1,3),
554      &      gloc_compon(iene,i),(gcomponx(iene,j,i),j=1,3),
555      &      gloc_compon(iene,nres+i)
556         enddo
557
558       enddo
559 #endif
560 #ifdef MPI
561 C FG slaves call the following matching MPI_Bcast in ERGASTULUM
562       if (nfgtasks.gt.1 .and. fg_rank.eq.0) 
563      &    call MPI_Bcast(1,1,MPI_INTEGER,king,FG_COMM,IERROR)
564 #endif
565
566       do iene=1,n_ene
567
568       gradbufc = 0.0d0
569       gradbufc_sum=0.0d0
570
571       do i=0,nct
572         do j=1,3
573           gradbufc_sum(j,i)=gcompon(iene,j,i)
574         enddo
575       enddo 
576       do j=1,3
577         gradbufc(j,nres-1)=gradbufc_sum(j,nres)
578       enddo
579       do i=nres-2,-1,-1
580         do j=1,3
581           gradbufc(j,i)=gradbufc(j,i+1)+gradbufc_sum(j,i+1)
582         enddo
583       enddo
584 #ifdef DEBUG
585       write (iout,*) "gradbufc after summing"
586       do i=1,nres
587         write (iout,'(i3,3f10.5)') i,(gradbufc(j,i),j=1,3)
588       enddo
589       call flush(iout)
590 #endif
591       do k=1,3
592         gradbufc(k,nres)=0.0d0
593       enddo
594       do i=-1,nct
595         do j=1,3
596           gcompon(iene,j,i)=gradbufc(j,i)+gcompon_loc(iene,j,i)
597         enddo
598       enddo 
599
600       enddo ! iprot
601 #ifdef MPI
602       if (nfgtasks.gt.1) then
603         gcompon_loc = gcompon
604         time00=MPI_Wtime()
605         call MPI_Reduce(gcompon_loc(1,1,1),gcompon(1,1,1),
606      &    3*nres*max_ene,
607      &    MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR)
608         gcompon_loc = gcomponx
609         call MPI_Reduce(gcompon_loc(1,1,1),gcomponx(1,1,1),
610      &    3*nres*max_ene,
611      &    MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR)
612         glocbuf = gloc_compon
613         call MPI_Reduce(glocbuf(1,1),gloc_compon(1,1),4*nres*max_ene,
614      &    MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR)
615         time_reduce=time_reduce+MPI_Wtime()-time00
616       endif
617 c#define DEBUG
618 #endif
619 #ifdef DEBUG
620       write (iout,*) "Gradient components in sum_gradient_compon"
621       do iene=1,n_ene
622
623         write (iout,'(a,i3,1x,a)') "Component",iene,ename(iene)
624
625         do i=1,nres
626           write (iout,'(a4,i5,3e15.5,5x,e15.5,5x,3e15.5,5x,e15.5)') 
627      &      restyp(itype(i)),i,(gcompon(iene,j,i),j=1,3),
628      &      gloc_compon(iene,i),(gcomponx(iene,j,i),j=1,3),
629      &      gloc_compon(iene,nres+i)
630         enddo
631
632       enddo
633 #endif
634 #ifdef TIMING
635       time_sumgradient=time_sumgradient+MPI_Wtime()-time01
636 #endif
637       return
638       end