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