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