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