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