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