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