update new files
[unres.git] / source / maxlik / src_MD_T_maxlik-NEWCORR-PMF-PDB / make_list_sc.F
1       subroutine make_list(lprn,first_call,nvarr,x)
2       implicit none
3       include "DIMENSIONS"
4       include "DIMENSIONS.ZSCOPT"
5 #ifdef MPI
6       include "mpif.h"
7       integer IERROR,ErrCode,Status(MPI_STATUS_SIZE,10)
8       integer req(10),msg_in(5),msg_out(5),address,size
9       character*1 buffer(8*(2*maxstr_proc*nntyp+8000))
10       include "COMMON.MPI"
11 #endif
12       include "COMMON.WEIGHTS"
13       include "COMMON.WEIGHTDER"
14       include "COMMON.ENERGIES"
15       include "COMMON.IOUNITS"
16       include "COMMON.VMCPAR"
17       include "COMMON.NAMES"
18       include "COMMON.INTERACT"
19       include "COMMON.TIME1"
20       include "COMMON.CHAIN"
21       include "COMMON.PROTFILES"
22       include "COMMON.VAR"
23       include "COMMON.GEO"
24       include "COMMON.OPTIM"
25       include "COMMON.CLASSES"
26       include "COMMON.COMPAR"
27       include "COMMON.TORSION"
28 C Define local variables
29       integer nvarr
30       integer i,ii,iii,ibatch,kkk,jj,j,k,kk,l,iprot,ib,nn
31       integer ipass_conf,istart_conf,iend_conf,Previous,Next
32       double precision energia(0:max_ene)
33       double precision etoti,elowesti,dene
34       double precision tcpu_ini,tcpu_fin,tcpu
35       double precision elowest_t(2,maxT),etot_aux,enepsjk,
36      &  emini,elowest_aux(2,maxT)
37       integer iroof,icant
38       external iroof,icant
39       logical lprn,first_call
40       logical*1 lflag(maxstr)
41       integer ipermut(maxstr)
42       integer list_conf_(maxstr,maxprot)
43       double precision x(nvarr)
44       double precision ftune_eps
45       external ftune_eps
46 #ifdef MPI
47       double precision e_total_(maxstr_proc)
48 #endif
49       write (iout,*) "Making the worklist of conformations."
50       write (iout,*) "enecut",(enecut(i),i=1,nprot)
51       write (iout,*) "first_call ",first_call," nvarr",nvarr
52       tcpu_ini = tcpu()
53       do iprot=1,nprot
54         do i=1,ntot(iprot)
55           list_conf(i,iprot)=i
56         enddo
57         ntot_work(iprot)=ntot(iprot)
58       enddo
59 #ifdef MPI
60 c
61 c Divide the whole database between processors
62 c
63       Previous = me1-1
64       Next = me1+1
65       if (Previous.lt.0) Previous = Nprocs-1
66       if (Next.ge.Nprocs) Next = 0
67       call work_partition(lprn)
68 #endif
69 c
70 c Loop over proteins
71 c
72       DO iprot=1,nprot
73
74        call restore_molinfo(iprot)
75 c      write (iout,*) "Processor",me," iprot",iprot,
76 c     & " indstart",indstart(me1,iprot)," indend",indend(me1,iprot),
77 c     & " init_ene",init_ene," mod_fourier",mod_fourier(0),
78 c     & " mod_elec",mod_elec," mod_scp",mod_scp
79       call restore_molinfo(iprot)
80 c
81 c Loop over the conformations of protein IPROT assigned to the current processor. 
82 c The conformations are read off a DA scratchfile and processed in passes, each
83 c of which requires no more than MAXSTR_PROC conformations to be stored in memory
84 c simultaneously.
85 c
86       if (.not.init_ene .and. mod_fourier(nloctyp).eq.0 
87      &  .and. mod_fouriertor(nloctyp).eq.0
88      &  .and. .not. mod_elec .and. .not. mod_scp) then
89         open (ientin,file=benefiles(iprot),status="old",
90      &    form="unformatted",access="direct",recl=lenrec_ene(iprot))
91       else
92         open (icbase,file=bprotfiles(iprot),status="old",
93      &    form="unformatted",access="direct",recl=lenrec(iprot))
94 c Change AL 12/30/2017
95         if (.not. mod_other_params) 
96      &  open (ientout,file=benefiles(iprot),status="old",
97      &    form="unformatted",access="direct",recl=lenrec_ene(iprot))
98       endif
99
100 #ifdef MPI
101       nchunk_ene(iprot) = iroof(scount(me1,iprot),maxstr)
102       nchunk_conf(iprot) = iroof(scount(me1,iprot),maxstr_proc)
103       write (iout,*)"Protein",iprot," energy evaluation in",
104      &   nchunk_conf(iprot)," passes."
105       ipass_conf=0
106       jj=0
107       do i=indstart(me1,iprot),indend(me1,iprot),maxstr_proc
108       ipass_conf=ipass_conf+1
109       write (iout,*) "MAKE_LIST: Pass",ipass_conf
110       istart_conf=i
111       iend_conf=min0(i+maxstr_proc-1,indend(me1,iprot))
112 #else
113       nchunk_ene(iprot) = iroof(ntot(iprot),maxstr)
114       nchunk_conf(iprot) = iroof(ntot(iprot),maxstr_proc)
115       write (iout,*)"Protein",iprot," energy evaluation in",
116      &   nchunk_conf," passes."
117       ipass_conf=0
118       do i=1,ntot(iprot),maxstr_proc
119       ipass_conf=ipass_conf+1
120       write (iout,*) "MAKE_LIST: Pass",ipass_conf
121       istart_conf=i
122       iend_conf=min0(i+maxstr_proc-1,ntot(iprot))
123 #endif
124       ii=0
125 c
126 c Read the chunk of conformations off a DA scratchfile.
127 c
128       if (.not.init_ene .and. mod_fourier(nloctyp).eq.0
129      &  .and. mod_fouriertor(nloctyp).eq.0 .and. .not. mod_elec 
130      &  .and. .not. mod_scp) then
131 c
132 c If energy components have been pre-computed read them off a DA file.
133 c
134         call daread_ene(iprot,istart_conf,iend_conf)
135         do iii=istart_conf,iend_conf
136           ii=ii+1
137           jj=jj+1
138           if (mod_side) then
139             enetb(ii,1,iprot)=0.0d0
140             do j=1,ntyp
141               do k=1,j
142                 enetb(ii,1,iprot)=enetb(ii,1,iprot)+ftune_eps(eps(j,k))*
143      &                 eneps(1,icant(j,k),ii,iprot)+
144      &                 eps(j,k)*eneps(2,icant(j,k),ii,iprot)
145               enddo
146             enddo
147           endif
148           etot_aux=0.0d0
149           do j=1,n_ene
150             etot_aux=etot_aux+ww(j)*enetb(ii,j,iprot)
151           enddo
152           e_total(iii,iprot)=etot_aux
153 #ifdef DEBUG
154           write (iout,'(i5,16(1pe12.4))') iii,
155      &    (enetb(ii,j,iprot),j=1,n_ene),e_total(iii,iprot)
156           call flush(iout)
157 #endif
158         enddo
159         if (first_call .and. mod_side) then
160           write (iout,*) "Callig x2w"
161           call flush(iout)
162           call x2w(nvarr,x_orig)
163           write (iout,*) "After x2w"
164           call flush(iout)
165           ii=0
166           do iii=istart_conf,iend_conf
167             ii=ii+1
168             enetb_oorig(ii,1,iprot)=0.0d0
169             do j=1,ntyp
170               do k=1,j
171                 enetb_oorig(ii,1,iprot)=enetb_oorig(ii,1,iprot)+
172      &              ftune_eps(eps(j,k))*eneps(1,icant(j,k),ii,iprot)+
173      &              eps(j,k)*eneps(2,icant(j,k),ii,iprot)  
174               enddo
175             enddo
176             do j=2,n_ene
177               enetb_oorig(ii,j,iprot)=enetb(ii,j,iprot)
178             enddo
179           enddo
180           call x2w(nvarr,x)
181         endif
182       else
183         call daread_ccoords(iprot,istart_conf,iend_conf)
184 c
185 c Compute the energies of the conformations currently in memory and compute
186 c the lowest energies.
187 c
188         do iii=istart_conf,iend_conf
189           ii=ii+1
190           jj=jj+1
191           call restore_ccoords(iprot,ii)
192           call int_from_cart1(.false.)
193 #ifdef DEBUG
194           write (iout,*) "Before etotal",iii,i
195           call flush(iout)
196 #endif
197           call etotal(energia(0))
198 #ifdef DEBUG
199           write (iout,*) "After etotal",i
200           call flush(iout)
201           call enerprint(energia(0))
202 c          write (iout,'(i5,16(1pe12.4))') i,
203 c     &    (energia(j),j=1,n_ene),energia(0)
204           call flush(iout)
205 #endif
206 #ifdef DEBUG
207           write (iout,*) "Conformation:",i
208           write (iout,'(8f10.5)') ((c(l,k),l=1,3),k=1,nres)
209           write (iout,'(8f10.5)') ((c(l,k+nres),l=1,3),k=nnt,nct)
210           write (iout,'(8f10.4)') (vbld(k),k=2,nres)
211           write (iout,'(8f10.4)') (vbld(k+nres),k=nnt,nct)
212           write (iout,'(8f10.4)') (rad2deg*theta(k),k=3,nres)
213           write (iout,'(8f10.4)') (rad2deg*phi(k),k=4,nres)
214           write (iout,'(8f10.4)') (rad2deg*alph(k),k=2,nres-1)
215           write (iout,'(8f10.4)') (rad2deg*omeg(k),k=2,nres-1)
216           call enerprint(energia(0))
217 #endif
218           e_total(iii,iprot)=energia(0)
219           do j=1,n_ene
220             enetb(ii,j,iprot)=energia(j)
221           enddo
222           do j=1,ntyp
223             do k=1,j
224              eneps(1,icant(j,k),ii,iprot)=eneps_temp(1,icant(j,k))
225              eneps(2,icant(j,k),ii,iprot)=eneps_temp(2,icant(j,k))
226             enddo
227           enddo
228 #ifdef DEBUG
229           write (iout,'(i5,20(1pe12.4))') iii,
230      &    (energia(j),j=1,n_ene),energia(0),eini(iii,iprot),
231      &    entfac(iii,iprot)
232           call flush(iout)
233 #endif
234           if (energia(0).ge.1.0d99) then
235             write (iout,*) 
236      &    "MAKE_LIST:CHUJ NASTAPIL in energy evaluation for",
237      &    " point",i,". Probably NaNs in some of the energy components."
238             write (iout,*) "The components of the energy are:"
239             call enerprint(energia(0))
240             write (iout,*) "Conformation:",i
241             write (iout,'(8f10.5)') ((c(l,k),l=1,3),k=1,nres)
242             write (iout,'(8f10.5)') ((c(l,k+nres),l=1,3),k=nnt,nct)
243             write (iout,'(8f10.4)') (vbld(k),k=2,nres)
244             write (iout,'(8f10.4)') (vbld(k+nres),k=nnt,nct)
245             write (iout,'(8f10.4)') (rad2deg*theta(k),k=3,nres)
246             write (iout,'(8f10.4)') (rad2deg*phi(k),k=4,nres)
247             write (iout,'(8f10.4)') (rad2deg*alph(k),k=2,nres-1)
248             write (iout,'(8f10.4)') (rad2deg*omeg(k),k=2,nres-1)
249             write (iout,*) "Calculation terminated at this point.",
250      &       " Check the database of conformations"
251 #ifdef MPI
252             call MPI_Abort(MPI_COMM_WORLD,ErrCode,IERROR)
253 #endif
254             stop "SEVERE error in energy calculation" 
255           endif
256         enddo ! iii
257         if (first_call) then
258           do ii=1,iend_conf-istart_conf+1
259             do j=1,n_ene
260               enetb_orig(ii,j,iprot)=enetb(ii,j,iprot)
261             enddo
262           enddo
263           ii=0
264           call x2w(nvarr,x_orig)
265 #ifdef DEBUG
266           write (iout,*) "x,xorig"
267           do k=1,nvarr
268             write (iout,'(i5,2f10.5)') k,x(k),x_orig(k)
269           enddo
270 #endif
271           do iii=istart_conf,iend_conf
272             ii=ii+1
273             call restore_ccoords(iprot,ii)
274             call int_from_cart1(.false.)
275             call etotal(energia(0))
276 #ifdef DEBUG
277             write (iout,*) "Conformation:",iii,ii
278             write (iout,'(8f10.5)') ((c(l,k),l=1,3),k=1,nres)
279             write (iout,'(8f10.5)') ((c(l,k+nres),l=1,3),k=nnt,nct)
280             write (iout,'(8f10.4)') (vbld(k),k=2,nres)
281             write (iout,'(8f10.4)') (vbld(k+nres),k=nnt,nct)
282             write (iout,'(8f10.4)') (rad2deg*theta(k),k=3,nres)
283             write (iout,'(8f10.4)') (rad2deg*phi(k),k=4,nres)
284             write (iout,'(8f10.4)') (rad2deg*alph(k),k=2,nres-1)
285             write (iout,'(8f10.4)') (rad2deg*omeg(k),k=2,nres-1)
286             call enerprint(energia(0))
287 #endif
288             do j=1,n_ene
289               enetb_oorig(ii,j,iprot)=energia(j)
290             enddo
291 #ifdef DEBUG
292             write (iout,'(2i5,20(1pe12.4))') iii,ii,
293      &      (energia(j),j=1,n_ene),energia(0)
294             write (iout,'(2i5,20(1pe12.4))') iii,ii,
295      &      (enetb(ii,j,iprot),j=1,n_ene)
296             call flush(iout)
297 #endif
298             if (energia(0).ge.1.0d99) then
299               write (iout,*) "CHUJ NASTAPIL in energy evaluation for",
300      &   " point",i,". Probably NaNs in some of the energy components."
301               write (iout,*) "The components of the energy are:"
302               call enerprint(energia(0))
303               write (iout,*) "Calculation terminated at this point.",
304      &         " Check the database of conformations"
305 #ifdef MPI
306               call MPI_Abort(MPI_COMM_WORLD,ErrCode,IERROR)
307 #endif
308               stop "SEVERE error in energy calculation" 
309             endif
310           enddo ! iii
311           write (iout,*) "MAKE_LIST Callig x2w"
312           call flush(iout)
313           call x2w(nvarr,x)
314           write (iout,*) "After x2w"
315           call flush(iout)
316         endif
317         write (iout,*) "make_list: calling dawrite_ene"
318         write (iout,*) "istart_conf",istart_conf,
319      &    " iend_conf",iend_conf
320 C Change AL 12/30/2017
321         if (.not. mod_other_params) 
322      &  call dawrite_ene(iprot,istart_conf,iend_conf,ientout)
323 #ifdef MPI
324 c Distribute energy components through ring
325         call MPI_Barrier(WHAM_COMM,IERROR)
326         write (iout,*) "Processes synchronized in make_list"
327         call flush(iout)
328         msg_out(1)= 5*me1+1000*ipass_conf+1
329         msg_out(2)= 5*me1+1000*ipass_conf+2
330         msg_out(3)= 5*me1+1000*ipass_conf+3
331         msg_out(4)= 5*me1+1000*ipass_conf+4
332         msg_out(5)= 5*me1+1000*ipass_conf+5
333         do iii=1,Nprocs-1
334 c Send the current energy tables to the right neighbor
335 c Receive the energy tables produced by processor kkk from the left neighbor
336           kkk = mod(me1-iii+NProcs,Nprocs)
337           msg_in(1)= 5*kkk+1000*ipass_conf+1
338           msg_in(2)= 5*kkk+1000*ipass_conf+2
339           msg_in(3)= 5*kkk+1000*ipass_conf+3
340           msg_in(4)= 5*kkk+1000*ipass_conf+4
341           msg_in(5)= 5*kkk+1000*ipass_conf+5
342           write (iout,*) "me1",me1," iii",iii," Previous",Previous,
343      &     " Next",Next," kkk",kkk
344           write (iout,*) "msg_in",msg_in
345           write (iout,*) "msg_out",msg_out
346           call flush(iout)
347           write (iout,*) "Processor",me1," Start Send and receive"
348           call flush(iout)
349           call MPI_Send(istart_conf,1,MPI_INTEGER,Next,msg_out(1),
350      &     WHAM_COMM,IERROR)
351           write (iout,*) "Send",msg_out(1)," complete"
352           call flush(iout)
353           call MPI_Recv(istart_conf,1,MPI_INTEGER,Previous,
354      &     msg_in(1),WHAM_COMM,STATUS(1,6),IERROR)
355           write (iout,*) "Recv",msg_in(1)," complete"
356           call flush(iout)
357           call MPI_Send(iend_conf,1,MPI_INTEGER,Next,msg_out(2),
358      &     WHAM_COMM,IERROR)
359           write (iout,*) "Send",msg_out(2)," complete"
360           call flush(iout)
361           call MPI_Recv(iend_conf,1,MPI_INTEGER,Previous,
362      &     msg_in(2),WHAM_COMM,STATUS(1,7),IERROR)
363           write (iout,*) "Recv",msg_in(2)," complete"
364           call flush(iout)
365           call MPI_Buffer_Attach(buffer(1),8*(2*maxstr_proc*n_ene+800),
366      &     IERROR)
367           call MPI_BSend(enetb(1,1,iprot),maxstr_proc*n_ene,
368      &    MPI_DOUBLE_PRECISION,Next,msg_out(3),
369      &     WHAM_COMM,IERROR)
370           write (iout,*) "Send",msg_out(3)," complete (enetb)"
371           call flush(iout)
372           call MPI_Recv(enetb(1,1,iprot),maxstr_proc*n_ene,
373      &    MPI_DOUBLE_PRECISION,Previous,msg_in(3),WHAM_COMM,
374      &     STATUS(1,8),IERROR)
375           write (iout,*) "Recv",msg_in(3)," complete (enetb)"
376           call flush(iout)
377           call MPI_Buffer_Detach(address,size,IERROR)
378 c          write (iout,*) "MPI_Buffer_Detach complete (enetb)"
379 c          call flush(iout)
380           call MPI_Buffer_Attach(buffer(1),8*(2*maxstr_proc*nntyp+800),
381      &     IERROR)
382 c          write (iout,*) "MPI_Buffer_Attach complete (eneps)"
383 c          call flush(iout)
384           call MPI_BSend(eneps(1,1,1,iprot),2*maxstr_proc*nntyp,
385      &     MPI_DOUBLE_PRECISION,Next,msg_out(4),
386      &     WHAM_COMM,IERROR)
387           write (iout,*) "Send",msg_out(4)," complete (eneps)"
388           call flush(iout)
389           call MPI_Recv(eneps(1,1,1,iprot),2*maxstr_proc*nntyp,
390      &    MPI_DOUBLE_PRECISION,Previous,msg_in(4),WHAM_COMM,
391      &     STATUS(1,9),IERROR)
392           write (iout,*) "Recv",msg_in(4)," complete (eneps)"
393           call flush(iout)
394           call MPI_Buffer_Detach(address,size,IERROR)
395           call MPI_Buffer_Attach(buffer(1),8*(2*maxstr_proc*maxnatlike
396      &     +800),IERROR)
397           call MPI_BSend(nu(1,1,1,iprot),
398      &    maxstr_proc*maxnatlike*maxdimnat,
399      &    MPI_DOUBLE_PRECISION,Next,msg_out(5),
400      &     WHAM_COMM,IERROR)
401           write (iout,*) "Send",msg_out(5)," complete (nu)"
402           call flush(iout)
403           call MPI_Recv(nu(1,1,1,iprot),
404      &    maxstr_proc*maxnatlike*maxdimnat,
405      &    MPI_DOUBLE_PRECISION,Previous,msg_in(5),WHAM_COMM,
406      &     STATUS(1,10),IERROR)
407           write (iout,*) "Recv",msg_in(5)," complete (nu)"
408           call flush(iout)
409           call MPI_Buffer_Detach(address,size,IERROR)
410           write (iout,*) "Send and receive complete"
411           call flush(iout)
412           if (.not. mod_other_params) then
413           write (iout,*) "Processor",me1," calling dawrite_ene",
414      &     " istart_conf",istart_conf," iend_conf",iend_conf
415           call flush(iout)
416           endif
417           if (first_call) then
418 c            write (iout,*) "assignment of enetb_orig"
419             do ii=1,iend_conf-istart_conf+1
420               do j=1,n_ene
421                 enetb_orig(ii,j,iprot)=enetb(ii,j,iprot)
422               enddo
423 c              write (iout,'(i5,20f8.2)') 
424 c     &         ii,(enetb_orig(ii,j,iprot),j=1,n_ene)
425             enddo
426           endif
427           if (.not. mod_other_params) 
428      &      call dawrite_ene(iprot,istart_conf,iend_conf,ientout)
429           do k=1,5
430             msg_out(k)=msg_in(k)
431           enddo
432         enddo
433 #endif
434       endif
435       enddo ! i
436       if (.not.init_ene .and. mod_fourier(nloctyp).eq.0
437      &  .and. mod_fouriertor(nloctyp).eq.0 
438      &  .and. .not. mod_elec .and. .not. mod_scp) then
439         close (ientin)
440       else
441         close (icbase)
442 c Change AL 12/30/2017
443         if (.not. mod_other_params) close (ientout)
444       endif
445
446       ENDDO ! iprot
447       init_ene=.false.
448 c Lowest free energies of structural classes
449       DO IPROT=1,NPROT
450
451
452       IF (NCHUNK_CONF(IPROT).EQ.1) THEN
453
454 #ifdef MPI
455       do i=1,ntot_work(iprot)
456         i2ii(i,iprot)=0
457       enddo
458       ii=0
459       do i=indstart(me1,iprot),indend(me1,iprot)
460         ii=ii+1
461         i2ii(i,iprot)=ii
462       enddo
463       istart_conf=indstart(me1,iprot)
464       iend_conf=indend(me1,iprot)
465 #else
466       do i=1,ntot_work(iprot)
467         i2ii(i,iprot)=i
468       endif
469       istart_conf=1
470       iend_conf=ntot_work(iprot)
471 #endif
472 #ifdef DEBUG
473       write (iout,*) "i2ii at make_list"
474       do i=1,ntot_work(iprot)
475         write (iout,*) "i",i," i2ii",i2ii(i,iprot)
476       enddo
477       call flush(iout)
478 #endif
479       if (.not. mod_other_params) then
480       open (ientin,file=benefiles(iprot),status="old",
481      &    form="unformatted",access="direct",recl=lenrec_ene(iprot))
482 c Change AL 12/30/2017
483       call daread_ene(iprot,istart_conf,iend_conf)
484       endif
485       call emin_search(iprot)
486
487       ELSE
488
489       if (.not. mod_other_params) 
490      &open (ientin,file=benefiles(iprot),status="old",
491      &    form="unformatted",access="direct",recl=lenrec_ene(iprot))
492       ipass_conf=0
493 #ifdef MPI
494       do istart_conf=indstart(me1,iprot),indend(me1,iprot),maxstr_proc
495         iend_conf=min0(istart_conf+maxstr_proc-1,indend(me1,iprot))
496 #else
497       do istart_conf=1,ntot_work(iprot),maxstr_proc
498         iend_conf=min0(istart_conf+maxstr_proc-1,ntot_work(iprot))
499 #endif
500 c
501 c Read the chunk of energies and derivatives off a DA scratchfile.
502 c
503         ipass_conf=ipass_conf+1
504         do i=1,ntot_work(iprot)
505           i2ii(i,iprot)=0
506         enddo
507         ii=0
508         do i=istart_conf,iend_conf
509           ii=ii+1
510           i2ii(i,iprot)=ii
511         enddo
512 #ifdef DEBUG
513         write (iout,*) "ipass_conf",ipass_conf,
514      &    " istart_conf",istart_conf," iend_conf",iend_conf
515         do i=1,ntot_work(iprot)
516           write (iout,*) "i",i," i2ii",i2ii(i,iprot)
517         enddo
518         call flush(iout)
519 #endif
520 c Change AL 12/30/2017
521         if (.not. mod_other_params)
522      &  call daread_ene(iprot,istart_conf,iend_conf)
523         call emin_search(iprot)
524       enddo
525
526       close(ientin)
527       ENDIF
528
529       ENDDO
530 #ifdef MPI
531 c Complete the calculation of the lowest energies over all classes and
532 c distribute the values to all procs
533       do iprot=1,nprot
534         do ibatch=1,natlike(iprot)+2
535 #ifdef DEBUG
536           do ib=1,nbeta(ibatch,iprot)
537             write (iout,'(7hELOWEST,3i3,f15.3,i12)') iprot,ibatch,ib,
538      &        elowest(ib,ibatch,iprot),ind_lowest(ib,ibatch,iprot)
539           enddo
540 #endif
541           do ib=1,nbeta(ibatch,iprot)
542             elowest_aux(1,ib)=elowest(ib,ibatch,iprot)
543             elowest_aux(2,ib)=ind_lowest(ib,ibatch,iprot)
544           enddo
545           call MPI_Allreduce(elowest_aux(1,1),elowest_t(1,1),
546      &      nbeta(ibatch,iprot), 
547      &      MPI_2DOUBLE_PRECISION, MPI_MINLOC, Comm1, IERROR)
548 #ifdef DEBUG
549           do ib=1,nbeta(ibatch,iprot)
550             write (iout,*) "beta=",betaT(ib,ibatch,iprot)
551             write (iout,'(9helowest_t,10f15.3)')
552      &        elowest_t(1,ib),elowest_t(2,ib)
553           enddo
554           write (iout,*) "Processor",me,me1," finished MPI_Reduce: 2"
555 #endif
556           do ib=1,nbeta(ibatch,iprot)
557             elowest(ib,ibatch,iprot)=elowest_t(1,ib)
558             ind_lowest(ib,ibatch,iprot)=elowest_t(2,ib)
559           enddo
560         enddo ! ibatc
561       enddo ! iprot
562       do iprot=1,nprot
563         do ibatch=1,natlike(iprot)+2
564           do ib=1,nbeta(ibatch,iprot)
565             write (iout,*) "iprot",iprot," ibatch",ibatch," ib",ib,
566      &         " elowest",elowest(ib,ibatch,iprot)
567           enddo
568         enddo
569       enddo
570 c
571 c Allgather to provide all energies to all processors
572 c
573       do iprot=1,nprot
574         do i=1,scount(me1,iprot)
575           e_total_(i)=e_total(indstart(me1,iprot)+i-1,iprot)
576         enddo
577         call MPI_Allgatherv(e_total_(1),
578      &     scount(me1,iprot),MPI_DOUBLE_PRECISION,e_total(1,iprot),
579      &     scount(0,iprot),idispl(0,iprot),MPI_DOUBLE_PRECISION,
580      &     Comm1, IERROR)
581 c        call MPI_Allgatherv(e_total(indstart(me1,iprot),iprot),
582 c     &     scount(me1,iprot),MPI_DOUBLE_PRECISION,e_total(1,iprot),
583 c     &     scount(0,iprot),idispl(0,iprot),MPI_DOUBLE_PRECISION,
584 c     &     Comm1, IERROR)
585       enddo
586 #endif
587 c
588 c Now determine which conformations will enter the database.
589 c
590       do iprot=1,nprot
591         call restore_molinfo(iprot)
592 c Clear the list of conformations
593         do i=1,min0(ntot(iprot),maxstr)
594           list_conf(i,iprot)=0
595         enddo
596         do i=1,ntot(iprot)
597           lflag(i)=.false.
598         enddo
599 c Make the list of conformations based on energy cut-off.
600         nn=0
601 #ifdef DEBUG
602         write (iout,*) "iprot",iprot," ibatch",ibatch," betmin",
603      &    betmin(ibatch,iprot)
604 #endif
605 #ifdef DEBUG
606         write (iout,*) "e_lowb",e_lowb(iprot)
607         write (iout,*) "t_lowb",t_lowb(iprot)
608 #endif
609         do ibatch=1,natlike(iprot)+2
610         do i=1,ntot(iprot)
611           jj = i2ii(i,iprot)
612 #ifdef MPI
613 c              write (iout,*) "i",i," ii",ii," indstart",
614 c     &         indstart(me1,iprot)," indend",indend(me1,iprot)
615           if (i.ge.indstart(me1,iprot).and.i.le.indend(me1,iprot)) 
616      &    then
617 #endif
618 c              write (iout,*) "i",i," ii",ii," kbatch",kbatch(i,iprot),
619 c     &          " flag",lflag(i)
620           if (.not.lflag(i)) then
621 c            if (eini(i,iprot).lt.e_lowb(iprot) .or.
622 c     &          entfac(i,iprot).lt.t_lowb(iprot)) then
623 c#ifdef DEBUG
624 c              write (iout,*) "Conformation",i," eini",eini(i,iprot),
625 c     &           " entfac",entfac(i,iprot)," e_lowb",
626 c     &         e_lowb(iprot),
627 c     &           " t_lowb",t_lowb(iprot)
628 c#endif
629 c              lflag(i)=.true.
630 c              goto 122
631 c            endif
632             do ib=1,nbeta(ibatch,iprot)
633               dene=betaT(ib,ibatch,iprot)*(e_total(i,iprot)
634      &          -elowest(ib,ibatch,iprot))+entfac(i,iprot)
635 #ifdef DEBUG
636               write (iout,*) "beta",betaT(ib,ibatch,iprot),
637      &           " i",i," e_total",
638      &        e_total(i,iprot),
639      &         " elowest",elowest(ib,ibatch,iprot)," dene",dene,
640      &         " enecut",enecut(iprot)
641 #endif
642               if (dene.lt.enecut(iprot)) then
643                 nn=nn+1     
644                 list_conf(nn,iprot)=i
645                 lflag(i)=.true.
646                 goto 122
647               endif
648             enddo
649   122       continue
650           endif
651 #ifdef MPI
652           endif
653 #endif
654         enddo ! j
655         enddo ! ibatch
656         ntot_work(iprot)=nn
657         if (nn.gt.maxstr) then
658           write (iout,*) "Error - after applying cutoff the number",
659      &     " of conformations for protein ",i," exceeds MAXSTR:",
660      &     nn,maxstr 
661           write (iout,*) "The calculation is terminating."
662           call flush(iout)
663 #ifdef MPI
664           call MPI_Finalize(ierror)
665 #endif
666           stop
667         endif
668         call imysort(ntot_work(iprot),list_conf(1,iprot),ipermut)
669         write (iout,*) "Protein",iprot,ntot_work(iprot),
670      &   " conformations within scaled energy cut-off=",enecut(iprot),
671      &   " found at processor",me
672 #ifdef MPI
673 c
674 c All workers get the complete list of conformations.
675 c
676         call MPI_Allgather(ntot_work(iprot),1,MPI_INTEGER,
677      &    scount(0,iprot),1,MPI_INTEGER,Comm1,IERROR)
678         idispl(0,iprot)=0
679         do i=1,nprocs1-1
680           idispl(i,iprot)=idispl(i-1,iprot)+scount(i-1,iprot)
681         enddo
682 #ifdef DEBUG
683         write (iout,*) "Protein",iprot," Scount and Idispl"
684         do i=0,nprocs1-1
685           write (iout,*) i,scount(i,iprot),idispl(i,iprot)
686         enddo
687         write (iout,*) "Protein",i,
688      &    " local list of conformations of processor",me
689         do i=1,ntot_work(iprot)
690           write(iout,*) i,list_conf(i,iprot) 
691         enddo
692         write (iout,*) "Before REDUCE: ntot_work",ntot_work(iprot)
693         call flush(iout)
694 #endif
695         call MPI_Allreduce(ntot_work(iprot),nn,1,MPI_INTEGER,MPI_SUM,
696      &    Comm1,IERROR)
697         ntot_work(iprot)=nn
698 #ifdef DEBUG
699         write (iout,*) "After REDUCE: ntot_work",ntot_work(iprot)
700         call flush(iout)
701 #endif
702         call MPI_Allgatherv(list_conf(1,iprot),
703      &    scount(me1,iprot),MPI_INTEGER,list_conf_(1,iprot),
704      &    scount(0,iprot),idispl(0,iprot),MPI_INTEGER,Comm1,IERROR)
705         do i=1,ntot_work(iprot)
706           list_conf(i,iprot)=list_conf_(i,iprot)
707         enddo
708 #ifdef DEBUG
709         write (iout,*) "Protein",i,
710      &    " global list of conformations of processor",me
711         do i=1,ntot_work(iprot)
712           write(iout,'(2i5,e15.5,33f7.3)')i,list_conf(i,iprot),
713      &      e_total(list_conf(i,iprot),iprot)
714 c     &      ,(nu(k,i,iprot),k=1,natlike(iprot))
715         enddo
716         call flush(iout)
717 #endif
718 #endif
719 c
720 c Construct the mapping of the new list to the original numbers of 
721 c conformations.
722 c
723         do i=1,ntot(iprot)
724           tsil(i,iprot)=0
725         enddo
726         do i=1,ntot_work(iprot)
727           tsil(list_conf(i,iprot),iprot)=i
728         enddo
729 #ifdef DEBUG
730         write (iout,*) "Protein",i," List-to-conformation mapping"
731         do i=1,ntot(iprot)
732           write(iout,*) i,tsil(i,iprot)
733         enddo
734 #endif
735       enddo       ! iprot
736 c
737 c Divide the work again based on the current lists of conformations
738 c
739       call work_partition(.true.)
740 c
741 c If the conformations fit into memory, read them off a DA scratchfile.
742 c
743       do iprot=1,nprot
744         call restore_molinfo(iprot)
745 #ifdef MPI
746         nchunk_conf(iprot)=iroof(scount(me1,iprot),maxstr_proc)
747 #else
748         nchunk_conf(iprot)=iroof(ntot_work(iprot),maxstr_proc)
749 #endif
750         if (nchunk_conf(iprot).eq.1) then
751           write (iout,*) "Protein",iprot,
752      &     " in-memory storage of conformations."
753           if (init_ene .or. mod_fourier(nloctyp).gt.0
754      &    .or.mod_fouriertor(nloctyp).gt.0.or.mod_elec.or.mod_scp) then
755             open (icbase,file=bprotfiles(iprot),status="old",
756      &      form="unformatted",access="direct",recl=lenrec(iprot))
757 #ifdef MPI
758             call daread_ccoords(iprot,indstart(me1,iprot),
759      &        indend(me1,iprot))
760 #else
761             call daread_ccoords(iprot,1,ntot_work(iprot))
762 #endif
763             close(icbase)
764           else 
765             open (ientin,file=benefiles(iprot),status="old",
766      &        form="unformatted",access="direct",recl=lenrec_ene(iprot))
767 #ifdef MPI
768             call daread_ene(iprot,indstart(me1,iprot),
769      &        indend(me1,iprot))
770 #else
771             call daread_ene(iprot,1,ntot_work(iprot))
772 #endif
773             close(ientin)
774             open (icbase,file=bprotfiles(iprot),status="old",
775      &      form="unformatted",access="direct",recl=lenrec(iprot))
776 #ifdef MPI
777             call daread_ccoords(iprot,indstart(me1,iprot),
778      &        indend(me1,iprot))
779 #else
780             call daread_ccoords(iprot,1,ntot_work(iprot))
781 #endif
782             close(icbase)
783           endif
784 #ifdef DEBUG
785           write (iout,*) "Protein",i,
786      &    " global list of conformations of processor",me
787           do i=1,ntot_work(iprot)
788             write(iout,'(2i5,e15.5,33f7.3)')i,list_conf(i,iprot),
789      &      e_total(list_conf(i,iprot),iprot)
790           enddo
791 #endif
792         else
793             write (iout,*) "Protein",iprot,
794      &     " off-memory storage of conformations; ",
795      &     "energy will be evaluated in",nchunk_conf(iprot)," passes."
796         endif
797       enddo
798       do i=1,n_ene
799         ww_orig(i)=ww(i)
800       enddo 
801       do iprot=1,nprot
802 #ifdef DEBUG
803         write (iout,*) "E_TOTAL and ETOT_ORIG of protein",iprot
804 #endif
805         do i=1,ntot_work(iprot)
806           etot_orig(i,iprot)=e_total(list_conf(i,iprot),iprot)
807 #ifdef DEBUG
808           write (iout,*) i,list_conf(i,iprot),
809      &     e_total(list_conf(i,iprot),iprot),etot_orig(i,iprot)
810 #endif
811         enddo
812         do ibatch=1,natlike(iprot)+2
813         do ib=1,nbeta(ibatch,iprot)
814           emini = elowest(ib,ibatch,iprot)
815           if (elowest(ib,ibatch,iprot) .lt. emini) 
816      &        emini = elowest(ib,ibatch,iprot)
817           emin_orig(ib,ibatch,iprot)=emini
818         enddo
819         enddo
820       enddo
821 #ifdef MPI
822       if (me1 .eq. Master) call write_conf_count
823 #else
824       call write_conf_count
825 #endif
826 #ifdef DEBUG
827       write (iout,*) "ELOWEST at the end of MAKE_LIST"
828       do iprot=1,nprot
829         do ibatch=1,natlike(iprot)+2
830           do ib=1,nbeta(ibatch,iprot)
831             write (iout,*) "iprot",iprot," ibatch",ibatch," ib",ib,
832      &         " elowest",elowest(ib,ibatch,iprot)
833           enddo
834         enddo
835       enddo
836 #endif
837       tcpu_fin = tcpu() - tcpu_ini
838       write (iout,*) "Time for creating list of conformations",tcpu_fin
839       call flush(iout)
840       t_func = t_func + tcpu_fin
841       return
842       end
843 c----------------------------------------------------------------------------
844       subroutine emin_search(iprot)
845       implicit none
846       include "DIMENSIONS"
847       include "DIMENSIONS.ZSCOPT"
848 #ifdef MPI
849       include "mpif.h"
850       integer IERROR,ErrCode,Status(MPI_STATUS_SIZE,10)
851       integer req(10),msg_in(5),msg_out(5),address,size
852       character*1 buffer(8*(2*maxstr_proc*nntyp+8000))
853       include "COMMON.MPI"
854 #endif
855       include "COMMON.WEIGHTS"
856       include "COMMON.WEIGHTDER"
857       include "COMMON.COMPAR"
858       include "COMMON.ENERGIES"
859       include "COMMON.IOUNITS"
860       include "COMMON.VMCPAR"
861       include "COMMON.NAMES"
862       include "COMMON.INTERACT"
863       include "COMMON.TIME1"
864       include "COMMON.CHAIN"
865       include "COMMON.PROTFILES"
866       include "COMMON.VAR"
867       include "COMMON.GEO"
868       include "COMMON.CLASSES"
869 C Define local variables
870       integer i,ii,iii,kkk,jj,j,k,kk,l,iprot,ib,ibatch,nn
871       integer ipass_conf,istart_conf,iend_conf,Previous,Next
872       double precision energia(0:max_ene)
873       double precision etoti,elowesti,dene
874       double precision tcpu_ini,tcpu_fin,tcpu
875       double precision etot_aux,enepsjk,
876      &  emini,elowest_aux(2,maxT)
877       integer iroof,icant
878       external iroof,icant
879       logical lprn
880
881       do ibatch=1,natlike(iprot)+2
882         do ib=1,nbeta(ibatch,iprot)
883           elowest(ib,ibatch,iprot)=1.0d20
884         enddo
885       enddo
886       do ibatch=1,natlike(iprot)+2
887         do k=1,ntot(iprot)
888         jj=i2ii(k,iprot)
889 #ifdef MPI
890         if (jj.gt.0) then
891 #endif
892           do ib=1,nbeta(ibatch,iprot)
893             etoti=0.0d0
894             do kk=1,n_ene
895               etoti=etoti+ww(kk)*enetb(jj,kk,iprot)
896      &          *escal(kk,ib,ibatch,iprot)
897             enddo
898             if (ib.eq.1 .and. etoti.lt.elowest_ent(1,ibatch,iprot)) 
899      &      then
900               elowest_ent(1,ibatch,iprot)=etoti
901               elowest_ent(2,ibatch,iprot)=entfac(k,iprot)
902             endif
903 c            if (ib.gt.1) 
904 c     &      etoti=etoti+entfac(k,iprot)/betaT(ib,ibatch,iprot)
905             etoti=etoti+entfac(k,iprot)/betaT(ib,ibatch,iprot)
906             if (etoti.lt.elowest(ib,ibatch,iprot)) then
907               elowest(ib,ibatch,iprot)=etoti
908               ind_lowest(ib,ibatch,iprot)=k
909             endif
910 c              write (iout,*) ib,betaT(ib,ibatch,iprot),etoti,
911 c     &          entfac(k,iprot)
912 #ifdef DEBUG
913             write (iout,'(2i5,2e15.5,f8.3)') k,jj,
914      &      etoti,elowest(ib,ibatch,iprot),betaT(ib,ibatch,iprot)
915 #endif
916           enddo ! ib
917 #ifdef DEBUG
918           write (iout,'(2i5,20f8.2)') j,jj,(enetb(jj,kk,iprot),
919      &         kk=1,n_ene)
920 #endif
921 #ifdef MPI
922           endif
923 #endif
924         enddo ! k
925       enddo ! ibatch
926       return
927       end