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