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