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