update new files
[unres.git] / source / maxlik / src_FPy.org / zscorez_new.F
1       subroutine maxlik_init(nvarr,xrange,comm,mepython)
2 c Optimize the UNRES energy function by minimization of a quartic target
3 c function or by the VMC method.
4       implicit none
5 #ifndef ISNAN
6       external proc_proc
7 #endif
8 #ifdef WINPGI
9 cMS$ATTRIBUTES C ::  proc_proc
10 #endif
11       include "DIMENSIONS"
12       include "DIMENSIONS.ZSCOPT"
13 #ifdef MPI
14       include "mpif.h"
15       integer IERROR,ERRCODE,kolor,key,comm,mepython
16       include "COMMON.MPI"
17 #endif
18       include "COMMON.IOUNITS"
19       include "COMMON.OPTIM"
20       include "COMMON.XBOUND"
21       integer nvarr,iparm
22       double precision rr,x(max_paropt)
23       integer idumm
24       integer i
25       double precision xrange(maxvar,2)
26       integer number_of_variables
27       common /patch/ number_of_variables
28 #ifdef PYTHON
29 Cf2py intent(in) comm
30 Cf2py intent(out) nvarr
31 Cf2Py intent(out) xrange
32 #endif
33       print *,"Starting..."
34 #ifdef MPI
35 #ifndef PYTHON
36 c      print *,"Initializing MPI..."
37       call MPI_Init( IERROR )
38       ALL_COMM = MPI_COMM_WORLD 
39 #else
40       ALL_COMM = comm
41       print *,"comm",comm
42 #endif
43       print *,"Before MPI_Comm_rank"
44       call MPI_Comm_rank( ALL_COMM, me, IERROR )
45       print *,"After MPI_Comm_rank"
46       call MPI_Comm_size( ALL_COMM, nprocs, IERROR )
47       print *,"Finished initializing MPI... nprocs",nprocs
48       Master = 0
49       Master1 = 1
50 c      print *,"Me",me," Master",master," Ierror",ierror
51 #ifndef PYTHON
52       if (me.eq.Master .and. ierror.gt.0) then
53         write(iout,*) "SEVERE ERROR - Can't initialize MPI."
54         call mpi_finalize(ierror)
55         stop
56       endif
57 #endif
58 #endif
59 #ifndef ISNAN
60 c NaNQ initialization
61       i=-1
62       rr=dacos(100.0d0)
63 #ifdef WINPGI
64       idumm=proc_proc(rr,i)
65 #else
66       call proc_proc(rr,i)
67 #endif
68 #endif
69       call initialize
70 c      print *,"calling openunits"
71       call openunits(mepython)
72 c      print *,"openunits called"
73       call read_general_data(*10)
74       if (me.eq.Master) then
75       write (iout,'(80(1h-)/10x,
76      & "Maximum likelihood optimization of UNRES energy function",
77      & " v. 05/10/16"/80(1h-))')
78       call flush(iout)
79       call cinfo
80       endif
81 c      call promienie(*10)
82       if (me.eq.Master)write (iout,*) "Finished READ_GENERAL_DATA"
83       call flush(iout)
84       do iparm=1,nparmset
85         call parmread(iparm,*10)
86       enddo
87       if (me.eq.Master) write (iout,*) "Finished parmread"
88       call flush(iout)
89       call read_optim_parm(*10)
90       call print_general_data(*10)
91       call read_protein_data(*10)
92       if (me.eq.Master) write (iout,*) "Finished READ_PROTEIN_DATA"
93       call flush(iout)
94       call read_database(*10)
95       if (me.eq.Master) write (iout,*) "Finished READ_DATABASE"
96       call flush(iout)
97 #ifdef MPI 
98 c      write (iout,*) Me,' calling PROC_GROUPS'
99       call proc_groups
100 c      write (iout,*) Me,' calling WORK_PARTITION_MAP'
101 c      call work_partition_map(nvarr)
102 #endif
103       call proc_data(nvarr,x,*10)
104 #ifdef PYTHON
105       number_of_variables=nvarr
106       do i=1,nvarr
107         xrange(i,1)=x_low(i)
108         xrange(i,2)=x_up(i)
109       enddo
110       if (me.eq.Master) write (iout,*) "xrange from MAXLIK_INIT"
111       do i=1,nvarr
112         if (me.eq.Master) write (iout,*) i,xrange(i,1),xrange(i,2)
113       enddo
114       if (me.eq.Master) 
115      & write (iout,*) "================ maxlik intiialization completed"
116       call flush (iout)
117 #endif
118       call read_thermal
119       return
120    10 if (me.eq.Master)write (iout,*) "Error termination of the program"
121       call MPI_Finalize( IERROR )
122       return
123       end
124 c------------------------------------------------------------------------
125       subroutine maxlik_optim(x,xmin,fmin)
126 c Optimize the UNRES energy function by minimization of a quartic target
127 c function or by the VMC method.
128       implicit none
129 #ifndef ISNAN
130       external proc_proc
131 #endif
132 #ifdef WINPGI
133 cMS$ATTRIBUTES C ::  proc_proc
134 #endif
135       include "DIMENSIONS"
136       include "DIMENSIONS.ZSCOPT"
137 #ifdef MPI
138       include "mpif.h"
139       integer IERROR,ERRCODE,kolor,key
140       include "COMMON.MPI"
141 #endif
142       include "COMMON.IOUNITS"
143       include "COMMON.OPTIM"
144       integer nvarr,iparm,i
145       double precision rr,x(max_paropt),xmin(max_paropt),fmin
146       double precision tcpu,t1,t1w,t1_ini,t1w_ini
147       integer number_of_variables
148       common /patch/ number_of_variables
149       nvarr=number_of_variables
150 #ifdef PYTHON
151 Cf2py intent(in) x
152 Cf2py intent(out) xmin
153 Cf2py intent(out) fmin
154       if (me.eq.Master) then
155         write (*,*) "MAXLIK_OPTIM: Variables from PYTHON:",nvarr
156         write (iout,*) "MAXLIK_OPTIM: Variables from PYTHON:",nvarr
157         do i=1,nvarr
158           write (iout,*) i,x(i)
159         enddo 
160       endif
161 #else
162       if (me.eq.Master) then
163         write (*,*) "MAXLIK_OPTIM: Variables from maxlik_init:",nvarr
164         write (iout,*) "MAXLIK_OPTIM: Variables from maxlik_init:",nvarr
165         do i=1,nvarr
166           write (iout,*) i,x(i)
167         enddo 
168       endif
169 #endif
170 #ifdef MPI
171       if (me.eq.Master) then
172         t1w_ini = MPI_WTIME()
173         t1_ini = tcpu()
174 #endif
175         call maxlikopt(nvarr,x,xmin,fmin)
176         write (iout,*) "fmin from MAXLIK_OPTIM:",fmin
177 #ifdef MPI
178         call jebadelko(nvarr)
179       else
180         call jebadelko(nvarr)
181       endif
182       call bilans
183       if (me.eq.Master) then
184         t1w = mpi_wtime() - t1w_ini
185         t1 = tcpu() - t1_ini
186         write (iout,*)
187         write (iout,*) "CPU time",t1," wall clock time",t1w
188       call flush(iout)
189       endif
190 #ifdef PYTHON
191       if (me.eq.Master) 
192      & write (iout,'(30(1h-),"Minimization completed",30(1h-))')
193       call flush(iout)
194 #else
195       call MPI_Finalize( IERROR )
196 #endif
197 #else
198       call bilans
199 #endif
200       return
201       end