correct in readpdb and in fcopy allocate
[unres4.git] / source / unres / REMD.f90
1       module REMD
2 !-----------------------------------------------------------------------------
3       use io_units
4       use MD_data
5       use REMD_data
6       use muca_md
7
8       implicit none
9 !-----------------------------------------------------------------------------
10 !
11 !
12 !-----------------------------------------------------------------------------
13       contains
14 !-----------------------------------------------------------------------------
15 ! lagrangian_lesyng.F
16 !-----------------------------------------------------------------------------
17       subroutine lagrangian
18 !-------------------------------------------------------------------------       
19 !  This subroutine contains the total lagrangain from which the accelerations
20 !  are obtained.  For numerical gradient checking, the derivetive of the     
21 !  lagrangian in the velocities and coordinates are calculated seperately      
22 !-------------------------------------------------------------------------
23 !       implicit real*8 (a-h,o-z)
24 !       include 'DIMENSIONS'
25       use comm_cipiszcze
26       use energy_data
27       use geometry_data, only: nres
28       use control_data    !el, only: mucadyn,lmuca
29 #ifdef MPI
30        include 'mpif.h'
31       real(kind=8) :: time00
32 #endif
33 !       include 'COMMON.VAR'
34 !       include 'COMMON.CHAIN'
35 !       include 'COMMON.DERIV'
36 !       include 'COMMON.GEO'
37 !       include 'COMMON.LOCAL'
38 !       include 'COMMON.INTERACT'
39 !       include 'COMMON.MD'
40 !       include 'COMMON.IOUNITS'
41 !       include 'COMMON.CONTROL'
42 !       include 'COMMON.MUCA'
43 !      include 'COMMON.TIME1'
44        integer ::i,j,ind,itime,mnum
45        real(kind=8) :: zapas(6*nres) !,muca_factor      !maxres6=6*maxres
46        real(kind=8) :: rs(dimen),xsolv(dimen)
47 #ifdef CHECK5DSOL
48        real(kind=8) :: rscheck(dimen),rsold(dimen)
49 #endif
50        logical :: lprn = .false.
51 !el       common /cipiszcze/ itime
52        itime = itt_comm
53
54 #ifdef TIMING
55        time00=MPI_Wtime()
56 #endif
57 #ifdef FIVEDIAG
58       if (lprn) then
59         write (iout,*) "Potential forces backbone"
60         do i=nnt,nct
61           write (iout,'(i5,3e15.5,5x,3e15.5)')i,(-gcart(j,i),j=1,3)
62         enddo
63         write (iout,*) "Potential forces sidechain"
64         do i=nnt,nct
65           mnum=molnum(i)
66           if (itype(i,1).ne.10 .and. itype(i,1).ne.ntyp1.and.mnum.eq.5)&
67             write (iout,'(i5,3e15.5,5x,3e15.5)') i,(-gxcart(j,i),j=1,3)
68         enddo
69       endif
70        do j=1,3
71          ind=1
72          do i=nnt,nct
73            mnum=molnum(i)     
74            if (itype(i,1).eq.10 .or. itype(i,1).eq.ntyp1 .or. mnum.eq.5)then
75             rs(ind)=-gcart(j,i)-gxcart(j,i)
76             ind=ind+1
77            else
78             rs(ind)=-gcart(j,i)
79             rs(ind+1)=-gxcart(j,i)
80             ind=ind+2
81            end if 
82          enddo
83 #ifdef CHECK5DSOL
84          rsold=rs 
85 #endif
86          if (lprn) then
87            write(iout,*) "RHS of the 5-diag equations system"
88            do i=1,dimen
89              write(iout,*) i,rs(i)
90            enddo
91          endif
92          
93          call FDISYS (dimen,DM,DU1,DU2,rs,xsolv)
94          if (lprn) then
95            write (iout,*) "Solution of the 5-diagonal equations system"
96            do i=1,dimen
97              write (iout,'(i5,f10.5)') i,xsolv(i)
98            enddo
99          endif
100 #ifdef CHECK5DSOL
101 ! Check the solution
102          rscheck(1)=DMorig(1)*xsolv(1)+DU1orig(1)*xsolv(2)+&
103            DU2orig(1)*xsolv(3) 
104          rscheck(2)=DU1orig(1)*xsolv(1)+DMorig(2)*xsolv(2)+&
105            DU1orig(2)*xsolv(3)+DU2orig(2)*xsolv(4)
106          
107          do i=3,dimen-2
108           rscheck(i)=DU2orig(i-2)*xsolv(i-2)+DU1orig(i-1)*&
109             xsolv(i-1)+DMorig(i)*xsolv(i)+DU1orig(i)*&
110              xsolv(i+1)+DU2orig(i)*xsolv(i+2)
111          enddo
112        rscheck(dimen-1)=DU2orig(dimen-3)*xsolv(dimen-3)+&
113          DU1orig(dimen-2)*xsolv(dimen-2)+DMorig(dimen-1)*&
114           xsolv(dimen-1)+DU1orig(dimen-1)*&
115            xsolv(dimen)
116        rscheck(dimen)=DU2orig(dimen-2)*xsolv(dimen-2)+DU1orig(dimen-1)*&
117           xsolv(dimen-1)+DMorig(dimen)*xsolv(dimen)
118          
119          do i=1,dimen
120             write(iout,*) "i",i,"rsold",rsold(i),"rscheck",rscheck(i),&
121             "ratio",rscheck(i)/rsold(i)
122            enddo
123 ! end check
124 #endif
125          ind=1
126          do i=nnt,nct
127           mnum=molnum(i)
128           if (itype(i,1).eq.10 .or. itype(i,mnum).eq.ntyp1_molec(mnum))then 
129            d_a(j,i)=xsolv(ind)
130            ind=ind+1
131           else
132            d_a(j,i)=xsolv(ind)
133            d_a(j,i+nres)=xsolv(ind+1)
134           ind=ind+2
135           end if 
136          enddo
137        enddo
138        if (lprn) then
139        do i=nnt,nct
140         write (iout,'(i3,3f10.5)') i,(d_a(j,i),j=1,3)
141        enddo
142        do i=nnt,nct
143         write (iout,'(i3,3f10.5)') i,(d_a(j,i+nres),j=1,3)
144        enddo
145        endif
146        do j=1,3
147          d_a(j,0)=d_a(j,nnt)
148        enddo
149        do i=nnt,nct
150           mnum=molnum(i)
151 !         if (itype(i,1).eq.10 .or. itype(i,1).eq.ntyp1) then
152           if (itype(i,1).eq.10 .or. itype(i,mnum).eq.ntyp1_molec(mnum))then
153            do j=1,3
154              d_a(j,i)=d_a(j,i+1)-d_a(j,i)
155            enddo
156          else 
157            do j=1,3
158              d_a(j,i+nres)=d_a(j,i+nres)-d_a(j,i)
159              d_a(j,i)=d_a(j,i+1)-d_a(j,i)
160            enddo
161        
162
163          end if
164        enddo
165
166       if (lprn) then
167         write(iout,*) 'acceleration 3D FIVEDIAG'
168         write (iout,'(i3,3f10.5,3x,3f10.5)') 0,(d_a(j,0),j=1,3)
169         do i=nnt,nct-1
170          write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3)
171         enddo
172         do i=nnt,nct
173          write (iout,'(i3,3f10.5,3x,3f10.5)') &
174            i+nres,(d_a(j,i+nres),j=1,3)
175         enddo
176       endif
177 #else
178 ! Old procedure
179        do j=1,3
180          zapas(j)=-gcart(j,0)
181        enddo
182       ind=3      
183       if (lprn) then
184         write (iout,*) "Potential forces backbone"
185       endif
186       do i=nnt,nct-1
187         if (lprn) write (iout,'(i5,3e15.5,5x,3e15.5)') &
188           i,(-gcart(j,i),j=1,3)
189         do j=1,3
190           ind=ind+1
191           zapas(ind)=-gcart(j,i)
192         enddo
193       enddo
194       if (lprn) write (iout,*) "Potential forces sidechain"
195       do i=nnt,nct
196           mnum=molnum(i)
197         if (itype(i,1).ne.10 .and. itype(i,mnum).ne.ntyp1_molec(mnum)&
198           .and.mnum.ne.5) then
199           if (lprn) write (iout,'(i5,3e15.5,5x,3e15.5)') &
200              i,(-gxcart(j,i),j=1,3)
201           do j=1,3
202             ind=ind+1
203             zapas(ind)=-gxcart(j,i)
204           enddo
205         endif
206       enddo
207
208       call ginv_mult(zapas,d_a_work)
209        
210       do j=1,3
211         d_a(j,0)=d_a_work(j)
212       enddo
213       ind=3
214       do i=nnt,nct-1
215         do j=1,3
216           ind=ind+1
217           d_a(j,i)=d_a_work(ind)
218         enddo
219       enddo
220       do i=nnt,nct
221        mnum=molnum(i)
222 !        if (itype(i,1).ne.10 .and. itype(i,1).ne.ntyp1) then
223          if (itype(i,1).ne.10 .and. itype(i,mnum).ne.ntyp1_molec(mnum)&
224           .and.mnum.ne.5) then
225           do j=1,3
226             ind=ind+1
227             d_a(j,i+nres)=d_a_work(ind)
228           enddo
229         endif
230       enddo
231
232 #endif      
233       if(lmuca) then
234        imtime=imtime+1
235        if(mucadyn.gt.0) call muca_update(potE)       
236        factor=muca_factor(potE)*t_bath*Rb
237
238 !d       print *,'lmuca ',factor,potE
239        do j=1,3
240           d_a(j,0)=d_a(j,0)*factor
241        enddo
242        do i=nnt,nct-1
243          do j=1,3
244           d_a(j,i)=d_a(j,i)*factor              
245          enddo
246        enddo
247        do i=nnt,nct
248          do j=1,3
249           d_a(j,i+nres)=d_a(j,i+nres)*factor              
250          enddo
251        enddo       
252        
253       endif
254       
255       if (lprn) then
256         write(iout,*) 'acceleration 3D'
257         write (iout,'(i3,3f10.5,3x,3f10.5)') 0,(d_a(j,0),j=1,3)
258         do i=nnt,nct-1
259          write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3)
260         enddo
261         do i=nnt,nct
262          write (iout,'(i3,3f10.5,3x,3f10.5)') &
263            i+nres,(d_a(j,i+nres),j=1,3)
264         enddo
265       endif
266 #ifdef TIMING
267       time_lagrangian=time_lagrangian+MPI_Wtime()-time00
268 #endif
269       return
270       end subroutine lagrangian
271 !-----------------------------------------------------------------------------
272       subroutine setup_MD_matrices
273
274       use geometry_data, only: nres,nside
275       use control_data
276       use MPI_data
277       use energy_data
278       use geometry, only:int_bounds
279       use md_calc
280 !      implicit real*8 (a-h,o-z)
281 !      include 'DIMENSIONS'
282 #ifdef MPI
283       include 'mpif.h'
284       integer :: ierror
285       real(kind=8) :: time00
286       real(kind=8) ,allocatable, dimension(:)  :: DDM,DDU1,DDU2
287 #endif
288 !      include 'COMMON.SETUP'
289 !      include 'COMMON.VAR'
290 !      include 'COMMON.CHAIN'
291 !      include 'COMMON.DERIV'
292 !      include 'COMMON.GEO'
293 !      include 'COMMON.LOCAL'
294 !      include 'COMMON.INTERACT'
295 !      include 'COMMON.MD'
296 !#ifndef LANG0
297 !      include 'COMMON.LANGEVIN'
298 !#else
299 !      include 'COMMON.LANGEVIN.lang0'
300 !#endif
301 !      include 'COMMON.IOUNITS'
302 !      include 'COMMON.TIME1'
303       logical :: lprn = .false.
304       logical :: osob
305 #ifndef FIVEDIAG
306       real(kind=8),allocatable,dimension(:,:) :: Bmat,matmult
307 #endif
308       real(kind=8) :: dtdi
309       real(kind=8),dimension(2*nres) :: massvec,sqreig  !(maxres2) maxres2=2*maxres
310       real(kind=8) :: relfeh,eps1,eps2
311 !el      real(kind=8),dimension(:),allocatable :: Ghalf
312 !el      real(kind=8),dimension(2*nres*(2*nres+1)/2) :: Ghalf   !(mmaxres2) (mmaxres2=(maxres2*(maxres2+1)/2))
313 !el      real(kind=8),dimension(2*nres,2*nres) :: Gcopy !(maxres2,maxres2)
314 !el      real(kind=8),dimension(:,:),allocatable :: Gcopy
315       real(kind=8),dimension(8*6*nres) :: work  !(8*maxres6)
316       integer,dimension(6*nres) :: iwork        !(maxres6) maxres6=6*maxres
317 !el      common /przechowalnia/ Gcopy,Ghalf
318       real(kind=8) :: coeff,mscab
319       integer :: i,j,ind,ind1,k,l,ii,jj,m,m1,ii1,iti,nres2,ierr,nind,mark
320       real(kind=8) :: ip4
321       integer :: iz,mnum
322       print *,"just entered"
323       relfeh=1.0d-14
324       nres2=2*nres
325       print *,"FIVEDIAG"
326       write (iout,*) "before FIVEDIAG"
327 #ifndef FIVEDIAG
328       write (iout,*) "ALLOCATE"
329       print *,"ALLOCATE"
330       if(.not.allocated(Gcopy)) allocate(Gcopy(nres2,nres2)) !(maxres2,maxres2)
331       if(.not.allocated(Ghalf)) allocate(Ghalf(nres2*(nres2+1)/2)) !mmaxres2=(maxres2*(maxres+1)/2)
332       if(.not.allocated(Bmat))  allocate(Bmat(nres2,nres2))
333       if(.not.allocated(matmult)) allocate(matmult(nres2,nres2))
334 #endif
335 !
336 ! Set up the matrix of the (dC,dX)-->(C,X) transformation (A), the
337 ! inertia matrix (Gmat) and the inverse of the inertia matrix (Ginv)
338 !
339 ! Determine the number of degrees of freedom (dimen) and the number of 
340 ! sites (dimen1)
341       dimen=(nct-nnt+1)+nside
342       dimen1=(nct-nnt)+(nct-nnt+1)
343       dimen3=dimen*3
344       write (iout,*) "nnt",nnt," nct",nct," nside",nside
345 #ifdef FIVEDIAG
346        ip4=ip/4
347        DM(1)=mp/4+ip4
348        if (iabs(itype(nnt)).eq.10) then
349          DM(1)=DM(1)+msc(10)
350          ind=2
351          nind=1
352        else
353          DM(1)=DM(1)+isc(iabs(itype(nnt)))
354          DM(2)=msc(iabs(itype(nnt)))+isc(iabs(itype(nnt)))
355          ind=3
356          nind=2
357        endif      
358        do i=nnt+1,nct-1
359 !         if (iabs(itype(i,1)).eq.ntyp1) cycle
360          DM(ind)=2*ip4+mp/2
361          if (iabs(itype(i,1)).eq.10 .or. iabs(itype(i,1)).eq.ntyp1) then
362            if (iabs(itype(i,1)).eq.10) DM(ind)=DM(ind)+msc(10)
363            ind=ind+1
364          else
365            DM(ind)=DM(ind)+isc(iabs(itype(i,1)))
366            DM(ind+1)=msc(iabs(itype(i,1)))+isc(iabs(itype(i,1)))
367            ind=ind+2
368          endif 
369        enddo
370        if (nct.gt.nnt) then
371        DM(ind)=ip4+mp/4
372        if (iabs(itype(nct)).eq.10) then
373          DM(ind)=DM(ind)+msc(10)
374          nind=ind
375        else
376          DM(ind)=DM(ind)+isc(iabs(itype(nct)))
377          DM(ind+1)=msc(iabs(itype(nct)))+isc(iabs(itype(nct)))
378          nind=ind+1
379         endif
380         endif
381         
382         
383         ind=1
384         do i=nnt,nct
385         mnum=molnum(i)
386         if (iabs(itype(i,1)).ne.10 .and.iabs(itype((i))).ne.ntyp1 &
387           .and.mnum.eq.5) then
388           DU1(ind)=-isc(iabs(itype(i,1)))
389           DU1(ind+1)=0.0d0
390           ind=ind+2
391         else
392           DU1(ind)=mp/4-ip4
393           ind=ind+1
394         endif
395         enddo
396
397        ind=1
398        do i=nnt,nct-1
399        mnum=molnum(i)
400 !        if (iabs(itype(i,1)).eq.ntyp1) cycle
401         write (iout,*) "i",i," itype",itype(i,1),ntyp1
402         if (iabs(itype(i,1)).ne.10 .and. &
403           iabs(itype(i,mnum)).ne.ntyp1_molec(mnum) .and. mnum.ne.5) then
404          DU2(ind)=mp/4-ip4
405          DU2(ind+1)=0.0d0
406          ind=ind+2
407         else
408          DU2(ind)=0.0d0
409          DU2(ind+1)=0.0d0
410          ind=ind+1
411         endif
412        enddo
413        DMorig=DM
414        DU1orig=DU1
415        DU2orig=DU2
416        write (iout,*) "nind",nind," dimen",dimen
417        if (nind.ne.dimen) then
418          write (iout,*) "ERROR, dimensions differ"
419 #ifdef MPI
420          call MPI_Finalize(ierr)
421 #endif
422          stop
423        endif
424      write (iout,*) "The upper part of the five-diagonal inertia matrix"
425        do i=1,dimen
426          if (i.lt.dimen-1) then
427            write (iout,'(i3,3f10.5)') i,DM(i),DU1(i),DU2(i)
428          else if (i.eq.dimen-1) then  
429            write (iout,'(i3,3f10.5)') i,DM(i),DU1(i)
430          else
431            write (iout,'(i3,3f10.5)') i,DM(i)
432          endif 
433        enddo
434
435         call FDISYP (dimen, DM, DU1, DU2, MARK)
436
437         if (mark.eq.-1) then
438      write (iout,*) "ERROR: the inertia matrix is not positive definite"
439 #ifdef MPI
440          call MPI_Finalize(ierr)
441 #endif
442          stop
443         else if (mark.eq.0) then
444      write (iout,*) "ERROR: the inertia matrix is singular"
445 #ifdef MPI
446          call MPI_Finalize(ierr)
447 #endif
448         else if (mark.eq.1) then
449           write (iout,*) "The transformed inertia matrix"
450           do i=1,dimen
451             if (i.lt.dimen-1) then
452               write (iout,'(i3,3f10.5)') i,DM(i),DU1(i),DU2(i)
453             else if (i.eq.dimen-1) then  
454               write (iout,'(i3,3f10.5)') i,DM(i),DU1(i)
455             else
456               write (iout,'(i3,3f10.5)') i,DM(i)
457             endif 
458           enddo
459         endif
460 ! Diagonalization of the pentadiagonal matrix
461         allocate(DDU1(2*nres))
462         allocate(DDU2(2*nres))
463         allocate(DDM(2*nres))
464         do i=1,dimen-1
465           DDU1(i+1)=DU1orig(i)
466         enddo
467         do i=1,dimen-2
468           DDU2(i+2)=DU2orig(i)
469         enddo
470         DDM=DMorig
471         call quindibisect2(DDM,DDU1,DDU2,dimen,1,dimen,eps1,relfeh,eps2,iz,geigen)
472         if (lprn) then
473           write (iout,*) &
474            "Eigenvalues of the five-diagonal inertia matrix"
475           do i=1,dimen
476             write (iout,'(i5,f10.5)') i,geigen(i)
477           enddo
478         endif
479         deallocate(DDU1)
480         deallocate(DDU2)
481         deallocate(DDM)
482 #else 
483 ! Old Gmatrix
484 #ifdef MPI
485       if (nfgtasks.gt.1) then
486       time00=MPI_Wtime()
487       call MPI_Bcast(5,1,MPI_INTEGER,king,FG_COMM,IERROR)
488       time_Bcast=time_Bcast+MPI_Wtime()-time00
489       call int_bounds(dimen,igmult_start,igmult_end)
490       igmult_start=igmult_start-1
491       call MPI_Allgather(3*igmult_start,1,MPI_INTEGER,&
492           ng_start(0),1,MPI_INTEGER,FG_COMM,IERROR)
493       my_ng_count=igmult_end-igmult_start
494       call MPI_Allgather(3*my_ng_count,1,MPI_INTEGER,ng_counts(0),1,&
495           MPI_INTEGER,FG_COMM,IERROR)
496       write (iout,*) 'Processor:',fg_rank,' CG group',kolor,&
497        ' absolute rank',myrank,' igmult_start',igmult_start,&
498        ' igmult_end',igmult_end,' count',my_ng_count
499       write (iout,*) "ng_start",(ng_start(i),i=0,nfgtasks-1)
500       write (iout,*) "ng_counts",(ng_counts(i),i=0,nfgtasks-1)
501       call flush(iout)
502       else
503 #endif
504       igmult_start=1
505       igmult_end=dimen
506       my_ng_count=dimen
507 #ifdef MPI
508       endif
509 #endif
510 !      write (iout,*) "dimen",dimen," dimen1",dimen1," dimen3",dimen3
511 !  Zeroing out A and fricmat
512       do i=1,dimen
513         do j=1,dimen
514           A(i,j)=0.0D0     
515         enddo    
516       enddo
517 !  Diagonal elements of the dC part of A and the respective friction coefficients
518       ind=1
519       ind1=0
520       print *,"TUTUTUT?!",nnt,nct-1
521       do i=nnt,nct-1
522         mnum=molnum(i)
523         ind=ind+1
524         ind1=ind1+1
525         coeff=0.25d0*IP(mnum)
526         if (mnum.eq.5) mp(mnum)=msc(itype(i,mnum),mnum)
527         print *,i,coeff,mp(mnum)
528         massvec(ind1)=mp(mnum)
529         Gmat(ind,ind)=coeff
530         print *,"i",mp(mnum)
531         A(ind1,ind)=0.5d0
532       enddo
533       
534 !  Off-diagonal elements of the dC part of A 
535       k=3
536       do i=1,nct-nnt
537         do j=1,i
538           A(i,j)=1.0d0
539         enddo
540       enddo
541 !  Diagonal elements of the dX part of A and the respective friction coefficients
542       m=nct-nnt
543       m1=nct-nnt+1
544       ind=0
545       ind1=0
546       do i=1,2
547       msc(ntyp1_molec(i),i)=1.0d0
548       enddo
549       do i=nnt,nct
550         ind=ind+1
551         ii = ind+m
552         mnum=molnum(i)
553         iti=itype(i,mnum)
554         if (mnum.eq.5) then
555         mscab=0.0
556         else
557         mscab=msc(iabs(iti),mnum)
558         endif
559         massvec(ii)=mscab
560
561         if (iti.ne.10 .and. iti.ne.ntyp1_molec(mnum) .and. mnum.ne.5) then
562           ind1=ind1+1
563           ii1= ind1+m1
564           A(ii,ii1)=1.0d0
565           Gmat(ii1,ii1)=ISC(iabs(iti),mnum)
566         endif
567       enddo
568 !  Off-diagonal elements of the dX part of A
569       ind=0
570       k=nct-nnt
571       do i=nnt,nct
572         iti=itype(i,1)
573         ind=ind+1
574         do j=nnt,i
575           ii = ind
576           jj = j-nnt+1
577           A(k+ii,jj)=1.0d0
578         enddo
579       enddo
580       if (lprn) then
581         write (iout,*)
582         write (iout,*) "Vector massvec"
583         do i=1,dimen1
584           write (iout,*) i,massvec(i)
585         enddo
586         write (iout,'(//a)') "A"
587         call matout(dimen,dimen1,nres2,nres2,A)
588       endif
589
590 ! Calculate the G matrix (store in Gmat)
591       do k=1,dimen
592        do i=1,dimen
593          dtdi=0.0d0
594          do j=1,dimen1
595            dtdi=dtdi+A(j,k)*A(j,i)*massvec(j)
596          enddo
597          Gmat(k,i)=Gmat(k,i)+dtdi
598        enddo
599       enddo 
600       
601       if (lprn) then
602         write (iout,'(//a)') "Gmat"
603         call matout(dimen,dimen,nres2,nres2,Gmat)
604       endif
605       do i=1,dimen
606         do j=1,dimen
607           Ginv(i,j)=0.0d0
608           Gcopy(i,j)=Gmat(i,j)
609         enddo
610         Ginv(i,i)=1.0d0
611       enddo
612 ! Invert the G matrix
613       call MATINVERT(dimen,nres2,Gcopy,Ginv,osob)
614       if (lprn) then
615         write (iout,'(//a)') "Ginv"
616         call matout(dimen,dimen,nres2,nres2,Ginv)
617       endif
618
619       Bmat=0.0d0
620       Bmat(1,1)=1.0d0
621       j=2
622       do i=nnt,nct
623         mnum=molnum(i)
624         if ((itype(i,1).eq.10).or.(itype(i,mnum).eq.ntyp1_molec(mnum))&
625           .or.(mnum.eq.5)) then
626           Bmat(i-nnt+2,j-1)=-1
627           Bmat(i-nnt+2,j)=1
628           j=j+1
629         else
630           if (i.lt.nct) then
631             Bmat(i-nnt+2,j-1)=-1
632             Bmat(i-nnt+2,j+1)=1
633             Bmat(i-nnt+1+(nct-nnt)+1,j-1)=-1
634             Bmat(i-nnt+1+(nct-nnt)+1,j)=1
635             j=j+2
636           else
637             Bmat(i-nnt+1+(nct-nnt)+1,j-1)=-1
638             Bmat(i-nnt+1+(nct-nnt)+1,j)=1
639           endif
640         endif
641       enddo
642       write (iout,*) "j",j," dimen",dimen 
643       if (lprn) then
644         write (iout,'(//a)') "Bmat"
645         call matout(dimen,dimen,nres2,nres2,Bmat)
646       endif
647       Gcopy=0.0d0
648       write(iout,*) "before Gcopy",dimen,nres*2
649 #ifdef TESTFIVEDIAG
650       do i=1,dimen
651         do j=1,dimen
652 !            write (iout,*) "myij",i,j
653             do k=1,dimen
654                do l=1,dimen
655 !                 write(iout,*) "i,j,k,l",i,j,k,l
656                  Gcopy(i,j)=Gcopy(i,j)+Bmat(k,i)*Gmat(k,l)*Bmat(l,j)
657                enddo
658             enddo
659         enddo
660       enddo
661 #endif
662       if (lprn) then
663         write (iout,'(//a)') "Gmat-transformed"
664         call matout(dimen,dimen,nres2,nres2,Gcopy)
665       endif
666 #ifdef MPI
667       if (nfgtasks.gt.1) then
668         myginv_ng_count=nres2*my_ng_count
669         call MPI_Allgather(nres2*igmult_start,1,MPI_INTEGER,&
670           nginv_start(0),1,MPI_INTEGER,FG_COMM,IERROR)
671         call MPI_Allgather(myginv_ng_count,1,MPI_INTEGER,&
672           nginv_counts(0),1,MPI_INTEGER,FG_COMM,IERROR)
673         write (iout,*) "nginv_start",(nginv_start(i),i=0,nfgtasks-1)
674         write (iout,*) "nginv_counts",(nginv_counts(i),i=0,nfgtasks-1)
675         call flush(iout)
676 !        call MPI_Scatterv(ginv(1,1),nginv_counts(0),
677 !     &    nginv_start(0),MPI_DOUBLE_PRECISION,ginv,
678 !     &    myginv_ng_count,MPI_DOUBLE_PRECISION,king,FG_COMM,IERR)
679 !        call MPI_Barrier(FG_COMM,IERR)
680         time00=MPI_Wtime()
681         call MPI_Scatterv(ginv(1,1),nginv_counts(0),&
682           nginv_start(0),MPI_DOUBLE_PRECISION,gcopy(1,1),&
683           myginv_ng_count,MPI_DOUBLE_PRECISION,king,FG_COMM,IERR)
684 #ifdef TIMING
685         time_scatter_ginv=time_scatter_ginv+MPI_Wtime()-time00
686 #endif
687         do i=1,dimen
688           do j=1,2*my_ng_count
689             ginv(j,i)=gcopy(i,j)
690           enddo
691         enddo
692 !        write (iout,*) "Master's chunk of ginv"
693 !        call MATOUT2(my_ng_count,dimen,maxres2,maxres2,ginv)
694       endif
695 #endif
696       if (osob) then
697         write (iout,*) "The G matrix is singular."
698         stop
699       endif
700 ! Compute G**(-1/2) and G**(1/2) 
701       ind=0
702       do i=1,dimen
703         do j=1,i
704           ind=ind+1
705           Ghalf(ind)=Gmat(i,j)
706         enddo
707       enddo
708       call gldiag(nres2,dimen,dimen,Ghalf,work,Geigen,Gvec,&
709         ierr,iwork)
710       if (lprn) then
711         write (iout,'(//a)') &
712          "Eigenvectors and eigenvalues of the G matrix"
713         call eigout(dimen,dimen,nres2,nres2,Gvec,Geigen)
714       endif
715       do i=1,dimen
716         sqreig(i)=dsqrt(Geigen(i))
717       enddo
718       do i=1,dimen
719         do j=1,dimen
720           Gsqrp(i,j)=0.0d0
721           Gsqrm(i,j)=0.0d0
722           Gcopy(i,j)=0.0d0
723           do k=1,dimen
724             Gsqrp(i,j)=Gsqrp(i,j)+Gvec(i,k)*Gvec(j,k)*sqreig(k)
725             Gsqrm(i,j)=Gsqrm(i,j)+Gvec(i,k)*Gvec(j,k)/sqreig(k)
726             Gcopy(i,j)=Gcopy(i,j)+Gvec(i,k)*Gvec(j,k)*Geigen(k)
727           enddo
728         enddo
729       enddo
730       if (lprn) then
731         write (iout,*) "Comparison of original and restored G"
732         do i=1,dimen
733           do j=1,dimen
734             write (iout,'(2i5,5f10.5)') i,j,Gmat(i,j),Gcopy(i,j),&
735               Gmat(i,j)-Gcopy(i,j),Gsqrp(i,j),Gsqrm(i,j)
736           enddo
737         enddo
738       endif
739       deallocate(Gcopy)
740 #endif
741       return
742       end subroutine setup_MD_matrices
743 !-----------------------------------------------------------------------------
744       subroutine EIGOUT(NC,NR,LM2,LM3,A,B)
745
746 !      implicit real*8 (a-h,o-z)
747 !      include 'DIMENSIONS'
748 !      include 'COMMON.IOUNITS'
749       integer :: LM2,LM3,NC,NR,KA,KC,KB,I,J,N
750       real(kind=8) :: A(LM2,LM3),B(LM2)
751       KA=1
752       KC=6
753     1 KB=MIN0(KC,NC)
754       WRITE(IOUT,600) (I,I=KA,KB)
755       WRITE(IOUT,601) (B(I),I=KA,KB)
756       WRITE(IOUT,602)
757     2 N=0
758       DO 3  I=1,NR
759       WRITE(IOUT,603) I,(A(I,J),J=KA,KB)
760       N=N+1
761       IF(N.LT.10) GO TO 3
762       WRITE(IOUT,602)
763       N=0
764     3 CONTINUE
765     4 IF (KB.EQ.NC) RETURN
766       KA=KC+1
767       KC=KC+6
768       GO TO 1
769   600 FORMAT (// 9H ROOT NO.,I4,9I11)
770   601 FORMAT (/5X,10(1PE11.4))
771   602 FORMAT (2H  )
772   603 FORMAT (I5,10F11.5)
773   604 FORMAT (1H1)
774       end subroutine EIGOUT
775 !-----------------------------------------------------------------------------
776       subroutine MATOUT(NC,NR,LM2,LM3,A)
777
778 !      implicit real*8 (a-h,o-z)
779 !      include 'DIMENSIONS'
780 !      include 'COMMON.IOUNITS'
781       integer :: LM2,LM3,NC,NR,KA,KC,KB,I,J,N
782       real(kind=8) :: A(LM2,LM3)
783       KA=1
784       KC=6
785     1 KB=MIN0(KC,NC)
786       WRITE(IOUT,600) (I,I=KA,KB)
787       WRITE(IOUT,602)
788     2 N=0
789       DO 3  I=1,NR
790       WRITE(IOUT,603) I,(A(I,J),J=KA,KB)
791       N=N+1
792       IF(N.LT.10) GO TO 3
793       WRITE(IOUT,602)
794       N=0
795     3 CONTINUE
796     4 IF (KB.EQ.NC) RETURN
797       KA=KC+1
798       KC=KC+6
799       GO TO 1
800   600 FORMAT (//5x,9I11)
801   602 FORMAT (2H  )
802   603 FORMAT (I5,10F11.3)
803   604 FORMAT (1H1)
804       end subroutine MATOUT
805 !-----------------------------------------------------------------------------
806       subroutine MATOUT1(NC,NR,LM2,LM3,A)
807
808 !      implicit real*8 (a-h,o-z)
809 !      include 'DIMENSIONS'
810 !      include 'COMMON.IOUNITS'
811       integer :: LM2,LM3,NC,NR,KA,KC,KB,I,J,N
812       real(kind=8) :: A(LM2,LM3)
813       KA=1
814       KC=21
815     1 KB=MIN0(KC,NC)
816       WRITE(IOUT,600) (I,I=KA,KB)
817       WRITE(IOUT,602)
818     2 N=0
819       DO 3  I=1,NR
820       WRITE(IOUT,603) I,(A(I,J),J=KA,KB)
821       N=N+1
822       IF(N.LT.3) GO TO 3
823       WRITE(IOUT,602)
824       N=0
825     3 CONTINUE
826     4 IF (KB.EQ.NC) RETURN
827       KA=KC+1
828       KC=KC+21
829       GO TO 1
830   600 FORMAT (//5x,7(3I5,2x))
831   602 FORMAT (2H  )
832   603 FORMAT (I5,7(3F5.1,2x))
833   604 FORMAT (1H1)
834       end subroutine MATOUT1
835 !-----------------------------------------------------------------------------
836       subroutine MATOUT2(NC,NR,LM2,LM3,A)
837
838 !      implicit real*8 (a-h,o-z)
839 !      include 'DIMENSIONS'
840 !      include 'COMMON.IOUNITS'
841       integer :: I,J,KA,KC,KB,N
842       integer :: LM2,LM3,NC,NR
843       real(kind=8) :: A(LM2,LM3)
844       KA=1
845       KC=12
846     1 KB=MIN0(KC,NC)
847       WRITE(IOUT,600) (I,I=KA,KB)
848       WRITE(IOUT,602)
849     2 N=0
850       DO 3  I=1,NR
851       WRITE(IOUT,603) I,(A(I,J),J=KA,KB)
852       N=N+1
853       IF(N.LT.3) GO TO 3
854       WRITE(IOUT,602)
855       N=0
856     3 CONTINUE
857     4 IF (KB.EQ.NC) RETURN
858       KA=KC+1
859       KC=KC+12
860       GO TO 1
861   600 FORMAT (//5x,4(3I9,2x))
862   602 FORMAT (2H  )
863   603 FORMAT (I5,4(3F9.3,2x))
864   604 FORMAT (1H1)
865       end subroutine MATOUT2
866 !-----------------------------------------------------------------------------
867       subroutine ginv_mult(z,d_a_tmp)
868
869       use geometry_data, only: nres
870       use control_data
871       use MPI_data
872 !      implicit real*8 (a-h,o-z)
873 !      include 'DIMENSIONS'
874 #ifdef MPI
875       include 'mpif.h'
876       integer :: ierr,ierror
877 #endif
878 !      include 'COMMON.SETUP'
879 !      include 'COMMON.TIME1'
880 !      include 'COMMON.MD'
881       real(kind=8),dimension(dimen3) :: z,z1,d_a_tmp
882       real(kind=8),dimension(6*nres) :: temp    !(maxres6) maxres6=6*maxres
883       real(kind=8) :: time00,time01
884       integer :: i,j,k,ind
885 #ifdef MPI
886       if (nfgtasks.gt.1) then
887         if (fg_rank.eq.0) then
888 ! The matching BROADCAST for fg processors is called in ERGASTULUM
889           time00=MPI_Wtime()
890           call MPI_Bcast(4,1,MPI_INTEGER,king,FG_COMM,IERROR)
891           time_Bcast=time_Bcast+MPI_Wtime()-time00
892 !          print *,"Processor",myrank," BROADCAST iorder in GINV_MULT"
893         endif
894 !        write (2,*) "time00",time00
895 !        write (2,*) "Before Scatterv"
896 !        call flush(2)
897 !        write (2,*) "Whole z (for FG master)"
898 !        do i=1,dimen
899 !          write (2,*) i,z(i)
900 !        enddo
901 !        call MPI_Barrier(FG_COMM,IERROR)
902         time00=MPI_Wtime()
903 !elwrite(iout,*) "do tej pory jest OK, MPI_Scatterv w ginv_mult"
904         call MPI_Scatterv(z,ng_counts(0),ng_start(0),&
905           MPI_DOUBLE_PRECISION,&
906           z1,3*my_ng_count,MPI_DOUBLE_PRECISION,king,FG_COMM,IERR)
907 !        write (2,*) "My chunk of z"
908 !        do i=1,3*my_ng_count
909 !          write (2,*) i,z(i)
910 !        enddo
911 !        write (2,*) "After SCATTERV"
912 !        call flush(2)
913 !        write (2,*) "MPI_Wtime",MPI_Wtime()
914         time_scatter=time_scatter+MPI_Wtime()-time00
915 #ifdef TIMING
916         time_scatter_ginvmult=time_scatter_ginvmult+MPI_Wtime()-time00
917 #endif
918 !        write (2,*) "time_scatter",time_scatter
919 !        write (2,*) "dimen",dimen," dimen3",dimen3," my_ng_count",
920 !     &    my_ng_count
921 !        call flush(2)
922         time01=MPI_Wtime()
923         do k=0,2
924           do i=1,dimen
925             ind=(i-1)*3+k+1
926             temp(ind)=0.0d0
927             do j=1,my_ng_count
928 !              write (2,*) "k,i,j,ind",k,i,j,ind,(j-1)*3+k+1,
929 !     &         Ginv(i,j),z((j-1)*3+k+1),
930 !     &          Ginv(i,j)*z((j-1)*3+k+1)
931 !              temp(ind)=temp(ind)+Ginv(i,j)*z((j-1)*3+k+1)
932               temp(ind)=temp(ind)+Ginv(j,i)*z1((j-1)*3+k+1)
933             enddo
934           enddo 
935         enddo
936         time_ginvmult=time_ginvmult+MPI_Wtime()-time01
937 !        write (2,*) "Before REDUCE"
938 !        call flush(2)
939 !        write (2,*) "z before reduce"
940 !        do i=1,dimen
941 !          write (2,*) i,temp(i)
942 !        enddo
943         time00=MPI_Wtime()
944         call MPI_Reduce(temp(1),d_a_tmp(1),dimen3,MPI_DOUBLE_PRECISION,&
945             MPI_SUM,king,FG_COMM,IERR)
946         time_reduce=time_reduce+MPI_Wtime()-time00
947 !        write (2,*) "After REDUCE"
948 !        call flush(2)
949       else
950 #endif
951 #ifdef TIMING
952         time01=MPI_Wtime()
953 #endif
954         do k=0,2
955           do i=1,dimen
956             ind=(i-1)*3+k+1
957             d_a_tmp(ind)=0.0d0
958             do j=1,dimen
959 !              write (2,*) "k,i,j,ind",k,i,j,ind,(j-1)*3+k+1
960 !              call flush(2)
961 !     &         Ginv(i,j),z((j-1)*3+k+1),
962 !     &          Ginv(i,j)*z((j-1)*3+k+1)
963               d_a_tmp(ind)=d_a_tmp(ind) &
964                                +Ginv(j,i)*z((j-1)*3+k+1)
965 !              d_a_tmp(ind)=d_a_tmp(ind)
966 !     &                         +Ginv(i,j)*z((j-1)*3+k+1)
967             enddo
968           enddo 
969         enddo
970 #ifdef TIMING
971         time_ginvmult=time_ginvmult+MPI_Wtime()-time01
972 #endif
973 #ifdef MPI
974       endif
975 #endif
976       return
977       end subroutine ginv_mult
978 !-----------------------------------------------------------------------------
979 #ifdef GINV_MULT
980       subroutine ginv_mult_test(z,d_a_tmp)
981
982 !      include 'DIMENSIONS'
983 !el      integer :: dimen
984 !      include 'COMMON.MD'
985       real(kind=8),dimension(dimen) :: z,d_a_tmp
986       real(kind=8),dimension(dimen/3) :: ztmp,dtmp
987       integer :: i,j,k,ind
988 !      do i=1,dimen
989 !        d_a_tmp(i)=0.0d0
990 !        do j=1,dimen
991 !          d_a_tmp(i)=d_a_tmp(i)+Ginv(i,j)*z(j)
992 !        enddo
993 !      enddo
994 !
995 !      return
996
997 !ibm* unroll(3)
998       do k=0,2
999        do j=1,dimen/3
1000         ztmp(j)=z((j-1)*3+k+1)
1001        enddo
1002
1003        call alignx(16,ztmp(1))
1004        call alignx(16,dtmp(1))
1005        call alignx(16,Ginv(1,1)) 
1006
1007        do i=1,dimen/3
1008         dtmp(i)=0.0d0
1009         do j=1,dimen/3
1010            dtmp(i)=dtmp(i)+Ginv(i,j)*ztmp(j)
1011         enddo
1012        enddo
1013        do i=1,dimen/3
1014         ind=(i-1)*3+k+1
1015         d_a_tmp(ind)=dtmp(i)
1016        enddo 
1017       enddo
1018       return
1019       end subroutine ginv_mult_test
1020 #endif
1021 !-----------------------------------------------------------------------------
1022       subroutine fricmat_mult(z,d_a_tmp)
1023
1024       use geometry_data, only: nres
1025       use control_data
1026       use MPI_data
1027 !      include 'DIMENSIONS'
1028 #ifdef MPI
1029       include 'mpif.h'
1030       integer :: IERROR,ierr
1031 #endif
1032 !      include 'COMMON.MD'
1033 !      include 'COMMON.IOUNITS'
1034 !      include 'COMMON.SETUP'
1035 !      include 'COMMON.TIME1'
1036 !#ifndef LANG0
1037 !      include 'COMMON.LANGEVIN'
1038 !#else
1039 !      include 'COMMON.LANGEVIN.lang0'
1040 !#endif
1041       real(kind=8),dimension(dimen3) :: z,z1,d_a_tmp
1042       real(kind=8),dimension(6*nres) :: temp    !(maxres6) maxres6=6*maxres
1043       real(kind=8) :: time00,time01
1044       integer :: i,j,k,ind,nres2
1045       nres2=2*nres
1046 !el      if(.not.allocated(fricmat)) allocate(fricmat(nres2,nres2))
1047
1048 #ifdef MPI
1049       if (nfgtasks.gt.1) then
1050         if (fg_rank.eq.0) then
1051 ! The matching BROADCAST for fg processors is called in ERGASTULUM
1052           time00=MPI_Wtime()
1053           call MPI_Bcast(9,1,MPI_INTEGER,king,FG_COMM,IERROR)
1054           time_Bcast=time_Bcast+MPI_Wtime()-time00
1055 !          print *,"Processor",myrank," BROADCAST iorder in FRICMAT_MULT"
1056         endif
1057 !        call MPI_Barrier(FG_COMM,IERROR)
1058         time00=MPI_Wtime()
1059         call MPI_Scatterv(z,ng_counts(0),ng_start(0),&
1060           MPI_DOUBLE_PRECISION,&
1061           z1,3*my_ng_count,MPI_DOUBLE_PRECISION,king,FG_COMM,IERR)
1062 !        write (2,*) "My chunk of z"
1063 !        do i=1,3*my_ng_count
1064 !          write (2,*) i,z(i)
1065 !        enddo
1066         time_scatter=time_scatter+MPI_Wtime()-time00
1067 #ifdef TIMING
1068         time_scatter_fmatmult=time_scatter_fmatmult+MPI_Wtime()-time00
1069 #endif
1070         time01=MPI_Wtime()
1071         do k=0,2
1072           do i=1,dimen
1073             ind=(i-1)*3+k+1
1074             temp(ind)=0.0d0
1075             do j=1,my_ng_count
1076               temp(ind)=temp(ind)-fricmat(j,i)*z1((j-1)*3+k+1)
1077             enddo
1078           enddo 
1079         enddo
1080         time_fricmatmult=time_fricmatmult+MPI_Wtime()-time01
1081 !        write (2,*) "Before REDUCE"
1082 !        write (2,*) "d_a_tmp before reduce"
1083 !        do i=1,dimen3
1084 !          write (2,*) i,temp(i)
1085 !        enddo
1086 !        call flush(2)
1087         time00=MPI_Wtime()
1088         call MPI_Reduce(temp(1),d_a_tmp(1),dimen3,MPI_DOUBLE_PRECISION,&
1089             MPI_SUM,king,FG_COMM,IERR)
1090         time_reduce=time_reduce+MPI_Wtime()-time00
1091 !        write (2,*) "After REDUCE"
1092 !        call flush(2)
1093       else
1094 #endif
1095 #ifdef TIMING
1096         time01=MPI_Wtime()
1097 #endif
1098         do k=0,2
1099          do i=1,dimen
1100           ind=(i-1)*3+k+1
1101           d_a_tmp(ind)=0.0d0
1102           do j=1,dimen
1103              d_a_tmp(ind)=d_a_tmp(ind) &
1104                                  -fricmat(j,i)*z((j-1)*3+k+1)
1105           enddo
1106          enddo 
1107         enddo
1108 #ifdef TIMING
1109         time_fricmatmult=time_fricmatmult+MPI_Wtime()-time01
1110 #endif
1111 #ifdef MPI
1112       endif
1113 #endif
1114 !      write (iout,*) "Vector d_a"
1115 !      do i=1,dimen3
1116 !        write (2,*) i,d_a_tmp(i)
1117 !      enddo
1118       return
1119       end subroutine fricmat_mult
1120 !-----------------------------------------------------------------------------
1121 !-----------------------------------------------------------------------------
1122       end module REMD