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