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