86f23f9ecf696bf8f1c78d2664cce17c23406f06
[unres4.git] / control.F90
1       module control
2 !-----------------------------------------------------------------------------
3       use io_units
4       use names
5       use MPI_data
6       use geometry_data
7       use energy_data
8       use control_data
9       use minim_data
10       use geometry, only:int_bounds
11 #ifndef CLUSTER
12       use csa_data
13 #ifdef WHAM_RUN
14       use wham_data
15 #endif
16 #endif
17       implicit none
18 !-----------------------------------------------------------------------------
19 ! commom.control
20 !      common /cntrl/
21 !      integer :: modecalc,iscode,indpdb,indback,indphi,iranconf,&
22 !       icheckgrad,iprint,i2ndstr,mucadyn,constr_dist,symetr
23 !      logical :: minim,refstr,pdbref,outpdb,outmol2,overlapsc,&
24 !       energy_dec,sideadd,lsecondary,read_cart,unres_pdb,&
25 !       vdisulf,searchsc,lmuca,dccart,extconf,out1file,&
26 !       gnorm_check,gradout,split_ene
27 !... minim = .true. means DO minimization.
28 !... energy_dec = .true. means print energy decomposition matrix
29 !-----------------------------------------------------------------------------
30 ! common.time1
31 !     FOUND_NAN - set by calcf to stop sumsl via stopx
32 !      COMMON/TIME1/
33       real(kind=8) :: STIME,BATIME,PREVTIM,RSTIME
34 !el      real(kind=8) :: TIMLIM,SAFETY
35 !el      real(kind=8) :: WALLTIME
36 !      COMMON/STOPTIM/
37       integer :: ISTOP
38 !      common /sumsl_flag/
39       logical :: FOUND_NAN
40 !      common /timing/
41       real(kind=8) :: t_init
42 !       time_bcast,time_reduce,time_gather,&
43 !       time_sendrecv,time_barrier_e,time_barrier_g,time_scatter,&
44        !t_eelecij,
45 !       time_allreduce,&
46 !       time_lagrangian,time_cartgrad,&
47 !       time_sumgradient,time_intcartderiv,time_inttocart,time_intfcart,&
48 !       time_mat,time_fricmatmult,&
49 !       time_scatter_fmat,time_scatter_ginv,&
50 !       time_scatter_fmatmult,time_scatter_ginvmult,&
51 !       t_eshort,t_elong,t_etotal
52 !-----------------------------------------------------------------------------
53 ! initialize_p.F
54 !-----------------------------------------------------------------------------
55 !      block data
56 !      integer,parameter :: MaxMoveType = 4
57 !      character(len=14),dimension(-1:MaxMoveType+1) :: MovTypID=(/'pool','chain regrow',&
58 !      character :: MovTypID(-1:MaxMoveType+1)=(/'pool','chain regrow',&
59 !       'multi-bond','phi','theta','side chain','total'/)
60 ! Conversion from poises to molecular unit and the gas constant
61 !el      real(kind=8) :: cPoise=2.9361d0, Rb=0.001986d0
62 !-----------------------------------------------------------------------------
63 !      common /przechowalnia/ subroutines: init_int_table,add_int,add_int_from
64       integer,dimension(:),allocatable :: iturn3_start_all,&
65         iturn3_end_all,iturn4_start_all,iturn4_end_all,iatel_s_all,&
66         iatel_e_all !(0:max_fg_procs)
67       integer,dimension(:,:),allocatable :: ielstart_all,&
68         ielend_all !(maxres,0:max_fg_procs-1)
69
70 !      common /przechowalnia/ subroutine: init_int_table
71       integer,dimension(:),allocatable :: ntask_cont_from_all,&
72         ntask_cont_to_all !(0:max_fg_procs-1)
73       integer,dimension(:,:),allocatable :: itask_cont_from_all,&
74         itask_cont_to_all !(0:max_fg_procs-1,0:max_fg_procs-1)
75 !-----------------------------------------------------------------------------
76 !
77 !
78 !-----------------------------------------------------------------------------
79       contains
80 !-----------------------------------------------------------------------------
81 ! initialize_p.F
82 !-----------------------------------------------------------------------------
83       subroutine initialize
84 !
85 ! Define constants and zero out tables.
86 !
87       use comm_iofile
88       use comm_machsw
89       use MCM_data, only: MovTypID
90 !      implicit real*8 (a-h,o-z)
91 !      include 'DIMENSIONS'
92 #ifdef MPI
93       include 'mpif.h'
94 #endif
95 #ifndef ISNAN
96       external proc_proc
97 #ifdef WINPGI
98 !MS$ATTRIBUTES C ::  proc_proc
99 #endif
100 #endif
101 !      include 'COMMON.IOUNITS'
102 !      include 'COMMON.CHAIN'
103 !      include 'COMMON.INTERACT'
104 !      include 'COMMON.GEO'
105 !      include 'COMMON.LOCAL'
106 !      include 'COMMON.TORSION'
107 !      include 'COMMON.FFIELD'
108 !      include 'COMMON.SBRIDGE'
109 !      include 'COMMON.MCM'
110 !      include 'COMMON.MINIM' 
111 !      include 'COMMON.DERIV'
112 !      include 'COMMON.SPLITELE'
113 !      implicit none
114 ! Common blocks from the diagonalization routines
115 !el      integer :: IR,IW,IP,IJK,IPK,IDAF,NAV,IODA(400)
116 !el      integer :: KDIAG,ICORFL,IXDR
117 !el      COMMON /IOFILE/ IR,IW,IP,IJK,IPK,IDAF,NAV,IODA
118 !el      COMMON /MACHSW/ KDIAG,ICORFL,IXDR
119       logical :: mask_r
120 !      real*8 text1 /'initial_i'/
121       real(kind=4) :: rr
122
123 !local variables el
124       integer :: i,j,k,l,ichir1,ichir2,iblock,m,maxit
125
126 #if .not. defined(WHAM_RUN) && .not. defined(CLUSTER)
127       mask_r=.false.
128 #ifndef ISNAN
129 ! NaNQ initialization
130       i=-1
131       rr=dacos(100.0d0)
132 #ifdef WINPGI
133       idumm=proc_proc(rr,i)
134 #elif defined(WHAM_RUN)
135       call proc_proc(rr,i)
136 #endif
137 #endif
138
139       kdiag=0
140       icorfl=0
141       iw=2
142       
143       allocate(MovTypID(-1:MaxMoveType+1))
144       MovTypID=(/'pool          ','chain regrow  ',&
145        'multi-bond    ','phi           ','theta         ',&
146        'side chain    ','total         '/)
147 #endif
148 !
149 ! The following is just to define auxiliary variables used in angle conversion
150 !
151       pi=4.0D0*datan(1.0D0)
152       dwapi=2.0D0*pi
153       dwapi3=dwapi/3.0D0
154       pipol=0.5D0*pi
155       deg2rad=pi/180.0D0
156       rad2deg=1.0D0/deg2rad
157       angmin=10.0D0*deg2rad
158 !el#ifdef CLUSTER
159 !el      Rgas = 1.987D-3
160 !el#endif
161 !
162 ! Define I/O units.
163 !
164       inp=    1
165       iout=   2
166       ipdbin= 3
167       ipdb=   7
168 #ifdef CLUSTER
169       imol2= 18
170       jplot= 19
171 !el      jstatin=10
172       imol2=  4
173       jrms=30
174 #else
175       icart = 30
176       imol2=  4
177       ithep_pdb=51
178       irotam_pdb=52
179       irest1=55
180       irest2=56
181       iifrag=57
182       ientin=18
183       ientout=19
184 !rc for write_rmsbank1  
185       izs1=21
186 !dr  include secondary structure prediction bias
187       isecpred=27
188 #endif
189       igeom=  8
190       intin=  9
191       ithep= 11
192       irotam=12
193       itorp= 13
194       itordp= 23
195       ielep= 14
196       isidep=15
197 #if defined(WHAM_RUN) || defined(CLUSTER)
198       isidep1=22 !wham
199 #else
200 !
201 ! CSA I/O units (separated from others especially for Jooyoung)
202 !
203       icsa_rbank=30
204       icsa_seed=31
205       icsa_history=32
206       icsa_bank=33
207       icsa_bank1=34
208       icsa_alpha=35
209       icsa_alpha1=36
210       icsa_bankt=37
211       icsa_int=39
212       icsa_bank_reminimized=38
213       icsa_native_int=41
214       icsa_in=40
215 !rc for ifc error 118
216       icsa_pdb=42
217 #endif
218       iscpp=25
219       icbase=16
220       ifourier=20
221       istat= 17
222       ibond = 28
223       isccor = 29
224 #ifdef WHAM_RUN
225 !
226 ! WHAM files
227 !
228       ihist=30
229       iweight=31
230       izsc=32
231 #endif
232       ibond_nucl=126
233       ithep_nucl=127
234       irotam_nucl=128
235       itorp_nucl= 129
236       itordp_nucl= 130
237 !      ielep_nucl= 131
238       isidep_nucl=132
239
240
241
242       iliptranpar=60
243       itube=61
244 #if defined(WHAM_RUN) || defined(CLUSTER)
245 !
246 ! setting the mpi variables for WHAM
247 !
248       fgprocs=1
249       nfgtasks=1
250       nfgtasks1=1
251 #endif
252 !
253 ! Set default weights of the energy terms.
254 !
255       wsc=1.0D0 ! in wham:  wlong=1.0D0
256       welec=1.0D0
257       wtor =1.0D0
258       wang =1.0D0
259       wscloc=1.0D0
260       wstrain=1.0D0
261 !
262 ! Zero out tables.
263 !
264 !      print '(a,$)','Inside initialize'
265 !      call memmon_print_usage()
266       
267 !      do i=1,maxres2
268 !       do j=1,3
269 !         c(j,i)=0.0D0
270 !         dc(j,i)=0.0D0
271 !       enddo
272 !      enddo
273 !      do i=1,maxres
274 !       do j=1,3
275 !         xloc(j,i)=0.0D0
276 !        enddo
277 !      enddo
278 !      do i=1,ntyp
279 !       do j=1,ntyp
280 !         aa(i,j)=0.0D0
281 !         bb(i,j)=0.0D0
282 !         augm(i,j)=0.0D0
283 !         sigma(i,j)=0.0D0
284 !         r0(i,j)=0.0D0
285 !         chi(i,j)=0.0D0
286 !        enddo
287 !       do j=1,2
288 !         bad(i,j)=0.0D0
289 !        enddo
290 !       chip(i)=0.0D0
291 !       alp(i)=0.0D0
292 !       sigma0(i)=0.0D0
293 !       sigii(i)=0.0D0
294 !       rr0(i)=0.0D0
295 !       a0thet(i)=0.0D0
296 !       do j=1,2
297 !         do ichir1=-1,1
298 !          do ichir2=-1,1
299 !          athet(j,i,ichir1,ichir2)=0.0D0
300 !          bthet(j,i,ichir1,ichir2)=0.0D0
301 !          enddo
302 !         enddo
303 !        enddo
304 !        do j=0,3
305 !         polthet(j,i)=0.0D0
306 !        enddo
307 !       do j=1,3
308 !         gthet(j,i)=0.0D0
309 !        enddo
310 !       theta0(i)=0.0D0
311 !       sig0(i)=0.0D0
312 !       sigc0(i)=0.0D0
313 !       do j=1,maxlob
314 !         bsc(j,i)=0.0D0
315 !         do k=1,3
316 !           censc(k,j,i)=0.0D0
317 !          enddo
318 !          do k=1,3
319 !           do l=1,3
320 !             gaussc(l,k,j,i)=0.0D0
321 !            enddo
322 !          enddo
323 !         nlob(i)=0
324 !        enddo
325 !      enddo
326 !      nlob(ntyp1)=0
327 !      dsc(ntyp1)=0.0D0
328 !      do i=-maxtor,maxtor
329 !        itortyp(i)=0
330 !c      write (iout,*) "TU DOCHODZE",i,itortyp(i)
331 !       do iblock=1,2
332 !        do j=-maxtor,maxtor
333 !          do k=1,maxterm
334 !            v1(k,j,i,iblock)=0.0D0
335 !            v2(k,j,i,iblock)=0.0D0
336 !          enddo
337 !        enddo
338 !        enddo
339 !      enddo
340 !      do iblock=1,2
341 !       do i=-maxtor,maxtor
342 !        do j=-maxtor,maxtor
343 !         do k=-maxtor,maxtor
344 !          do l=1,maxtermd_1
345 !            v1c(1,l,i,j,k,iblock)=0.0D0
346 !            v1s(1,l,i,j,k,iblock)=0.0D0
347 !            v1c(2,l,i,j,k,iblock)=0.0D0
348 !            v1s(2,l,i,j,k,iblock)=0.0D0
349 !          enddo !l
350 !          do l=1,maxtermd_2
351 !           do m=1,maxtermd_2
352 !            v2c(m,l,i,j,k,iblock)=0.0D0
353 !            v2s(m,l,i,j,k,iblock)=0.0D0
354 !           enddo !m
355 !          enddo !l
356 !        enddo !k
357 !       enddo !j
358 !      enddo !i
359 !      enddo !iblock
360
361 !      do i=1,maxres
362 !       itype(i,1)=0
363 !       itel(i)=0
364 !      enddo
365 ! Initialize the bridge arrays
366       ns=0
367       nss=0 
368       nhpb=0
369 !      do i=1,maxss
370 !       iss(i)=0
371 !      enddo
372 !      do i=1,maxdim
373 !       dhpb(i)=0.0D0
374 !      enddo
375 !      do i=1,maxres
376 !       ihpb(i)=0
377 !       jhpb(i)=0
378 !      enddo
379 !
380 ! Initialize timing.
381 !
382       call set_timers
383 !
384 ! Initialize variables used in minimization.
385 !   
386 !c     maxfun=5000
387 !c     maxit=2000
388       maxfun=500
389       maxit=200
390       tolf=1.0D-2
391       rtolf=5.0D-4
392
393 ! Initialize the variables responsible for the mode of gradient storage.
394 !
395       nfl=0
396       icg=1
397       
398 #ifdef WHAM_RUN
399       allocate(iww(max_eneW))
400       do i=1,14
401         do j=1,14
402           if (print_order(i).eq.j) then
403             iww(print_order(i))=j
404             goto 1121
405           endif
406         enddo
407 1121    continue
408       enddo
409 #endif
410  
411 #if defined(WHAM_RUN) || defined(CLUSTER)
412       ndih_constr=0
413
414 !      allocate(ww0(max_eneW))
415 !      ww0 = reshape((/1.0d0,1.0d0,1.0d0,1.0d0,1.0d0,1.0d0,1.0d0,1.0d0,&
416 !          1.0d0,1.0d0,1.0d0,1.0d0,1.0d0,1.0d0,1.0d0,1.0d0,0.4d0,1.0d0,&
417 !          1.0d0,0.0d0,0.0/), shape(ww0))
418 !
419       calc_grad=.false.
420 ! Set timers and counters for the respective routines
421       t_func = 0.0d0
422       t_grad = 0.0d0
423       t_fhel = 0.0d0
424       t_fbet = 0.0d0
425       t_ghel = 0.0d0
426       t_gbet = 0.0d0
427       t_viol = 0.0d0
428       t_gviol = 0.0d0
429       n_func = 0
430       n_grad = 0
431       n_fhel = 0
432       n_fbet = 0
433       n_ghel = 0
434       n_gbet = 0
435       n_viol = 0
436       n_gviol = 0
437       n_map = 0
438 #endif
439 !
440 ! Initialize constants used to split the energy into long- and short-range
441 ! components
442 !
443       r_cut=2.0d0
444       rlamb=0.3d0
445 #ifndef SPLITELE
446       nprint_ene=nprint_ene-1
447 #endif
448       return
449       end subroutine initialize
450 !-----------------------------------------------------------------------------
451       subroutine init_int_table
452
453       use geometry, only:int_bounds1
454 !el      use MPI_data
455 !el      implicit none
456 !      implicit real*8 (a-h,o-z)
457 !      include 'DIMENSIONS'
458 #ifdef MPI
459       include 'mpif.h'
460       integer,dimension(15) :: blocklengths,displs
461 #endif
462 !      include 'COMMON.CONTROL'
463 !      include 'COMMON.SETUP'
464 !      include 'COMMON.CHAIN'
465 !      include 'COMMON.INTERACT'
466 !      include 'COMMON.LOCAL'
467 !      include 'COMMON.SBRIDGE'
468 !      include 'COMMON.TORCNSTR'
469 !      include 'COMMON.IOUNITS'
470 !      include 'COMMON.DERIV'
471 !      include 'COMMON.CONTACTS'
472 !el      integer,dimension(0:nfgtasks) :: iturn3_start_all,iturn3_end_all,&
473 !el        iturn4_start_all,iturn4_end_all,iatel_s_all,iatel_e_all  !(0:max_fg_procs)
474 !el      integer,dimension(nres,0:nfgtasks) :: ielstart_all,&
475 !el        ielend_all !(maxres,0:max_fg_procs-1)
476 !el      integer,dimension(0:nfgtasks-1) :: ntask_cont_from_all,&
477 !el        ntask_cont_to_all !(0:max_fg_procs-1),
478 !el      integer,dimension(0:nfgtasks-1,0:nfgtasks-1) :: itask_cont_from_all,&
479 !el        itask_cont_to_all !(0:max_fg_procs-1,0:max_fg_procs-1)
480
481 !el      common /przechowalnia/ iturn3_start_all,iturn3_end_all,&
482 !el        iturn4_start_all,iturn4_end_all,iatel_s_all,iatel_e_all,&
483 !el        ielstart_all,ielend_all,ntask_cont_from_all,itask_cont_from_all,&
484 !el        ntask_cont_to_all,itask_cont_to_all
485
486       integer :: FG_GROUP,CONT_FROM_GROUP,CONT_TO_GROUP
487       logical :: scheck,lprint,flag
488
489 !el local variables
490       integer :: ind_scint=0,ind_scint_old,ii,jj,i,j,iint,itmp
491
492 #ifdef MPI
493       integer :: my_sc_int(0:nfgtasks-1),my_ele_int(0:nfgtasks-1)
494       integer :: my_sc_intt(0:nfgtasks),my_ele_intt(0:nfgtasks)
495       integer :: n_sc_int_tot,my_sc_inde,my_sc_inds,ind_sctint,npept
496       integer :: nele_int_tot,my_ele_inds,my_ele_inde,ind_eleint_old,&
497             ind_eleint,ijunk,nele_int_tot_vdw,my_ele_inds_vdw,&
498             my_ele_inde_vdw,ind_eleint_vdw,ind_eleint_vdw_old,&
499             nscp_int_tot,my_scp_inds,my_scp_inde,ind_scpint,&
500             ind_scpint_old,nsumgrad,nlen,ngrad_start,ngrad_end,&
501             ierror,k,ierr,iaux,ncheck_to,ncheck_from,ind_typ,&
502             ichunk,int_index_old
503       integer,dimension(5) :: nct_molec,nnt_molec
504 !el      allocate(itask_cont_from(0:nfgtasks-1)) !(0:max_fg_procs-1)
505 !el      allocate(itask_cont_to(0:nfgtasks-1)) !(0:max_fg_procs-1)
506
507 !... Determine the numbers of start and end SC-SC interaction
508 !... to deal with by current processor.
509 !write (iout,*) '******INIT_INT_TABLE nres=',nres,' nnt=',nnt,' nct=',nct
510       do i=0,nfgtasks-1
511         itask_cont_from(i)=fg_rank
512         itask_cont_to(i)=fg_rank
513       enddo
514       lprint=energy_dec
515       itmp=0
516       do i=1,5
517        if (nres_molec(i).eq.0) cycle
518       itmp=itmp+nres_molec(i)
519       if (itype(itmp,i).eq.ntyp1_molec(i)) then
520       nct_molec(i)=itmp-1
521       else
522       nct_molec(i)=itmp
523       endif
524       enddo
525 !      nct_molec(1)=nres_molec(1)-1
526       itmp=0
527       do i=2,5
528        itmp=itmp+nres_molec(i-1)
529       if (itype(itmp+1,i).eq.ntyp1_molec(i)) then
530       nnt_molec(i)=itmp+2
531       else
532       nnt_molec(i)=itmp+1
533       endif
534       enddo
535       print *,"nres_molec",nres_molec(:)
536       print *,"nnt_molec",nnt_molec(:)
537       print *,"nct_molec",nct_molec(:)
538 !      lprint=.true.
539       if (lprint) &
540        write (iout,*)'INIT_INT_TABLE nres=',nres,' nnt=',nnt,' nct=',nct
541       n_sc_int_tot=(nct_molec(1)-nnt+1)*(nct_molec(1)-nnt)/2-nss
542       call int_bounds(n_sc_int_tot,my_sc_inds,my_sc_inde)
543 !write (iout,*) 'INIT_INT_TABLE nres=',nres,' nnt=',nnt,' nct=',nct
544       if (lprint) &
545         write (iout,*) 'Processor',fg_rank,' CG group',kolor,&
546         ' absolute rank',MyRank,&
547         ' n_sc_int_tot',n_sc_int_tot,' my_sc_inds=',my_sc_inds,&
548         ' my_sc_inde',my_sc_inde
549       ind_sctint=0
550       iatsc_s=0
551       iatsc_e=0
552 #endif
553 !el       common /przechowalnia/
554       allocate(iturn3_start_all(0:nfgtasks))
555       allocate(iturn3_end_all(0:nfgtasks))
556       allocate(iturn4_start_all(0:nfgtasks))
557       allocate(iturn4_end_all(0:nfgtasks))
558       allocate(iatel_s_all(0:nfgtasks))
559       allocate(iatel_e_all(0:nfgtasks))
560       allocate(ielstart_all(nres,0:nfgtasks-1))
561       allocate(ielend_all(nres,0:nfgtasks-1))
562
563       allocate(ntask_cont_from_all(0:nfgtasks-1))
564       allocate(ntask_cont_to_all(0:nfgtasks-1))
565       allocate(itask_cont_from_all(0:nfgtasks-1,0:nfgtasks-1))
566       allocate(itask_cont_to_all(0:nfgtasks-1,0:nfgtasks-1))
567 !el----------
568 !      lprint=.false.
569         print *,"NCT",nct_molec(1),nct
570       do i=1,nres !el   !maxres
571         nint_gr(i)=0
572         nscp_gr(i)=0
573         ielstart(i)=0
574         ielend(i)=0
575         do j=1,maxint_gr
576           istart(i,j)=0
577           iend(i,j)=0
578           iscpstart(i,j)=0
579           iscpend(i,j)=0    
580         enddo
581       enddo
582       ind_scint=0
583       ind_scint_old=0
584 !d    write (iout,*) 'ns=',ns,' nss=',nss,' ihpb,jhpb',
585 !d   &   (ihpb(i),jhpb(i),i=1,nss)
586 !       print *,nnt,nct_molec(1)
587       do i=nnt,nct_molec(1)-1
588 !        print*, "inloop",i
589         scheck=.false.
590         if (dyn_ss) goto 10
591         do ii=1,nss
592           if (ihpb(ii).eq.i+nres) then
593             scheck=.true.
594             jj=jhpb(ii)-nres
595             goto 10
596           endif
597         enddo
598    10   continue
599 !        print *,'i=',i,' scheck=',scheck,' jj=',jj
600 !d      write (iout,*) 'i=',i,' scheck=',scheck,' jj=',jj
601         if (scheck) then
602           if (jj.eq.i+1) then
603 #ifdef MPI
604 !            write (iout,*) 'jj=i+1'
605             call int_partition(ind_scint,my_sc_inds,my_sc_inde,i,&
606        iatsc_s,iatsc_e,i+2,nct_molec(1),nint_gr(i),istart(i,1),iend(i,1),*12)
607 #else
608             nint_gr(i)=1
609             istart(i,1)=i+2
610             iend(i,1)=nct
611 #endif
612           else if (jj.eq.nct_molec(1)) then
613 #ifdef MPI
614 !            write (iout,*) 'jj=nct'
615             call int_partition(ind_scint,my_sc_inds,my_sc_inde,i,&
616         iatsc_s,iatsc_e,i+1,nct_molec(1)-1,nint_gr(i),istart(i,1),iend(i,1),*12)
617 #else
618             nint_gr(i)=1
619             istart(i,1)=i+1
620             iend(i,1)=nct_molecule(1)-1
621 #endif
622           else
623 #ifdef MPI
624             call int_partition(ind_scint,my_sc_inds,my_sc_inde,i,&
625        iatsc_s,iatsc_e,i+1,jj-1,nint_gr(i),istart(i,1),iend(i,1),*12)
626             ii=nint_gr(i)+1
627             call int_partition(ind_scint,my_sc_inds,my_sc_inde,i,&
628        iatsc_s,iatsc_e,jj+1,nct_molec(1),nint_gr(i),istart(i,ii),iend(i,ii),*12)
629          
630 #else
631             nint_gr(i)=2
632             istart(i,1)=i+1
633             iend(i,1)=jj-1
634             istart(i,2)=jj+1
635             iend(i,2)=nct_molec(1)
636 #endif
637           endif
638         else
639 #ifdef MPI
640 !          print *,"i for EVDW",iatsc_s,iatsc_e,istart(i,1),iend(i,1),&
641 !          i+1,nct_molec(1),nint_gr(i),ind_scint,my_sc_inds,my_sc_inde,i
642           call int_partition(ind_scint,my_sc_inds,my_sc_inde,i,&
643           iatsc_s,iatsc_e,i+1,nct_molec(1),nint_gr(i), &
644           istart(i,1),iend(i,1),*12)
645 !          print *,"i for EVDW",iatsc_s,iatsc_e,istart(i,1),iend(i,1)
646 #else
647           nint_gr(i)=1
648           istart(i,1)=i+1
649           iend(i,1)=nct_molec(1)
650           ind_scint=ind_scint+nct_molec(1)-i
651 #endif
652         endif
653 #ifdef MPI
654         ind_scint_old=ind_scint
655 #endif
656       enddo
657    12 continue
658 !      print *,"i for EVDW",iatsc_s,iatsc_e,istart(i,1),iend(i,1)
659
660 #ifndef MPI
661       iatsc_s=nnt
662       iatsc_e=nct-1
663 #endif
664       if (iatsc_s.eq.0) iatsc_s=1
665 #ifdef MPI
666       if (lprint) write (*,*) 'Processor',fg_rank,' CG Group',kolor,&
667          ' absolute rank',myrank,' iatsc_s=',iatsc_s,' iatsc_e=',iatsc_e
668 #endif
669 !      lprint=.true.
670       if (lprint) then
671       write (iout,'(a)') 'Interaction array:'
672       do i=iatsc_s,iatsc_e
673         write (iout,'(i3,2(2x,2i3))') &
674        i,(istart(i,iint),iend(i,iint),iint=1,nint_gr(i))
675       enddo
676       endif
677 !      lprint=.false.
678       ispp=4 !?? wham ispp=2
679 #ifdef MPI
680 ! Now partition the electrostatic-interaction array
681       if (itype(nres_molec(1),1).eq.ntyp1_molec(1)) then
682       npept=nres_molec(1)-nnt-1
683       else
684       npept=nres_molec(1)-nnt
685       endif
686       nele_int_tot=(npept-ispp)*(npept-ispp+1)/2
687       call int_bounds(nele_int_tot,my_ele_inds,my_ele_inde)
688       if (lprint) &
689        write (*,*) 'Processor',fg_rank,' CG group',kolor,&
690         ' absolute rank',MyRank,&
691         ' nele_int_tot',nele_int_tot,' my_ele_inds=',my_ele_inds,&
692                     ' my_ele_inde',my_ele_inde
693       iatel_s=0
694       iatel_e=0
695       ind_eleint=0
696       ind_eleint_old=0
697 !      if (itype(nres_molec(1),1).eq.ntyp1_molec(1)) then
698 !      nct_molec(1)=nres_molec(1)-1
699 !      else
700 !      nct_molec(1)=nres_molec(1)
701 !      endif
702 !       print *,"nct",nct,nct_molec(1),itype(nres_molec(1),1),ntyp_molec(1)
703       do i=nnt,nct_molec(1)-3
704         ijunk=0
705         call int_partition(ind_eleint,my_ele_inds,my_ele_inde,i,&
706           iatel_s,iatel_e,i+ispp,nct_molec(1)-1,ijunk,ielstart(i),ielend(i),*13)
707       enddo ! i 
708    13 continue
709       if (iatel_s.eq.0) iatel_s=1
710       nele_int_tot_vdw=(npept-2)*(npept-2+1)/2
711 !      write (iout,*) "nele_int_tot_vdw",nele_int_tot_vdw
712       call int_bounds(nele_int_tot_vdw,my_ele_inds_vdw,my_ele_inde_vdw)
713 !      write (iout,*) "my_ele_inds_vdw",my_ele_inds_vdw,
714 !     & " my_ele_inde_vdw",my_ele_inde_vdw
715       ind_eleint_vdw=0
716       ind_eleint_vdw_old=0
717       iatel_s_vdw=0
718       iatel_e_vdw=0
719       do i=nnt,nct_molec(1)-3
720         ijunk=0
721         call int_partition(ind_eleint_vdw,my_ele_inds_vdw,&
722           my_ele_inde_vdw,i,&
723           iatel_s_vdw,iatel_e_vdw,i+2,nct_molec(1)-1,ijunk,ielstart_vdw(i),&
724           ielend_vdw(i),*15)
725 !        write (iout,*) i," ielstart_vdw",ielstart_vdw(i),
726 !     &   " ielend_vdw",ielend_vdw(i)
727       enddo ! i 
728       if (iatel_s_vdw.eq.0) iatel_s_vdw=1
729    15 continue
730 #else
731       iatel_s=nnt
732       iatel_e=nct_molec(1)-5 ! ?? wham iatel_e=nct-3
733       do i=iatel_s,iatel_e
734         ielstart(i)=i+4 ! ?? wham +2
735         ielend(i)=nct_molec(1)-1
736       enddo
737       iatel_s_vdw=nnt
738       iatel_e_vdw=nct_molec(1)-3
739       do i=iatel_s_vdw,iatel_e_vdw
740         ielstart_vdw(i)=i+2
741         ielend_vdw(i)=nct_molec(1)-1
742       enddo
743 #endif
744       if (lprint) then
745         write (*,'(a)') 'Processor',fg_rank,' CG group',kolor,&
746         ' absolute rank',MyRank
747         write (iout,*) 'Electrostatic interaction array:'
748         do i=iatel_s,iatel_e
749           write (iout,'(i3,2(2x,2i3))') i,ielstart(i),ielend(i)
750         enddo
751       endif ! lprint
752 !     iscp=3
753       iscp=2
754 ! Partition the SC-p interaction array
755 #ifdef MPI
756       nscp_int_tot=(npept-iscp+1)*(npept-iscp+1)
757       call int_bounds(nscp_int_tot,my_scp_inds,my_scp_inde)
758       if (lprint) write (iout,*) 'Processor',fg_rank,' CG group',kolor,&
759         ' absolute rank',myrank,&
760         ' nscp_int_tot',nscp_int_tot,' my_scp_inds=',my_scp_inds,&
761                     ' my_scp_inde',my_scp_inde
762       iatscp_s=0
763       iatscp_e=0
764       ind_scpint=0
765       ind_scpint_old=0
766       do i=nnt,nct_molec(1)-1
767         if (i.lt.nnt+iscp) then
768 !d        write (iout,*) 'i.le.nnt+iscp'
769           call int_partition(ind_scpint,my_scp_inds,my_scp_inde,i,&
770             iatscp_s,iatscp_e,i+iscp,nct_molec(1),nscp_gr(i),iscpstart(i,1),&
771             iscpend(i,1),*14)
772         else if (i.gt.nct-iscp) then
773 !d        write (iout,*) 'i.gt.nct-iscp'
774           call int_partition(ind_scpint,my_scp_inds,my_scp_inde,i,&
775             iatscp_s,iatscp_e,nnt,i-iscp,nscp_gr(i),iscpstart(i,1),&
776             iscpend(i,1),*14)
777         else
778           call int_partition(ind_scpint,my_scp_inds,my_scp_inde,i,&
779             iatscp_s,iatscp_e,nnt,i-iscp,nscp_gr(i),iscpstart(i,1),&
780            iscpend(i,1),*14)
781           ii=nscp_gr(i)+1
782           call int_partition(ind_scpint,my_scp_inds,my_scp_inde,i,&
783             iatscp_s,iatscp_e,i+iscp,nct_molec(1),nscp_gr(i),iscpstart(i,ii),&
784             iscpend(i,ii),*14)
785         endif
786       enddo ! i
787    14 continue
788 #else
789       iatscp_s=nnt
790       iatscp_e=nct_molec(1)-1
791       do i=nnt,nct_molec(1)-1
792         if (i.lt.nnt+iscp) then
793           nscp_gr(i)=1
794           iscpstart(i,1)=i+iscp
795           iscpend(i,1)=nct_molec(1)
796         elseif (i.gt.nct-iscp) then
797           nscp_gr(i)=1
798           iscpstart(i,1)=nnt
799           iscpend(i,1)=i-iscp
800         else
801           nscp_gr(i)=2
802           iscpstart(i,1)=nnt
803           iscpend(i,1)=i-iscp
804           iscpstart(i,2)=i+iscp
805           iscpend(i,2)=nct_molec(1)
806         endif 
807       enddo ! i
808 #endif
809       if (iatscp_s.eq.0) iatscp_s=1
810       if (lprint) then
811         write (iout,'(a)') 'SC-p interaction array:'
812         do i=iatscp_s,iatscp_e
813           write (iout,'(i3,2(2x,2i3))') &
814               i,(iscpstart(i,j),iscpend(i,j),j=1,nscp_gr(i))
815         enddo
816       endif ! lprint
817 ! Partition local interactions
818 #ifdef MPI
819       call int_bounds(nres_molec(1)-2,loc_start,loc_end)
820       loc_start=loc_start+1
821       loc_end=loc_end+1
822       call int_bounds(nres_molec(1)-2,ithet_start,ithet_end)
823       ithet_start=ithet_start+2
824       ithet_end=ithet_end+2
825       call int_bounds(nres_molec(2)-2,ithet_nucl_start,ithet_nucl_end)
826       ithet_nucl_start=ithet_nucl_start+2+nres_molec(1)
827       ithet_nucl_end=ithet_nucl_end+2+nres_molec(1)
828       call int_bounds(nct_molec(1)-nnt-2,iturn3_start,iturn3_end) 
829       iturn3_start=iturn3_start+nnt
830       iphi_start=iturn3_start+2
831       iturn3_end=iturn3_end+nnt
832       iphi_end=iturn3_end+2
833       iturn3_start=iturn3_start-1
834       iturn3_end=iturn3_end-1
835       call int_bounds(nct_molec(2)-nnt_molec(2)-2,iphi_nucl_start,iphi_nucl_end)
836       iphi_nucl_start=iphi_nucl_start+nnt_molec(2)+2
837       iphi_nucl_end=iphi_nucl_end+nnt_molec(2)+2
838       print *,"KURDE",iphi_nucl_start,iphi_nucl_end
839       call int_bounds(nres_molec(1)-3,itau_start,itau_end)
840       itau_start=itau_start+3
841       itau_end=itau_end+3
842       call int_bounds(nres_molec(1)-3,iphi1_start,iphi1_end)
843       iphi1_start=iphi1_start+3
844       iphi1_end=iphi1_end+3
845       call int_bounds(nct_molec(1)-nnt-3,iturn4_start,iturn4_end) 
846       iturn4_start=iturn4_start+nnt
847       iphid_start=iturn4_start+2
848       iturn4_end=iturn4_end+nnt
849       iphid_end=iturn4_end+2
850       iturn4_start=iturn4_start-1
851       iturn4_end=iturn4_end-1
852 !      print *,"TUTUTU",nres_molec(1),nres
853       call int_bounds(nres_molec(1)-2,ibond_start,ibond_end) 
854       ibond_start=ibond_start+1
855       ibond_end=ibond_end+1
856 !      print *,ibond_start,ibond_end
857       call int_bounds(nct_molec(1)-nnt,ibondp_start,ibondp_end) 
858       ibondp_start=ibondp_start+nnt
859       ibondp_end=ibondp_end+nnt
860      call int_bounds(nres_molec(2)-2,ibond_nucl_start,ibond_nucl_end)
861       ibond_nucl_start=ibond_nucl_start+nnt_molec(2)-1
862       ibond_nucl_end=ibond_nucl_end+nnt_molec(2)-1
863       print *,"NUCLibond",ibond_nucl_start,ibond_nucl_end
864       print *, "before devision",nnt_molec(2),nct_molec(2)-nnt_molec(2)
865       call int_bounds(nct_molec(2)-nnt_molec(2),ibondp_nucl_start,ibondp_nucl_end)
866       ibondp_nucl_start=ibondp_nucl_start+nnt_molec(2)
867       ibondp_nucl_end=ibondp_nucl_end+nnt_molec(2)
868       print *,"NUCLibond2",ibondp_nucl_start,ibondp_nucl_end
869
870
871       call int_bounds1(nres_molec(1)-1,ivec_start,ivec_end) 
872 !      print *,"Processor",myrank,fg_rank,fg_rank1,
873 !     &  " ivec_start",ivec_start," ivec_end",ivec_end
874       iset_start=loc_start+2
875       iset_end=loc_end+2
876       call int_bounds(nres_molec(1),ilip_start,ilip_end)
877       ilip_start=ilip_start
878       ilip_end=ilip_end
879       call int_bounds(nres_molec(1)-1,itube_start,itube_end)
880       itube_start=itube_start
881       itube_end=itube_end
882       if (ndih_constr.eq.0) then
883         idihconstr_start=1
884         idihconstr_end=0
885       else
886         call int_bounds(ndih_constr,idihconstr_start,idihconstr_end)
887       endif
888       if (ntheta_constr.eq.0) then
889         ithetaconstr_start=1
890         ithetaconstr_end=0
891       else
892         call int_bounds &
893        (ntheta_constr,ithetaconstr_start,ithetaconstr_end)
894       endif
895
896 !      nsumgrad=(nres-nnt)*(nres-nnt+1)/2
897 !      nlen=nres-nnt+1
898       nsumgrad=(nres-nnt)*(nres-nnt+1)/2
899       nlen=nres-nnt+1
900       call int_bounds(nsumgrad,ngrad_start,ngrad_end)
901       igrad_start=((2*nlen+1) &
902          -sqrt(float((2*nlen-1)**2-8*(ngrad_start-1))))/2
903       igrad_end=((2*nlen+1) &
904          -sqrt(float((2*nlen-1)**2-8*(ngrad_end-1))))/2
905 !el      allocate(jgrad_start(igrad_start:igrad_end))
906 !el      allocate(jgrad_end(igrad_start:igrad_end)) !(maxres)
907       jgrad_start(igrad_start)= &
908          ngrad_start-(2*nlen-igrad_start)*(igrad_start-1)/2 &
909          +igrad_start
910       jgrad_end(igrad_start)=nres
911       if (igrad_end.gt.igrad_start) jgrad_start(igrad_end)=igrad_end+1
912       jgrad_end(igrad_end)=ngrad_end-(2*nlen-igrad_end)*(igrad_end-1)/2 &
913           +igrad_end
914       do i=igrad_start+1,igrad_end-1
915         jgrad_start(i)=i+1
916         jgrad_end(i)=nres
917       enddo
918       if (lprint) then 
919         write (*,*) 'Processor:',fg_rank,' CG group',kolor,&
920        ' absolute rank',myrank,&
921        ' loc_start',loc_start,' loc_end',loc_end,&
922        ' ithet_start',ithet_start,' ithet_end',ithet_end,&
923        ' iphi_start',iphi_start,' iphi_end',iphi_end,&
924        ' iphid_start',iphid_start,' iphid_end',iphid_end,&
925        ' ibond_start',ibond_start,' ibond_end',ibond_end,&
926        ' ibondp_start',ibondp_start,' ibondp_end',ibondp_end,&
927        ' iturn3_start',iturn3_start,' iturn3_end',iturn3_end,&
928        ' iturn4_start',iturn4_start,' iturn4_end',iturn4_end,&
929        ' ivec_start',ivec_start,' ivec_end',ivec_end,&
930        ' iset_start',iset_start,' iset_end',iset_end,&
931        ' idihconstr_start',idihconstr_start,' idihconstr_end',&
932          idihconstr_end
933        write (*,*) 'Processor:',fg_rank,myrank,' igrad_start',&
934          igrad_start,' igrad_end',igrad_end,' ngrad_start',ngrad_start,&
935          ' ngrad_end',ngrad_end
936        do i=igrad_start,igrad_end
937          write(*,*) 'Processor:',fg_rank,myrank,i,&
938           jgrad_start(i),jgrad_end(i)
939        enddo
940       endif
941       if (nfgtasks.gt.1) then
942         call MPI_Allgather(ivec_start,1,MPI_INTEGER,ivec_displ(0),1,&
943           MPI_INTEGER,FG_COMM1,IERROR)
944         iaux=ivec_end-ivec_start+1
945         call MPI_Allgather(iaux,1,MPI_INTEGER,ivec_count(0),1,&
946           MPI_INTEGER,FG_COMM1,IERROR)
947         call MPI_Allgather(iset_start-2,1,MPI_INTEGER,iset_displ(0),1,&
948           MPI_INTEGER,FG_COMM,IERROR)
949         iaux=iset_end-iset_start+1
950         call MPI_Allgather(iaux,1,MPI_INTEGER,iset_count(0),1,&
951           MPI_INTEGER,FG_COMM,IERROR)
952         call MPI_Allgather(ibond_start,1,MPI_INTEGER,ibond_displ(0),1,&
953           MPI_INTEGER,FG_COMM,IERROR)
954         iaux=ibond_end-ibond_start+1
955         call MPI_Allgather(iaux,1,MPI_INTEGER,ibond_count(0),1,&
956           MPI_INTEGER,FG_COMM,IERROR)
957         call MPI_Allgather(ithet_start,1,MPI_INTEGER,ithet_displ(0),1,&
958           MPI_INTEGER,FG_COMM,IERROR)
959         iaux=ithet_end-ithet_start+1
960         call MPI_Allgather(iaux,1,MPI_INTEGER,ithet_count(0),1,&
961           MPI_INTEGER,FG_COMM,IERROR)
962         call MPI_Allgather(iphi_start,1,MPI_INTEGER,iphi_displ(0),1,&
963           MPI_INTEGER,FG_COMM,IERROR)
964         iaux=iphi_end-iphi_start+1
965         call MPI_Allgather(iaux,1,MPI_INTEGER,iphi_count(0),1,&
966           MPI_INTEGER,FG_COMM,IERROR)
967         call MPI_Allgather(iphi1_start,1,MPI_INTEGER,iphi1_displ(0),1,&
968           MPI_INTEGER,FG_COMM,IERROR)
969         iaux=iphi1_end-iphi1_start+1
970         call MPI_Allgather(iaux,1,MPI_INTEGER,iphi1_count(0),1,&
971           MPI_INTEGER,FG_COMM,IERROR)
972         do i=0,nfgtasks-1
973           do j=1,nres
974             ielstart_all(j,i)=0
975             ielend_all(j,i)=0
976           enddo
977         enddo
978         call MPI_Allgather(iturn3_start,1,MPI_INTEGER,&
979           iturn3_start_all(0),1,MPI_INTEGER,FG_COMM,IERROR)
980         call MPI_Allgather(iturn4_start,1,MPI_INTEGER,&
981           iturn4_start_all(0),1,MPI_INTEGER,FG_COMM,IERROR)
982         call MPI_Allgather(iturn3_end,1,MPI_INTEGER,&
983           iturn3_end_all(0),1,MPI_INTEGER,FG_COMM,IERROR)
984         call MPI_Allgather(iturn4_end,1,MPI_INTEGER,&
985           iturn4_end_all(0),1,MPI_INTEGER,FG_COMM,IERROR)
986         call MPI_Allgather(iatel_s,1,MPI_INTEGER,&
987           iatel_s_all(0),1,MPI_INTEGER,FG_COMM,IERROR)
988         call MPI_Allgather(iatel_e,1,MPI_INTEGER,&
989           iatel_e_all(0),1,MPI_INTEGER,FG_COMM,IERROR)
990         call MPI_Allgather(ielstart(1),nres,MPI_INTEGER,&
991           ielstart_all(1,0),nres,MPI_INTEGER,FG_COMM,IERROR)
992         call MPI_Allgather(ielend(1),nres,MPI_INTEGER,&
993           ielend_all(1,0),nres,MPI_INTEGER,FG_COMM,IERROR)
994         if (lprint) then
995         write (iout,*) "iatel_s_all",(iatel_s_all(i),i=0,nfgtasks)
996         write (iout,*) "iatel_e_all",(iatel_e_all(i),i=0,nfgtasks)
997         write (iout,*) "iturn3_start_all",&
998           (iturn3_start_all(i),i=0,nfgtasks-1)
999         write (iout,*) "iturn3_end_all",&
1000           (iturn3_end_all(i),i=0,nfgtasks-1)
1001         write (iout,*) "iturn4_start_all",&
1002           (iturn4_start_all(i),i=0,nfgtasks-1)
1003         write (iout,*) "iturn4_end_all",&
1004           (iturn4_end_all(i),i=0,nfgtasks-1)
1005         write (iout,*) "The ielstart_all array"
1006         do i=nnt,nct
1007           write (iout,'(20i4)') i,(ielstart_all(i,j),j=0,nfgtasks-1)
1008         enddo
1009         write (iout,*) "The ielend_all array"
1010         do i=nnt,nct
1011           write (iout,'(20i4)') i,(ielend_all(i,j),j=0,nfgtasks-1)
1012         enddo
1013         call flush(iout)
1014         endif
1015         ntask_cont_from=0
1016         ntask_cont_to=0
1017         itask_cont_from(0)=fg_rank
1018         itask_cont_to(0)=fg_rank
1019         flag=.false.
1020 !el        allocate(iturn3_sent(4,iturn3_start:iturn3_end))
1021 !el        allocate(iturn4_sent(4,iturn4_start:iturn4_end)) !(4,maxres)
1022         do ii=iturn3_start,iturn3_end
1023           call add_int(ii,ii+2,iturn3_sent(1,ii),&
1024                       ntask_cont_to,itask_cont_to,flag)
1025         enddo
1026         do ii=iturn4_start,iturn4_end
1027           call add_int(ii,ii+3,iturn4_sent(1,ii),&
1028                       ntask_cont_to,itask_cont_to,flag)
1029         enddo
1030         do ii=iturn3_start,iturn3_end
1031           call add_int_from(ii,ii+2,ntask_cont_from,itask_cont_from)
1032         enddo
1033         do ii=iturn4_start,iturn4_end
1034           call add_int_from(ii,ii+3,ntask_cont_from,itask_cont_from)
1035         enddo
1036         if (lprint) then
1037         write (iout,*) "After turn3 ntask_cont_from",ntask_cont_from,&
1038          " ntask_cont_to",ntask_cont_to
1039         write (iout,*) "itask_cont_from",&
1040           (itask_cont_from(i),i=1,ntask_cont_from)
1041         write (iout,*) "itask_cont_to",&
1042           (itask_cont_to(i),i=1,ntask_cont_to)
1043         call flush(iout)
1044         endif
1045 !        write (iout,*) "Loop forward"
1046 !        call flush(iout)
1047         do i=iatel_s,iatel_e
1048 !          write (iout,*) "from loop i=",i
1049 !          call flush(iout)
1050           do j=ielstart(i),ielend(i)
1051             call add_int_from(i,j,ntask_cont_from,itask_cont_from)
1052           enddo
1053         enddo
1054 !        write (iout,*) "Loop backward iatel_e-1",iatel_e-1,
1055 !     &     " iatel_e",iatel_e
1056 !        call flush(iout)
1057         nat_sent=0
1058         do i=iatel_s,iatel_e
1059 !          write (iout,*) "i",i," ielstart",ielstart(i),
1060 !     &      " ielend",ielend(i)
1061 !          call flush(iout)
1062           flag=.false.
1063           do j=ielstart(i),ielend(i)
1064             call add_int(i,j,iint_sent(1,j,nat_sent+1),ntask_cont_to,&
1065                         itask_cont_to,flag)
1066           enddo
1067           if (flag) then
1068             nat_sent=nat_sent+1
1069             iat_sent(nat_sent)=i
1070           endif
1071         enddo
1072         if (lprint) then
1073         write (iout,*)"After longrange ntask_cont_from",ntask_cont_from,&
1074          " ntask_cont_to",ntask_cont_to
1075         write (iout,*) "itask_cont_from",&
1076           (itask_cont_from(i),i=1,ntask_cont_from)
1077         write (iout,*) "itask_cont_to",&
1078           (itask_cont_to(i),i=1,ntask_cont_to)
1079         call flush(iout)
1080         write (iout,*) "iint_sent"
1081         do i=1,nat_sent
1082           ii=iat_sent(i)
1083           write (iout,'(20i4)') ii,(j,(iint_sent(k,j,i),k=1,4),&
1084             j=ielstart(ii),ielend(ii))
1085         enddo
1086         write (iout,*) "iturn3_sent iturn3_start",iturn3_start,&
1087           " iturn3_end",iturn3_end
1088         write (iout,'(20i4)') (i,(iturn3_sent(j,i),j=1,4),&
1089            i=iturn3_start,iturn3_end)
1090         write (iout,*) "iturn4_sent iturn4_start",iturn4_start,&
1091           " iturn4_end",iturn4_end
1092         write (iout,'(20i4)') (i,(iturn4_sent(j,i),j=1,4),&
1093            i=iturn4_start,iturn4_end)
1094         call flush(iout)
1095         endif
1096         call MPI_Gather(ntask_cont_from,1,MPI_INTEGER,&
1097          ntask_cont_from_all,1,MPI_INTEGER,king,FG_COMM,IERR)
1098 !        write (iout,*) "Gather ntask_cont_from ended"
1099 !        call flush(iout)
1100         call MPI_Gather(itask_cont_from(0),nfgtasks,MPI_INTEGER,&
1101          itask_cont_from_all(0,0),nfgtasks,MPI_INTEGER,king,&
1102          FG_COMM,IERR)
1103 !        write (iout,*) "Gather itask_cont_from ended"
1104 !        call flush(iout)
1105         call MPI_Gather(ntask_cont_to,1,MPI_INTEGER,ntask_cont_to_all,&
1106          1,MPI_INTEGER,king,FG_COMM,IERR)
1107 !        write (iout,*) "Gather ntask_cont_to ended"
1108 !        call flush(iout)
1109         call MPI_Gather(itask_cont_to,nfgtasks,MPI_INTEGER,&
1110          itask_cont_to_all,nfgtasks,MPI_INTEGER,king,FG_COMM,IERR)
1111 !        write (iout,*) "Gather itask_cont_to ended"
1112 !        call flush(iout)
1113         if (fg_rank.eq.king) then
1114           write (iout,*)"Contact receive task map (proc, #tasks, tasks)"
1115           do i=0,nfgtasks-1
1116             write (iout,'(20i4)') i,ntask_cont_from_all(i),&
1117               (itask_cont_from_all(j,i),j=1,ntask_cont_from_all(i)) 
1118           enddo
1119           write (iout,*)
1120           call flush(iout)
1121           write (iout,*) "Contact send task map (proc, #tasks, tasks)"
1122           do i=0,nfgtasks-1
1123             write (iout,'(20i4)') i,ntask_cont_to_all(i),&
1124              (itask_cont_to_all(j,i),j=1,ntask_cont_to_all(i)) 
1125           enddo
1126           write (iout,*)
1127           call flush(iout)
1128 ! Check if every send will have a matching receive
1129           ncheck_to=0
1130           ncheck_from=0
1131           do i=0,nfgtasks-1
1132             ncheck_to=ncheck_to+ntask_cont_to_all(i)
1133             ncheck_from=ncheck_from+ntask_cont_from_all(i)
1134           enddo
1135           write (iout,*) "Control sums",ncheck_from,ncheck_to
1136           if (ncheck_from.ne.ncheck_to) then
1137             write (iout,*) "Error: #receive differs from #send."
1138             write (iout,*) "Terminating program...!"
1139             call flush(iout)
1140             flag=.false.
1141           else
1142             flag=.true.
1143             do i=0,nfgtasks-1
1144               do j=1,ntask_cont_to_all(i)
1145                 ii=itask_cont_to_all(j,i)
1146                 do k=1,ntask_cont_from_all(ii)
1147                   if (itask_cont_from_all(k,ii).eq.i) then
1148                     if(lprint)write(iout,*)"Matching send/receive",i,ii
1149                     exit
1150                   endif
1151                 enddo
1152                 if (k.eq.ntask_cont_from_all(ii)+1) then
1153                   flag=.false.
1154                   write (iout,*) "Error: send by",j," to",ii,&
1155                     " would have no matching receive"
1156                 endif
1157               enddo
1158             enddo
1159           endif
1160           if (.not.flag) then
1161             write (iout,*) "Unmatched sends; terminating program"
1162             call flush(iout)
1163           endif
1164         endif
1165         call MPI_Bcast(flag,1,MPI_LOGICAL,king,FG_COMM,IERROR)
1166 !        write (iout,*) "flag broadcast ended flag=",flag
1167 !        call flush(iout)
1168         if (.not.flag) then
1169           call MPI_Finalize(IERROR)
1170           stop "Error in INIT_INT_TABLE: unmatched send/receive."
1171         endif
1172         call MPI_Comm_group(FG_COMM,fg_group,IERR)
1173 !        write (iout,*) "MPI_Comm_group ended"
1174 !        call flush(iout)
1175         call MPI_Group_incl(fg_group,ntask_cont_from+1,&
1176           itask_cont_from(0),CONT_FROM_GROUP,IERR)
1177         call MPI_Group_incl(fg_group,ntask_cont_to+1,itask_cont_to(0),&
1178           CONT_TO_GROUP,IERR)
1179         do i=1,nat_sent
1180           ii=iat_sent(i)
1181           iaux=4*(ielend(ii)-ielstart(ii)+1)
1182           call MPI_Group_translate_ranks(fg_group,iaux,&
1183             iint_sent(1,ielstart(ii),i),CONT_TO_GROUP,&
1184             iint_sent_local(1,ielstart(ii),i),IERR )
1185 !          write (iout,*) "Ranks translated i=",i
1186 !          call flush(iout)
1187         enddo
1188         iaux=4*(iturn3_end-iturn3_start+1)
1189         call MPI_Group_translate_ranks(fg_group,iaux,&
1190            iturn3_sent(1,iturn3_start),CONT_TO_GROUP,&
1191            iturn3_sent_local(1,iturn3_start),IERR)
1192         iaux=4*(iturn4_end-iturn4_start+1)
1193         call MPI_Group_translate_ranks(fg_group,iaux,&
1194            iturn4_sent(1,iturn4_start),CONT_TO_GROUP,&
1195            iturn4_sent_local(1,iturn4_start),IERR)
1196         if (lprint) then
1197         write (iout,*) "iint_sent_local"
1198         do i=1,nat_sent
1199           ii=iat_sent(i)
1200           write (iout,'(20i4)') ii,(j,(iint_sent_local(k,j,i),k=1,4),&
1201             j=ielstart(ii),ielend(ii))
1202           call flush(iout)
1203         enddo
1204         write (iout,*) "iturn3_sent_local iturn3_start",iturn3_start,&
1205           " iturn3_end",iturn3_end
1206         write (iout,'(20i4)') (i,(iturn3_sent_local(j,i),j=1,4),&
1207            i=iturn3_start,iturn3_end)
1208         write (iout,*) "iturn4_sent_local iturn4_start",iturn4_start,&
1209           " iturn4_end",iturn4_end
1210         write (iout,'(20i4)') (i,(iturn4_sent_local(j,i),j=1,4),&
1211            i=iturn4_start,iturn4_end)
1212         call flush(iout)
1213         endif
1214         call MPI_Group_free(fg_group,ierr)
1215         call MPI_Group_free(cont_from_group,ierr)
1216         call MPI_Group_free(cont_to_group,ierr)
1217         call MPI_Type_contiguous(3,MPI_DOUBLE_PRECISION,MPI_UYZ,IERROR)
1218         call MPI_Type_commit(MPI_UYZ,IERROR)
1219         call MPI_Type_contiguous(18,MPI_DOUBLE_PRECISION,MPI_UYZGRAD,&
1220           IERROR)
1221         call MPI_Type_commit(MPI_UYZGRAD,IERROR)
1222         call MPI_Type_contiguous(2,MPI_DOUBLE_PRECISION,MPI_MU,IERROR)
1223         call MPI_Type_commit(MPI_MU,IERROR)
1224         call MPI_Type_contiguous(4,MPI_DOUBLE_PRECISION,MPI_MAT1,IERROR)
1225         call MPI_Type_commit(MPI_MAT1,IERROR)
1226         call MPI_Type_contiguous(8,MPI_DOUBLE_PRECISION,MPI_MAT2,IERROR)
1227         call MPI_Type_commit(MPI_MAT2,IERROR)
1228         call MPI_Type_contiguous(6,MPI_DOUBLE_PRECISION,MPI_THET,IERROR)
1229         call MPI_Type_commit(MPI_THET,IERROR)
1230         call MPI_Type_contiguous(9,MPI_DOUBLE_PRECISION,MPI_GAM,IERROR)
1231         call MPI_Type_commit(MPI_GAM,IERROR)
1232
1233 !el        allocate(lentyp(0:nfgtasks-1))
1234 #ifndef MATGATHER
1235 ! 9/22/08 Derived types to send matrices which appear in correlation terms
1236         do i=0,nfgtasks-1
1237           if (ivec_count(i).eq.ivec_count(0)) then
1238             lentyp(i)=0
1239           else
1240             lentyp(i)=1
1241           endif
1242         enddo
1243         do ind_typ=lentyp(0),lentyp(nfgtasks-1)
1244         if (ind_typ.eq.0) then
1245           ichunk=ivec_count(0)
1246         else
1247           ichunk=ivec_count(1)
1248         endif
1249 !        do i=1,4
1250 !          blocklengths(i)=4
1251 !        enddo
1252 !        displs(1)=0
1253 !        do i=2,4
1254 !          displs(i)=displs(i-1)+blocklengths(i-1)*maxres
1255 !        enddo
1256 !        do i=1,4
1257 !          blocklengths(i)=blocklengths(i)*ichunk
1258 !        enddo
1259 !        write (iout,*) "blocklengths and displs"
1260 !        do i=1,4
1261 !          write (iout,*) i,blocklengths(i),displs(i)
1262 !        enddo
1263 !        call flush(iout)
1264 !        call MPI_Type_indexed(4,blocklengths(1),displs(1),
1265 !     &    MPI_DOUBLE_PRECISION,MPI_ROTAT1(ind_typ),IERROR)
1266 !        call MPI_Type_commit(MPI_ROTAT1(ind_typ),IERROR)
1267 !        write (iout,*) "MPI_ROTAT1",MPI_ROTAT1 
1268 !        do i=1,4
1269 !          blocklengths(i)=2
1270 !        enddo
1271 !        displs(1)=0
1272 !        do i=2,4
1273 !          displs(i)=displs(i-1)+blocklengths(i-1)*maxres
1274 !        enddo
1275 !        do i=1,4
1276 !          blocklengths(i)=blocklengths(i)*ichunk
1277 !        enddo
1278 !        write (iout,*) "blocklengths and displs"
1279 !        do i=1,4
1280 !          write (iout,*) i,blocklengths(i),displs(i)
1281 !        enddo
1282 !        call flush(iout)
1283 !        call MPI_Type_indexed(4,blocklengths(1),displs(1),
1284 !     &    MPI_DOUBLE_PRECISION,MPI_ROTAT2(ind_typ),IERROR)
1285 !        call MPI_Type_commit(MPI_ROTAT2(ind_typ),IERROR)
1286 !        write (iout,*) "MPI_ROTAT2",MPI_ROTAT2 
1287         do i=1,8
1288           blocklengths(i)=2
1289         enddo
1290         displs(1)=0
1291         do i=2,8
1292           displs(i)=displs(i-1)+blocklengths(i-1)*nres !maxres
1293         enddo
1294         do i=1,15
1295           blocklengths(i)=blocklengths(i)*ichunk
1296         enddo
1297         call MPI_Type_indexed(8,blocklengths,displs,&
1298           MPI_DOUBLE_PRECISION,MPI_PRECOMP11(ind_typ),IERROR)
1299         call MPI_Type_commit(MPI_PRECOMP11(ind_typ),IERROR)
1300         do i=1,8
1301           blocklengths(i)=4
1302         enddo
1303         displs(1)=0
1304         do i=2,8
1305           displs(i)=displs(i-1)+blocklengths(i-1)*nres !maxres
1306         enddo
1307         do i=1,15
1308           blocklengths(i)=blocklengths(i)*ichunk
1309         enddo
1310         call MPI_Type_indexed(8,blocklengths,displs,&
1311           MPI_DOUBLE_PRECISION,MPI_PRECOMP12(ind_typ),IERROR)
1312         call MPI_Type_commit(MPI_PRECOMP12(ind_typ),IERROR)
1313         do i=1,6
1314           blocklengths(i)=4
1315         enddo
1316         displs(1)=0
1317         do i=2,6
1318           displs(i)=displs(i-1)+blocklengths(i-1)*nres !maxres
1319         enddo
1320         do i=1,6
1321           blocklengths(i)=blocklengths(i)*ichunk
1322         enddo
1323         call MPI_Type_indexed(6,blocklengths,displs,&
1324           MPI_DOUBLE_PRECISION,MPI_PRECOMP22(ind_typ),IERROR)
1325         call MPI_Type_commit(MPI_PRECOMP22(ind_typ),IERROR)
1326         do i=1,2
1327           blocklengths(i)=8
1328         enddo
1329         displs(1)=0
1330         do i=2,2
1331           displs(i)=displs(i-1)+blocklengths(i-1)*nres !maxres
1332         enddo
1333         do i=1,2
1334           blocklengths(i)=blocklengths(i)*ichunk
1335         enddo
1336         call MPI_Type_indexed(2,blocklengths,displs,&
1337           MPI_DOUBLE_PRECISION,MPI_PRECOMP23(ind_typ),IERROR)
1338         call MPI_Type_commit(MPI_PRECOMP23(ind_typ),IERROR)
1339         do i=1,4
1340           blocklengths(i)=1
1341         enddo
1342         displs(1)=0
1343         do i=2,4
1344           displs(i)=displs(i-1)+blocklengths(i-1)*nres !maxres
1345         enddo
1346         do i=1,4
1347           blocklengths(i)=blocklengths(i)*ichunk
1348         enddo
1349         call MPI_Type_indexed(4,blocklengths,displs,&
1350           MPI_DOUBLE_PRECISION,MPI_ROTAT_OLD(ind_typ),IERROR)
1351         call MPI_Type_commit(MPI_ROTAT_OLD(ind_typ),IERROR)
1352         enddo
1353 #endif
1354       endif
1355       iint_start=ivec_start+1
1356       iint_end=ivec_end+1
1357       do i=0,nfgtasks-1
1358           iint_count(i)=ivec_count(i)
1359           iint_displ(i)=ivec_displ(i)
1360           ivec_displ(i)=ivec_displ(i)-1
1361           iset_displ(i)=iset_displ(i)-1
1362           ithet_displ(i)=ithet_displ(i)-1
1363           iphi_displ(i)=iphi_displ(i)-1
1364           iphi1_displ(i)=iphi1_displ(i)-1
1365           ibond_displ(i)=ibond_displ(i)-1
1366       enddo
1367       if (nfgtasks.gt.1 .and. fg_rank.eq.king &
1368           .and. (me.eq.0 .or. .not. out1file)) then
1369         write (iout,*) "IVEC_DISPL, IVEC_COUNT, ISET_START, ISET_COUNT"
1370         do i=0,nfgtasks-1
1371           write (iout,*) i,ivec_displ(i),ivec_count(i),iset_displ(i),&
1372             iset_count(i)
1373         enddo
1374         write (iout,*) "iphi_start",iphi_start," iphi_end",iphi_end,&
1375           " iphi1_start",iphi1_start," iphi1_end",iphi1_end
1376         write (iout,*)"IPHI_COUNT, IPHI_DISPL, IPHI1_COUNT, IPHI1_DISPL"
1377         do i=0,nfgtasks-1
1378           write (iout,*) i,iphi_count(i),iphi_displ(i),iphi1_count(i),&
1379             iphi1_displ(i)
1380         enddo
1381         write(iout,'(i10,a,i10,a,i10,a/a,i3,a)') n_sc_int_tot,' SC-SC ',&
1382           nele_int_tot,' electrostatic and ',nscp_int_tot,&
1383           ' SC-p interactions','were distributed among',nfgtasks,&
1384           ' fine-grain processors.'
1385       endif
1386 #else
1387       loc_start=2
1388       loc_end=nres_molec(1)-1
1389       ithet_start=3 
1390       ithet_end=nres_molec(1)
1391       ithet_nucl_start=3+nres_molec(1)
1392       ithet_nucl_end=nres_molec(1)+nres_molec(2)
1393       iturn3_start=nnt
1394       iturn3_end=nct_molec(1)-3
1395       iturn4_start=nnt
1396       iturn4_end=nct_molec(1)-4
1397       iphi_start=nnt+3
1398       iphi_end=nct_molec(1)
1399       iphi1_start=4
1400       iphi1_end=nres_molec(1)
1401       iphi_nucl_start=4+nres_molec(1)
1402       iphi_nucl_end=nres_molec(1)+nres_molec(2)
1403       idihconstr_start=1
1404       idihconstr_end=ndih_constr
1405       ithetaconstr_start=1
1406       ithetaconstr_end=ntheta_constr
1407       iphid_start=iphi_start
1408       iphid_end=iphi_end-1
1409       itau_start=4
1410       itau_end=nres_molec(1)
1411       ibond_start=2
1412       ibond_end=nres_molec(1)-1
1413       ibond_nucl_start=2+nres_molec(1)
1414       ibond_nucl_end=nres_molec(2)-1
1415       ibondp_start=nnt
1416       ibondp_end=nct_molec(1)-1
1417       ibondp_nucl_start=nnt_molec(2)
1418       ibondp_nucl_end=nct_molec(2)
1419       ivec_start=1
1420       ivec_end=nres_molec(1)-1
1421       iset_start=3
1422       iset_end=nres_molec(1)+1
1423       iint_start=2
1424       iint_end=nres_molec(1)-1
1425       ilip_start=1
1426       ilip_end=nres_molec(1)
1427       itube_start=1
1428       itube_end=nres_molec(1)
1429 #endif
1430 !el       common /przechowalnia/
1431 !      deallocate(iturn3_start_all)
1432 !      deallocate(iturn3_end_all)
1433 !      deallocate(iturn4_start_all)
1434 !      deallocate(iturn4_end_all)
1435 !      deallocate(iatel_s_all)
1436 !      deallocate(iatel_e_all)
1437 !      deallocate(ielstart_all)
1438 !      deallocate(ielend_all)
1439
1440 !      deallocate(ntask_cont_from_all)
1441 !      deallocate(ntask_cont_to_all)
1442 !      deallocate(itask_cont_from_all)
1443 !      deallocate(itask_cont_to_all)
1444 !el----------
1445       return
1446       end subroutine init_int_table
1447 #ifdef MPI
1448 !-----------------------------------------------------------------------------
1449       subroutine add_int(ii,jj,itask,ntask_cont_to,itask_cont_to,flag)
1450
1451 !el      implicit none
1452 !      include "DIMENSIONS"
1453 !      include "COMMON.INTERACT"
1454 !      include "COMMON.SETUP"
1455 !      include "COMMON.IOUNITS"
1456       integer :: ii,jj,ntask_cont_to
1457       integer,dimension(4) :: itask
1458       integer :: itask_cont_to(0:nfgtasks-1)    !(0:max_fg_procs-1)
1459       logical :: flag
1460 !el      integer,dimension(0:nfgtasks) :: iturn3_start_all,iturn3_end_all,iturn4_start_all,&
1461 !el       iturn4_end_all,iatel_s_all,iatel_e_all        !(0:max_fg_procs)
1462 !el      integer,dimension(nres,0:nfgtasks-1) :: ielstart_all,ielend_all        !(maxres,0:max_fg_procs-1)
1463 !el      common /przechowalnia/ iturn3_start_all,iturn3_end_all,iturn4_start_all,&
1464 !el       iturn4_end_all,iatel_s_all,iatel_e_all,ielstart_all,ielend_all
1465       integer :: iproc,isent,k,l
1466 ! Determines whether to send interaction ii,jj to other processors; a given
1467 ! interaction can be sent to at most 2 processors.
1468 ! Sets flag=.true. if interaction ii,jj needs to be sent to at least 
1469 ! one processor, otherwise flag is unchanged from the input value.
1470       isent=0
1471       itask(1)=fg_rank
1472       itask(2)=fg_rank
1473       itask(3)=fg_rank
1474       itask(4)=fg_rank
1475 !      write (iout,*) "ii",ii," jj",jj
1476 ! Loop over processors to check if anybody could need interaction ii,jj
1477       do iproc=0,fg_rank-1
1478 ! Check if the interaction matches any turn3 at iproc
1479         do k=iturn3_start_all(iproc),iturn3_end_all(iproc)
1480           l=k+2
1481           if (k.eq.ii-1 .and. l.eq.jj-1 .or. k.eq.ii-1 .and. l.eq.jj+1 &
1482          .or. k.eq.ii+1 .and. l.eq.jj+1 .or. k.eq.ii+1 .and. l.eq.jj-1) &
1483           then 
1484 !            write (iout,*) "turn3 to iproc",iproc," ij",ii,jj,"kl",k,l
1485 !            call flush(iout)
1486             flag=.true.
1487             if (iproc.ne.itask(1).and.iproc.ne.itask(2) &
1488               .and.iproc.ne.itask(3).and.iproc.ne.itask(4)) then
1489               isent=isent+1
1490               itask(isent)=iproc
1491               call add_task(iproc,ntask_cont_to,itask_cont_to)
1492             endif
1493           endif
1494         enddo
1495 ! Check if the interaction matches any turn4 at iproc
1496         do k=iturn4_start_all(iproc),iturn4_end_all(iproc)
1497           l=k+3
1498           if (k.eq.ii-1 .and. l.eq.jj-1 .or. k.eq.ii-1 .and. l.eq.jj+1 &
1499          .or. k.eq.ii+1 .and. l.eq.jj+1 .or. k.eq.ii+1 .and. l.eq.jj-1) &
1500           then 
1501 !            write (iout,*) "turn3 to iproc",iproc," ij",ii,jj," kl",k,l
1502 !            call flush(iout)
1503             flag=.true.
1504             if (iproc.ne.itask(1).and.iproc.ne.itask(2) &
1505               .and.iproc.ne.itask(3).and.iproc.ne.itask(4)) then
1506               isent=isent+1
1507               itask(isent)=iproc
1508               call add_task(iproc,ntask_cont_to,itask_cont_to)
1509             endif
1510           endif
1511         enddo
1512         if (iatel_s_all(iproc).gt.0 .and. iatel_e_all(iproc).gt.0 .and. &
1513         iatel_s_all(iproc).le.ii-1 .and. iatel_e_all(iproc).ge.ii-1)then
1514           if (ielstart_all(ii-1,iproc).le.jj-1.and. &
1515               ielend_all(ii-1,iproc).ge.jj-1) then
1516             flag=.true.
1517             if (iproc.ne.itask(1).and.iproc.ne.itask(2) &
1518               .and.iproc.ne.itask(3).and.iproc.ne.itask(4)) then
1519               isent=isent+1
1520               itask(isent)=iproc
1521               call add_task(iproc,ntask_cont_to,itask_cont_to)
1522             endif
1523           endif
1524           if (ielstart_all(ii-1,iproc).le.jj+1.and. &
1525               ielend_all(ii-1,iproc).ge.jj+1) then
1526             flag=.true.
1527             if (iproc.ne.itask(1).and.iproc.ne.itask(2) &
1528               .and.iproc.ne.itask(3).and.iproc.ne.itask(4)) then
1529               isent=isent+1
1530               itask(isent)=iproc
1531               call add_task(iproc,ntask_cont_to,itask_cont_to)
1532             endif
1533           endif
1534         endif
1535       enddo
1536       return
1537       end subroutine add_int
1538 !-----------------------------------------------------------------------------
1539       subroutine add_int_from(ii,jj,ntask_cont_from,itask_cont_from)
1540
1541 !el      use MPI_data
1542 !el      implicit none
1543 !      include "DIMENSIONS"
1544 !      include "COMMON.INTERACT"
1545 !      include "COMMON.SETUP"
1546 !      include "COMMON.IOUNITS"
1547       integer :: ii,jj,itask(2),ntask_cont_from,&
1548        itask_cont_from(0:nfgtasks-1)    !(0:max_fg_procs)
1549       logical :: flag
1550 !el      integer,dimension(0:nfgtasks) :: iturn3_start_all,iturn3_end_all,&
1551 !el       iturn4_start_all,iturn4_end_all,iatel_s_all,iatel_e_all       !(0:max_fg_procs)
1552 !el      integer,dimension(nres,0:nfgtasks-1) :: ielstart_all,ielend_all        !(maxres,0:max_fg_procs-1)
1553 !el      common /przechowalnia/ iturn3_start_all,iturn3_end_all,iturn4_start_all,&
1554 !el       iturn4_end_all,iatel_s_all,iatel_e_all,ielstart_all,ielend_all
1555       integer :: iproc,k,l
1556       do iproc=fg_rank+1,nfgtasks-1
1557         do k=iturn3_start_all(iproc),iturn3_end_all(iproc)
1558           l=k+2
1559           if (k.eq.ii+1 .and. l.eq.jj+1 .or. k.eq.ii+1.and.l.eq.jj-1 &
1560          .or. k.eq.ii-1 .and. l.eq.jj-1 .or. k.eq.ii-1 .and. l.eq.jj+1) &
1561           then
1562 !            write (iout,*)"turn3 from iproc",iproc," ij",ii,jj," kl",k,l
1563             call add_task(iproc,ntask_cont_from,itask_cont_from)
1564           endif
1565         enddo 
1566         do k=iturn4_start_all(iproc),iturn4_end_all(iproc)
1567           l=k+3
1568           if (k.eq.ii+1 .and. l.eq.jj+1 .or. k.eq.ii+1.and.l.eq.jj-1 &
1569          .or. k.eq.ii-1 .and. l.eq.jj-1 .or. k.eq.ii-1 .and. l.eq.jj+1) &
1570           then
1571 !            write (iout,*)"turn4 from iproc",iproc," ij",ii,jj," kl",k,l
1572             call add_task(iproc,ntask_cont_from,itask_cont_from)
1573           endif
1574         enddo 
1575         if (iatel_s_all(iproc).gt.0 .and. iatel_e_all(iproc).gt.0) then
1576           if (ii+1.ge.iatel_s_all(iproc).and.ii+1.le.iatel_e_all(iproc)) &
1577           then
1578             if (jj+1.ge.ielstart_all(ii+1,iproc).and. &
1579                 jj+1.le.ielend_all(ii+1,iproc)) then
1580               call add_task(iproc,ntask_cont_from,itask_cont_from)
1581             endif            
1582             if (jj-1.ge.ielstart_all(ii+1,iproc).and. &
1583                 jj-1.le.ielend_all(ii+1,iproc)) then
1584               call add_task(iproc,ntask_cont_from,itask_cont_from)
1585             endif
1586           endif
1587           if (ii-1.ge.iatel_s_all(iproc).and.ii-1.le.iatel_e_all(iproc)) &
1588           then
1589             if (jj-1.ge.ielstart_all(ii-1,iproc).and. &
1590                 jj-1.le.ielend_all(ii-1,iproc)) then
1591               call add_task(iproc,ntask_cont_from,itask_cont_from)
1592             endif
1593             if (jj+1.ge.ielstart_all(ii-1,iproc).and. &
1594                 jj+1.le.ielend_all(ii-1,iproc)) then
1595                call add_task(iproc,ntask_cont_from,itask_cont_from)
1596             endif
1597           endif
1598         endif
1599       enddo
1600       return
1601       end subroutine add_int_from
1602 !-----------------------------------------------------------------------------
1603       subroutine add_task(iproc,ntask_cont,itask_cont)
1604
1605 !el      use MPI_data
1606 !el      implicit none
1607 !      include "DIMENSIONS"
1608       integer :: iproc,ntask_cont,itask_cont(0:nfgtasks-1)      !(0:max_fg_procs-1)
1609       integer :: ii
1610       do ii=1,ntask_cont
1611         if (itask_cont(ii).eq.iproc) return
1612       enddo
1613       ntask_cont=ntask_cont+1
1614       itask_cont(ntask_cont)=iproc
1615       return
1616       end subroutine add_task
1617 #endif
1618 !-----------------------------------------------------------------------------
1619 #if defined MPI || defined WHAM_RUN
1620       subroutine int_partition(int_index,lower_index,upper_index,atom,&
1621        at_start,at_end,first_atom,last_atom,int_gr,jat_start,jat_end,*)
1622
1623 !      implicit real*8 (a-h,o-z)
1624 !      include 'DIMENSIONS'
1625 !      include 'COMMON.IOUNITS'
1626       integer :: int_index,lower_index,upper_index,atom,at_start,at_end,&
1627        first_atom,last_atom,int_gr,jat_start,jat_end,int_index_old
1628       logical :: lprn
1629       lprn=.false.
1630       if (lprn) write (iout,*) 'int_index=',int_index
1631       int_index_old=int_index
1632       int_index=int_index+last_atom-first_atom+1
1633       if (lprn) &
1634          write (iout,*) 'int_index=',int_index,&
1635                      ' int_index_old',int_index_old,&
1636                      ' lower_index=',lower_index,&
1637                      ' upper_index=',upper_index,&
1638                      ' atom=',atom,' first_atom=',first_atom,&
1639                      ' last_atom=',last_atom
1640       if (int_index.ge.lower_index) then
1641         int_gr=int_gr+1
1642         if (at_start.eq.0) then
1643           at_start=atom
1644           jat_start=first_atom-1+lower_index-int_index_old
1645         else
1646           jat_start=first_atom
1647         endif
1648         if (lprn) write (iout,*) 'jat_start',jat_start
1649         if (int_index.ge.upper_index) then
1650           at_end=atom
1651           jat_end=first_atom-1+upper_index-int_index_old
1652           return 1
1653         else
1654           jat_end=last_atom
1655         endif
1656         if (lprn) write (iout,*) 'jat_end',jat_end
1657       endif
1658       return
1659       end subroutine int_partition
1660 #endif
1661 !-----------------------------------------------------------------------------
1662 #ifndef CLUSTER
1663       subroutine hpb_partition
1664
1665 !      implicit real*8 (a-h,o-z)
1666 !      include 'DIMENSIONS'
1667 #ifdef MPI
1668       include 'mpif.h'
1669 #endif
1670 !      include 'COMMON.SBRIDGE'
1671 !      include 'COMMON.IOUNITS'
1672 !      include 'COMMON.SETUP'
1673 #ifdef MPI
1674       call int_bounds(nhpb,link_start,link_end)
1675       write (iout,*) 'Processor',fg_rank,' CG group',kolor,&
1676         ' absolute rank',MyRank,&
1677         ' nhpb',nhpb,' link_start=',link_start,&
1678         ' link_end',link_end
1679 #else
1680       link_start=1
1681       link_end=nhpb
1682 #endif
1683       return
1684       end subroutine hpb_partition
1685 #endif
1686 !-----------------------------------------------------------------------------
1687 ! misc.f in module io_base
1688 !-----------------------------------------------------------------------------
1689 !-----------------------------------------------------------------------------
1690 ! parmread.F
1691 !-----------------------------------------------------------------------------
1692       subroutine getenv_loc(var, val)
1693
1694       character(*) :: var, val
1695
1696 #ifdef WINIFL
1697       character(len=2000) :: line
1698 !el      external ilen
1699
1700       open (196,file='env',status='old',readonly,shared)
1701       iread=0
1702 !      write(*,*)'looking for ',var
1703 10    read(196,*,err=11,end=11)line
1704       iread=index(line,var)
1705 !      write(*,*)iread,' ',var,' ',line
1706       if (iread.eq.0) go to 10 
1707 !      write(*,*)'---> ',line
1708 11    continue
1709       if(iread.eq.0) then
1710 !       write(*,*)'CHUJ'
1711        val=''
1712       else
1713        iread=iread+ilen(var)+1
1714        read (line(iread:),*,err=12,end=12) val
1715 !       write(*,*)'OK: ',var,' = ',val
1716       endif
1717       close(196)
1718       return
1719 12    val=''
1720       close(196)
1721 #elif (defined CRAY)
1722       integer :: lennam,lenval,ierror
1723 !
1724 !        getenv using a POSIX call, useful on the T3D
1725 !        Sept 1996, comment out error check on advice of H. Pritchard
1726 !
1727       lennam = len(var)
1728       if(lennam.le.0) stop '--error calling getenv--'
1729       call pxfgetenv(var,lennam,val,lenval,ierror)
1730 !-HP- if(ierror.ne.0) stop '--error returned by pxfgetenv--'
1731 #else
1732       call getenv(var,val)
1733 #endif
1734
1735       return
1736       end subroutine getenv_loc
1737 !-----------------------------------------------------------------------------
1738 ! readrtns_CSA.F
1739 !-----------------------------------------------------------------------------
1740       subroutine setup_var
1741
1742       integer :: i
1743 !      implicit real*8 (a-h,o-z)
1744 !      include 'DIMENSIONS'
1745 !      include 'COMMON.IOUNITS'
1746 !      include 'COMMON.GEO'
1747 !      include 'COMMON.VAR'
1748 !      include 'COMMON.INTERACT'
1749 !      include 'COMMON.LOCAL'
1750 !      include 'COMMON.NAMES'
1751 !      include 'COMMON.CHAIN'
1752 !      include 'COMMON.FFIELD'
1753 !      include 'COMMON.SBRIDGE'
1754 !      include 'COMMON.HEADER'
1755 !      include 'COMMON.CONTROL'
1756 !      include 'COMMON.DBASE'
1757 !      include 'COMMON.THREAD'
1758 !      include 'COMMON.TIME1'
1759 ! Set up variable list.
1760       ntheta=nres-2
1761       nphi=nres-3
1762       nvar=ntheta+nphi
1763       nside=0
1764       do i=2,nres-1
1765 #ifdef WHAM_RUN
1766         if (itype(i,1).ne.10) then
1767 #else
1768         if (itype(i,1).ne.10 .and. itype(i,1).ne.ntyp1) then
1769 #endif
1770           nside=nside+1
1771           ialph(i,1)=nvar+nside
1772           ialph(nside,2)=i
1773         endif
1774       enddo
1775       if (indphi.gt.0) then
1776         nvar=nphi
1777       else if (indback.gt.0) then
1778         nvar=nphi+ntheta
1779       else
1780         nvar=nvar+2*nside
1781       endif
1782 !d    write (iout,'(3i4)') (i,ialph(i,1),ialph(i,2),i=2,nres-1)
1783       return
1784       end subroutine setup_var
1785 !-----------------------------------------------------------------------------
1786 ! rescode.f
1787 !-----------------------------------------------------------------------------
1788       integer function rescode(iseq,nam,itype,molecule)
1789
1790       use io_base, only: ucase
1791 !      implicit real*8 (a-h,o-z)
1792 !      include 'DIMENSIONS'
1793 !      include 'COMMON.NAMES'
1794 !      include 'COMMON.IOUNITS'
1795       character(len=3) :: nam   !,ucase
1796       integer :: iseq,itype,i
1797       integer :: molecule
1798       print *,molecule,nam
1799       if (molecule.eq.1) then 
1800       if (itype.eq.0) then
1801
1802       do i=-ntyp1_molec(molecule),ntyp1_molec(molecule)
1803         if (ucase(nam).eq.restyp(i,molecule)) then
1804           rescode=i
1805           return
1806         endif
1807       enddo
1808
1809       else
1810
1811       do i=-ntyp1_molec(molecule),ntyp1_molec(molecule)
1812         if (nam(1:1).eq.onelet(i)) then
1813           rescode=i
1814           return  
1815         endif  
1816       enddo
1817
1818       endif
1819       else if (molecule.eq.2) then
1820       do i=1,ntyp1_molec(molecule)
1821          print *,nam(1:1),restyp(i,molecule)(1:1) 
1822         if (nam(2:2).eq.restyp(i,molecule)(1:1)) then
1823           rescode=i
1824           return
1825         endif
1826       enddo
1827       else if (molecule.eq.3) then
1828        write(iout,*) "SUGAR not yet implemented"
1829        stop
1830       else if (molecule.eq.4) then
1831        write(iout,*) "Explicit LIPID not yet implemented"
1832        stop
1833       else if (molecule.eq.5) then
1834       do i=1,ntyp1_molec(molecule)
1835         print *,i,restyp(i,molecule)(1:2)
1836         if (ucase(nam(1:2)).eq.restyp(i,molecule)(1:2)) then
1837           rescode=i
1838           return
1839         endif
1840       enddo
1841       else   
1842        write(iout,*) "molecule not defined"
1843       endif
1844       write (iout,10) iseq,nam
1845       stop
1846    10 format ('**** Error - residue',i4,' has an unresolved name ',a3)
1847       end function rescode
1848       integer function sugarcode(sugar,ires)
1849       character sugar
1850       integer ires
1851       if (sugar.eq.'D') then
1852         sugarcode=1
1853       else if (sugar.eq.' ') then
1854         sugarcode=2
1855       else
1856         write (iout,*) 'UNKNOWN sugar type for residue',ires,' ',sugar
1857         stop
1858       endif
1859       return
1860       end function sugarcode
1861
1862 !-----------------------------------------------------------------------------
1863 ! timing.F
1864 !-----------------------------------------------------------------------------
1865 ! $Date: 1994/10/05 16:41:52 $
1866 ! $Revision: 2.2 $
1867 !
1868       subroutine set_timers
1869 !
1870 !el      implicit none
1871 !el      real(kind=8) :: tcpu
1872 !      include 'COMMON.TIME1'
1873 !#ifdef MP
1874 #ifdef MPI
1875       include 'mpif.h'
1876 #endif
1877 ! Diminish the assigned time limit a little so that there is some time to
1878 ! end a batch job
1879 !     timlim=batime-150.0
1880 ! Calculate the initial time, if it is not zero (e.g. for the SUN).
1881       stime=tcpu()
1882 #if .not. defined(WHAM_RUN) && .not. defined(CLUSTER)
1883 #ifdef MPI
1884       walltime=MPI_WTIME()
1885       time_reduce=0.0d0
1886       time_allreduce=0.0d0
1887       time_bcast=0.0d0
1888       time_gather=0.0d0
1889       time_sendrecv=0.0d0
1890       time_scatter=0.0d0
1891       time_scatter_fmat=0.0d0
1892       time_scatter_ginv=0.0d0
1893       time_scatter_fmatmult=0.0d0
1894       time_scatter_ginvmult=0.0d0
1895       time_barrier_e=0.0d0
1896       time_barrier_g=0.0d0
1897       time_enecalc=0.0d0
1898       time_sumene=0.0d0
1899       time_lagrangian=0.0d0
1900       time_sumgradient=0.0d0
1901       time_intcartderiv=0.0d0
1902       time_inttocart=0.0d0
1903       time_ginvmult=0.0d0
1904       time_fricmatmult=0.0d0
1905       time_cartgrad=0.0d0
1906       time_bcastc=0.0d0
1907       time_bcast7=0.0d0
1908       time_bcastw=0.0d0
1909       time_intfcart=0.0d0
1910       time_vec=0.0d0
1911       time_mat=0.0d0
1912       time_fric=0.0d0
1913       time_stoch=0.0d0
1914       time_fricmatmult=0.0d0
1915       time_fsample=0.0d0
1916 #endif
1917 #endif
1918 !d    print *,' in SET_TIMERS stime=',stime
1919       return
1920       end subroutine set_timers
1921 !-----------------------------------------------------------------------------
1922 #ifndef CLUSTER
1923       logical function stopx(nf)
1924 ! This function returns .true. if one of the following reasons to exit SUMSL
1925 ! occurs. The "reason" code is stored in WHATSUP passed thru a COMMON block:
1926 !
1927 !... WHATSUP = 0 - go on, no reason to stop. Stopx will return .false.
1928 !...           1 - Time up in current node;
1929 !...           2 - STOP signal was received from another node because the
1930 !...               node's task was accomplished (parallel only);
1931 !...          -1 - STOP signal was received from another node because of error;
1932 !...          -2 - STOP signal was received from another node, because 
1933 !...               the node's time was up.
1934 !      implicit real*8 (a-h,o-z)
1935 !      include 'DIMENSIONS'
1936 !el#ifdef WHAM_RUN
1937 !el      use control_data, only:WhatsUp
1938 !el#endif
1939 #ifdef MP
1940 !el      use MPI_data   !include 'COMMON.INFO'
1941       include 'mpif.h'
1942 #endif
1943       integer :: nf
1944 !el      logical :: ovrtim
1945
1946 !      include 'COMMON.IOUNITS'
1947 !      include 'COMMON.TIME1'
1948       integer :: Kwita
1949
1950 !d    print *,'Processor',MyID,' NF=',nf
1951 !d      write (iout,*) "stopx: ",nf
1952 #ifndef WHAM_RUN
1953 #ifndef MPI
1954       if (ovrtim()) then
1955 ! Finish if time is up.
1956          stopx = .true.
1957          WhatsUp=1
1958 #ifdef MPL
1959       else if (mod(nf,100).eq.0) then
1960 ! Other processors might have finished. Check this every 100th function 
1961 ! evaluation.
1962 ! Master checks if any other processor has sent accepted conformation(s) to it. 
1963          if (MyID.ne.MasterID) call receive_mcm_info
1964          if (MyID.eq.MasterID) call receive_conf
1965 !d       print *,'Processor ',MyID,' is checking STOP: nf=',nf
1966          call recv_stop_sig(Kwita)
1967          if (Kwita.eq.-1) then
1968            write (iout,'(a,i4,a,i5)') 'Processor',&
1969            MyID,' has received STOP signal in STOPX; NF=',nf
1970            write (*,'(a,i4,a,i5)') 'Processor',&
1971            MyID,' has received STOP signal in STOPX; NF=',nf
1972            stopx=.true.
1973            WhatsUp=2
1974          elseif (Kwita.eq.-2) then
1975            write (iout,*) &
1976           'Processor',MyID,' received TIMEUP-STOP signal in SUMSL.'
1977            write (*,*) &
1978           'Processor',MyID,' received TIMEUP-STOP signal in SUMSL.'
1979            WhatsUp=-2
1980            stopx=.true.  
1981          else if (Kwita.eq.-3) then
1982            write (iout,*) &
1983           'Processor',MyID,' received ERROR-STOP signal in SUMSL.'
1984            write (*,*) &
1985           'Processor',MyID,' received ERROR-STOP signal in SUMSL.'
1986            WhatsUp=-1
1987            stopx=.true.
1988          else
1989            stopx=.false.
1990            WhatsUp=0
1991          endif
1992 #endif
1993       else
1994          stopx = .false.
1995          WhatsUp=0
1996       endif
1997 #else
1998       stopx=.false.
1999 !d      write (iout,*) "stopx set at .false."
2000 #endif
2001
2002 #ifdef OSF
2003 ! Check for FOUND_NAN flag
2004       if (FOUND_NAN) then
2005         write(iout,*)"   ***   stopx : Found a NaN"
2006         stopx=.true.
2007       endif
2008 #endif
2009 #else
2010       if (ovrtim()) then
2011 ! Finish if time is up.
2012          stopx = .true.
2013          WhatsUp=1
2014       else if (cutoffviol) then
2015         stopx = .true.
2016         WhatsUp=2
2017       else
2018         stopx=.false.
2019       endif
2020 #endif
2021       return
2022       end function stopx
2023 !-----------------------------------------------------------------------------
2024 #else
2025       logical function stopx(nf)
2026 !
2027 !     ..................................................................
2028 !
2029 !     *****PURPOSE...
2030 !     THIS FUNCTION MAY SERVE AS THE STOPX (ASYNCHRONOUS INTERRUPTION)
2031 !     FUNCTION FOR THE NL2SOL (NONLINEAR LEAST-SQUARES) PACKAGE AT
2032 !     THOSE INSTALLATIONS WHICH DO NOT WISH TO IMPLEMENT A
2033 !     DYNAMIC STOPX.
2034 !
2035 !     *****ALGORITHM NOTES...
2036 !     AT INSTALLATIONS WHERE THE NL2SOL SYSTEM IS USED
2037 !     INTERACTIVELY, THIS DUMMY STOPX SHOULD BE REPLACED BY A
2038 !     FUNCTION THAT RETURNS .TRUE. IF AND ONLY IF THE INTERRUPT
2039 !     (BREAK) KEY HAS BEEN PRESSED SINCE THE LAST CALL ON STOPX.
2040 !
2041 !     $$$ MODIFIED FOR USE AS  THE TIMER ROUTINE.
2042 !     $$$                              WHEN THE TIME LIMIT HAS BEEN
2043 !     $$$ REACHED     STOPX IS SET TO .TRUE  AND INITIATES (IN ITSUM)
2044 !     $$$ AND ORDERLY EXIT OUT OF SUMSL.  IF ARRAYS IV AND V ARE
2045 !     $$$ SAVED, THE SUMSL ROUTINES CAN BE RESTARTED AT THE SAME
2046 !     $$$ POINT AT WHICH THEY WERE INTERRUPTED.
2047 !
2048 !     ..................................................................
2049 !
2050 !      include 'DIMENSIONS'
2051       integer :: nf
2052 !      logical ovrtim
2053 !      include 'COMMON.IOUNITS'
2054 !      include 'COMMON.TIME1'
2055 #ifdef MPL
2056 !     include 'COMMON.INFO'
2057       integer :: Kwita
2058
2059 !d    print *,'Processor',MyID,' NF=',nf
2060 #endif
2061       if (ovrtim()) then
2062 ! Finish if time is up.
2063          stopx = .true.
2064 #ifdef MPL
2065       else if (mod(nf,100).eq.0) then
2066 ! Other processors might have finished. Check this every 100th function 
2067 ! evaluation.
2068 !d       print *,'Processor ',MyID,' is checking STOP: nf=',nf
2069          call recv_stop_sig(Kwita)
2070          if (Kwita.eq.-1) then
2071            write (iout,'(a,i4,a,i5)') 'Processor',&
2072            MyID,' has received STOP signal in STOPX; NF=',nf
2073            write (*,'(a,i4,a,i5)') 'Processor',&
2074            MyID,' has received STOP signal in STOPX; NF=',nf
2075            stopx=.true.
2076          else
2077            stopx=.false.
2078          endif
2079 #endif
2080       else
2081          stopx = .false.
2082       endif
2083       return
2084       end function stopx
2085 #endif
2086 !-----------------------------------------------------------------------------
2087       logical function ovrtim()
2088
2089 !      include 'DIMENSIONS'
2090 !      include 'COMMON.IOUNITS'
2091 !      include 'COMMON.TIME1'
2092 !el      real(kind=8) :: tcpu
2093       real(kind=8) :: curtim
2094 #ifdef MPI
2095       include "mpif.h"
2096       curtim = MPI_Wtime()-walltime
2097 #else
2098       curtim= tcpu()
2099 #endif
2100 !  curtim is the current time in seconds.
2101 !      write (iout,*) "curtim",curtim," timlim",timlim," safety",safety
2102 #ifndef WHAM_RUN
2103       if (curtim .ge. timlim - safety) then
2104         write (iout,'(a,f10.2,a,f10.2,a,f10.2,a)') &
2105         "***************** Elapsed time (",curtim,&
2106         " s) is within the safety limit (",safety,&
2107         " s) of the allocated time (",timlim," s). Terminating."
2108         ovrtim=.true.
2109       else
2110         ovrtim=.false.
2111       endif
2112 #else
2113       ovrtim=.false.
2114 #endif
2115 !elwrite (iout,*) "ovrtim",ovrtim
2116       return
2117       end function ovrtim
2118 !-----------------------------------------------------------------------------
2119       real(kind=8) function tcpu()
2120
2121 !      include 'COMMON.TIME1'
2122       real(kind=8) :: seconds
2123 #ifdef ES9000
2124 !***************************
2125 ! Next definition for EAGLE (ibm-es9000)
2126       real(kind=8) :: micseconds
2127       integer :: rcode
2128       tcpu=cputime(micseconds,rcode)
2129       tcpu=(micseconds/1.0E6) - stime
2130 !***************************
2131 #endif
2132 #ifdef SUN
2133 !***************************
2134 ! Next definitions for sun
2135       REAL(kind=8) ::  ECPU,ETIME,ETCPU
2136       real(kind=8),dimension(2) :: tarray
2137       tcpu=etime(tarray)
2138       tcpu=tarray(1)
2139 !***************************
2140 #endif
2141 #ifdef KSR
2142 !***************************
2143 ! Next definitions for ksr
2144 ! this function uses the ksr timer ALL_SECONDS from the PMON library to
2145 ! return the elapsed time in seconds
2146       tcpu= all_seconds() - stime
2147 !***************************
2148 #endif
2149 #ifdef SGI
2150 !***************************
2151 ! Next definitions for sgi
2152       real(kind=4) :: timar(2), etime
2153       seconds = etime(timar)
2154 !d    print *,'seconds=',seconds,' stime=',stime
2155 !      usrsec = timar(1)
2156 !      syssec = timar(2)
2157       tcpu=seconds - stime
2158 !***************************
2159 #endif
2160
2161 #ifdef LINUX
2162 !***************************
2163 ! Next definitions for sgi
2164       real(kind=4) :: timar(2), etime
2165       seconds = etime(timar)
2166 !d    print *,'seconds=',seconds,' stime=',stime
2167 !      usrsec = timar(1)
2168 !      syssec = timar(2)
2169       tcpu=seconds - stime
2170 !***************************
2171 #endif
2172
2173
2174 #ifdef CRAY
2175 !***************************
2176 ! Next definitions for Cray
2177 !     call date(curdat)
2178 !     curdat=curdat(1:9)
2179 !     call clock(curtim)
2180 !     curtim=curtim(1:8)
2181       cpusec = second()
2182       tcpu=cpusec - stime
2183 !***************************
2184 #endif
2185 #ifdef AIX
2186 !***************************
2187 ! Next definitions for RS6000
2188        integer(kind=4) :: i1,mclock
2189        i1 = mclock()
2190        tcpu = (i1+0.0D0)/100.0D0
2191 #endif
2192 #ifdef WINPGI
2193 !***************************
2194 ! next definitions for windows NT Digital fortran
2195        real(kind=4) :: time_real
2196        call cpu_time(time_real)
2197        tcpu = time_real
2198 #endif
2199 #ifdef WINIFL
2200 !***************************
2201 ! next definitions for windows NT Digital fortran
2202        real(kind=4) :: time_real
2203        call cpu_time(time_real)
2204        tcpu = time_real
2205 #endif
2206       tcpu = 0d0 !el
2207       return
2208       end function tcpu
2209 !-----------------------------------------------------------------------------
2210 #ifndef CLUSTER
2211       subroutine dajczas(rntime,hrtime,mintime,sectime)
2212
2213 !      include 'COMMON.IOUNITS'
2214       integer :: ihr,imn,isc
2215       real(kind=8) :: rntime,hrtime,mintime,sectime 
2216       hrtime=rntime/3600.0D0 
2217       hrtime=aint(hrtime)
2218       mintime=aint((rntime-3600.0D0*hrtime)/60.0D0)
2219       sectime=aint((rntime-3600.0D0*hrtime-60.0D0*mintime)+0.5D0)
2220       if (sectime.eq.60.0D0) then
2221         sectime=0.0D0
2222         mintime=mintime+1.0D0
2223       endif
2224       ihr=hrtime
2225       imn=mintime
2226       isc=sectime
2227       write (iout,328) ihr,imn,isc
2228   328 FORMAT(//'***** Computation time: ',I4  ,' hours ',I2  ,&
2229                ' minutes ', I2  ,' seconds *****')       
2230       return
2231       end subroutine dajczas
2232 !-----------------------------------------------------------------------------
2233       subroutine print_detailed_timing
2234
2235 !el      use MPI_data
2236 !      implicit real*8 (a-h,o-z)
2237 !      include 'DIMENSIONS'
2238 #ifdef MPI
2239       include 'mpif.h'
2240 #endif
2241 !      include 'COMMON.IOUNITS'
2242 !      include 'COMMON.TIME1'
2243 !      include 'COMMON.SETUP'
2244       real(kind=8) :: time1,time_barrier
2245       time_barrier = 0.0d0
2246 #ifdef MPI !el
2247       time1=MPI_WTIME()
2248 #endif !el
2249          write (iout,'(80(1h=)/a/(80(1h=)))') &
2250           "Details of FG communication time"
2251          write (*,'(7(a40,1pe15.5/),40(1h-)/a40,1pe15.5/80(1h=))') &
2252           "BROADCAST:",time_bcast,"REDUCE:",time_reduce,&
2253           "GATHER:",time_gather,&
2254           "SCATTER:",time_scatter,"SENDRECV:",time_sendrecv,&
2255           "BARRIER ene",time_barrier_e,&
2256           "BARRIER grad",time_barrier_g,&
2257           "TOTAL:",&
2258           time_bcast+time_reduce+time_gather+time_scatter+time_sendrecv
2259          write (*,*) fg_rank,myrank,&
2260            ': Total wall clock time',time1-walltime,' sec'
2261          write (*,*) "Processor",fg_rank,myrank,&
2262            ": BROADCAST time",time_bcast," REDUCE time",&
2263             time_reduce," GATHER time",time_gather," SCATTER time",&
2264             time_scatter,&
2265            " SCATTER fmatmult",time_scatter_fmatmult,&
2266            " SCATTER ginvmult",time_scatter_ginvmult,&
2267            " SCATTER fmat",time_scatter_fmat,&
2268            " SCATTER ginv",time_scatter_ginv,&
2269             " SENDRECV",time_sendrecv,&
2270             " BARRIER ene",time_barrier_e,&
2271             " BARRIER GRAD",time_barrier_g,&
2272             " BCAST7",time_bcast7," BCASTC",time_bcastc,&
2273             " BCASTW",time_bcastw," ALLREDUCE",time_allreduce,&
2274             " TOTAL",&
2275             time_bcast+time_reduce+time_gather+time_scatter+ &
2276             time_sendrecv+time_barrier+time_bcastc
2277 !el#endif
2278          write (*,*) "Processor",fg_rank,myrank," enecalc",time_enecalc
2279          write (*,*) "Processor",fg_rank,myrank," sumene",time_sumene
2280          write (*,*) "Processor",fg_rank,myrank," intfromcart",&
2281            time_intfcart
2282          write (*,*) "Processor",fg_rank,myrank," vecandderiv",&
2283            time_vec
2284          write (*,*) "Processor",fg_rank,myrank," setmatrices",&
2285            time_mat
2286          write (*,*) "Processor",fg_rank,myrank," ginvmult",&
2287            time_ginvmult
2288          write (*,*) "Processor",fg_rank,myrank," fricmatmult",&
2289            time_fricmatmult
2290          write (*,*) "Processor",fg_rank,myrank," inttocart",&
2291            time_inttocart
2292          write (*,*) "Processor",fg_rank,myrank," sumgradient",&
2293            time_sumgradient
2294          write (*,*) "Processor",fg_rank,myrank," intcartderiv",&
2295            time_intcartderiv
2296          if (fg_rank.eq.0) then
2297            write (*,*) "Processor",fg_rank,myrank," lagrangian",&
2298              time_lagrangian
2299            write (*,*) "Processor",fg_rank,myrank," cartgrad",&
2300              time_cartgrad
2301          endif
2302       return
2303       end subroutine print_detailed_timing
2304 #endif
2305 !-----------------------------------------------------------------------------
2306 !-----------------------------------------------------------------------------
2307       end module control