76bd2808cdc71755c87ba8246389b1fb335af0c0
[unres.git] / source / unres / src-HCD-5D / unres.F
1 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
2 C                                                                              C
3 C                                U N R E S                                     C
4 C                                                                              C
5 C Program to carry out conformational search of proteins in an united-residue  C
6 C approximation.                                                               C
7 C                                                                              C
8 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
9       implicit none
10       include 'DIMENSIONS'
11
12
13 #ifdef MPI
14       include 'mpif.h'
15       include 'COMMON.SETUP'
16 #endif
17       include 'COMMON.TIME1'
18       include 'COMMON.INTERACT'
19       include 'COMMON.NAMES'
20       include 'COMMON.GEO'
21       include 'COMMON.HEADER'
22       include 'COMMON.CONTROL'
23       include 'COMMON.CONTACTS'
24       include 'COMMON.CHAIN'
25       include 'COMMON.VAR'
26       include 'COMMON.IOUNITS'
27       include 'COMMON.FFIELD'
28       include 'COMMON.REMD'
29       include 'COMMON.MD'
30       include 'COMMON.SBRIDGE'
31       double precision hrtime,mintime,sectime
32       character*64 text_mode_calc(-2:14) /'test',
33      & 'SC rotamer distribution',
34      & 'Energy evaluation or minimization',
35      & 'Regularization of PDB structure',
36      & 'Threading of a sequence on PDB structures',
37      & 'Monte Carlo (with minimization) ',
38      & 'Energy minimization of multiple conformations',
39      & 'Checking energy gradient',
40      & 'Entropic sampling Monte Carlo (with minimization)',
41      & 'Energy map',
42      & 'CSA calculations',
43      & 'Not used 9',
44      & 'Not used 10',
45      & 'Soft regularization of PDB structure',
46      & 'Mesoscopic molecular dynamics (MD) ',
47      & 'Not used 13',
48      & 'Replica exchange molecular dynamics (REMD)'/
49       integer ilen
50       external ilen
51       integer ierr
52
53 c      call memmon_print_usage()
54
55       call init_task
56       if (me.eq.king)
57      & write(iout,*)'### LAST MODIFIED  4/25/08 7:29PM by adam'  
58       if (me.eq.king) call cinfo
59 C Read force field parameters and job setup data
60       call readrtns
61 C
62 c      write (iout,*) "After readrtns"
63       call cartprint
64       call intout
65       if (me.eq.king .or. .not. out1file) then
66        write (iout,'(/2a/)') 
67      & text_mode_calc(modecalc)(:ilen(text_mode_calc(modecalc))),
68      & ' calculation.' 
69        if (minim) write (iout,'(a)') 
70      &  'Conformations will be energy-minimized.'
71        write (iout,'(80(1h*)/)') 
72       endif
73       call flush(iout)
74 C
75       if (modecalc.eq.-2) then
76         call test
77         stop
78       else if (modecalc.eq.-1) then
79         write(iout,*) "call check_sc_map next"
80         call check_bond
81         stop
82       endif
83 #ifdef MPI
84       if (fg_rank.gt.0) then
85 C Fine-grain slaves just do energy and gradient components.
86         call ergastulum ! slave workhouse in Latin
87       else
88 #endif
89       if (modecalc.eq.0) then
90         write (iout,*) "Calling exec_eeval_or_minim"
91         call cartprint
92         call exec_eeval_or_minim
93       else if (modecalc.eq.1) then
94         call exec_regularize
95       else if (modecalc.eq.2) then
96         call exec_thread
97       else if (modecalc.eq.3 .or. modecalc .eq.6) then
98         call exec_MC
99       else if (modecalc.eq.4) then
100         call exec_mult_eeval_or_minim
101       else if (modecalc.eq.5) then
102         call exec_checkgrad
103       else if (ModeCalc.eq.7) then
104         call exec_map
105       else if (ModeCalc.eq.8) then
106         call exec_CSA
107       else if (modecalc.eq.11) then
108         call exec_softreg
109       else if (modecalc.eq.12) then
110         call exec_MD
111       else if (modecalc.eq.14) then
112 #ifdef MPI
113         call exec_MREMD
114 #else
115         write (iout,*) "Need a parallel version to run MREMD."
116         stop
117 #endif
118       else
119         write (iout,'(a)') 'This calculation type is not supported',
120      &   ModeCalc
121       endif
122 #ifdef MPI
123       endif
124 C Finish task.
125       if (fg_rank.eq.0) call finish_task
126 c      call memmon_print_usage()
127 #ifdef TIMING
128        call print_detailed_timing
129 #endif
130       call MPI_Finalize(ierr)
131       stop 'Bye Bye...'
132 #else
133       call dajczas(tcpu(),hrtime,mintime,sectime)
134       stop '********** Program terminated normally.'
135 #endif
136       end
137 c--------------------------------------------------------------------------
138       subroutine exec_MD
139       implicit none
140       include 'DIMENSIONS'
141 #ifdef MPI
142       include "mpif.h"
143 #endif
144       include 'COMMON.SETUP'
145       include 'COMMON.CONTROL'
146       include 'COMMON.IOUNITS'
147 c      if (me.eq.king .or. .not. out1file) then
148 c        write (iout,*) "Calling chainbuild"
149 c        call flush(iout)
150 c      endif
151       call chainbuild
152 c      if (me.eq.king .or. .not. out1file) then
153 c        write (iout,*) "Calling MD"
154 c        call flush(iout)
155 c      endif
156       call MD
157       return
158       end
159 c---------------------------------------------------------------------------
160 #ifdef MPI
161       subroutine exec_MREMD
162       implicit none
163       include 'DIMENSIONS'
164 #ifdef MPI
165       include "mpif.h"
166 #endif
167       include 'COMMON.SETUP'
168       include 'COMMON.CONTROL'
169       include 'COMMON.IOUNITS'
170       include 'COMMON.REMD'
171       integer i
172       if (me.eq.king .or. .not. out1file)
173      &   write (iout,*) "Calling chainbuild"
174       call chainbuild
175       if (me.eq.king .or. .not. out1file)
176      &   write (iout,*) "Calling REMD"
177       if (remd_mlist) then 
178         call MREMD
179       else
180         do i=1,nrep
181           remd_m(i)=1
182         enddo
183         call MREMD
184       endif
185       return
186       end
187 #endif
188 c---------------------------------------------------------------------------
189       subroutine exec_eeval_or_minim
190       implicit none
191       include 'DIMENSIONS'
192 #ifdef MPI
193       include 'mpif.h'
194 #endif
195       include 'COMMON.SETUP'
196       include 'COMMON.TIME1'
197       include 'COMMON.INTERACT'
198       include 'COMMON.NAMES'
199       include 'COMMON.GEO'
200       include 'COMMON.HEADER'
201       include 'COMMON.CONTROL'
202       include 'COMMON.CONTACTS'
203       include 'COMMON.CHAIN'
204       include 'COMMON.VAR'
205       include 'COMMON.IOUNITS'
206       include 'COMMON.FFIELD'
207       include 'COMMON.REMD'
208       include 'COMMON.MD'
209       include 'COMMON.SBRIDGE'
210       integer i,icall,iretcode,nfun
211       common /srutu/ icall
212       integer nharp,iharp(4,maxres/3)
213       integer nft_sc
214       logical fail
215       double precision energy(0:n_ene),etot,etota
216       double precision energy_long(0:n_ene),energy_short(0:n_ene)
217       double precision rms,frac,frac_nn,co
218       double precision varia(maxvar)
219       double precision time00,time1,time_ene,evals
220       integer ilen
221       if (indpdb.eq.0)     call chainbuild
222       if (indpdb.ne.0) then
223       dc(1,0)=c(1,1)
224       dc(2,0)=c(2,1)
225       dc(3,0)=c(3,1)
226       endif
227 #ifdef MPI
228       time00=MPI_Wtime()
229 #else
230       time00=tcpu()
231 #endif
232       write (iout,*) "Energy evaluation/minimization"
233       call chainbuild_cart
234 c      print *,'dc',dc(1,0),dc(2,0),dc(3,0)
235       if (split_ene) then
236        print *,"Processor",myrank," after chainbuild"
237        icall=1
238        call etotal_long(energy_long(0))
239        write (iout,*) "Printing long range energy"
240        call enerprint(energy_long(0))
241        call etotal_short(energy_short(0))
242        write (iout,*) "Printing short range energy"
243        call enerprint(energy_short(0))
244        do i=0,n_ene
245          energy(i)=energy_long(i)+energy_short(i)
246 c         write (iout,*) i,energy_long(i),energy_short(i),energy(i)
247        enddo
248        write (iout,*) "Printing long+short range energy"
249        call enerprint(energy(0))
250       endif
251 c      write(iout,*)"before etotal"
252 c      call flush(iout)
253       call etotal(energy(0))
254 c      write(iout,*)"after etotal"
255 c      call flush(iout)
256 #ifdef MPI
257       time_ene=MPI_Wtime()-time00
258 #else 
259       time_ene=tcpu()-time00
260 #endif
261       write (iout,*) "Time for energy evaluation",time_ene
262 c      print *,"after etotal"
263       etota = energy(0)
264       etot =etota
265       call enerprint(energy(0))
266       call hairpin(.true.,nharp,iharp)
267 c        print *,'after hairpin'
268       call secondary2(.true.)
269 c        print *,'after secondary'
270       if (minim) then
271 crc overlap test
272         if (overlapsc) then 
273           print *, 'Calling OVERLAP_SC'
274           call overlap_sc(fail)
275           print *,"After overlap_sc"
276         endif 
277
278         if (searchsc) then 
279           call sc_move(2,nres-1,10,1d10,nft_sc,etot)
280           print *,'SC_move',nft_sc,etot
281           write(iout,*) 'SC_move',nft_sc,etot
282         endif 
283
284         if (dccart) then
285           print *, 'Calling MINIM_DC'
286 #ifdef MPI
287           time1=MPI_WTIME()
288 #else
289           time1=tcpu()
290 #endif
291           call minim_dc(etot,iretcode,nfun)
292         else
293           if (indpdb.ne.0) then 
294             call bond_regular
295             call chainbuild_extconf
296           endif
297           call geom_to_var(nvar,varia)
298           print *,'Calling MINIMIZE.'
299 #ifdef MPI
300           time1=MPI_WTIME()
301 #else
302           time1=tcpu()
303 #endif
304           call minimize(etot,varia,iretcode,nfun)
305         endif
306         print *,'SUMSL return code is',iretcode,' eval ',nfun
307 #ifdef MPI
308         evals=nfun/(MPI_WTIME()-time1)
309 #else
310         evals=nfun/(tcpu()-time1)
311 #endif
312         print *,'# eval/s',evals
313         print *,'refstr=',refstr
314         call hairpin(.false.,nharp,iharp)
315         print *,'after hairpin'
316         call secondary2(.true.)
317         print *,'after secondary'
318         call etotal(energy(0))
319         etot = energy(0)
320         call enerprint(energy(0))
321
322         call intout
323         if (out_int) call briefout(0,etot)
324         if (out_cart) then
325           cartname=prefix(:ilen(prefix))//'.x'
326           potE=etot
327           call cartoutx(0.0d0)
328         endif
329         if (refstr) call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
330           write (iout,'(a,i3)') 'SUMSL return code:',iretcode
331           write (iout,'(a,i20)') '# of energy evaluations:',nfun+1
332           write (iout,'(a,f16.3)')'# of energy evaluations/sec:',evals
333       else
334         print *,'refstr=',refstr
335         if (refstr) call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
336         call briefout(0,etot)
337       endif
338       if (outpdb) call pdbout(etot,titel(:50),ipdb)
339       if (outmol2) call mol2out(etot,titel(:32))
340       return
341       end
342 c---------------------------------------------------------------------------
343       subroutine exec_regularize
344       implicit none
345       include 'DIMENSIONS'
346 #ifdef MPI
347       include 'mpif.h'
348 #endif
349       include 'COMMON.SETUP'
350       include 'COMMON.TIME1'
351       include 'COMMON.INTERACT'
352       include 'COMMON.NAMES'
353       include 'COMMON.GEO'
354       include 'COMMON.HEADER'
355       include 'COMMON.CONTROL'
356       include 'COMMON.CONTACTS'
357       include 'COMMON.CHAIN'
358       include 'COMMON.VAR'
359       include 'COMMON.IOUNITS'
360       include 'COMMON.FFIELD'
361       include 'COMMON.REMD'
362       include 'COMMON.MD'
363       include 'COMMON.SBRIDGE'
364       double precision energy(0:n_ene)
365       double precision etot,rms,frac,frac_nn,co
366       integer iretcode
367
368       call gen_dist_constr
369       call sc_conf
370       call intout
371       call regularize(nct-nnt+1,etot,rms,cref(1,nnt),iretcode)
372       call etotal(energy(0))
373       energy(0)=energy(0)-energy(14)
374       etot=energy(0)
375       call enerprint(energy(0))
376       call intout
377       call briefout(0,etot)
378       if (outpdb) call pdbout(etot,titel(:50),ipdb)
379       if (outmol2) call mol2out(etot,titel(:32))
380       if (refstr) call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
381       write (iout,'(a,i3)') 'SUMSL return code:',iretcode
382       return
383       end
384 c---------------------------------------------------------------------------
385       subroutine exec_thread
386       implicit none
387       include 'DIMENSIONS'
388 #ifdef MP
389       include "mpif.h"
390 #endif
391       include "COMMON.SETUP"
392       call thread_seq
393       return
394       end
395 c---------------------------------------------------------------------------
396       subroutine exec_MC
397       implicit none
398       include 'DIMENSIONS'
399       character*10 nodeinfo
400       integer ipar
401       double precision varia(maxvar)
402 #ifdef MPI
403       include "mpif.h"
404 #endif
405       include "COMMON.SETUP"
406       include 'COMMON.CONTROL'
407       call mcm_setup
408       if (minim) then
409 #ifdef MPI
410         if (modecalc.eq.3) then
411           call do_mcm(ipar)
412         else
413           call entmcm
414         endif
415 #else
416         if (modecalc.eq.3) then
417           call do_mcm(ipar)
418         else
419           call entmcm
420         endif
421 #endif
422       else
423         call monte_carlo
424       endif
425       return
426       end
427 c---------------------------------------------------------------------------
428       subroutine exec_mult_eeval_or_minim
429       implicit none
430       include 'DIMENSIONS'
431 #ifdef MPI
432       include 'mpif.h'
433       integer muster(mpi_status_size)
434       integer ierr,ierror
435 #endif
436       include 'COMMON.SETUP'
437       include 'COMMON.TIME1'
438       include 'COMMON.INTERACT'
439       include 'COMMON.NAMES'
440       include 'COMMON.GEO'
441       include 'COMMON.HEADER'
442       include 'COMMON.CONTROL'
443       include 'COMMON.CONTACTS'
444       include 'COMMON.CHAIN'
445       include 'COMMON.VAR'
446       include 'COMMON.IOUNITS'
447       include 'COMMON.FFIELD'
448       include 'COMMON.REMD'
449       include 'COMMON.MD'
450       include 'COMMON.SBRIDGE'
451       double precision varia(maxvar)
452       integer i,j,iconf,ind(6)
453       integer n,it,man,nf_mcmf,nmin,imm,mm,nft
454       double precision energy(0:max_ene),ene,etot,ene0
455       double precision rms,frac,frac_nn,co
456       double precision time
457       logical eof
458       eof=.false.
459 #ifdef MPI
460       if(me.ne.king) then
461         call minim_mcmf
462         return
463       endif
464
465       close (intin)
466       open(intin,file=intinname,status='old')
467       write (istat,'(a5,20a12)')"#    ",
468      &  (wname(print_order(i)),i=1,nprint_ene)
469       if (refstr) then
470         write (istat,'(a5,20a12)')"#    ",
471      &   (ename(print_order(i)),i=1,nprint_ene),
472      &   "ETOT total","RMSD","nat.contact","nnt.contact"        
473       else
474         write (istat,'(a5,20a12)')"#    ",
475      &    (ename(print_order(i)),i=1,nprint_ene),"ETOT total"
476       endif
477
478       if (.not.minim) then
479         do while (.not. eof)
480           if (read_cart) then
481             read (intin,'(e15.10,e15.5)',end=1100,err=1100) time,ene
482             call read_x(intin,*11)
483 #ifdef MPI
484 c Broadcast the order to compute internal coordinates to the slaves.
485             if (nfgtasks.gt.1)
486      &        call MPI_Bcast(6,1,MPI_INTEGER,king,FG_COMM,IERROR)
487 #endif
488             call int_from_cart1(.false.)
489           else
490             read (intin,'(i5)',end=1100,err=1100) iconf
491             call read_angles(intin,*11)
492             call geom_to_var(nvar,varia)
493             call chainbuild
494           endif
495           write (iout,'(a,i7)') 'Conformation #',iconf
496           call etotal(energy(0))
497           call briefout(iconf,energy(0))
498           call enerprint(energy(0))
499           etot=energy(0)
500           if (refstr) then 
501             call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
502             write (istat,'(i5,20(f12.3))') iconf,
503      &      (energy(print_order(i)),i=1,nprint_ene),etot,
504      &       rms,frac,frac_nn,co
505 cjlee end
506           else
507             write (istat,'(i5,16(f12.3))') iconf,
508      &     (energy(print_order(i)),i=1,nprint_ene),etot
509           endif
510         enddo
511 1100    continue
512         goto 1101
513       endif
514
515       mm=0
516       imm=0
517       nft=0
518       ene0=0.0d0
519       n=0
520       iconf=0
521 c      do n=1,nzsc
522       do while (.not. eof)
523         mm=mm+1
524         if (mm.lt.nodes) then
525           if (read_cart) then
526             read (intin,'(e15.10,e15.5)',end=11,err=11) time,ene
527             call read_x(intin,*11)
528 #ifdef MPI
529 c Broadcast the order to compute internal coordinates to the slaves.
530             if (nfgtasks.gt.1) 
531      &        call MPI_Bcast(6,1,MPI_INTEGER,king,FG_COMM,IERROR)
532 #endif
533             call int_from_cart1(.false.)
534           else
535             read (intin,'(i5)',end=11,err=11) iconf
536             call read_angles(intin,*11)
537             call geom_to_var(nvar,varia)
538             call chainbuild
539           endif
540           write (iout,'(a,i7)') 'Conformation #',iconf
541           n=n+1
542          imm=imm+1
543          ind(1)=1
544          ind(2)=n
545          ind(3)=0
546          ind(4)=0
547          ind(5)=0
548          ind(6)=0
549          ene0=0.0d0
550          call mpi_send(ind,6,mpi_integer,mm,idint,CG_COMM,
551      *                  ierr)
552          call mpi_send(varia,nvar,mpi_double_precision,mm,
553      *                  idreal,CG_COMM,ierr)
554          call mpi_send(ene0,1,mpi_double_precision,mm,
555      *                  idreal,CG_COMM,ierr)
556 c         print *,'task ',n,' sent to worker ',mm,nvar
557         else
558          call mpi_recv(ind,6,mpi_integer,mpi_any_source,idint,
559      *                 CG_COMM,muster,ierr)
560          man=muster(mpi_source)
561 c         print *,'receiving result from worker ',man,' (',iii1,iii,')'
562          call mpi_recv(varia,nvar,mpi_double_precision, 
563      *               man,idreal,CG_COMM,muster,ierr)
564          call mpi_recv(ene,1,
565      *               mpi_double_precision,man,idreal,
566      *               CG_COMM,muster,ierr)
567          call mpi_recv(ene0,1,
568      *               mpi_double_precision,man,idreal,
569      *               CG_COMM,muster,ierr)
570 c         print *,'result received from worker ',man,' sending now'
571
572           call var_to_geom(nvar,varia)
573           call chainbuild
574           call etotal(energy(0))
575           iconf=ind(2)
576           write (iout,*)
577           write (iout,*)
578           write (iout,'(a,2i7)') 'Conformation #',iconf,ind(5)
579
580           etot=energy(0)
581           call enerprint(energy(0))
582           call briefout(it,etot)
583 c          if (minim) call briefout(it,etot)
584           if (refstr) then 
585             call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
586             write (istat,'(i5,19(f12.3))') iconf,
587      &     (energy(print_order(i)),i=1,nprint_ene),etot,
588      &     rms,frac,frac_nn,co
589           else
590             write (istat,'(i5,15(f12.3))') iconf,
591      &     (energy(print_order(i)),i=1,nprint_ene),etot
592           endif
593
594           imm=imm-1
595           if (read_cart) then
596             read (intin,'(e15.10,e15.5)',end=1101,err=1101) time,ene
597             call read_x(intin,*11)
598 #ifdef MPI
599 c Broadcast the order to compute internal coordinates to the slaves.
600             if (nfgtasks.gt.1)
601      &        call MPI_Bcast(6,1,MPI_INTEGER,king,FG_COMM,IERROR)
602 #endif
603             call int_from_cart1(.false.)
604           else
605             read (intin,'(i5)',end=1101,err=1101) iconf
606             call read_angles(intin,*11)
607             call geom_to_var(nvar,varia)
608             call chainbuild
609           endif
610           n=n+1
611           imm=imm+1
612           ind(1)=1
613           ind(2)=n
614           ind(3)=0
615           ind(4)=0
616           ind(5)=0
617           ind(6)=0
618           call mpi_send(ind,6,mpi_integer,man,idint,CG_COMM,
619      *                  ierr)
620           call mpi_send(varia,nvar,mpi_double_precision,man, 
621      *                  idreal,CG_COMM,ierr)
622           call mpi_send(ene0,1,mpi_double_precision,man,
623      *                  idreal,CG_COMM,ierr)
624           nf_mcmf=nf_mcmf+ind(4)
625           nmin=nmin+1
626         endif
627       enddo
628 11    continue
629       do j=1,imm
630         call mpi_recv(ind,6,mpi_integer,mpi_any_source,idint,
631      *               CG_COMM,muster,ierr)
632         man=muster(mpi_source)
633         call mpi_recv(varia,nvar,mpi_double_precision, 
634      *               man,idreal,CG_COMM,muster,ierr)
635         call mpi_recv(ene,1,
636      *               mpi_double_precision,man,idreal,
637      *               CG_COMM,muster,ierr)
638         call mpi_recv(ene0,1,
639      *               mpi_double_precision,man,idreal,
640      *               CG_COMM,muster,ierr)
641
642         call var_to_geom(nvar,varia)
643         call chainbuild
644         call etotal(energy(0))
645         iconf=ind(2)
646         write (iout,*)
647         write (iout,*)
648         write (iout,'(a,2i7)') 'Conformation #',iconf,ind(5)
649
650         etot=energy(0)
651         call enerprint(energy(0))
652         call briefout(it,etot)
653         if (refstr) then 
654           call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
655           write (istat,'(i5,19(f12.3))') iconf,
656      &   (energy(print_order(i)),i=1,nprint_ene),etot,
657      &   rms,frac,frac_nn,co
658         else
659           write (istat,'(i5,15(f12.3))') iconf,
660      &    (energy(print_order(i)),i=1,nprint_ene),etot
661         endif
662         nmin=nmin+1
663       enddo
664 1101  continue
665       do i=1, nodes-1
666          ind(1)=0
667          ind(2)=0
668          ind(3)=0
669          ind(4)=0
670          ind(5)=0
671          ind(6)=0
672          call mpi_send(ind,6,mpi_integer,i,idint,CG_COMM,
673      *                  ierr)
674       enddo
675 #else
676       close (intin)
677       open(intin,file=intinname,status='old')
678       write (istat,'(a5,20a12)')"#    ",
679      &   (wname(print_order(i)),i=1,nprint_ene)
680       write (istat,'("#    ",20(1pe12.4))')
681      &   (weights(print_order(i)),i=1,nprint_ene)
682       if (refstr) then
683         write (istat,'(a5,20a12)')"#    ",
684      &   (ename(print_order(i)),i=1,nprint_ene),
685      &   "ETOT total","RMSD","nat.contact","nnt.contact"
686       else
687         write (istat,'(a5,14a12)')"#    ",
688      &   (ename(print_order(i)),i=1,nprint_ene),"ETOT total"
689       endif
690       do while (.not. eof)
691           if (read_cart) then
692             read (intin,'(e15.10,e15.5)',end=11,err=11) time,ene
693             call read_x(intin,*11)
694 #ifdef MPI
695 c Broadcast the order to compute internal coordinates to the slaves.
696             if (nfgtasks.gt.1)
697      &        call MPI_Bcast(6,1,MPI_INTEGER,king,FG_COMM,IERROR)
698 #endif
699             call int_from_cart1(.false.)
700           else
701             read (intin,'(i5)',end=11,err=11) iconf
702             call read_angles(intin,*11)
703             call geom_to_var(nvar,varia)
704             call chainbuild
705           endif
706         write (iout,'(a,i7)') 'Conformation #',iconf
707         if (minim) call minimize(etot,varia,iretcode,nfun)
708         call etotal(energy(0))
709
710         etot=energy(0)
711         call enerprint(energy(0))
712         if (minim) call briefout(it,etot) 
713         if (refstr) then 
714           call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
715           write (istat,'(i5,18(f12.3))') iconf,
716      &   (energy(print_order(i)),i=1,nprint_ene),
717      &   etot,rms,frac,frac_nn,co
718 cjlee end
719         else
720           write (istat,'(i5,14(f12.3))') iconf,
721      &   (energy(print_order(i)),i=1,nprint_ene),etot
722         endif
723       enddo
724    11 continue
725 #endif
726       return
727       end
728 c---------------------------------------------------------------------------
729       subroutine exec_checkgrad
730       implicit none
731       include 'DIMENSIONS'
732 #ifdef MPI
733       include 'mpif.h'
734 #endif
735       include 'COMMON.SETUP'
736       include 'COMMON.TIME1'
737       include 'COMMON.INTERACT'
738       include 'COMMON.NAMES'
739       include 'COMMON.GEO'
740       include 'COMMON.HEADER'
741       include 'COMMON.CONTROL'
742       include 'COMMON.CONTACTS'
743       include 'COMMON.CHAIN'
744       include 'COMMON.VAR'
745       include 'COMMON.IOUNITS'
746       include 'COMMON.FFIELD'
747       include 'COMMON.REMD'
748       include 'COMMON.MD'
749       include 'COMMON.QRESTR'
750       include 'COMMON.SBRIDGE'
751       integer icall
752       common /srutu/ icall
753       double precision energy(0:max_ene)
754 c      print *,"A TU?"
755 c      do i=2,nres
756 c        vbld(i)=vbld(i)+ran_number(-0.1d0,0.1d0)
757 c        if (itype(i).ne.10) 
758 c     &      vbld(i+nres)=vbld(i+nres)+ran_number(-0.001d0,0.001d0)
759 c      enddo
760       if (indpdb.eq.0) call chainbuild
761 c      do i=0,nres
762 c        do j=1,3
763 c          dc(j,i)=dc(j,i)+ran_number(-0.2d0,0.2d0)
764 c        enddo
765 c      enddo
766 c      do i=1,nres-1
767 c        if (itype(i).ne.10) then
768 c          do j=1,3
769 c            dc(j,i+nres)=dc(j,i+nres)+ran_number(-0.2d0,0.2d0)
770 c          enddo
771 c        endif
772 c      enddo
773 c      do j=1,3
774 c        dc(j,0)=ran_number(-0.2d0,0.2d0)
775 c      enddo
776 #ifdef UMB
777       usampl=.true.
778       scale_umb=.false.
779 #ifdef PMF
780       adaptive=.true.
781 #endif
782       totT=1.d0
783       eq_time=0.0d0
784       call read_fragments
785       iset=1
786       nperm=1
787       print *, "AFTER read fragments"
788            write (iout,*) "iset",iset
789            if (loc_qlike) then
790            write(iout,*) "fragment, weights, q0:"
791            do i=1,nfrag_back
792             write(iout,'(2i5,3(f8.1,f8.2))') ifrag_back(1,i,iset),
793      &         ifrag_back(2,i,iset),
794      &         wfrag_back(1,i,iset),qin_back(1,i,iset),
795      &         wfrag_back(2,i,iset),qin_back(2,i,iset),
796      &         wfrag_back(3,i,iset),qin_back(3,i,iset)
797            enddo
798            else
799            write(iout,*) "fragment, weights:"
800            do i=1,nfrag_back
801             write(iout,'(2i5,3f8.1)') ifrag_back(1,i,iset),
802      &         ifrag_back(2,i,iset),wfrag_back(1,i,iset),
803      &         wfrag_back(2,i,iset),wfrag_back(3,i,iset)
804            enddo
805            endif
806 #ifdef PMF
807       call read_REMDpar
808       call PMFread
809 #endif
810       call rescale_weights(t_bath)
811       call chainbuild_cart
812       print *,"chainbuild_cart"
813       call cartprint
814       print *,"After cartprint"
815       call intout
816       icall=1
817       print *,"before ETOT"
818       write (iout,*) "usampl",usampl
819       call etotal(energy(0))
820       etot = energy(0)
821       call enerprint(energy(0))
822       write (iout,*) "Uconst",Uconst," Uconst_back",uconst_back
823       print *,'icheckgrad=',icheckgrad
824 #endif
825       goto (10,20,30) icheckgrad
826   10  call check_ecartint
827       return
828   20  call check_cartgrad
829       return
830   30  call check_eint
831       return
832       end
833 c---------------------------------------------------------------------------
834       subroutine exec_map
835 C Energy maps
836       call map_read
837       call map
838       return
839       end
840 c---------------------------------------------------------------------------
841       subroutine exec_CSA
842       implicit none
843 #ifdef MPI
844       include "mpif.h"
845 #endif
846       include 'DIMENSIONS'
847       include 'COMMON.IOUNITS'
848 C Conformational Space Annealling programmed by Jooyoung Lee.
849 C This method works only with parallel machines!
850 #ifdef MPI
851       call together
852 #else
853       write (iout,*) "CSA works on parallel machines only"
854 #endif
855       return
856       end
857 c---------------------------------------------------------------------------
858       subroutine exec_softreg
859       implicit none
860       include 'DIMENSIONS'
861       include 'COMMON.IOUNITS'
862       include 'COMMON.CONTROL'
863       double precision energy(0:max_ene),etot
864       double precision rms,frac,frac_nn,co
865       call chainbuild
866       call etotal(energy(0))
867       call enerprint(energy(0))
868       if (.not.lsecondary) then
869         write(iout,*) 'Calling secondary structure recognition'
870         call secondary2(.true.)
871       else
872         write(iout,*) 'Using secondary structure supplied in pdb'
873       endif
874
875       call softreg
876
877       call etotal(energy(0))
878       etot=energy(0)
879       call enerprint(energy(0))
880       call intout
881       call briefout(0,etot)
882       call secondary2(.true.)
883       if (refstr) call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
884       return
885       end