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