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