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