working energy for shield and lipid wrong gradient
[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 #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       call int_bounds(nres,ilip_start,ilip_end)
800       ilip_start=ilip_start
801       ilip_end=ilip_end
802       if (ndih_constr.eq.0) then
803         idihconstr_start=1
804         idihconstr_end=0
805       else
806         call int_bounds(ndih_constr,idihconstr_start,idihconstr_end)
807       endif
808 !      nsumgrad=(nres-nnt)*(nres-nnt+1)/2
809 !      nlen=nres-nnt+1
810       nsumgrad=(nres-nnt)*(nres-nnt+1)/2
811       nlen=nres-nnt+1
812       call int_bounds(nsumgrad,ngrad_start,ngrad_end)
813       igrad_start=((2*nlen+1) &
814          -sqrt(float((2*nlen-1)**2-8*(ngrad_start-1))))/2
815       igrad_end=((2*nlen+1) &
816          -sqrt(float((2*nlen-1)**2-8*(ngrad_end-1))))/2
817 !el      allocate(jgrad_start(igrad_start:igrad_end))
818 !el      allocate(jgrad_end(igrad_start:igrad_end)) !(maxres)
819       jgrad_start(igrad_start)= &
820          ngrad_start-(2*nlen-igrad_start)*(igrad_start-1)/2 &
821          +igrad_start
822       jgrad_end(igrad_start)=nres
823       if (igrad_end.gt.igrad_start) jgrad_start(igrad_end)=igrad_end+1
824       jgrad_end(igrad_end)=ngrad_end-(2*nlen-igrad_end)*(igrad_end-1)/2 &
825           +igrad_end
826       do i=igrad_start+1,igrad_end-1
827         jgrad_start(i)=i+1
828         jgrad_end(i)=nres
829       enddo
830       if (lprint) then 
831         write (*,*) 'Processor:',fg_rank,' CG group',kolor,&
832        ' absolute rank',myrank,&
833        ' loc_start',loc_start,' loc_end',loc_end,&
834        ' ithet_start',ithet_start,' ithet_end',ithet_end,&
835        ' iphi_start',iphi_start,' iphi_end',iphi_end,&
836        ' iphid_start',iphid_start,' iphid_end',iphid_end,&
837        ' ibond_start',ibond_start,' ibond_end',ibond_end,&
838        ' ibondp_start',ibondp_start,' ibondp_end',ibondp_end,&
839        ' iturn3_start',iturn3_start,' iturn3_end',iturn3_end,&
840        ' iturn4_start',iturn4_start,' iturn4_end',iturn4_end,&
841        ' ivec_start',ivec_start,' ivec_end',ivec_end,&
842        ' iset_start',iset_start,' iset_end',iset_end,&
843        ' idihconstr_start',idihconstr_start,' idihconstr_end',&
844          idihconstr_end
845        write (*,*) 'Processor:',fg_rank,myrank,' igrad_start',&
846          igrad_start,' igrad_end',igrad_end,' ngrad_start',ngrad_start,&
847          ' ngrad_end',ngrad_end
848        do i=igrad_start,igrad_end
849          write(*,*) 'Processor:',fg_rank,myrank,i,&
850           jgrad_start(i),jgrad_end(i)
851        enddo
852       endif
853       if (nfgtasks.gt.1) then
854         call MPI_Allgather(ivec_start,1,MPI_INTEGER,ivec_displ(0),1,&
855           MPI_INTEGER,FG_COMM1,IERROR)
856         iaux=ivec_end-ivec_start+1
857         call MPI_Allgather(iaux,1,MPI_INTEGER,ivec_count(0),1,&
858           MPI_INTEGER,FG_COMM1,IERROR)
859         call MPI_Allgather(iset_start-2,1,MPI_INTEGER,iset_displ(0),1,&
860           MPI_INTEGER,FG_COMM,IERROR)
861         iaux=iset_end-iset_start+1
862         call MPI_Allgather(iaux,1,MPI_INTEGER,iset_count(0),1,&
863           MPI_INTEGER,FG_COMM,IERROR)
864         call MPI_Allgather(ibond_start,1,MPI_INTEGER,ibond_displ(0),1,&
865           MPI_INTEGER,FG_COMM,IERROR)
866         iaux=ibond_end-ibond_start+1
867         call MPI_Allgather(iaux,1,MPI_INTEGER,ibond_count(0),1,&
868           MPI_INTEGER,FG_COMM,IERROR)
869         call MPI_Allgather(ithet_start,1,MPI_INTEGER,ithet_displ(0),1,&
870           MPI_INTEGER,FG_COMM,IERROR)
871         iaux=ithet_end-ithet_start+1
872         call MPI_Allgather(iaux,1,MPI_INTEGER,ithet_count(0),1,&
873           MPI_INTEGER,FG_COMM,IERROR)
874         call MPI_Allgather(iphi_start,1,MPI_INTEGER,iphi_displ(0),1,&
875           MPI_INTEGER,FG_COMM,IERROR)
876         iaux=iphi_end-iphi_start+1
877         call MPI_Allgather(iaux,1,MPI_INTEGER,iphi_count(0),1,&
878           MPI_INTEGER,FG_COMM,IERROR)
879         call MPI_Allgather(iphi1_start,1,MPI_INTEGER,iphi1_displ(0),1,&
880           MPI_INTEGER,FG_COMM,IERROR)
881         iaux=iphi1_end-iphi1_start+1
882         call MPI_Allgather(iaux,1,MPI_INTEGER,iphi1_count(0),1,&
883           MPI_INTEGER,FG_COMM,IERROR)
884         do i=0,nfgtasks-1
885           do j=1,nres
886             ielstart_all(j,i)=0
887             ielend_all(j,i)=0
888           enddo
889         enddo
890         call MPI_Allgather(iturn3_start,1,MPI_INTEGER,&
891           iturn3_start_all(0),1,MPI_INTEGER,FG_COMM,IERROR)
892         call MPI_Allgather(iturn4_start,1,MPI_INTEGER,&
893           iturn4_start_all(0),1,MPI_INTEGER,FG_COMM,IERROR)
894         call MPI_Allgather(iturn3_end,1,MPI_INTEGER,&
895           iturn3_end_all(0),1,MPI_INTEGER,FG_COMM,IERROR)
896         call MPI_Allgather(iturn4_end,1,MPI_INTEGER,&
897           iturn4_end_all(0),1,MPI_INTEGER,FG_COMM,IERROR)
898         call MPI_Allgather(iatel_s,1,MPI_INTEGER,&
899           iatel_s_all(0),1,MPI_INTEGER,FG_COMM,IERROR)
900         call MPI_Allgather(iatel_e,1,MPI_INTEGER,&
901           iatel_e_all(0),1,MPI_INTEGER,FG_COMM,IERROR)
902         call MPI_Allgather(ielstart(1),nres,MPI_INTEGER,&
903           ielstart_all(1,0),nres,MPI_INTEGER,FG_COMM,IERROR)
904         call MPI_Allgather(ielend(1),nres,MPI_INTEGER,&
905           ielend_all(1,0),nres,MPI_INTEGER,FG_COMM,IERROR)
906         if (lprint) then
907         write (iout,*) "iatel_s_all",(iatel_s_all(i),i=0,nfgtasks)
908         write (iout,*) "iatel_e_all",(iatel_e_all(i),i=0,nfgtasks)
909         write (iout,*) "iturn3_start_all",&
910           (iturn3_start_all(i),i=0,nfgtasks-1)
911         write (iout,*) "iturn3_end_all",&
912           (iturn3_end_all(i),i=0,nfgtasks-1)
913         write (iout,*) "iturn4_start_all",&
914           (iturn4_start_all(i),i=0,nfgtasks-1)
915         write (iout,*) "iturn4_end_all",&
916           (iturn4_end_all(i),i=0,nfgtasks-1)
917         write (iout,*) "The ielstart_all array"
918         do i=nnt,nct
919           write (iout,'(20i4)') i,(ielstart_all(i,j),j=0,nfgtasks-1)
920         enddo
921         write (iout,*) "The ielend_all array"
922         do i=nnt,nct
923           write (iout,'(20i4)') i,(ielend_all(i,j),j=0,nfgtasks-1)
924         enddo
925         call flush(iout)
926         endif
927         ntask_cont_from=0
928         ntask_cont_to=0
929         itask_cont_from(0)=fg_rank
930         itask_cont_to(0)=fg_rank
931         flag=.false.
932 !el        allocate(iturn3_sent(4,iturn3_start:iturn3_end))
933 !el        allocate(iturn4_sent(4,iturn4_start:iturn4_end)) !(4,maxres)
934         do ii=iturn3_start,iturn3_end
935           call add_int(ii,ii+2,iturn3_sent(1,ii),&
936                       ntask_cont_to,itask_cont_to,flag)
937         enddo
938         do ii=iturn4_start,iturn4_end
939           call add_int(ii,ii+3,iturn4_sent(1,ii),&
940                       ntask_cont_to,itask_cont_to,flag)
941         enddo
942         do ii=iturn3_start,iturn3_end
943           call add_int_from(ii,ii+2,ntask_cont_from,itask_cont_from)
944         enddo
945         do ii=iturn4_start,iturn4_end
946           call add_int_from(ii,ii+3,ntask_cont_from,itask_cont_from)
947         enddo
948         if (lprint) then
949         write (iout,*) "After turn3 ntask_cont_from",ntask_cont_from,&
950          " ntask_cont_to",ntask_cont_to
951         write (iout,*) "itask_cont_from",&
952           (itask_cont_from(i),i=1,ntask_cont_from)
953         write (iout,*) "itask_cont_to",&
954           (itask_cont_to(i),i=1,ntask_cont_to)
955         call flush(iout)
956         endif
957 !        write (iout,*) "Loop forward"
958 !        call flush(iout)
959         do i=iatel_s,iatel_e
960 !          write (iout,*) "from loop i=",i
961 !          call flush(iout)
962           do j=ielstart(i),ielend(i)
963             call add_int_from(i,j,ntask_cont_from,itask_cont_from)
964           enddo
965         enddo
966 !        write (iout,*) "Loop backward iatel_e-1",iatel_e-1,
967 !     &     " iatel_e",iatel_e
968 !        call flush(iout)
969         nat_sent=0
970         do i=iatel_s,iatel_e
971 !          write (iout,*) "i",i," ielstart",ielstart(i),
972 !     &      " ielend",ielend(i)
973 !          call flush(iout)
974           flag=.false.
975           do j=ielstart(i),ielend(i)
976             call add_int(i,j,iint_sent(1,j,nat_sent+1),ntask_cont_to,&
977                         itask_cont_to,flag)
978           enddo
979           if (flag) then
980             nat_sent=nat_sent+1
981             iat_sent(nat_sent)=i
982           endif
983         enddo
984         if (lprint) then
985         write (iout,*)"After longrange ntask_cont_from",ntask_cont_from,&
986          " ntask_cont_to",ntask_cont_to
987         write (iout,*) "itask_cont_from",&
988           (itask_cont_from(i),i=1,ntask_cont_from)
989         write (iout,*) "itask_cont_to",&
990           (itask_cont_to(i),i=1,ntask_cont_to)
991         call flush(iout)
992         write (iout,*) "iint_sent"
993         do i=1,nat_sent
994           ii=iat_sent(i)
995           write (iout,'(20i4)') ii,(j,(iint_sent(k,j,i),k=1,4),&
996             j=ielstart(ii),ielend(ii))
997         enddo
998         write (iout,*) "iturn3_sent iturn3_start",iturn3_start,&
999           " iturn3_end",iturn3_end
1000         write (iout,'(20i4)') (i,(iturn3_sent(j,i),j=1,4),&
1001            i=iturn3_start,iturn3_end)
1002         write (iout,*) "iturn4_sent iturn4_start",iturn4_start,&
1003           " iturn4_end",iturn4_end
1004         write (iout,'(20i4)') (i,(iturn4_sent(j,i),j=1,4),&
1005            i=iturn4_start,iturn4_end)
1006         call flush(iout)
1007         endif
1008         call MPI_Gather(ntask_cont_from,1,MPI_INTEGER,&
1009          ntask_cont_from_all,1,MPI_INTEGER,king,FG_COMM,IERR)
1010 !        write (iout,*) "Gather ntask_cont_from ended"
1011 !        call flush(iout)
1012         call MPI_Gather(itask_cont_from(0),nfgtasks,MPI_INTEGER,&
1013          itask_cont_from_all(0,0),nfgtasks,MPI_INTEGER,king,&
1014          FG_COMM,IERR)
1015 !        write (iout,*) "Gather itask_cont_from ended"
1016 !        call flush(iout)
1017         call MPI_Gather(ntask_cont_to,1,MPI_INTEGER,ntask_cont_to_all,&
1018          1,MPI_INTEGER,king,FG_COMM,IERR)
1019 !        write (iout,*) "Gather ntask_cont_to ended"
1020 !        call flush(iout)
1021         call MPI_Gather(itask_cont_to,nfgtasks,MPI_INTEGER,&
1022          itask_cont_to_all,nfgtasks,MPI_INTEGER,king,FG_COMM,IERR)
1023 !        write (iout,*) "Gather itask_cont_to ended"
1024 !        call flush(iout)
1025         if (fg_rank.eq.king) then
1026           write (iout,*)"Contact receive task map (proc, #tasks, tasks)"
1027           do i=0,nfgtasks-1
1028             write (iout,'(20i4)') i,ntask_cont_from_all(i),&
1029               (itask_cont_from_all(j,i),j=1,ntask_cont_from_all(i)) 
1030           enddo
1031           write (iout,*)
1032           call flush(iout)
1033           write (iout,*) "Contact send task map (proc, #tasks, tasks)"
1034           do i=0,nfgtasks-1
1035             write (iout,'(20i4)') i,ntask_cont_to_all(i),&
1036              (itask_cont_to_all(j,i),j=1,ntask_cont_to_all(i)) 
1037           enddo
1038           write (iout,*)
1039           call flush(iout)
1040 ! Check if every send will have a matching receive
1041           ncheck_to=0
1042           ncheck_from=0
1043           do i=0,nfgtasks-1
1044             ncheck_to=ncheck_to+ntask_cont_to_all(i)
1045             ncheck_from=ncheck_from+ntask_cont_from_all(i)
1046           enddo
1047           write (iout,*) "Control sums",ncheck_from,ncheck_to
1048           if (ncheck_from.ne.ncheck_to) then
1049             write (iout,*) "Error: #receive differs from #send."
1050             write (iout,*) "Terminating program...!"
1051             call flush(iout)
1052             flag=.false.
1053           else
1054             flag=.true.
1055             do i=0,nfgtasks-1
1056               do j=1,ntask_cont_to_all(i)
1057                 ii=itask_cont_to_all(j,i)
1058                 do k=1,ntask_cont_from_all(ii)
1059                   if (itask_cont_from_all(k,ii).eq.i) then
1060                     if(lprint)write(iout,*)"Matching send/receive",i,ii
1061                     exit
1062                   endif
1063                 enddo
1064                 if (k.eq.ntask_cont_from_all(ii)+1) then
1065                   flag=.false.
1066                   write (iout,*) "Error: send by",j," to",ii,&
1067                     " would have no matching receive"
1068                 endif
1069               enddo
1070             enddo
1071           endif
1072           if (.not.flag) then
1073             write (iout,*) "Unmatched sends; terminating program"
1074             call flush(iout)
1075           endif
1076         endif
1077         call MPI_Bcast(flag,1,MPI_LOGICAL,king,FG_COMM,IERROR)
1078 !        write (iout,*) "flag broadcast ended flag=",flag
1079 !        call flush(iout)
1080         if (.not.flag) then
1081           call MPI_Finalize(IERROR)
1082           stop "Error in INIT_INT_TABLE: unmatched send/receive."
1083         endif
1084         call MPI_Comm_group(FG_COMM,fg_group,IERR)
1085 !        write (iout,*) "MPI_Comm_group ended"
1086 !        call flush(iout)
1087         call MPI_Group_incl(fg_group,ntask_cont_from+1,&
1088           itask_cont_from(0),CONT_FROM_GROUP,IERR)
1089         call MPI_Group_incl(fg_group,ntask_cont_to+1,itask_cont_to(0),&
1090           CONT_TO_GROUP,IERR)
1091         do i=1,nat_sent
1092           ii=iat_sent(i)
1093           iaux=4*(ielend(ii)-ielstart(ii)+1)
1094           call MPI_Group_translate_ranks(fg_group,iaux,&
1095             iint_sent(1,ielstart(ii),i),CONT_TO_GROUP,&
1096             iint_sent_local(1,ielstart(ii),i),IERR )
1097 !          write (iout,*) "Ranks translated i=",i
1098 !          call flush(iout)
1099         enddo
1100         iaux=4*(iturn3_end-iturn3_start+1)
1101         call MPI_Group_translate_ranks(fg_group,iaux,&
1102            iturn3_sent(1,iturn3_start),CONT_TO_GROUP,&
1103            iturn3_sent_local(1,iturn3_start),IERR)
1104         iaux=4*(iturn4_end-iturn4_start+1)
1105         call MPI_Group_translate_ranks(fg_group,iaux,&
1106            iturn4_sent(1,iturn4_start),CONT_TO_GROUP,&
1107            iturn4_sent_local(1,iturn4_start),IERR)
1108         if (lprint) then
1109         write (iout,*) "iint_sent_local"
1110         do i=1,nat_sent
1111           ii=iat_sent(i)
1112           write (iout,'(20i4)') ii,(j,(iint_sent_local(k,j,i),k=1,4),&
1113             j=ielstart(ii),ielend(ii))
1114           call flush(iout)
1115         enddo
1116         write (iout,*) "iturn3_sent_local iturn3_start",iturn3_start,&
1117           " iturn3_end",iturn3_end
1118         write (iout,'(20i4)') (i,(iturn3_sent_local(j,i),j=1,4),&
1119            i=iturn3_start,iturn3_end)
1120         write (iout,*) "iturn4_sent_local iturn4_start",iturn4_start,&
1121           " iturn4_end",iturn4_end
1122         write (iout,'(20i4)') (i,(iturn4_sent_local(j,i),j=1,4),&
1123            i=iturn4_start,iturn4_end)
1124         call flush(iout)
1125         endif
1126         call MPI_Group_free(fg_group,ierr)
1127         call MPI_Group_free(cont_from_group,ierr)
1128         call MPI_Group_free(cont_to_group,ierr)
1129         call MPI_Type_contiguous(3,MPI_DOUBLE_PRECISION,MPI_UYZ,IERROR)
1130         call MPI_Type_commit(MPI_UYZ,IERROR)
1131         call MPI_Type_contiguous(18,MPI_DOUBLE_PRECISION,MPI_UYZGRAD,&
1132           IERROR)
1133         call MPI_Type_commit(MPI_UYZGRAD,IERROR)
1134         call MPI_Type_contiguous(2,MPI_DOUBLE_PRECISION,MPI_MU,IERROR)
1135         call MPI_Type_commit(MPI_MU,IERROR)
1136         call MPI_Type_contiguous(4,MPI_DOUBLE_PRECISION,MPI_MAT1,IERROR)
1137         call MPI_Type_commit(MPI_MAT1,IERROR)
1138         call MPI_Type_contiguous(8,MPI_DOUBLE_PRECISION,MPI_MAT2,IERROR)
1139         call MPI_Type_commit(MPI_MAT2,IERROR)
1140         call MPI_Type_contiguous(6,MPI_DOUBLE_PRECISION,MPI_THET,IERROR)
1141         call MPI_Type_commit(MPI_THET,IERROR)
1142         call MPI_Type_contiguous(9,MPI_DOUBLE_PRECISION,MPI_GAM,IERROR)
1143         call MPI_Type_commit(MPI_GAM,IERROR)
1144
1145 !el        allocate(lentyp(0:nfgtasks-1))
1146 #ifndef MATGATHER
1147 ! 9/22/08 Derived types to send matrices which appear in correlation terms
1148         do i=0,nfgtasks-1
1149           if (ivec_count(i).eq.ivec_count(0)) then
1150             lentyp(i)=0
1151           else
1152             lentyp(i)=1
1153           endif
1154         enddo
1155         do ind_typ=lentyp(0),lentyp(nfgtasks-1)
1156         if (ind_typ.eq.0) then
1157           ichunk=ivec_count(0)
1158         else
1159           ichunk=ivec_count(1)
1160         endif
1161 !        do i=1,4
1162 !          blocklengths(i)=4
1163 !        enddo
1164 !        displs(1)=0
1165 !        do i=2,4
1166 !          displs(i)=displs(i-1)+blocklengths(i-1)*maxres
1167 !        enddo
1168 !        do i=1,4
1169 !          blocklengths(i)=blocklengths(i)*ichunk
1170 !        enddo
1171 !        write (iout,*) "blocklengths and displs"
1172 !        do i=1,4
1173 !          write (iout,*) i,blocklengths(i),displs(i)
1174 !        enddo
1175 !        call flush(iout)
1176 !        call MPI_Type_indexed(4,blocklengths(1),displs(1),
1177 !     &    MPI_DOUBLE_PRECISION,MPI_ROTAT1(ind_typ),IERROR)
1178 !        call MPI_Type_commit(MPI_ROTAT1(ind_typ),IERROR)
1179 !        write (iout,*) "MPI_ROTAT1",MPI_ROTAT1 
1180 !        do i=1,4
1181 !          blocklengths(i)=2
1182 !        enddo
1183 !        displs(1)=0
1184 !        do i=2,4
1185 !          displs(i)=displs(i-1)+blocklengths(i-1)*maxres
1186 !        enddo
1187 !        do i=1,4
1188 !          blocklengths(i)=blocklengths(i)*ichunk
1189 !        enddo
1190 !        write (iout,*) "blocklengths and displs"
1191 !        do i=1,4
1192 !          write (iout,*) i,blocklengths(i),displs(i)
1193 !        enddo
1194 !        call flush(iout)
1195 !        call MPI_Type_indexed(4,blocklengths(1),displs(1),
1196 !     &    MPI_DOUBLE_PRECISION,MPI_ROTAT2(ind_typ),IERROR)
1197 !        call MPI_Type_commit(MPI_ROTAT2(ind_typ),IERROR)
1198 !        write (iout,*) "MPI_ROTAT2",MPI_ROTAT2 
1199         do i=1,8
1200           blocklengths(i)=2
1201         enddo
1202         displs(1)=0
1203         do i=2,8
1204           displs(i)=displs(i-1)+blocklengths(i-1)*nres !maxres
1205         enddo
1206         do i=1,15
1207           blocklengths(i)=blocklengths(i)*ichunk
1208         enddo
1209         call MPI_Type_indexed(8,blocklengths,displs,&
1210           MPI_DOUBLE_PRECISION,MPI_PRECOMP11(ind_typ),IERROR)
1211         call MPI_Type_commit(MPI_PRECOMP11(ind_typ),IERROR)
1212         do i=1,8
1213           blocklengths(i)=4
1214         enddo
1215         displs(1)=0
1216         do i=2,8
1217           displs(i)=displs(i-1)+blocklengths(i-1)*nres !maxres
1218         enddo
1219         do i=1,15
1220           blocklengths(i)=blocklengths(i)*ichunk
1221         enddo
1222         call MPI_Type_indexed(8,blocklengths,displs,&
1223           MPI_DOUBLE_PRECISION,MPI_PRECOMP12(ind_typ),IERROR)
1224         call MPI_Type_commit(MPI_PRECOMP12(ind_typ),IERROR)
1225         do i=1,6
1226           blocklengths(i)=4
1227         enddo
1228         displs(1)=0
1229         do i=2,6
1230           displs(i)=displs(i-1)+blocklengths(i-1)*nres !maxres
1231         enddo
1232         do i=1,6
1233           blocklengths(i)=blocklengths(i)*ichunk
1234         enddo
1235         call MPI_Type_indexed(6,blocklengths,displs,&
1236           MPI_DOUBLE_PRECISION,MPI_PRECOMP22(ind_typ),IERROR)
1237         call MPI_Type_commit(MPI_PRECOMP22(ind_typ),IERROR)
1238         do i=1,2
1239           blocklengths(i)=8
1240         enddo
1241         displs(1)=0
1242         do i=2,2
1243           displs(i)=displs(i-1)+blocklengths(i-1)*nres !maxres
1244         enddo
1245         do i=1,2
1246           blocklengths(i)=blocklengths(i)*ichunk
1247         enddo
1248         call MPI_Type_indexed(2,blocklengths,displs,&
1249           MPI_DOUBLE_PRECISION,MPI_PRECOMP23(ind_typ),IERROR)
1250         call MPI_Type_commit(MPI_PRECOMP23(ind_typ),IERROR)
1251         do i=1,4
1252           blocklengths(i)=1
1253         enddo
1254         displs(1)=0
1255         do i=2,4
1256           displs(i)=displs(i-1)+blocklengths(i-1)*nres !maxres
1257         enddo
1258         do i=1,4
1259           blocklengths(i)=blocklengths(i)*ichunk
1260         enddo
1261         call MPI_Type_indexed(4,blocklengths,displs,&
1262           MPI_DOUBLE_PRECISION,MPI_ROTAT_OLD(ind_typ),IERROR)
1263         call MPI_Type_commit(MPI_ROTAT_OLD(ind_typ),IERROR)
1264         enddo
1265 #endif
1266       endif
1267       iint_start=ivec_start+1
1268       iint_end=ivec_end+1
1269       do i=0,nfgtasks-1
1270           iint_count(i)=ivec_count(i)
1271           iint_displ(i)=ivec_displ(i)
1272           ivec_displ(i)=ivec_displ(i)-1
1273           iset_displ(i)=iset_displ(i)-1
1274           ithet_displ(i)=ithet_displ(i)-1
1275           iphi_displ(i)=iphi_displ(i)-1
1276           iphi1_displ(i)=iphi1_displ(i)-1
1277           ibond_displ(i)=ibond_displ(i)-1
1278       enddo
1279       if (nfgtasks.gt.1 .and. fg_rank.eq.king &
1280           .and. (me.eq.0 .or. .not. out1file)) then
1281         write (iout,*) "IVEC_DISPL, IVEC_COUNT, ISET_START, ISET_COUNT"
1282         do i=0,nfgtasks-1
1283           write (iout,*) i,ivec_displ(i),ivec_count(i),iset_displ(i),&
1284             iset_count(i)
1285         enddo
1286         write (iout,*) "iphi_start",iphi_start," iphi_end",iphi_end,&
1287           " iphi1_start",iphi1_start," iphi1_end",iphi1_end
1288         write (iout,*)"IPHI_COUNT, IPHI_DISPL, IPHI1_COUNT, IPHI1_DISPL"
1289         do i=0,nfgtasks-1
1290           write (iout,*) i,iphi_count(i),iphi_displ(i),iphi1_count(i),&
1291             iphi1_displ(i)
1292         enddo
1293         write(iout,'(i10,a,i10,a,i10,a/a,i3,a)') n_sc_int_tot,' SC-SC ',&
1294           nele_int_tot,' electrostatic and ',nscp_int_tot,&
1295           ' SC-p interactions','were distributed among',nfgtasks,&
1296           ' fine-grain processors.'
1297       endif
1298 #else
1299       loc_start=2
1300       loc_end=nres-1
1301       ithet_start=3 
1302       ithet_end=nres
1303       iturn3_start=nnt
1304       iturn3_end=nct-3
1305       iturn4_start=nnt
1306       iturn4_end=nct-4
1307       iphi_start=nnt+3
1308       iphi_end=nct
1309       iphi1_start=4
1310       iphi1_end=nres
1311       idihconstr_start=1
1312       idihconstr_end=ndih_constr
1313       iphid_start=iphi_start
1314       iphid_end=iphi_end-1
1315       itau_start=4
1316       itau_end=nres
1317       ibond_start=2
1318       ibond_end=nres-1
1319       ibondp_start=nnt
1320       ibondp_end=nct-1
1321       ivec_start=1
1322       ivec_end=nres-1
1323       iset_start=3
1324       iset_end=nres+1
1325       iint_start=2
1326       iint_end=nres-1
1327       ilip_start=1
1328       ilip_end=nres
1329 #endif
1330 !el       common /przechowalnia/
1331 !      deallocate(iturn3_start_all)
1332 !      deallocate(iturn3_end_all)
1333 !      deallocate(iturn4_start_all)
1334 !      deallocate(iturn4_end_all)
1335 !      deallocate(iatel_s_all)
1336 !      deallocate(iatel_e_all)
1337 !      deallocate(ielstart_all)
1338 !      deallocate(ielend_all)
1339
1340 !      deallocate(ntask_cont_from_all)
1341 !      deallocate(ntask_cont_to_all)
1342 !      deallocate(itask_cont_from_all)
1343 !      deallocate(itask_cont_to_all)
1344 !el----------
1345       return
1346       end subroutine init_int_table
1347 #ifdef MPI
1348 !-----------------------------------------------------------------------------
1349       subroutine add_int(ii,jj,itask,ntask_cont_to,itask_cont_to,flag)
1350
1351 !el      implicit none
1352 !      include "DIMENSIONS"
1353 !      include "COMMON.INTERACT"
1354 !      include "COMMON.SETUP"
1355 !      include "COMMON.IOUNITS"
1356       integer :: ii,jj,ntask_cont_to
1357       integer,dimension(4) :: itask
1358       integer :: itask_cont_to(0:nfgtasks-1)    !(0:max_fg_procs-1)
1359       logical :: flag
1360 !el      integer,dimension(0:nfgtasks) :: iturn3_start_all,iturn3_end_all,iturn4_start_all,&
1361 !el       iturn4_end_all,iatel_s_all,iatel_e_all        !(0:max_fg_procs)
1362 !el      integer,dimension(nres,0:nfgtasks-1) :: ielstart_all,ielend_all        !(maxres,0:max_fg_procs-1)
1363 !el      common /przechowalnia/ iturn3_start_all,iturn3_end_all,iturn4_start_all,&
1364 !el       iturn4_end_all,iatel_s_all,iatel_e_all,ielstart_all,ielend_all
1365       integer :: iproc,isent,k,l
1366 ! Determines whether to send interaction ii,jj to other processors; a given
1367 ! interaction can be sent to at most 2 processors.
1368 ! Sets flag=.true. if interaction ii,jj needs to be sent to at least 
1369 ! one processor, otherwise flag is unchanged from the input value.
1370       isent=0
1371       itask(1)=fg_rank
1372       itask(2)=fg_rank
1373       itask(3)=fg_rank
1374       itask(4)=fg_rank
1375 !      write (iout,*) "ii",ii," jj",jj
1376 ! Loop over processors to check if anybody could need interaction ii,jj
1377       do iproc=0,fg_rank-1
1378 ! Check if the interaction matches any turn3 at iproc
1379         do k=iturn3_start_all(iproc),iturn3_end_all(iproc)
1380           l=k+2
1381           if (k.eq.ii-1 .and. l.eq.jj-1 .or. k.eq.ii-1 .and. l.eq.jj+1 &
1382          .or. k.eq.ii+1 .and. l.eq.jj+1 .or. k.eq.ii+1 .and. l.eq.jj-1) &
1383           then 
1384 !            write (iout,*) "turn3 to iproc",iproc," ij",ii,jj,"kl",k,l
1385 !            call flush(iout)
1386             flag=.true.
1387             if (iproc.ne.itask(1).and.iproc.ne.itask(2) &
1388               .and.iproc.ne.itask(3).and.iproc.ne.itask(4)) then
1389               isent=isent+1
1390               itask(isent)=iproc
1391               call add_task(iproc,ntask_cont_to,itask_cont_to)
1392             endif
1393           endif
1394         enddo
1395 ! Check if the interaction matches any turn4 at iproc
1396         do k=iturn4_start_all(iproc),iturn4_end_all(iproc)
1397           l=k+3
1398           if (k.eq.ii-1 .and. l.eq.jj-1 .or. k.eq.ii-1 .and. l.eq.jj+1 &
1399          .or. k.eq.ii+1 .and. l.eq.jj+1 .or. k.eq.ii+1 .and. l.eq.jj-1) &
1400           then 
1401 !            write (iout,*) "turn3 to iproc",iproc," ij",ii,jj," kl",k,l
1402 !            call flush(iout)
1403             flag=.true.
1404             if (iproc.ne.itask(1).and.iproc.ne.itask(2) &
1405               .and.iproc.ne.itask(3).and.iproc.ne.itask(4)) then
1406               isent=isent+1
1407               itask(isent)=iproc
1408               call add_task(iproc,ntask_cont_to,itask_cont_to)
1409             endif
1410           endif
1411         enddo
1412         if (iatel_s_all(iproc).gt.0 .and. iatel_e_all(iproc).gt.0 .and. &
1413         iatel_s_all(iproc).le.ii-1 .and. iatel_e_all(iproc).ge.ii-1)then
1414           if (ielstart_all(ii-1,iproc).le.jj-1.and. &
1415               ielend_all(ii-1,iproc).ge.jj-1) then
1416             flag=.true.
1417             if (iproc.ne.itask(1).and.iproc.ne.itask(2) &
1418               .and.iproc.ne.itask(3).and.iproc.ne.itask(4)) then
1419               isent=isent+1
1420               itask(isent)=iproc
1421               call add_task(iproc,ntask_cont_to,itask_cont_to)
1422             endif
1423           endif
1424           if (ielstart_all(ii-1,iproc).le.jj+1.and. &
1425               ielend_all(ii-1,iproc).ge.jj+1) then
1426             flag=.true.
1427             if (iproc.ne.itask(1).and.iproc.ne.itask(2) &
1428               .and.iproc.ne.itask(3).and.iproc.ne.itask(4)) then
1429               isent=isent+1
1430               itask(isent)=iproc
1431               call add_task(iproc,ntask_cont_to,itask_cont_to)
1432             endif
1433           endif
1434         endif
1435       enddo
1436       return
1437       end subroutine add_int
1438 !-----------------------------------------------------------------------------
1439       subroutine add_int_from(ii,jj,ntask_cont_from,itask_cont_from)
1440
1441 !el      use MPI_data
1442 !el      implicit none
1443 !      include "DIMENSIONS"
1444 !      include "COMMON.INTERACT"
1445 !      include "COMMON.SETUP"
1446 !      include "COMMON.IOUNITS"
1447       integer :: ii,jj,itask(2),ntask_cont_from,&
1448        itask_cont_from(0:nfgtasks-1)    !(0:max_fg_procs)
1449       logical :: flag
1450 !el      integer,dimension(0:nfgtasks) :: iturn3_start_all,iturn3_end_all,&
1451 !el       iturn4_start_all,iturn4_end_all,iatel_s_all,iatel_e_all       !(0:max_fg_procs)
1452 !el      integer,dimension(nres,0:nfgtasks-1) :: ielstart_all,ielend_all        !(maxres,0:max_fg_procs-1)
1453 !el      common /przechowalnia/ iturn3_start_all,iturn3_end_all,iturn4_start_all,&
1454 !el       iturn4_end_all,iatel_s_all,iatel_e_all,ielstart_all,ielend_all
1455       integer :: iproc,k,l
1456       do iproc=fg_rank+1,nfgtasks-1
1457         do k=iturn3_start_all(iproc),iturn3_end_all(iproc)
1458           l=k+2
1459           if (k.eq.ii+1 .and. l.eq.jj+1 .or. k.eq.ii+1.and.l.eq.jj-1 &
1460          .or. k.eq.ii-1 .and. l.eq.jj-1 .or. k.eq.ii-1 .and. l.eq.jj+1) &
1461           then
1462 !            write (iout,*)"turn3 from iproc",iproc," ij",ii,jj," kl",k,l
1463             call add_task(iproc,ntask_cont_from,itask_cont_from)
1464           endif
1465         enddo 
1466         do k=iturn4_start_all(iproc),iturn4_end_all(iproc)
1467           l=k+3
1468           if (k.eq.ii+1 .and. l.eq.jj+1 .or. k.eq.ii+1.and.l.eq.jj-1 &
1469          .or. k.eq.ii-1 .and. l.eq.jj-1 .or. k.eq.ii-1 .and. l.eq.jj+1) &
1470           then
1471 !            write (iout,*)"turn4 from iproc",iproc," ij",ii,jj," kl",k,l
1472             call add_task(iproc,ntask_cont_from,itask_cont_from)
1473           endif
1474         enddo 
1475         if (iatel_s_all(iproc).gt.0 .and. iatel_e_all(iproc).gt.0) then
1476           if (ii+1.ge.iatel_s_all(iproc).and.ii+1.le.iatel_e_all(iproc)) &
1477           then
1478             if (jj+1.ge.ielstart_all(ii+1,iproc).and. &
1479                 jj+1.le.ielend_all(ii+1,iproc)) then
1480               call add_task(iproc,ntask_cont_from,itask_cont_from)
1481             endif            
1482             if (jj-1.ge.ielstart_all(ii+1,iproc).and. &
1483                 jj-1.le.ielend_all(ii+1,iproc)) then
1484               call add_task(iproc,ntask_cont_from,itask_cont_from)
1485             endif
1486           endif
1487           if (ii-1.ge.iatel_s_all(iproc).and.ii-1.le.iatel_e_all(iproc)) &
1488           then
1489             if (jj-1.ge.ielstart_all(ii-1,iproc).and. &
1490                 jj-1.le.ielend_all(ii-1,iproc)) then
1491               call add_task(iproc,ntask_cont_from,itask_cont_from)
1492             endif
1493             if (jj+1.ge.ielstart_all(ii-1,iproc).and. &
1494                 jj+1.le.ielend_all(ii-1,iproc)) then
1495                call add_task(iproc,ntask_cont_from,itask_cont_from)
1496             endif
1497           endif
1498         endif
1499       enddo
1500       return
1501       end subroutine add_int_from
1502 !-----------------------------------------------------------------------------
1503       subroutine add_task(iproc,ntask_cont,itask_cont)
1504
1505 !el      use MPI_data
1506 !el      implicit none
1507 !      include "DIMENSIONS"
1508       integer :: iproc,ntask_cont,itask_cont(0:nfgtasks-1)      !(0:max_fg_procs-1)
1509       integer :: ii
1510       do ii=1,ntask_cont
1511         if (itask_cont(ii).eq.iproc) return
1512       enddo
1513       ntask_cont=ntask_cont+1
1514       itask_cont(ntask_cont)=iproc
1515       return
1516       end subroutine add_task
1517 #endif
1518 !-----------------------------------------------------------------------------
1519 #if defined MPI || defined WHAM_RUN
1520       subroutine int_partition(int_index,lower_index,upper_index,atom,&
1521        at_start,at_end,first_atom,last_atom,int_gr,jat_start,jat_end,*)
1522
1523 !      implicit real*8 (a-h,o-z)
1524 !      include 'DIMENSIONS'
1525 !      include 'COMMON.IOUNITS'
1526       integer :: int_index,lower_index,upper_index,atom,at_start,at_end,&
1527        first_atom,last_atom,int_gr,jat_start,jat_end,int_index_old
1528       logical :: lprn
1529       lprn=.false.
1530       if (lprn) write (iout,*) 'int_index=',int_index
1531       int_index_old=int_index
1532       int_index=int_index+last_atom-first_atom+1
1533       if (lprn) &
1534          write (iout,*) 'int_index=',int_index,&
1535                      ' int_index_old',int_index_old,&
1536                      ' lower_index=',lower_index,&
1537                      ' upper_index=',upper_index,&
1538                      ' atom=',atom,' first_atom=',first_atom,&
1539                      ' last_atom=',last_atom
1540       if (int_index.ge.lower_index) then
1541         int_gr=int_gr+1
1542         if (at_start.eq.0) then
1543           at_start=atom
1544           jat_start=first_atom-1+lower_index-int_index_old
1545         else
1546           jat_start=first_atom
1547         endif
1548         if (lprn) write (iout,*) 'jat_start',jat_start
1549         if (int_index.ge.upper_index) then
1550           at_end=atom
1551           jat_end=first_atom-1+upper_index-int_index_old
1552           return 1
1553         else
1554           jat_end=last_atom
1555         endif
1556         if (lprn) write (iout,*) 'jat_end',jat_end
1557       endif
1558       return
1559       end subroutine int_partition
1560 #endif
1561 !-----------------------------------------------------------------------------
1562 #ifndef CLUSTER
1563       subroutine hpb_partition
1564
1565 !      implicit real*8 (a-h,o-z)
1566 !      include 'DIMENSIONS'
1567 #ifdef MPI
1568       include 'mpif.h'
1569 #endif
1570 !      include 'COMMON.SBRIDGE'
1571 !      include 'COMMON.IOUNITS'
1572 !      include 'COMMON.SETUP'
1573 #ifdef MPI
1574       call int_bounds(nhpb,link_start,link_end)
1575       write (iout,*) 'Processor',fg_rank,' CG group',kolor,&
1576         ' absolute rank',MyRank,&
1577         ' nhpb',nhpb,' link_start=',link_start,&
1578         ' link_end',link_end
1579 #else
1580       link_start=1
1581       link_end=nhpb
1582 #endif
1583       return
1584       end subroutine hpb_partition
1585 #endif
1586 !-----------------------------------------------------------------------------
1587 ! misc.f in module io_base
1588 !-----------------------------------------------------------------------------
1589 !-----------------------------------------------------------------------------
1590 ! parmread.F
1591 !-----------------------------------------------------------------------------
1592       subroutine getenv_loc(var, val)
1593
1594       character(*) :: var, val
1595
1596 #ifdef WINIFL
1597       character(len=2000) :: line
1598 !el      external ilen
1599
1600       open (196,file='env',status='old',readonly,shared)
1601       iread=0
1602 !      write(*,*)'looking for ',var
1603 10    read(196,*,err=11,end=11)line
1604       iread=index(line,var)
1605 !      write(*,*)iread,' ',var,' ',line
1606       if (iread.eq.0) go to 10 
1607 !      write(*,*)'---> ',line
1608 11    continue
1609       if(iread.eq.0) then
1610 !       write(*,*)'CHUJ'
1611        val=''
1612       else
1613        iread=iread+ilen(var)+1
1614        read (line(iread:),*,err=12,end=12) val
1615 !       write(*,*)'OK: ',var,' = ',val
1616       endif
1617       close(196)
1618       return
1619 12    val=''
1620       close(196)
1621 #elif (defined CRAY)
1622       integer :: lennam,lenval,ierror
1623 !
1624 !        getenv using a POSIX call, useful on the T3D
1625 !        Sept 1996, comment out error check on advice of H. Pritchard
1626 !
1627       lennam = len(var)
1628       if(lennam.le.0) stop '--error calling getenv--'
1629       call pxfgetenv(var,lennam,val,lenval,ierror)
1630 !-HP- if(ierror.ne.0) stop '--error returned by pxfgetenv--'
1631 #else
1632       call getenv(var,val)
1633 #endif
1634
1635       return
1636       end subroutine getenv_loc
1637 !-----------------------------------------------------------------------------
1638 ! readrtns_CSA.F
1639 !-----------------------------------------------------------------------------
1640       subroutine setup_var
1641
1642       integer :: i
1643 !      implicit real*8 (a-h,o-z)
1644 !      include 'DIMENSIONS'
1645 !      include 'COMMON.IOUNITS'
1646 !      include 'COMMON.GEO'
1647 !      include 'COMMON.VAR'
1648 !      include 'COMMON.INTERACT'
1649 !      include 'COMMON.LOCAL'
1650 !      include 'COMMON.NAMES'
1651 !      include 'COMMON.CHAIN'
1652 !      include 'COMMON.FFIELD'
1653 !      include 'COMMON.SBRIDGE'
1654 !      include 'COMMON.HEADER'
1655 !      include 'COMMON.CONTROL'
1656 !      include 'COMMON.DBASE'
1657 !      include 'COMMON.THREAD'
1658 !      include 'COMMON.TIME1'
1659 ! Set up variable list.
1660       ntheta=nres-2
1661       nphi=nres-3
1662       nvar=ntheta+nphi
1663       nside=0
1664       do i=2,nres-1
1665 #ifdef WHAM_RUN
1666         if (itype(i).ne.10) then
1667 #else
1668         if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then
1669 #endif
1670           nside=nside+1
1671           ialph(i,1)=nvar+nside
1672           ialph(nside,2)=i
1673         endif
1674       enddo
1675       if (indphi.gt.0) then
1676         nvar=nphi
1677       else if (indback.gt.0) then
1678         nvar=nphi+ntheta
1679       else
1680         nvar=nvar+2*nside
1681       endif
1682 !d    write (iout,'(3i4)') (i,ialph(i,1),ialph(i,2),i=2,nres-1)
1683       return
1684       end subroutine setup_var
1685 !-----------------------------------------------------------------------------
1686 ! rescode.f
1687 !-----------------------------------------------------------------------------
1688       integer function rescode(iseq,nam,itype)
1689
1690       use io_base, only: ucase
1691 !      implicit real*8 (a-h,o-z)
1692 !      include 'DIMENSIONS'
1693 !      include 'COMMON.NAMES'
1694 !      include 'COMMON.IOUNITS'
1695       character(len=3) :: nam   !,ucase
1696       integer :: iseq,itype,i
1697
1698       if (itype.eq.0) then
1699
1700       do i=-ntyp1,ntyp1
1701         if (ucase(nam).eq.restyp(i)) then
1702           rescode=i
1703           return
1704         endif
1705       enddo
1706
1707       else
1708
1709       do i=-ntyp1,ntyp1
1710         if (nam(1:1).eq.onelet(i)) then
1711           rescode=i
1712           return  
1713         endif  
1714       enddo
1715
1716       endif
1717       write (iout,10) iseq,nam
1718       stop
1719    10 format ('**** Error - residue',i4,' has an unresolved name ',a3)
1720       end function rescode
1721 !-----------------------------------------------------------------------------
1722 ! timing.F
1723 !-----------------------------------------------------------------------------
1724 ! $Date: 1994/10/05 16:41:52 $
1725 ! $Revision: 2.2 $
1726 !
1727       subroutine set_timers
1728 !
1729 !el      implicit none
1730 !el      real(kind=8) :: tcpu
1731 !      include 'COMMON.TIME1'
1732 !#ifdef MP
1733 #ifdef MPI
1734       include 'mpif.h'
1735 #endif
1736 ! Diminish the assigned time limit a little so that there is some time to
1737 ! end a batch job
1738 !     timlim=batime-150.0
1739 ! Calculate the initial time, if it is not zero (e.g. for the SUN).
1740       stime=tcpu()
1741 #if .not. defined(WHAM_RUN) && .not. defined(CLUSTER)
1742 #ifdef MPI
1743       walltime=MPI_WTIME()
1744       time_reduce=0.0d0
1745       time_allreduce=0.0d0
1746       time_bcast=0.0d0
1747       time_gather=0.0d0
1748       time_sendrecv=0.0d0
1749       time_scatter=0.0d0
1750       time_scatter_fmat=0.0d0
1751       time_scatter_ginv=0.0d0
1752       time_scatter_fmatmult=0.0d0
1753       time_scatter_ginvmult=0.0d0
1754       time_barrier_e=0.0d0
1755       time_barrier_g=0.0d0
1756       time_enecalc=0.0d0
1757       time_sumene=0.0d0
1758       time_lagrangian=0.0d0
1759       time_sumgradient=0.0d0
1760       time_intcartderiv=0.0d0
1761       time_inttocart=0.0d0
1762       time_ginvmult=0.0d0
1763       time_fricmatmult=0.0d0
1764       time_cartgrad=0.0d0
1765       time_bcastc=0.0d0
1766       time_bcast7=0.0d0
1767       time_bcastw=0.0d0
1768       time_intfcart=0.0d0
1769       time_vec=0.0d0
1770       time_mat=0.0d0
1771       time_fric=0.0d0
1772       time_stoch=0.0d0
1773       time_fricmatmult=0.0d0
1774       time_fsample=0.0d0
1775 #endif
1776 #endif
1777 !d    print *,' in SET_TIMERS stime=',stime
1778       return
1779       end subroutine set_timers
1780 !-----------------------------------------------------------------------------
1781 #ifndef CLUSTER
1782       logical function stopx(nf)
1783 ! This function returns .true. if one of the following reasons to exit SUMSL
1784 ! occurs. The "reason" code is stored in WHATSUP passed thru a COMMON block:
1785 !
1786 !... WHATSUP = 0 - go on, no reason to stop. Stopx will return .false.
1787 !...           1 - Time up in current node;
1788 !...           2 - STOP signal was received from another node because the
1789 !...               node's task was accomplished (parallel only);
1790 !...          -1 - STOP signal was received from another node because of error;
1791 !...          -2 - STOP signal was received from another node, because 
1792 !...               the node's time was up.
1793 !      implicit real*8 (a-h,o-z)
1794 !      include 'DIMENSIONS'
1795 !el#ifdef WHAM_RUN
1796 !el      use control_data, only:WhatsUp
1797 !el#endif
1798 #ifdef MP
1799 !el      use MPI_data   !include 'COMMON.INFO'
1800       include 'mpif.h'
1801 #endif
1802       integer :: nf
1803 !el      logical :: ovrtim
1804
1805 !      include 'COMMON.IOUNITS'
1806 !      include 'COMMON.TIME1'
1807       integer :: Kwita
1808
1809 !d    print *,'Processor',MyID,' NF=',nf
1810 !d      write (iout,*) "stopx: ",nf
1811 #ifndef WHAM_RUN
1812 #ifndef MPI
1813       if (ovrtim()) then
1814 ! Finish if time is up.
1815          stopx = .true.
1816          WhatsUp=1
1817 #ifdef MPL
1818       else if (mod(nf,100).eq.0) then
1819 ! Other processors might have finished. Check this every 100th function 
1820 ! evaluation.
1821 ! Master checks if any other processor has sent accepted conformation(s) to it. 
1822          if (MyID.ne.MasterID) call receive_mcm_info
1823          if (MyID.eq.MasterID) call receive_conf
1824 !d       print *,'Processor ',MyID,' is checking STOP: nf=',nf
1825          call recv_stop_sig(Kwita)
1826          if (Kwita.eq.-1) then
1827            write (iout,'(a,i4,a,i5)') 'Processor',&
1828            MyID,' has received STOP signal in STOPX; NF=',nf
1829            write (*,'(a,i4,a,i5)') 'Processor',&
1830            MyID,' has received STOP signal in STOPX; NF=',nf
1831            stopx=.true.
1832            WhatsUp=2
1833          elseif (Kwita.eq.-2) then
1834            write (iout,*) &
1835           'Processor',MyID,' received TIMEUP-STOP signal in SUMSL.'
1836            write (*,*) &
1837           'Processor',MyID,' received TIMEUP-STOP signal in SUMSL.'
1838            WhatsUp=-2
1839            stopx=.true.  
1840          else if (Kwita.eq.-3) then
1841            write (iout,*) &
1842           'Processor',MyID,' received ERROR-STOP signal in SUMSL.'
1843            write (*,*) &
1844           'Processor',MyID,' received ERROR-STOP signal in SUMSL.'
1845            WhatsUp=-1
1846            stopx=.true.
1847          else
1848            stopx=.false.
1849            WhatsUp=0
1850          endif
1851 #endif
1852       else
1853          stopx = .false.
1854          WhatsUp=0
1855       endif
1856 #else
1857       stopx=.false.
1858 !d      write (iout,*) "stopx set at .false."
1859 #endif
1860
1861 #ifdef OSF
1862 ! Check for FOUND_NAN flag
1863       if (FOUND_NAN) then
1864         write(iout,*)"   ***   stopx : Found a NaN"
1865         stopx=.true.
1866       endif
1867 #endif
1868 #else
1869       if (ovrtim()) then
1870 ! Finish if time is up.
1871          stopx = .true.
1872          WhatsUp=1
1873       else if (cutoffviol) then
1874         stopx = .true.
1875         WhatsUp=2
1876       else
1877         stopx=.false.
1878       endif
1879 #endif
1880       return
1881       end function stopx
1882 !-----------------------------------------------------------------------------
1883 #else
1884       logical function stopx(nf)
1885 !
1886 !     ..................................................................
1887 !
1888 !     *****PURPOSE...
1889 !     THIS FUNCTION MAY SERVE AS THE STOPX (ASYNCHRONOUS INTERRUPTION)
1890 !     FUNCTION FOR THE NL2SOL (NONLINEAR LEAST-SQUARES) PACKAGE AT
1891 !     THOSE INSTALLATIONS WHICH DO NOT WISH TO IMPLEMENT A
1892 !     DYNAMIC STOPX.
1893 !
1894 !     *****ALGORITHM NOTES...
1895 !     AT INSTALLATIONS WHERE THE NL2SOL SYSTEM IS USED
1896 !     INTERACTIVELY, THIS DUMMY STOPX SHOULD BE REPLACED BY A
1897 !     FUNCTION THAT RETURNS .TRUE. IF AND ONLY IF THE INTERRUPT
1898 !     (BREAK) KEY HAS BEEN PRESSED SINCE THE LAST CALL ON STOPX.
1899 !
1900 !     $$$ MODIFIED FOR USE AS  THE TIMER ROUTINE.
1901 !     $$$                              WHEN THE TIME LIMIT HAS BEEN
1902 !     $$$ REACHED     STOPX IS SET TO .TRUE  AND INITIATES (IN ITSUM)
1903 !     $$$ AND ORDERLY EXIT OUT OF SUMSL.  IF ARRAYS IV AND V ARE
1904 !     $$$ SAVED, THE SUMSL ROUTINES CAN BE RESTARTED AT THE SAME
1905 !     $$$ POINT AT WHICH THEY WERE INTERRUPTED.
1906 !
1907 !     ..................................................................
1908 !
1909 !      include 'DIMENSIONS'
1910       integer :: nf
1911 !      logical ovrtim
1912 !      include 'COMMON.IOUNITS'
1913 !      include 'COMMON.TIME1'
1914 #ifdef MPL
1915 !     include 'COMMON.INFO'
1916       integer :: Kwita
1917
1918 !d    print *,'Processor',MyID,' NF=',nf
1919 #endif
1920       if (ovrtim()) then
1921 ! Finish if time is up.
1922          stopx = .true.
1923 #ifdef MPL
1924       else if (mod(nf,100).eq.0) then
1925 ! Other processors might have finished. Check this every 100th function 
1926 ! evaluation.
1927 !d       print *,'Processor ',MyID,' is checking STOP: nf=',nf
1928          call recv_stop_sig(Kwita)
1929          if (Kwita.eq.-1) then
1930            write (iout,'(a,i4,a,i5)') 'Processor',&
1931            MyID,' has received STOP signal in STOPX; NF=',nf
1932            write (*,'(a,i4,a,i5)') 'Processor',&
1933            MyID,' has received STOP signal in STOPX; NF=',nf
1934            stopx=.true.
1935          else
1936            stopx=.false.
1937          endif
1938 #endif
1939       else
1940          stopx = .false.
1941       endif
1942       return
1943       end function stopx
1944 #endif
1945 !-----------------------------------------------------------------------------
1946       logical function ovrtim()
1947
1948 !      include 'DIMENSIONS'
1949 !      include 'COMMON.IOUNITS'
1950 !      include 'COMMON.TIME1'
1951 !el      real(kind=8) :: tcpu
1952       real(kind=8) :: curtim
1953 #ifdef MPI
1954       include "mpif.h"
1955       curtim = MPI_Wtime()-walltime
1956 #else
1957       curtim= tcpu()
1958 #endif
1959 !  curtim is the current time in seconds.
1960 !      write (iout,*) "curtim",curtim," timlim",timlim," safety",safety
1961 #ifndef WHAM_RUN
1962       if (curtim .ge. timlim - safety) then
1963         write (iout,'(a,f10.2,a,f10.2,a,f10.2,a)') &
1964         "***************** Elapsed time (",curtim,&
1965         " s) is within the safety limit (",safety,&
1966         " s) of the allocated time (",timlim," s). Terminating."
1967         ovrtim=.true.
1968       else
1969         ovrtim=.false.
1970       endif
1971 #else
1972       ovrtim=.false.
1973 #endif
1974 !elwrite (iout,*) "ovrtim",ovrtim
1975       return
1976       end function ovrtim
1977 !-----------------------------------------------------------------------------
1978       real(kind=8) function tcpu()
1979
1980 !      include 'COMMON.TIME1'
1981       real(kind=8) :: seconds
1982 #ifdef ES9000
1983 !***************************
1984 ! Next definition for EAGLE (ibm-es9000)
1985       real(kind=8) :: micseconds
1986       integer :: rcode
1987       tcpu=cputime(micseconds,rcode)
1988       tcpu=(micseconds/1.0E6) - stime
1989 !***************************
1990 #endif
1991 #ifdef SUN
1992 !***************************
1993 ! Next definitions for sun
1994       REAL(kind=8) ::  ECPU,ETIME,ETCPU
1995       real(kind=8),dimension(2) :: tarray
1996       tcpu=etime(tarray)
1997       tcpu=tarray(1)
1998 !***************************
1999 #endif
2000 #ifdef KSR
2001 !***************************
2002 ! Next definitions for ksr
2003 ! this function uses the ksr timer ALL_SECONDS from the PMON library to
2004 ! return the elapsed time in seconds
2005       tcpu= all_seconds() - stime
2006 !***************************
2007 #endif
2008 #ifdef SGI
2009 !***************************
2010 ! Next definitions for sgi
2011       real(kind=4) :: timar(2), etime
2012       seconds = etime(timar)
2013 !d    print *,'seconds=',seconds,' stime=',stime
2014 !      usrsec = timar(1)
2015 !      syssec = timar(2)
2016       tcpu=seconds - stime
2017 !***************************
2018 #endif
2019
2020 #ifdef LINUX
2021 !***************************
2022 ! Next definitions for sgi
2023       real(kind=4) :: timar(2), etime
2024       seconds = etime(timar)
2025 !d    print *,'seconds=',seconds,' stime=',stime
2026 !      usrsec = timar(1)
2027 !      syssec = timar(2)
2028       tcpu=seconds - stime
2029 !***************************
2030 #endif
2031
2032
2033 #ifdef CRAY
2034 !***************************
2035 ! Next definitions for Cray
2036 !     call date(curdat)
2037 !     curdat=curdat(1:9)
2038 !     call clock(curtim)
2039 !     curtim=curtim(1:8)
2040       cpusec = second()
2041       tcpu=cpusec - stime
2042 !***************************
2043 #endif
2044 #ifdef AIX
2045 !***************************
2046 ! Next definitions for RS6000
2047        integer(kind=4) :: i1,mclock
2048        i1 = mclock()
2049        tcpu = (i1+0.0D0)/100.0D0
2050 #endif
2051 #ifdef WINPGI
2052 !***************************
2053 ! next definitions for windows NT Digital fortran
2054        real(kind=4) :: time_real
2055        call cpu_time(time_real)
2056        tcpu = time_real
2057 #endif
2058 #ifdef WINIFL
2059 !***************************
2060 ! next definitions for windows NT Digital fortran
2061        real(kind=4) :: time_real
2062        call cpu_time(time_real)
2063        tcpu = time_real
2064 #endif
2065       tcpu = 0d0 !el
2066       return
2067       end function tcpu
2068 !-----------------------------------------------------------------------------
2069 #ifndef CLUSTER
2070       subroutine dajczas(rntime,hrtime,mintime,sectime)
2071
2072 !      include 'COMMON.IOUNITS'
2073       integer :: ihr,imn,isc
2074       real(kind=8) :: rntime,hrtime,mintime,sectime 
2075       hrtime=rntime/3600.0D0 
2076       hrtime=aint(hrtime)
2077       mintime=aint((rntime-3600.0D0*hrtime)/60.0D0)
2078       sectime=aint((rntime-3600.0D0*hrtime-60.0D0*mintime)+0.5D0)
2079       if (sectime.eq.60.0D0) then
2080         sectime=0.0D0
2081         mintime=mintime+1.0D0
2082       endif
2083       ihr=hrtime
2084       imn=mintime
2085       isc=sectime
2086       write (iout,328) ihr,imn,isc
2087   328 FORMAT(//'***** Computation time: ',I4  ,' hours ',I2  ,&
2088                ' minutes ', I2  ,' seconds *****')       
2089       return
2090       end subroutine dajczas
2091 !-----------------------------------------------------------------------------
2092       subroutine print_detailed_timing
2093
2094 !el      use MPI_data
2095 !      implicit real*8 (a-h,o-z)
2096 !      include 'DIMENSIONS'
2097 #ifdef MPI
2098       include 'mpif.h'
2099 #endif
2100 !      include 'COMMON.IOUNITS'
2101 !      include 'COMMON.TIME1'
2102 !      include 'COMMON.SETUP'
2103       real(kind=8) :: time1,time_barrier
2104       time_barrier = 0.0d0
2105 #ifdef MPI !el
2106       time1=MPI_WTIME()
2107 #endif !el
2108          write (iout,'(80(1h=)/a/(80(1h=)))') &
2109           "Details of FG communication time"
2110          write (*,'(7(a40,1pe15.5/),40(1h-)/a40,1pe15.5/80(1h=))') &
2111           "BROADCAST:",time_bcast,"REDUCE:",time_reduce,&
2112           "GATHER:",time_gather,&
2113           "SCATTER:",time_scatter,"SENDRECV:",time_sendrecv,&
2114           "BARRIER ene",time_barrier_e,&
2115           "BARRIER grad",time_barrier_g,&
2116           "TOTAL:",&
2117           time_bcast+time_reduce+time_gather+time_scatter+time_sendrecv
2118          write (*,*) fg_rank,myrank,&
2119            ': Total wall clock time',time1-walltime,' sec'
2120          write (*,*) "Processor",fg_rank,myrank,&
2121            ": BROADCAST time",time_bcast," REDUCE time",&
2122             time_reduce," GATHER time",time_gather," SCATTER time",&
2123             time_scatter,&
2124            " SCATTER fmatmult",time_scatter_fmatmult,&
2125            " SCATTER ginvmult",time_scatter_ginvmult,&
2126            " SCATTER fmat",time_scatter_fmat,&
2127            " SCATTER ginv",time_scatter_ginv,&
2128             " SENDRECV",time_sendrecv,&
2129             " BARRIER ene",time_barrier_e,&
2130             " BARRIER GRAD",time_barrier_g,&
2131             " BCAST7",time_bcast7," BCASTC",time_bcastc,&
2132             " BCASTW",time_bcastw," ALLREDUCE",time_allreduce,&
2133             " TOTAL",&
2134             time_bcast+time_reduce+time_gather+time_scatter+ &
2135             time_sendrecv+time_barrier+time_bcastc
2136 !el#endif
2137          write (*,*) "Processor",fg_rank,myrank," enecalc",time_enecalc
2138          write (*,*) "Processor",fg_rank,myrank," sumene",time_sumene
2139          write (*,*) "Processor",fg_rank,myrank," intfromcart",&
2140            time_intfcart
2141          write (*,*) "Processor",fg_rank,myrank," vecandderiv",&
2142            time_vec
2143          write (*,*) "Processor",fg_rank,myrank," setmatrices",&
2144            time_mat
2145          write (*,*) "Processor",fg_rank,myrank," ginvmult",&
2146            time_ginvmult
2147          write (*,*) "Processor",fg_rank,myrank," fricmatmult",&
2148            time_fricmatmult
2149          write (*,*) "Processor",fg_rank,myrank," inttocart",&
2150            time_inttocart
2151          write (*,*) "Processor",fg_rank,myrank," sumgradient",&
2152            time_sumgradient
2153          write (*,*) "Processor",fg_rank,myrank," intcartderiv",&
2154            time_intcartderiv
2155          if (fg_rank.eq.0) then
2156            write (*,*) "Processor",fg_rank,myrank," lagrangian",&
2157              time_lagrangian
2158            write (*,*) "Processor",fg_rank,myrank," cartgrad",&
2159              time_cartgrad
2160          endif
2161       return
2162       end subroutine print_detailed_timing
2163 #endif
2164 !-----------------------------------------------------------------------------
2165 !-----------------------------------------------------------------------------
2166       end module control