ctest the first static disulfide test ene/min cart
[unres.git] / source / unres / src_MD-M / timing.F
1 C $Date: 1994/10/05 16:41:52 $
2 C $Revision: 2.2 $
3 C
4 C
5 C
6       subroutine set_timers
7 c
8       implicit none
9       double precision tcpu
10       include 'COMMON.TIME1'
11 #ifdef MP
12       include 'mpif.h'
13 #endif
14 C Diminish the assigned time limit a little so that there is some time to
15 C end a batch job
16 c     timlim=batime-150.0
17 C Calculate the initial time, if it is not zero (e.g. for the SUN).
18       stime=tcpu()
19 #ifdef MPI
20       walltime=MPI_WTIME()
21       time_reduce=0.0d0
22       time_allreduce=0.0d0
23       time_bcast=0.0d0
24       time_gather=0.0d0
25       time_sendrecv=0.0d0
26       time_scatter=0.0d0
27       time_scatter_fmat=0.0d0
28       time_scatter_ginv=0.0d0
29       time_scatter_fmatmult=0.0d0
30       time_scatter_ginvmult=0.0d0
31       time_barrier_e=0.0d0
32       time_barrier_g=0.0d0
33       time_enecalc=0.0d0
34       time_sumene=0.0d0
35       time_lagrangian=0.0d0
36       time_sumgradient=0.0d0
37       time_intcartderiv=0.0d0
38       time_inttocart=0.0d0
39       time_ginvmult=0.0d0
40       time_fricmatmult=0.0d0
41       time_cartgrad=0.0d0
42       time_bcastc=0.0d0
43       time_bcast7=0.0d0
44       time_bcastw=0.0d0
45       time_intfcart=0.0d0
46       time_vec=0.0d0
47       time_mat=0.0d0
48       time_fric=0.0d0
49       time_stoch=0.0d0
50       time_fricmatmult=0.0d0
51       time_fsample=0.0d0
52 #endif
53 cd    print *,' in SET_TIMERS stime=',stime
54       return 
55       end
56 C------------------------------------------------------------------------------
57       logical function stopx(nf)
58 C This function returns .true. if one of the following reasons to exit SUMSL
59 C occurs. The "reason" code is stored in WHATSUP passed thru a COMMON block:
60 C
61 C... WHATSUP = 0 - go on, no reason to stop. Stopx will return .false.
62 C...           1 - Time up in current node;
63 C...           2 - STOP signal was received from another node because the
64 C...               node's task was accomplished (parallel only);
65 C...          -1 - STOP signal was received from another node because of error;
66 C...          -2 - STOP signal was received from another node, because 
67 C...               the node's time was up.
68       implicit real*8 (a-h,o-z)
69       include 'DIMENSIONS'
70       integer nf
71       logical ovrtim
72 #ifdef MP
73       include 'mpif.h'
74       include 'COMMON.INFO'
75 #endif
76       include 'COMMON.IOUNITS'
77       include 'COMMON.TIME1'
78       integer Kwita
79
80 cd    print *,'Processor',MyID,' NF=',nf
81 #ifndef MPI
82       if (ovrtim()) then
83 C Finish if time is up.
84          stopx = .true.
85          WhatsUp=1
86 #ifdef MPL
87       else if (mod(nf,100).eq.0) then
88 C Other processors might have finished. Check this every 100th function 
89 C evaluation.
90 C Master checks if any other processor has sent accepted conformation(s) to it. 
91          if (MyID.ne.MasterID) call receive_mcm_info
92          if (MyID.eq.MasterID) call receive_conf
93 cd       print *,'Processor ',MyID,' is checking STOP: nf=',nf
94          call recv_stop_sig(Kwita)
95          if (Kwita.eq.-1) then
96            write (iout,'(a,i4,a,i5)') 'Processor',
97      &     MyID,' has received STOP signal in STOPX; NF=',nf
98            write (*,'(a,i4,a,i5)') 'Processor',
99      &     MyID,' has received STOP signal in STOPX; NF=',nf
100            stopx=.true.
101            WhatsUp=2
102          elseif (Kwita.eq.-2) then
103            write (iout,*)
104      &    'Processor',MyID,' received TIMEUP-STOP signal in SUMSL.'
105            write (*,*)
106      &    'Processor',MyID,' received TIMEUP-STOP signal in SUMSL.'
107            WhatsUp=-2
108            stopx=.true.  
109          else if (Kwita.eq.-3) then
110            write (iout,*)
111      &    'Processor',MyID,' received ERROR-STOP signal in SUMSL.'
112            write (*,*)
113      &    'Processor',MyID,' received ERROR-STOP signal in SUMSL.'
114            WhatsUp=-1
115            stopx=.true.
116          else
117            stopx=.false.
118            WhatsUp=0
119          endif
120 #endif
121       else
122          stopx = .false.
123          WhatsUp=0
124       endif
125 #else
126       stopx=.false.
127 #endif
128
129 #ifdef OSF
130 c Check for FOUND_NAN flag
131       if (FOUND_NAN) then
132         write(iout,*)"   ***   stopx : Found a NaN"
133         stopx=.true.
134       endif
135 #endif
136
137       return
138       end
139 C--------------------------------------------------------------------------
140       logical function ovrtim() 
141       include 'DIMENSIONS'
142       include 'COMMON.IOUNITS'
143       include 'COMMON.TIME1'
144       real*8 tcpu
145 #ifdef MPI
146       include "mpif.h"
147       curtim = MPI_Wtime()-walltime
148 #else
149       curtim= tcpu()
150 #endif
151 C  curtim is the current time in seconds.
152 c      write (iout,*) "curtim",curtim," timlim",timlim," safety",safety
153       if (curtim .ge. timlim - safety) then
154         write (iout,'(a,f10.2,a,f10.2,a,f10.2,a)') 
155      &  "***************** Elapsed time (",curtim,
156      &  " s) is within the safety limit (",safety,
157      &  " s) of the allocated time (",timlim," s). Terminating."
158         ovrtim=.true.
159       else
160         ovrtim=.false.
161       endif
162       return                                               
163       end
164 **************************************************************************      
165       double precision function tcpu()
166       include 'COMMON.TIME1'
167 #ifdef ES9000 
168 ****************************
169 C Next definition for EAGLE (ibm-es9000)
170       real*8 micseconds
171       integer rcode
172       tcpu=cputime(micseconds,rcode)
173       tcpu=(micseconds/1.0E6) - stime
174 ****************************
175 #endif
176 #ifdef SUN
177 ****************************
178 C Next definitions for sun
179       REAL*8  ECPU,ETIME,ETCPU
180       dimension tarray(2)
181       tcpu=etime(tarray)
182       tcpu=tarray(1)
183 ****************************
184 #endif
185 #ifdef KSR
186 ****************************
187 C Next definitions for ksr
188 C this function uses the ksr timer ALL_SECONDS from the PMON library to
189 C return the elapsed time in seconds
190       tcpu= all_seconds() - stime
191 ****************************
192 #endif
193 #ifdef SGI
194 ****************************
195 C Next definitions for sgi
196       real timar(2), etime
197       seconds = etime(timar)
198 Cd    print *,'seconds=',seconds,' stime=',stime
199 C      usrsec = timar(1)
200 C      syssec = timar(2)
201       tcpu=seconds - stime
202 ****************************
203 #endif
204
205 #ifdef LINUX
206 ****************************
207 C Next definitions for sgi
208       real timar(2), etime
209       seconds = etime(timar)
210 Cd    print *,'seconds=',seconds,' stime=',stime
211 C      usrsec = timar(1)
212 C      syssec = timar(2)
213       tcpu=seconds - stime
214 ****************************
215 #endif
216
217
218 #ifdef CRAY
219 ****************************
220 C Next definitions for Cray
221 C     call date(curdat)
222 C     curdat=curdat(1:9)
223 C     call clock(curtim)
224 C     curtim=curtim(1:8)
225       cpusec = second()
226       tcpu=cpusec - stime
227 ****************************
228 #endif
229 #ifdef AIX
230 ****************************
231 C Next definitions for RS6000
232        integer*4 i1,mclock
233        i1 = mclock()
234        tcpu = (i1+0.0D0)/100.0D0
235 #endif
236 #ifdef WINPGI
237 ****************************
238 c next definitions for windows NT Digital fortran
239        real time_real
240        call cpu_time(time_real)
241        tcpu = time_real
242 #endif
243 #ifdef WINIFL
244 ****************************
245 c next definitions for windows NT Digital fortran
246        real time_real
247        call cpu_time(time_real)
248        tcpu = time_real
249 #endif
250
251       return     
252       end  
253 C---------------------------------------------------------------------------
254       subroutine dajczas(rntime,hrtime,mintime,sectime)
255       include 'COMMON.IOUNITS'
256       real*8 rntime,hrtime,mintime,sectime 
257       hrtime=rntime/3600.0D0 
258       hrtime=aint(hrtime)
259       mintime=aint((rntime-3600.0D0*hrtime)/60.0D0)
260       sectime=aint((rntime-3600.0D0*hrtime-60.0D0*mintime)+0.5D0)
261       if (sectime.eq.60.0D0) then
262         sectime=0.0D0
263         mintime=mintime+1.0D0
264       endif
265       ihr=hrtime
266       imn=mintime
267       isc=sectime
268       write (iout,328) ihr,imn,isc
269   328 FORMAT(//'***** Computation time: ',I4  ,' hours ',I2  ,
270      1         ' minutes ', I2  ,' seconds *****')       
271       return
272       end
273 C---------------------------------------------------------------------------
274       subroutine print_detailed_timing
275       implicit real*8 (a-h,o-z)
276       include 'DIMENSIONS'
277 #ifdef MPI
278       include 'mpif.h'
279 #endif
280       include 'COMMON.IOUNITS'
281       include 'COMMON.TIME1'
282       include 'COMMON.SETUP'
283 #ifdef MPI
284       time1=MPI_WTIME()
285          write (iout,'(80(1h=)/a/(80(1h=)))') 
286      &    "Details of FG communication time"
287          write (*,'(7(a40,1pe15.5/),40(1h-)/a40,1pe15.5/80(1h=))') 
288      &    "BROADCAST:",time_bcast,"REDUCE:",time_reduce,
289      &    "GATHER:",time_gather,
290      &    "SCATTER:",time_scatter,"SENDRECV:",time_sendrecv,
291      &    "BARRIER ene",time_barrier_e,
292      &    "BARRIER grad",time_barrier_g,
293      &    "TOTAL:",
294      &    time_bcast+time_reduce+time_gather+time_scatter+time_sendrecv
295          write (*,*) fg_rank,myrank,
296      &     ': Total wall clock time',time1-walltime,' sec'
297          write (*,*) "Processor",fg_rank,myrank,
298      &     ": BROADCAST time",time_bcast," REDUCE time",
299      &      time_reduce," GATHER time",time_gather," SCATTER time",
300      &      time_scatter,
301      &     " SCATTER fmatmult",time_scatter_fmatmult,
302      &     " SCATTER ginvmult",time_scatter_ginvmult,
303      &     " SCATTER fmat",time_scatter_fmat,
304      &     " SCATTER ginv",time_scatter_ginv,
305      &      " SENDRECV",time_sendrecv,
306      &      " BARRIER ene",time_barrier_e,
307      &      " BARRIER GRAD",time_barrier_g,
308      &      " BCAST7",time_bcast7," BCASTC",time_bcastc,
309      &      " BCASTW",time_bcastw," ALLREDUCE",time_allreduce,
310      &      " TOTAL",
311      &      time_bcast+time_reduce+time_gather+time_scatter+
312      &      time_sendrecv+time_barrier+time_bcastc
313          write (*,*) "Processor",fg_rank,myrank," enecalc",time_enecalc
314          write (*,*) "Processor",fg_rank,myrank," sumene",time_sumene
315          write (*,*) "Processor",fg_rank,myrank," intfromcart",
316      &     time_intfcart
317          write (*,*) "Processor",fg_rank,myrank," vecandderiv",
318      &     time_vec
319          write (*,*) "Processor",fg_rank,myrank," setmatrices",
320      &     time_mat
321          write (*,*) "Processor",fg_rank,myrank," ginvmult",
322      &     time_ginvmult
323          write (*,*) "Processor",fg_rank,myrank," fricmatmult",
324      &     time_fricmatmult
325          write (*,*) "Processor",fg_rank,myrank," inttocart",
326      &     time_inttocart
327          write (*,*) "Processor",fg_rank,myrank," sumgradient",
328      &     time_sumgradient
329          write (*,*) "Processor",fg_rank,myrank," intcartderiv",
330      &     time_intcartderiv
331          if (fg_rank.eq.0) then
332            write (*,*) "Processor",fg_rank,myrank," lagrangian",
333      &       time_lagrangian
334            write (*,*) "Processor",fg_rank,myrank," cartgrad",
335      &       time_cartgrad
336          endif
337 #else
338          write (*,*) "enecalc",time_enecalc
339          write (*,*) "sumene",time_sumene
340          write (*,*) "intfromcart",time_intfcart
341          write (*,*) "vecandderiv",time_vec
342          write (*,*) "setmatrices",time_mat
343          write (*,*) "ginvmult",time_ginvmult
344          write (*,*) "fricmatmult",time_fricmatmult
345          write (*,*) "inttocart",time_inttocart
346          write (*,*) "sumgradient",time_sumgradient
347          write (*,*) "intcartderiv",time_intcartderiv
348          write (*,*) "lagrangian",time_lagrangian
349          write (*,*) "cartgrad",time_cartgrad
350 #endif
351       return
352       end