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