1 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
5 C Program to carry out conformational search of proteins in an united-residue C
8 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
9 implicit real*8 (a-h,o-z)
15 include 'COMMON.SETUP'
17 include 'COMMON.TIME1'
18 include 'COMMON.INTERACT'
19 include 'COMMON.NAMES'
21 include 'COMMON.HEADER'
22 include 'COMMON.CONTROL'
23 include 'COMMON.CONTACTS'
24 include 'COMMON.CHAIN'
26 include 'COMMON.IOUNITS'
27 include 'COMMON.FFIELD'
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)',
45 & 'Soft regularization of PDB structure',
46 & 'Mesoscopic molecular dynamics (MD) ',
48 & 'Replica exchange molecular dynamics (REMD)'/
51 c call memmon_print_usage()
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
61 if (me.eq.king .or. .not. out1file) then
63 & text_mode_calc(modecalc)(:ilen(text_mode_calc(modecalc))),
65 if (minim) write (iout,'(a)')
66 & 'Conformations will be energy-minimized.'
67 write (iout,'(80(1h*)/)')
71 if (modecalc.eq.-2) then
74 else if (modecalc.eq.-1) then
75 write(iout,*) "call check_sc_map next"
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
85 if (modecalc.eq.0) then
86 call exec_eeval_or_minim
87 else if (modecalc.eq.1) then
89 else if (modecalc.eq.2) then
91 else if (modecalc.eq.3 .or. modecalc .eq.6) then
93 else if (modecalc.eq.4) then
94 call exec_mult_eeval_or_minim
95 else if (modecalc.eq.5) then
97 else if (ModeCalc.eq.7) then
99 else if (ModeCalc.eq.8) then
101 else if (modecalc.eq.11) then
103 else if (modecalc.eq.12) then
105 else if (modecalc.eq.14) then
108 write (iout,'(a)') 'This calculation type is not supported',
114 if (fg_rank.eq.0) call finish_task
115 c call memmon_print_usage()
117 call print_detailed_timing
119 call MPI_Finalize(ierr)
122 call dajczas(tcpu(),hrtime,mintime,sectime)
123 stop '********** Program terminated normally.'
126 c--------------------------------------------------------------------------
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"
141 c---------------------------------------------------------------------------
142 subroutine exec_MREMD
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"
153 if (me.eq.king .or. .not. out1file)
154 & write (iout,*) "Calling REMD"
164 write (iout,*) "MREMD works on parallel machines only"
168 c---------------------------------------------------------------------------
169 subroutine exec_eeval_or_minim
170 implicit real*8 (a-h,o-z)
175 include 'COMMON.SETUP'
176 include 'COMMON.TIME1'
177 include 'COMMON.INTERACT'
178 include 'COMMON.NAMES'
180 include 'COMMON.HEADER'
181 include 'COMMON.CONTROL'
182 include 'COMMON.CONTACTS'
183 include 'COMMON.CHAIN'
185 include 'COMMON.IOUNITS'
186 include 'COMMON.FFIELD'
187 include 'COMMON.REMD'
189 include 'COMMON.SBRIDGE'
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
202 print *,"Processor",myrank," after chainbuild"
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))
211 energy(i)=energy_long(i)+energy_short(i)
212 write (iout,*) i,energy_long(i),energy_short(i),energy(i)
214 write (iout,*) "Printing long+short range energy"
215 call enerprint(energy(0))
217 call etotal(energy(0))
219 time_ene=MPI_Wtime()-time00
221 time_ene=tcpu()-time00
223 write (iout,*) "Time for energy evaluation",time_ene
224 print *,"after etotal"
227 call enerprint(energy(0))
228 call hairpin(.true.,nharp,iharp)
229 call secondary2(.true.)
233 print *, 'Calling OVERLAP_SC'
234 call overlap_sc(fail)
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
244 print *, 'Calling MINIM_DC'
250 call minim_dc(etot,iretcode,nfun)
252 if (indpdb.ne.0) then
256 call geom_to_var(nvar,varia)
257 print *,'Calling MINIMIZE.'
263 call minimize(etot,varia,iretcode,nfun)
265 print *,'SUMSL return code is',iretcode,' eval ',nfun
267 evals=nfun/(MPI_WTIME()-time1)
269 evals=nfun/(tcpu()-time1)
271 print *,'# eval/s',evals
272 print *,'refstr=',refstr
273 call hairpin(.true.,nharp,iharp)
274 call secondary2(.true.)
275 call etotal(energy(0))
277 call enerprint(energy(0))
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
286 print *,'refstr=',refstr
287 if (refstr) call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
288 call briefout(0,etot)
290 if (outpdb) call pdbout(etot,titel(:32),ipdb)
291 if (outmol2) call mol2out(etot,titel(:32))
294 c---------------------------------------------------------------------------
295 subroutine exec_regularize
296 implicit real*8 (a-h,o-z)
301 include 'COMMON.SETUP'
302 include 'COMMON.TIME1'
303 include 'COMMON.INTERACT'
304 include 'COMMON.NAMES'
306 include 'COMMON.HEADER'
307 include 'COMMON.CONTROL'
308 include 'COMMON.CONTACTS'
309 include 'COMMON.CHAIN'
311 include 'COMMON.IOUNITS'
312 include 'COMMON.FFIELD'
313 include 'COMMON.REMD'
315 include 'COMMON.SBRIDGE'
316 double precision energy(0:n_ene)
321 call regularize(nct-nnt+1,etot,rms,cref(1,nnt),iretcode)
322 call etotal(energy(0))
323 energy(0)=energy(0)-energy(14)
325 call enerprint(energy(0))
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
334 c---------------------------------------------------------------------------
335 subroutine exec_thread
340 include "COMMON.SETUP"
344 c---------------------------------------------------------------------------
346 implicit real*8 (a-h,o-z)
348 character*10 nodeinfo
349 double precision varia(maxvar)
353 include "COMMON.SETUP"
354 include 'COMMON.CONTROL'
358 if (modecalc.eq.3) then
364 if (modecalc.eq.3) then
375 c---------------------------------------------------------------------------
376 subroutine exec_mult_eeval_or_minim
377 implicit real*8 (a-h,o-z)
381 dimension muster(mpi_status_size)
383 include 'COMMON.SETUP'
384 include 'COMMON.TIME1'
385 include 'COMMON.INTERACT'
386 include 'COMMON.NAMES'
388 include 'COMMON.HEADER'
389 include 'COMMON.CONTROL'
390 include 'COMMON.CONTACTS'
391 include 'COMMON.CHAIN'
393 include 'COMMON.IOUNITS'
394 include 'COMMON.FFIELD'
395 include 'COMMON.REMD'
397 include 'COMMON.SBRIDGE'
398 double precision varia(maxvar)
400 double precision energy(0:n_ene)
410 open(intin,file=intinname,status='old')
411 write (istat,'(a5,30a12)')"# ",
412 & (wname(print_order(i)),i=1,nprint_ene)
414 write (istat,'(a5,30a12)')"# ",
415 & (ename(print_order(i)),i=1,nprint_ene),
416 & "ETOT total","RMSD","nat.contact","nnt.contact","cont.order"
418 write (istat,'(a5,30a12)')"# ",
419 & (ename(print_order(i)),i=1,nprint_ene),"ETOT total"
425 read (intin,'(e15.10,e15.5)',end=1100,err=1100) time,ene
426 call read_x(intin,*11)
428 c Broadcast the order to compute internal coordinates to the slaves.
430 & call MPI_Bcast(6,1,MPI_INTEGER,king,FG_COMM,IERROR)
432 call int_from_cart1(.false.)
434 read (intin,'(i5)',end=1100,err=1100) iconf
435 call read_angles(intin,*11)
436 call geom_to_var(nvar,varia)
439 write (iout,'(a,i7)') 'Conformation #',iconf
440 call etotal(energy(0))
441 call briefout(iconf,energy(0))
442 call enerprint(energy(0))
445 call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
446 write (istat,'(i5,30(f12.3))') iconf,
447 & (energy(print_order(i)),i=1,nprint_ene),etot,
448 & rms,frac,frac_nn,co
451 write (istat,'(i5,30(f12.3))') iconf,
452 & (energy(print_order(i)),i=1,nprint_ene),etot
468 if (mm.lt.nodes) then
470 read (intin,'(e15.10,e15.5)',end=11,err=11) time,ene
471 call read_x(intin,*11)
473 c Broadcast the order to compute internal coordinates to the slaves.
475 & call MPI_Bcast(6,1,MPI_INTEGER,king,FG_COMM,IERROR)
477 call int_from_cart1(.false.)
479 read (intin,'(i5)',end=11,err=11) iconf
480 call read_angles(intin,*11)
481 call geom_to_var(nvar,varia)
486 write (iout,*) 'Conformation #',iconf,' read'
495 call mpi_send(ind,6,mpi_integer,mm,idint,CG_COMM,
497 call mpi_send(varia,nvar,mpi_double_precision,mm,
498 * idreal,CG_COMM,ierr)
499 call mpi_send(ene0,1,mpi_double_precision,mm,
500 * idreal,CG_COMM,ierr)
501 c print *,'task ',n,' sent to worker ',mm,nvar
503 call mpi_recv(ind,6,mpi_integer,mpi_any_source,idint,
504 * CG_COMM,muster,ierr)
505 man=muster(mpi_source)
506 c print *,'receiving result from worker ',man,' (',iii1,iii,')'
507 call mpi_recv(varia,nvar,mpi_double_precision,
508 * man,idreal,CG_COMM,muster,ierr)
510 * mpi_double_precision,man,idreal,
511 * CG_COMM,muster,ierr)
512 call mpi_recv(ene0,1,
513 * mpi_double_precision,man,idreal,
514 * CG_COMM,muster,ierr)
515 c print *,'result received from worker ',man,' sending now'
517 call var_to_geom(nvar,varia)
519 call etotal(energy(0))
523 write (iout,*) 'Conformation #',iconf," sumsl return code ",
527 call enerprint(energy(0))
528 call briefout(it,etot)
529 c if (minim) call briefout(it,etot)
531 call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
532 write (istat,'(i5,30(f12.3))') iconf,
533 & (energy(print_order(i)),i=1,nprint_ene),etot,
534 & rms,frac,frac_nn,co
536 write (istat,'(i5,30(f12.3))') iconf,
537 & (energy(print_order(i)),i=1,nprint_ene),etot
542 read (intin,'(e15.10,e15.5)',end=11,err=11) time,ene
543 call read_x(intin,*11)
545 c Broadcast the order to compute internal coordinates to the slaves.
547 & call MPI_Bcast(6,1,MPI_INTEGER,king,FG_COMM,IERROR)
549 call int_from_cart1(.false.)
551 read (intin,'(i5)',end=11,err=11) iconf
552 call read_angles(intin,*11)
553 call geom_to_var(nvar,varia)
557 write (iout,*) 'Conformation #',iconf,' read'
565 call mpi_send(ind,6,mpi_integer,man,idint,CG_COMM,
567 call mpi_send(varia,nvar,mpi_double_precision,man,
568 * idreal,CG_COMM,ierr)
569 call mpi_send(ene0,1,mpi_double_precision,man,
570 * idreal,CG_COMM,ierr)
571 nf_mcmf=nf_mcmf+ind(4)
577 call mpi_recv(ind,6,mpi_integer,mpi_any_source,idint,
578 * CG_COMM,muster,ierr)
579 man=muster(mpi_source)
580 call mpi_recv(varia,nvar,mpi_double_precision,
581 * man,idreal,CG_COMM,muster,ierr)
583 * mpi_double_precision,man,idreal,
584 * CG_COMM,muster,ierr)
585 call mpi_recv(ene0,1,
586 * mpi_double_precision,man,idreal,
587 * CG_COMM,muster,ierr)
589 call var_to_geom(nvar,varia)
591 call etotal(energy(0))
595 write (iout,*) 'Conformation #',iconf," sumsl return code ",
599 call enerprint(energy(0))
600 call briefout(it,etot)
602 call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
603 write (istat,'(i5,30(f12.3))') iconf,
604 & (energy(print_order(i)),i=1,nprint_ene),etot,
605 & rms,frac,frac_nn,co
607 write (istat,'(i5,30(f12.3))') iconf,
608 & (energy(print_order(i)),i=1,nprint_ene),etot
620 call mpi_send(ind,6,mpi_integer,i,idint,CG_COMM,
625 open(intin,file=intinname,status='old')
626 write (istat,'(a5,20a12)')"# ",
627 & (wname(print_order(i)),i=1,nprint_ene)
628 write (istat,'("# ",20(1pe12.4))')
629 & (weights(print_order(i)),i=1,nprint_ene)
631 write (istat,'(a5,20a12)')"# ",
632 & (ename(print_order(i)),i=1,nprint_ene),
633 & "ETOT total","RMSD","nat.contact","nnt.contact"
635 write (istat,'(a5,14a12)')"# ",
636 & (ename(print_order(i)),i=1,nprint_ene),"ETOT total"
640 read (intin,'(e15.10,e15.5)',end=1100,err=1100) time,ene
641 call read_x(intin,*11)
643 c Broadcast the order to compute internal coordinates to the slaves.
645 & call MPI_Bcast(6,1,MPI_INTEGER,king,FG_COMM,IERROR)
647 call int_from_cart1(.false.)
649 read (intin,'(i5)',end=1100,err=1100) iconf
650 call read_angles(intin,*11)
651 call geom_to_var(nvar,varia)
654 write (iout,'(a,i7)') 'Conformation #',iconf
655 if (minim) call minimize(etot,varia,iretcode,nfun)
656 call etotal(energy(0))
659 call enerprint(energy(0))
660 if (minim) call briefout(it,etot)
662 call rms_nac_nnc(rms,frac,frac_nn,co,.true.)
663 write (istat,'(i5,18(f12.3))') iconf,
664 & (energy(print_order(i)),i=1,nprint_ene),
665 & etot,rms,frac,frac_nn,co
668 write (istat,'(i5,14(f12.3))') iconf,
669 & (energy(print_order(i)),i=1,nprint_ene),etot
677 c---------------------------------------------------------------------------
678 subroutine exec_checkgrad
679 implicit real*8 (a-h,o-z)
684 include 'COMMON.SETUP'
685 include 'COMMON.TIME1'
686 include 'COMMON.INTERACT'
687 include 'COMMON.NAMES'
689 include 'COMMON.HEADER'
690 include 'COMMON.CONTROL'
691 include 'COMMON.CONTACTS'
692 include 'COMMON.CHAIN'
694 include 'COMMON.IOUNITS'
695 include 'COMMON.FFIELD'
696 include 'COMMON.REMD'
698 include 'COMMON.SBRIDGE'
700 double precision energy(0:max_ene)
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)
706 if (indpdb.eq.0) call chainbuild
709 c dc(j,i)=dc(j,i)+ran_number(-0.2d0,0.2d0)
713 c if (itype(i).ne.10) then
715 c dc(j,i+nres)=dc(j,i+nres)+ran_number(-0.2d0,0.2d0)
720 c dc(j,0)=ran_number(-0.2d0,0.2d0)
727 call rescale_weights(t_bath)
732 call etotal(energy(0))
734 call enerprint(energy(0))
735 write (iout,*) "Uconst",Uconst," Uconst_back",uconst_back
736 print *,'icheckgrad=',icheckgrad
737 goto (10,20,30) icheckgrad
738 10 call check_ecartint
741 20 call check_cartgrad
746 c---------------------------------------------------------------------------
753 c---------------------------------------------------------------------------
759 include 'COMMON.IOUNITS'
760 C Conformational Space Annealling programmed by Jooyoung Lee.
761 C This method works only with parallel machines!
764 write (iout,*) "CSA is not supported in this version"
766 csa write (iout,*) "CSA works on parallel machines only"
767 write (iout,*) "CSA is not supported in this version"
771 c---------------------------------------------------------------------------
772 subroutine exec_softreg
773 implicit real*8 (a-h,o-z)
775 include 'COMMON.IOUNITS'
776 include 'COMMON.CONTROL'
777 double precision energy(0:max_ene)
778 logical debug /.false./
780 call etotal(energy(0))
781 call enerprint(energy(0))
782 if (.not.lsecondary) then
783 write(iout,*) 'Calling secondary structure recognition'
784 call secondary2(debug)
786 write(iout,*) 'Using secondary structure supplied in pdb'
791 call etotal(energy(0))
793 call enerprint(energy(0))
795 call briefout(0,etot)
796 call secondary2(.true.)
797 if (refstr) call rms_nac_nnc(rms,frac,frac_nn,co,.true.)