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