e2db47896619eb6a519ebd30fb2bccab65fb94e3
[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 #ifdef CRYST_BOND
347        ip4=ip/4
348        DM(1)=mp/4+ip4
349        if (iabs(itype(nnt)).eq.10) then
350          DM(1)=DM(1)+msc(10)
351          ind=2
352          nind=1
353        else
354          DM(1)=DM(1)+isc(iabs(itype(nnt)))
355          DM(2)=msc(iabs(itype(nnt)))+isc(iabs(itype(nnt)))
356          ind=3
357          nind=2
358        endif      
359        do i=nnt+1,nct-1
360 !         if (iabs(itype(i,1)).eq.ntyp1) cycle
361          DM(ind)=2*ip4+mp/2
362          if (iabs(itype(i,1)).eq.10 .or. iabs(itype(i,1)).eq.ntyp1) then
363            if (iabs(itype(i,1)).eq.10) DM(ind)=DM(ind)+msc(10)
364            ind=ind+1
365          else
366            DM(ind)=DM(ind)+isc(iabs(itype(i,1)))
367            DM(ind+1)=msc(iabs(itype(i,1)))+isc(iabs(itype(i,1)))
368            ind=ind+2
369          endif 
370        enddo
371        if (nct.gt.nnt) then
372        DM(ind)=ip4+mp/4
373        if (iabs(itype(nct)).eq.10) then
374          DM(ind)=DM(ind)+msc(10)
375          nind=ind
376        else
377          DM(ind)=DM(ind)+isc(iabs(itype(nct)))
378          DM(ind+1)=msc(iabs(itype(nct)))+isc(iabs(itype(nct)))
379          nind=ind+1
380         endif
381         endif
382         
383         
384         ind=1
385         do i=nnt,nct
386         mnum=molnum(i)
387         if (iabs(itype(i,1)).ne.10 .and.iabs(itype((i))).ne.ntyp1 &
388           .and.mnum.eq.5) then
389           DU1(ind)=-isc(iabs(itype(i,1)))
390           DU1(ind+1)=0.0d0
391           ind=ind+2
392         else
393           DU1(ind)=mp/4-ip4
394           ind=ind+1
395         endif
396         enddo
397
398        ind=1
399        do i=nnt,nct-1
400        mnum=molnum(i)
401 !        if (iabs(itype(i,1)).eq.ntyp1) cycle
402         write (iout,*) "i",i," itype",itype(i,1),ntyp1
403         if (iabs(itype(i,1)).ne.10 .and. &
404           iabs(itype(i,mnum)).ne.ntyp1_molec(mnum) .and. mnum.ne.5) then
405          DU2(ind)=mp/4-ip4
406          DU2(ind+1)=0.0d0
407          ind=ind+2
408         else
409          DU2(ind)=0.0d0
410          DU2(ind+1)=0.0d0
411          ind=ind+1
412         endif
413        enddo
414 #else
415        ip4=ip(1)/4
416        DM(1)=mp(1)/4+ip4
417        if (iabs(itype(nnt,1)).eq.10) then
418          DM(1)=DM(1)+msc(10,1)
419          ind=2
420          nind=1
421        else
422          DM(1)=DM(1)+isc(iabs(itype(nnt,1)),1)
423          DM(2)=msc(iabs(itype(nnt,1)),1)+isc(iabs(itype(nnt,1)),1)
424          ind=3
425          nind=2
426        endif     
427        do i=nnt+1,nct-1
428 !         if (iabs(itype(i,1)).eq.ntyp1) cycle
429          DM(ind)=2*ip4+mp(1)/2
430          if (iabs(itype(i,1)).eq.10 .or. iabs(itype(i,1)).eq.ntyp1) then
431            if (iabs(itype(i,1)).eq.10) DM(ind)=DM(ind)+msc(10,1)
432            ind=ind+1
433          else
434            DM(ind)=DM(ind)+isc(iabs(itype(i,1)),1)
435            DM(ind+1)=msc(iabs(itype(i,1)),1)+isc(iabs(itype(i,1)),1)
436            ind=ind+2
437          endif 
438        enddo
439        if (nct.gt.nnt) then
440        DM(ind)=ip4+mp(1)/4
441        if (iabs(itype(nct,1)).eq.10) then
442          DM(ind)=DM(ind)+msc(10,1)
443          nind=ind
444        else
445          DM(ind)=DM(ind)+isc(iabs(itype(nct,1)),1)
446          DM(ind+1)=msc(iabs(itype(nct,1)),1)+isc(iabs(itype(nct,1)),1)
447          nind=ind+1
448        endif
449         endif
450        
451        
452         ind=1
453         do i=nnt,nct
454         mnum=molnum(i)
455        if (iabs(itype(i,1)).ne.10 .and.iabs(itype(i,mnum)).ne.ntyp1 &
456           .and.mnum.eq.5) then
457           DU1(ind)=-isc(iabs(itype(i,1)),1)
458           DU1(ind+1)=0.0d0
459          ind=ind+2
460        else
461          DU1(ind)=mp(1)/4-ip4
462          ind=ind+1
463        endif
464        enddo
465
466        ind=1
467        do i=nnt,nct-1
468        mnum=molnum(i)
469 !        if (iabs(itype(i,1)).eq.ntyp1) cycle
470         write (iout,*) "i",i," itype",itype(i,1),ntyp1
471        if (iabs(itype(i,1)).ne.10 .and. &
472           iabs(itype(i,mnum)).ne.ntyp1_molec(mnum) .and. mnum.ne.5) then
473          DU2(ind)=mp(1)/4-ip4
474          DU2(ind+1)=0.0d0
475          ind=ind+2
476        else
477          DU2(ind)=0.0d0
478         DU2(ind+1)=0.0d0
479         ind=ind+1
480        endif
481        enddo
482 #endif
483        DMorig=DM
484        DU1orig=DU1
485        DU2orig=DU2
486        write (iout,*) "nind",nind," dimen",dimen
487        if (nind.ne.dimen) then
488          write (iout,*) "ERROR, dimensions differ"
489 #ifdef MPI
490          call MPI_Finalize(ierr)
491 #endif
492          stop
493        endif
494      write (iout,*) "The upper part of the five-diagonal inertia matrix"
495        do i=1,dimen
496          if (i.lt.dimen-1) then
497            write (iout,'(i3,3f10.5)') i,DM(i),DU1(i),DU2(i)
498          else if (i.eq.dimen-1) then  
499            write (iout,'(i3,3f10.5)') i,DM(i),DU1(i)
500          else
501            write (iout,'(i3,3f10.5)') i,DM(i)
502          endif 
503        enddo
504
505         call FDISYP (dimen, DM, DU1, DU2, MARK)
506
507         if (mark.eq.-1) then
508      write (iout,*) "ERROR: the inertia matrix is not positive definite"
509 #ifdef MPI
510          call MPI_Finalize(ierr)
511 #endif
512          stop
513         else if (mark.eq.0) then
514      write (iout,*) "ERROR: the inertia matrix is singular"
515 #ifdef MPI
516          call MPI_Finalize(ierr)
517 #endif
518         else if (mark.eq.1) then
519           write (iout,*) "The transformed inertia matrix"
520           do i=1,dimen
521             if (i.lt.dimen-1) then
522               write (iout,'(i3,3f10.5)') i,DM(i),DU1(i),DU2(i)
523             else if (i.eq.dimen-1) then  
524               write (iout,'(i3,3f10.5)') i,DM(i),DU1(i)
525             else
526               write (iout,'(i3,3f10.5)') i,DM(i)
527             endif 
528           enddo
529         endif
530 ! Diagonalization of the pentadiagonal matrix
531         allocate(DDU1(2*nres))
532         allocate(DDU2(2*nres))
533         allocate(DDM(2*nres))
534         do i=1,dimen-1
535           DDU1(i+1)=DU1orig(i)
536         enddo
537         do i=1,dimen-2
538           DDU2(i+2)=DU2orig(i)
539         enddo
540         DDM=DMorig
541         call quindibisect2(DDM,DDU1,DDU2,dimen,1,dimen,eps1,relfeh,eps2,iz,geigen)
542         if (lprn) then
543           write (iout,*) &
544            "Eigenvalues of the five-diagonal inertia matrix"
545           do i=1,dimen
546             write (iout,'(i5,f10.5)') i,geigen(i)
547           enddo
548         endif
549         deallocate(DDU1)
550         deallocate(DDU2)
551         deallocate(DDM)
552 #else 
553 ! Old Gmatrix
554 #ifdef MPI
555       if (nfgtasks.gt.1) then
556       time00=MPI_Wtime()
557       call MPI_Bcast(5,1,MPI_INTEGER,king,FG_COMM,IERROR)
558       time_Bcast=time_Bcast+MPI_Wtime()-time00
559       call int_bounds(dimen,igmult_start,igmult_end)
560       igmult_start=igmult_start-1
561       call MPI_Allgather(3*igmult_start,1,MPI_INTEGER,&
562           ng_start(0),1,MPI_INTEGER,FG_COMM,IERROR)
563       my_ng_count=igmult_end-igmult_start
564       call MPI_Allgather(3*my_ng_count,1,MPI_INTEGER,ng_counts(0),1,&
565           MPI_INTEGER,FG_COMM,IERROR)
566       write (iout,*) 'Processor:',fg_rank,' CG group',kolor,&
567        ' absolute rank',myrank,' igmult_start',igmult_start,&
568        ' igmult_end',igmult_end,' count',my_ng_count
569       write (iout,*) "ng_start",(ng_start(i),i=0,nfgtasks-1)
570       write (iout,*) "ng_counts",(ng_counts(i),i=0,nfgtasks-1)
571       call flush(iout)
572       else
573 #endif
574       igmult_start=1
575       igmult_end=dimen
576       my_ng_count=dimen
577 #ifdef MPI
578       endif
579 #endif
580 !      write (iout,*) "dimen",dimen," dimen1",dimen1," dimen3",dimen3
581 !  Zeroing out A and fricmat
582       do i=1,dimen
583         do j=1,dimen
584           A(i,j)=0.0D0     
585         enddo    
586       enddo
587 !  Diagonal elements of the dC part of A and the respective friction coefficients
588       ind=1
589       ind1=0
590       print *,"TUTUTUT?!",nnt,nct-1
591       do i=nnt,nct-1
592         mnum=molnum(i)
593         ind=ind+1
594         ind1=ind1+1
595         coeff=0.25d0*IP(mnum)
596         if (mnum.eq.5) mp(mnum)=msc(itype(i,mnum),mnum)
597         print *,i,coeff,mp(mnum)
598         massvec(ind1)=mp(mnum)
599         Gmat(ind,ind)=coeff
600         print *,"i",mp(mnum)
601         A(ind1,ind)=0.5d0
602       enddo
603       
604 !  Off-diagonal elements of the dC part of A 
605       k=3
606       do i=1,nct-nnt
607         do j=1,i
608           A(i,j)=1.0d0
609         enddo
610       enddo
611 !  Diagonal elements of the dX part of A and the respective friction coefficients
612       m=nct-nnt
613       m1=nct-nnt+1
614       ind=0
615       ind1=0
616       do i=1,2
617       msc(ntyp1_molec(i),i)=1.0d0
618       enddo
619       do i=nnt,nct
620         ind=ind+1
621         ii = ind+m
622         mnum=molnum(i)
623         iti=itype(i,mnum)
624         if (mnum.eq.5) then
625         mscab=0.0
626         else
627         mscab=msc(iabs(iti),mnum)
628         endif
629         massvec(ii)=mscab
630
631         if (iti.ne.10 .and. iti.ne.ntyp1_molec(mnum) .and. mnum.ne.5) then
632           ind1=ind1+1
633           ii1= ind1+m1
634           A(ii,ii1)=1.0d0
635           Gmat(ii1,ii1)=ISC(iabs(iti),mnum)
636         endif
637       enddo
638 !  Off-diagonal elements of the dX part of A
639       ind=0
640       k=nct-nnt
641       do i=nnt,nct
642         iti=itype(i,1)
643         ind=ind+1
644         do j=nnt,i
645           ii = ind
646           jj = j-nnt+1
647           A(k+ii,jj)=1.0d0
648         enddo
649       enddo
650       if (lprn) then
651         write (iout,*)
652         write (iout,*) "Vector massvec"
653         do i=1,dimen1
654           write (iout,*) i,massvec(i)
655         enddo
656         write (iout,'(//a)') "A"
657         call matout(dimen,dimen1,nres2,nres2,A)
658       endif
659
660 ! Calculate the G matrix (store in Gmat)
661       do k=1,dimen
662        do i=1,dimen
663          dtdi=0.0d0
664          do j=1,dimen1
665            dtdi=dtdi+A(j,k)*A(j,i)*massvec(j)
666          enddo
667          Gmat(k,i)=Gmat(k,i)+dtdi
668        enddo
669       enddo 
670       
671       if (lprn) then
672         write (iout,'(//a)') "Gmat"
673         call matout(dimen,dimen,nres2,nres2,Gmat)
674       endif
675       do i=1,dimen
676         do j=1,dimen
677           Ginv(i,j)=0.0d0
678           Gcopy(i,j)=Gmat(i,j)
679         enddo
680         Ginv(i,i)=1.0d0
681       enddo
682 ! Invert the G matrix
683       call MATINVERT(dimen,nres2,Gcopy,Ginv,osob)
684       if (lprn) then
685         write (iout,'(//a)') "Ginv"
686         call matout(dimen,dimen,nres2,nres2,Ginv)
687       endif
688
689       Bmat=0.0d0
690       Bmat(1,1)=1.0d0
691       j=2
692       do i=nnt,nct
693         mnum=molnum(i)
694         if ((itype(i,1).eq.10).or.(itype(i,mnum).eq.ntyp1_molec(mnum))&
695           .or.(mnum.eq.5)) then
696           Bmat(i-nnt+2,j-1)=-1
697           Bmat(i-nnt+2,j)=1
698           j=j+1
699         else
700           if (i.lt.nct) then
701             Bmat(i-nnt+2,j-1)=-1
702             Bmat(i-nnt+2,j+1)=1
703             Bmat(i-nnt+1+(nct-nnt)+1,j-1)=-1
704             Bmat(i-nnt+1+(nct-nnt)+1,j)=1
705             j=j+2
706           else
707             Bmat(i-nnt+1+(nct-nnt)+1,j-1)=-1
708             Bmat(i-nnt+1+(nct-nnt)+1,j)=1
709           endif
710         endif
711       enddo
712       write (iout,*) "j",j," dimen",dimen 
713       if (lprn) then
714         write (iout,'(//a)') "Bmat"
715         call matout(dimen,dimen,nres2,nres2,Bmat)
716       endif
717       Gcopy=0.0d0
718       write(iout,*) "before Gcopy",dimen,nres*2
719 #ifdef TESTFIVEDIAG
720       do i=1,dimen
721         do j=1,dimen
722 !            write (iout,*) "myij",i,j
723             do k=1,dimen
724                do l=1,dimen
725 !                 write(iout,*) "i,j,k,l",i,j,k,l
726                  Gcopy(i,j)=Gcopy(i,j)+Bmat(k,i)*Gmat(k,l)*Bmat(l,j)
727                enddo
728             enddo
729         enddo
730       enddo
731 #endif
732       if (lprn) then
733         write (iout,'(//a)') "Gmat-transformed"
734         call matout(dimen,dimen,nres2,nres2,Gcopy)
735       endif
736 #ifdef MPI
737       if (nfgtasks.gt.1) then
738         myginv_ng_count=nres2*my_ng_count
739         call MPI_Allgather(nres2*igmult_start,1,MPI_INTEGER,&
740           nginv_start(0),1,MPI_INTEGER,FG_COMM,IERROR)
741         call MPI_Allgather(myginv_ng_count,1,MPI_INTEGER,&
742           nginv_counts(0),1,MPI_INTEGER,FG_COMM,IERROR)
743         write (iout,*) "nginv_start",(nginv_start(i),i=0,nfgtasks-1)
744         write (iout,*) "nginv_counts",(nginv_counts(i),i=0,nfgtasks-1)
745         call flush(iout)
746 !        call MPI_Scatterv(ginv(1,1),nginv_counts(0),
747 !     &    nginv_start(0),MPI_DOUBLE_PRECISION,ginv,
748 !     &    myginv_ng_count,MPI_DOUBLE_PRECISION,king,FG_COMM,IERR)
749 !        call MPI_Barrier(FG_COMM,IERR)
750         time00=MPI_Wtime()
751         call MPI_Scatterv(ginv(1,1),nginv_counts(0),&
752           nginv_start(0),MPI_DOUBLE_PRECISION,gcopy(1,1),&
753           myginv_ng_count,MPI_DOUBLE_PRECISION,king,FG_COMM,IERR)
754 #ifdef TIMING
755         time_scatter_ginv=time_scatter_ginv+MPI_Wtime()-time00
756 #endif
757         do i=1,dimen
758           do j=1,2*my_ng_count
759             ginv(j,i)=gcopy(i,j)
760           enddo
761         enddo
762 !        write (iout,*) "Master's chunk of ginv"
763 !        call MATOUT2(my_ng_count,dimen,maxres2,maxres2,ginv)
764       endif
765 #endif
766       if (osob) then
767         write (iout,*) "The G matrix is singular."
768         stop
769       endif
770 ! Compute G**(-1/2) and G**(1/2) 
771       ind=0
772       do i=1,dimen
773         do j=1,i
774           ind=ind+1
775           Ghalf(ind)=Gmat(i,j)
776         enddo
777       enddo
778       call gldiag(nres2,dimen,dimen,Ghalf,work,Geigen,Gvec,&
779         ierr,iwork)
780       if (lprn) then
781         write (iout,'(//a)') &
782          "Eigenvectors and eigenvalues of the G matrix"
783         call eigout(dimen,dimen,nres2,nres2,Gvec,Geigen)
784       endif
785       do i=1,dimen
786         sqreig(i)=dsqrt(Geigen(i))
787       enddo
788       do i=1,dimen
789         do j=1,dimen
790           Gsqrp(i,j)=0.0d0
791           Gsqrm(i,j)=0.0d0
792           Gcopy(i,j)=0.0d0
793           do k=1,dimen
794             Gsqrp(i,j)=Gsqrp(i,j)+Gvec(i,k)*Gvec(j,k)*sqreig(k)
795             Gsqrm(i,j)=Gsqrm(i,j)+Gvec(i,k)*Gvec(j,k)/sqreig(k)
796             Gcopy(i,j)=Gcopy(i,j)+Gvec(i,k)*Gvec(j,k)*Geigen(k)
797           enddo
798         enddo
799       enddo
800       if (lprn) then
801         write (iout,*) "Comparison of original and restored G"
802         do i=1,dimen
803           do j=1,dimen
804             write (iout,'(2i5,5f10.5)') i,j,Gmat(i,j),Gcopy(i,j),&
805               Gmat(i,j)-Gcopy(i,j),Gsqrp(i,j),Gsqrm(i,j)
806           enddo
807         enddo
808       endif
809       deallocate(Gcopy)
810 #endif
811       return
812       end subroutine setup_MD_matrices
813 !-----------------------------------------------------------------------------
814       subroutine EIGOUT(NC,NR,LM2,LM3,A,B)
815
816 !      implicit real*8 (a-h,o-z)
817 !      include 'DIMENSIONS'
818 !      include 'COMMON.IOUNITS'
819       integer :: LM2,LM3,NC,NR,KA,KC,KB,I,J,N
820       real(kind=8) :: A(LM2,LM3),B(LM2)
821       KA=1
822       KC=6
823     1 KB=MIN0(KC,NC)
824       WRITE(IOUT,600) (I,I=KA,KB)
825       WRITE(IOUT,601) (B(I),I=KA,KB)
826       WRITE(IOUT,602)
827     2 N=0
828       DO 3  I=1,NR
829       WRITE(IOUT,603) I,(A(I,J),J=KA,KB)
830       N=N+1
831       IF(N.LT.10) GO TO 3
832       WRITE(IOUT,602)
833       N=0
834     3 CONTINUE
835     4 IF (KB.EQ.NC) RETURN
836       KA=KC+1
837       KC=KC+6
838       GO TO 1
839   600 FORMAT (// 9H ROOT NO.,I4,9I11)
840   601 FORMAT (/5X,10(1PE11.4))
841   602 FORMAT (2H  )
842   603 FORMAT (I5,10F11.5)
843   604 FORMAT (1H1)
844       end subroutine EIGOUT
845 !-----------------------------------------------------------------------------
846       subroutine MATOUT(NC,NR,LM2,LM3,A)
847
848 !      implicit real*8 (a-h,o-z)
849 !      include 'DIMENSIONS'
850 !      include 'COMMON.IOUNITS'
851       integer :: LM2,LM3,NC,NR,KA,KC,KB,I,J,N
852       real(kind=8) :: A(LM2,LM3)
853       KA=1
854       KC=6
855     1 KB=MIN0(KC,NC)
856       WRITE(IOUT,600) (I,I=KA,KB)
857       WRITE(IOUT,602)
858     2 N=0
859       DO 3  I=1,NR
860       WRITE(IOUT,603) I,(A(I,J),J=KA,KB)
861       N=N+1
862       IF(N.LT.10) GO TO 3
863       WRITE(IOUT,602)
864       N=0
865     3 CONTINUE
866     4 IF (KB.EQ.NC) RETURN
867       KA=KC+1
868       KC=KC+6
869       GO TO 1
870   600 FORMAT (//5x,9I11)
871   602 FORMAT (2H  )
872   603 FORMAT (I5,10F11.3)
873   604 FORMAT (1H1)
874       end subroutine MATOUT
875 !-----------------------------------------------------------------------------
876       subroutine MATOUT1(NC,NR,LM2,LM3,A)
877
878 !      implicit real*8 (a-h,o-z)
879 !      include 'DIMENSIONS'
880 !      include 'COMMON.IOUNITS'
881       integer :: LM2,LM3,NC,NR,KA,KC,KB,I,J,N
882       real(kind=8) :: A(LM2,LM3)
883       KA=1
884       KC=21
885     1 KB=MIN0(KC,NC)
886       WRITE(IOUT,600) (I,I=KA,KB)
887       WRITE(IOUT,602)
888     2 N=0
889       DO 3  I=1,NR
890       WRITE(IOUT,603) I,(A(I,J),J=KA,KB)
891       N=N+1
892       IF(N.LT.3) GO TO 3
893       WRITE(IOUT,602)
894       N=0
895     3 CONTINUE
896     4 IF (KB.EQ.NC) RETURN
897       KA=KC+1
898       KC=KC+21
899       GO TO 1
900   600 FORMAT (//5x,7(3I5,2x))
901   602 FORMAT (2H  )
902   603 FORMAT (I5,7(3F5.1,2x))
903   604 FORMAT (1H1)
904       end subroutine MATOUT1
905 !-----------------------------------------------------------------------------
906       subroutine MATOUT2(NC,NR,LM2,LM3,A)
907
908 !      implicit real*8 (a-h,o-z)
909 !      include 'DIMENSIONS'
910 !      include 'COMMON.IOUNITS'
911       integer :: I,J,KA,KC,KB,N
912       integer :: LM2,LM3,NC,NR
913       real(kind=8) :: A(LM2,LM3)
914       KA=1
915       KC=12
916     1 KB=MIN0(KC,NC)
917       WRITE(IOUT,600) (I,I=KA,KB)
918       WRITE(IOUT,602)
919     2 N=0
920       DO 3  I=1,NR
921       WRITE(IOUT,603) I,(A(I,J),J=KA,KB)
922       N=N+1
923       IF(N.LT.3) GO TO 3
924       WRITE(IOUT,602)
925       N=0
926     3 CONTINUE
927     4 IF (KB.EQ.NC) RETURN
928       KA=KC+1
929       KC=KC+12
930       GO TO 1
931   600 FORMAT (//5x,4(3I9,2x))
932   602 FORMAT (2H  )
933   603 FORMAT (I5,4(3F9.3,2x))
934   604 FORMAT (1H1)
935       end subroutine MATOUT2
936 !-----------------------------------------------------------------------------
937       subroutine ginv_mult(z,d_a_tmp)
938
939       use geometry_data, only: nres
940       use control_data
941       use MPI_data
942 !      implicit real*8 (a-h,o-z)
943 !      include 'DIMENSIONS'
944 #ifdef MPI
945       include 'mpif.h'
946       integer :: ierr,ierror
947 #endif
948 !      include 'COMMON.SETUP'
949 !      include 'COMMON.TIME1'
950 !      include 'COMMON.MD'
951       real(kind=8),dimension(dimen3) :: z,z1,d_a_tmp
952       real(kind=8),dimension(6*nres) :: temp    !(maxres6) maxres6=6*maxres
953       real(kind=8) :: time00,time01
954       integer :: i,j,k,ind
955 #ifdef MPI
956       if (nfgtasks.gt.1) then
957         if (fg_rank.eq.0) then
958 ! The matching BROADCAST for fg processors is called in ERGASTULUM
959           time00=MPI_Wtime()
960           call MPI_Bcast(4,1,MPI_INTEGER,king,FG_COMM,IERROR)
961           time_Bcast=time_Bcast+MPI_Wtime()-time00
962 !          print *,"Processor",myrank," BROADCAST iorder in GINV_MULT"
963         endif
964 !        write (2,*) "time00",time00
965 !        write (2,*) "Before Scatterv"
966 !        call flush(2)
967 !        write (2,*) "Whole z (for FG master)"
968 !        do i=1,dimen
969 !          write (2,*) i,z(i)
970 !        enddo
971 !        call MPI_Barrier(FG_COMM,IERROR)
972         time00=MPI_Wtime()
973 !elwrite(iout,*) "do tej pory jest OK, MPI_Scatterv w ginv_mult"
974         call MPI_Scatterv(z,ng_counts(0),ng_start(0),&
975           MPI_DOUBLE_PRECISION,&
976           z1,3*my_ng_count,MPI_DOUBLE_PRECISION,king,FG_COMM,IERR)
977 !        write (2,*) "My chunk of z"
978 !        do i=1,3*my_ng_count
979 !          write (2,*) i,z(i)
980 !        enddo
981 !        write (2,*) "After SCATTERV"
982 !        call flush(2)
983 !        write (2,*) "MPI_Wtime",MPI_Wtime()
984         time_scatter=time_scatter+MPI_Wtime()-time00
985 #ifdef TIMING
986         time_scatter_ginvmult=time_scatter_ginvmult+MPI_Wtime()-time00
987 #endif
988 !        write (2,*) "time_scatter",time_scatter
989 !        write (2,*) "dimen",dimen," dimen3",dimen3," my_ng_count",
990 !     &    my_ng_count
991 !        call flush(2)
992         time01=MPI_Wtime()
993         do k=0,2
994           do i=1,dimen
995             ind=(i-1)*3+k+1
996             temp(ind)=0.0d0
997             do j=1,my_ng_count
998 !              write (2,*) "k,i,j,ind",k,i,j,ind,(j-1)*3+k+1,
999 !     &         Ginv(i,j),z((j-1)*3+k+1),
1000 !     &          Ginv(i,j)*z((j-1)*3+k+1)
1001 !              temp(ind)=temp(ind)+Ginv(i,j)*z((j-1)*3+k+1)
1002               temp(ind)=temp(ind)+Ginv(j,i)*z1((j-1)*3+k+1)
1003             enddo
1004           enddo 
1005         enddo
1006         time_ginvmult=time_ginvmult+MPI_Wtime()-time01
1007 !        write (2,*) "Before REDUCE"
1008 !        call flush(2)
1009 !        write (2,*) "z before reduce"
1010 !        do i=1,dimen
1011 !          write (2,*) i,temp(i)
1012 !        enddo
1013         time00=MPI_Wtime()
1014         call MPI_Reduce(temp(1),d_a_tmp(1),dimen3,MPI_DOUBLE_PRECISION,&
1015             MPI_SUM,king,FG_COMM,IERR)
1016         time_reduce=time_reduce+MPI_Wtime()-time00
1017 !        write (2,*) "After REDUCE"
1018 !        call flush(2)
1019       else
1020 #endif
1021 #ifdef TIMING
1022         time01=MPI_Wtime()
1023 #endif
1024         do k=0,2
1025           do i=1,dimen
1026             ind=(i-1)*3+k+1
1027             d_a_tmp(ind)=0.0d0
1028             do j=1,dimen
1029 !              write (2,*) "k,i,j,ind",k,i,j,ind,(j-1)*3+k+1
1030 !              call flush(2)
1031 !     &         Ginv(i,j),z((j-1)*3+k+1),
1032 !     &          Ginv(i,j)*z((j-1)*3+k+1)
1033               d_a_tmp(ind)=d_a_tmp(ind) &
1034                                +Ginv(j,i)*z((j-1)*3+k+1)
1035 !              d_a_tmp(ind)=d_a_tmp(ind)
1036 !     &                         +Ginv(i,j)*z((j-1)*3+k+1)
1037             enddo
1038           enddo 
1039         enddo
1040 #ifdef TIMING
1041         time_ginvmult=time_ginvmult+MPI_Wtime()-time01
1042 #endif
1043 #ifdef MPI
1044       endif
1045 #endif
1046       return
1047       end subroutine ginv_mult
1048 !-----------------------------------------------------------------------------
1049 #ifdef GINV_MULT
1050       subroutine ginv_mult_test(z,d_a_tmp)
1051
1052 !      include 'DIMENSIONS'
1053 !el      integer :: dimen
1054 !      include 'COMMON.MD'
1055       real(kind=8),dimension(dimen) :: z,d_a_tmp
1056       real(kind=8),dimension(dimen/3) :: ztmp,dtmp
1057       integer :: i,j,k,ind
1058 !      do i=1,dimen
1059 !        d_a_tmp(i)=0.0d0
1060 !        do j=1,dimen
1061 !          d_a_tmp(i)=d_a_tmp(i)+Ginv(i,j)*z(j)
1062 !        enddo
1063 !      enddo
1064 !
1065 !      return
1066
1067 !ibm* unroll(3)
1068       do k=0,2
1069        do j=1,dimen/3
1070         ztmp(j)=z((j-1)*3+k+1)
1071        enddo
1072
1073        call alignx(16,ztmp(1))
1074        call alignx(16,dtmp(1))
1075        call alignx(16,Ginv(1,1)) 
1076
1077        do i=1,dimen/3
1078         dtmp(i)=0.0d0
1079         do j=1,dimen/3
1080            dtmp(i)=dtmp(i)+Ginv(i,j)*ztmp(j)
1081         enddo
1082        enddo
1083        do i=1,dimen/3
1084         ind=(i-1)*3+k+1
1085         d_a_tmp(ind)=dtmp(i)
1086        enddo 
1087       enddo
1088       return
1089       end subroutine ginv_mult_test
1090 #endif
1091 !-----------------------------------------------------------------------------
1092       subroutine fricmat_mult(z,d_a_tmp)
1093
1094       use geometry_data, only: nres
1095       use control_data
1096       use MPI_data
1097 !      include 'DIMENSIONS'
1098 #ifdef MPI
1099       include 'mpif.h'
1100       integer :: IERROR,ierr
1101 #endif
1102 !      include 'COMMON.MD'
1103 !      include 'COMMON.IOUNITS'
1104 !      include 'COMMON.SETUP'
1105 !      include 'COMMON.TIME1'
1106 !#ifndef LANG0
1107 !      include 'COMMON.LANGEVIN'
1108 !#else
1109 !      include 'COMMON.LANGEVIN.lang0'
1110 !#endif
1111       real(kind=8),dimension(dimen3) :: z,z1,d_a_tmp
1112       real(kind=8),dimension(6*nres) :: temp    !(maxres6) maxres6=6*maxres
1113       real(kind=8) :: time00,time01
1114       integer :: i,j,k,ind,nres2
1115       nres2=2*nres
1116 !el      if(.not.allocated(fricmat)) allocate(fricmat(nres2,nres2))
1117
1118 #ifdef MPI
1119       if (nfgtasks.gt.1) then
1120         if (fg_rank.eq.0) then
1121 ! The matching BROADCAST for fg processors is called in ERGASTULUM
1122           time00=MPI_Wtime()
1123           call MPI_Bcast(9,1,MPI_INTEGER,king,FG_COMM,IERROR)
1124           time_Bcast=time_Bcast+MPI_Wtime()-time00
1125 !          print *,"Processor",myrank," BROADCAST iorder in FRICMAT_MULT"
1126         endif
1127 !        call MPI_Barrier(FG_COMM,IERROR)
1128         time00=MPI_Wtime()
1129         call MPI_Scatterv(z,ng_counts(0),ng_start(0),&
1130           MPI_DOUBLE_PRECISION,&
1131           z1,3*my_ng_count,MPI_DOUBLE_PRECISION,king,FG_COMM,IERR)
1132 !        write (2,*) "My chunk of z"
1133 !        do i=1,3*my_ng_count
1134 !          write (2,*) i,z(i)
1135 !        enddo
1136         time_scatter=time_scatter+MPI_Wtime()-time00
1137 #ifdef TIMING
1138         time_scatter_fmatmult=time_scatter_fmatmult+MPI_Wtime()-time00
1139 #endif
1140         time01=MPI_Wtime()
1141         do k=0,2
1142           do i=1,dimen
1143             ind=(i-1)*3+k+1
1144             temp(ind)=0.0d0
1145             do j=1,my_ng_count
1146               temp(ind)=temp(ind)-fricmat(j,i)*z1((j-1)*3+k+1)
1147             enddo
1148           enddo 
1149         enddo
1150         time_fricmatmult=time_fricmatmult+MPI_Wtime()-time01
1151 !        write (2,*) "Before REDUCE"
1152 !        write (2,*) "d_a_tmp before reduce"
1153 !        do i=1,dimen3
1154 !          write (2,*) i,temp(i)
1155 !        enddo
1156 !        call flush(2)
1157         time00=MPI_Wtime()
1158         call MPI_Reduce(temp(1),d_a_tmp(1),dimen3,MPI_DOUBLE_PRECISION,&
1159             MPI_SUM,king,FG_COMM,IERR)
1160         time_reduce=time_reduce+MPI_Wtime()-time00
1161 !        write (2,*) "After REDUCE"
1162 !        call flush(2)
1163       else
1164 #endif
1165 #ifdef TIMING
1166         time01=MPI_Wtime()
1167 #endif
1168         do k=0,2
1169          do i=1,dimen
1170           ind=(i-1)*3+k+1
1171           d_a_tmp(ind)=0.0d0
1172           do j=1,dimen
1173              d_a_tmp(ind)=d_a_tmp(ind) &
1174                                  -fricmat(j,i)*z((j-1)*3+k+1)
1175           enddo
1176          enddo 
1177         enddo
1178 #ifdef TIMING
1179         time_fricmatmult=time_fricmatmult+MPI_Wtime()-time01
1180 #endif
1181 #ifdef MPI
1182       endif
1183 #endif
1184 !      write (iout,*) "Vector d_a"
1185 !      do i=1,dimen3
1186 !        write (2,*) i,d_a_tmp(i)
1187 !      enddo
1188       return
1189       end subroutine fricmat_mult
1190 !-----------------------------------------------------------------------------
1191 !-----------------------------------------------------------------------------
1192       end module REMD