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