update new files
[unres.git] / source / unres / src_MD-M-homology / map.f
1       subroutine map
2       implicit real*8 (a-h,o-z)
3       include 'DIMENSIONS'
4       include 'COMMON.MAP'
5       include 'COMMON.VAR'
6       include 'COMMON.GEO'
7       include 'COMMON.DERIV'
8       include 'COMMON.IOUNITS'
9       include 'COMMON.NAMES'
10       include 'COMMON.CONTROL'
11       include 'COMMON.TORCNSTR'
12       include 'COMMON.HEADER'
13       double precision energia(0:n_ene)
14       character*5 angid(4) /'PHI','THETA','ALPHA','OMEGA'/
15       double precision ang_list(10)
16       double precision g(maxvar),x(maxvar)
17       integer nn(10)
18       write (iout,'(a,i3,a)')'Energy map constructed in the following ',
19      &       nmap,' groups of variables:'
20       do i=1,nmap
21         write (iout,'(2a,i3,a,i3)') angid(kang(i)),' of residues ',
22      &   res1(i),' to ',res2(i)
23       enddo
24       write (iout,*) "outpdb",outpdb
25       nmax=nstep(1)
26       do i=2,nmap
27         if (nmax.lt.nstep(i)) nmax=nstep(i)
28       enddo
29       ntot=nmax**nmap
30       iii=0
31       write (istat,'(1h#,a14,29a15)') (" ",k=1,nmap),
32      &    (ename(print_order(k)),k=1,nprint_ene),"ETOT","GNORM"
33       do i=0,ntot-1
34         ii=i
35         do j=1,nmap
36           nn(j)=mod(ii,nmax)+1
37           ii=ii/nmax
38         enddo
39         do j=1,nmap
40           if (nn(j).gt.nstep(j)) goto 10
41         enddo
42         iii=iii+1
43 Cd      write (iout,*) i,iii,(nn(j),j=1,nmap)
44         do j=1,nmap
45           ang_list(j)=ang_from(j)
46      &       +(nn(j)-1)*(ang_to(j)-ang_from(j))/nstep(j)
47           do k=res1(j),res2(j)
48             goto (1,2,3,4), kang(j)
49     1       phi(k)=deg2rad*ang_list(j)
50             if (minim) phi0(k-res1(j)+1)=deg2rad*ang_list(j)
51             goto 5
52     2       theta(k)=deg2rad*ang_list(j)
53             goto 5
54     3       alph(k)=deg2rad*ang_list(j)
55             goto 5
56     4       omeg(k)=deg2rad*ang_list(j)
57     5       continue
58           enddo ! k
59         enddo ! j
60         call chainbuild_extconf
61         if (minim) then 
62          call geom_to_var(nvar,x)
63          call minimize(etot,x,iretcode,nfun)
64          print *,'SUMSL return code is',iretcode,' eval ',nfun
65 c         call intout
66         else
67          call zerograd
68          call geom_to_var(nvar,x)
69         endif
70          call etotal(energia(0))
71          etot = energia(0)
72          nf=1
73          nfl=3
74          call gradient(nvar,x,nf,g,uiparm,urparm,fdum)
75          gnorm=0.0d0
76          do k=1,nvar
77            gnorm=gnorm+g(k)**2
78          enddo
79         etot=energia(0)
80
81         gnorm=dsqrt(gnorm)
82 c        write (iout,'(6(1pe15.5))') (ang_list(k),k=1,nmap),etot,gnorm
83         write (istat,'(30e15.5)') (ang_list(k),k=1,nmap),
84      &   (energia(print_order(ii)),ii=1,nprint_ene),etot,gnorm
85         if (outpdb) call pdbout(etot,titel(:50),ipdb)
86 c        write (iout,*) 'POINT',I,' ANGLES:',(ang_list(k),k=1,nmap)
87 c        call intout
88 c        call enerprint(energia)
89    10   continue
90       enddo ! i
91       return
92       end