1 subroutine WHAM_CALC(islice,*)
2 ! Weighed Histogram Analysis Method (WHAM) code
3 ! Written by A. Liwo based on the work of Kumar et al.,
4 ! J.Comput.Chem., 13, 1011 (1992)
6 ! 2/1/05 Multiple temperatures allowed.
7 ! 2/2/05 Free energies calculated directly from data points
8 ! acc. to Eq. (21) of Kumar et al.; final histograms also
9 ! constructed based on this equation.
10 ! 2/12/05 Multiple parameter sets included
12 ! 2/2/05 Parallel version
15 include "DIMENSIONS.ZSCOPT"
16 include "DIMENSIONS.FREE"
17 integer MaxBinRms,MaxBinRgy
18 parameter (MaxBinRms=100,MaxBinRgy=100)
20 c parameter (MaxHdim=200000)
21 parameter (MaxHdim=200)
23 parameter (maxinde=200)
27 integer ierror,errcode,status(MPI_STATUS_SIZE)
29 include "COMMON.CONTROL"
30 include "COMMON.IOUNITS"
32 include "COMMON.ENERGIES"
33 include "COMMON.FFIELD"
34 include "COMMON.SBRIDGE"
36 include "COMMON.ENEPS"
37 integer MaxPoint,MaxPointProc
38 parameter (MaxPoint=MaxStr,
39 & MaxPointProc=MaxStr_Proc)
40 double precision finorm_max,potfac,entmin,entmax,expfac,vf
41 double precision entfac_min,entfac_min_t
42 parameter (finorm_max=1.0d0)
44 integer i,ii,j,jj,k,kk,l,m,ind,iter,t,tmax,ient,ientmax,iln
45 integer start,end,iharm,ib,iib,nbin1,nbin,nbin_rms,nbin_rgy,
46 & nbin_rmsrgy,liczba,iparm,nFi,indrgy,indrms
47 integer htot(0:MaxHdim),histent(0:2000)
48 double precision v(MaxPointProc,MaxR,MaxT_h,Max_Parm)
49 double precision energia(0:max_ene)
51 integer tmax_t,upindE_p
52 double precision fi_p(MaxR,MaxT_h,Max_Parm),
53 & fi_p_min(MaxR,MaxT_h,Max_Parm)
54 double precision sumW_p(0:Max_GridT,Max_Parm),
55 & sumE_p(0:Max_GridT,Max_Parm),sumEsq_p(0:Max_GridT,Max_Parm),
56 & sumQ_p(MaxQ1,0:Max_GridT,Max_Parm),
57 & sumQsq_p(MaxQ1,0:Max_GridT,Max_Parm),
58 & sumEQ_p(MaxQ1,0:Max_GridT,Max_Parm),
59 & sumEprim_p(MaxQ1,0:Max_GridT,Max_Parm),
60 & sumEbis_p(0:Max_GridT,Max_Parm)
61 double precision hfin_p(0:MaxHdim,maxT_h),
62 & hfin_ent_p(0:MaxHdim),histE_p(0:maxindE),sumH,
63 & hrmsrgy_p(0:MaxBinRgy,0:MaxBinRms,maxT_h)
64 double precision rgymin_t,rmsmin_t,rgymax_t,rmsmax_t
65 double precision potEmin_t_all(maxT_h,Max_Parm),entmin_p,entmax_p
66 integer histent_p(0:2000)
67 logical lprint /.true./
69 double precision rgymin,rmsmin,rgymax,rmsmax
70 double precision sumW(0:NGridT,Max_Parm),sumE(0:NGridT,Max_Parm),
71 & sumEsq(0:NGridT,Max_Parm),sumQ(MaxQ1,0:NGridT,Max_Parm),
72 & sumQsq(MaxQ1,0:NGridT,Max_Parm),sumEQ(MaxQ1,0:NGridT,Max_Parm),
73 & sumEprim(0:NGridT,Max_Parm),sumEbis(0:NGridT,Max_Parm),betaT,
75 double precision fi(MaxR,maxT_h,Max_Parm),
76 & fi_min(MaxR,maxT_h,Max_Parm),
77 & dd,dd1,dd2,hh,dmin,denom,finorm,avefi,pom,
78 & hfin(0:MaxHdim,maxT_h),histE(0:maxindE),
79 & hrmsrgy(0:MaxBinRgy,0:MaxBinRms,maxT_h),
80 & potEmin_all(maxT_h,Max_Parm),potEmin,potEmin_min,ent,
81 & hfin_ent(0:MaxHdim),vmax,aux
82 double precision fT(6),fTprim(6),fTbis(6),quot,quotl1,quotl,kfacl,
83 & eprim,ebis,temper,kfac/2.4d0/,T0/300.0d0/,
84 & eplus,eminus,logfac,tanhT,tt
85 double precision etot,evdw,evdw_t,evdw2,ees,evdw1,ebe,etors,
86 & escloc,ehpb,ecorr,ecorr5,ecorr6,eello_turn4,eello_turn3,
87 & eturn6,eel_loc,edihcnstr,etors_d,estr,evdw2_14,esccor
89 integer ind_point(maxpoint),upindE,indE
98 write(licz2,'(bz,i2.2)') islice
100 write (iout,'(//80(1h-)/"Solving WHAM equations for slice",
101 & i2/80(1h-)//)') islice
102 write (iout,*) "delta",delta," nbin1",nbin1
103 write (iout,*) "MaxN",MaxN," MaxQ",MaxQ," MaHdim",MaxHdim
109 potEmin_all(j,i)=1.0d10
124 if (q(nQ+1,i).lt.rmsmin) rmsmin=q(nQ+1,i)
125 if (q(nQ+1,i).gt.rmsmax) rmsmax=q(nQ+1,i)
126 if (q(nQ+2,i).lt.rgymin) rgymin=q(nQ+2,i)
127 if (q(nQ+2,i).gt.rgymax) rgymax=q(nQ+2,i)
130 ind=(q(j,i)-dmin+1.0d-8)/delta
132 ind_point(i)=ind_point(i)+ind
134 ind_point(i)=ind_point(i)+nbin1**(j-1)*ind
136 if (ind_point(i).lt.0 .or. ind_point(i).gt.MaxHdim) then
137 write (iout,*) "Error - index exceeds range for point",i,
138 & " q=",q(j,i)," ind",ind_point(i)
140 write (iout,*) "Processor",me1
142 call MPI_Abort(MPI_COMM_WORLD, Ierror, Errcode )
147 if (ind_point(i).gt.tmax) tmax=ind_point(i)
148 htot(ind_point(i))=htot(ind_point(i))+1
150 write (iout,*) "i",i,"q",(q(j,i),j=1,nQ)," ind",ind_point(i),
151 & " htot",htot(ind_point(i))
158 write (iout,'(a)') "Numbers of counts in Q bins"
160 if (htot(t).gt.0) then
161 write (iout,'(i15,$)') t
164 jj = mod(liczba,nbin1)
166 write (iout,'(i5,$)') jj
168 write (iout,'(i8)') htot(t)
172 write (iout,'(a,i3)') "Number of data points for parameter set",
174 write (iout,'(i7,$)') ((snk(m,ib,iparm,islice),m=1,nr(ib,iparm)),
176 write (iout,'(i8)') stot(islice)
182 call MPI_AllReduce(tmax,tmax_t,1,MPI_INTEGER,MPI_MAX,
185 call MPI_AllReduce(rmsmin,rmsmin_t,1,MPI_DOUBLE_PRECISION,
186 & MPI_MIN,WHAM_COMM,IERROR)
187 call MPI_AllReduce(rmsmax,rmsmax_t,1,MPI_DOUBLE_PRECISION,
188 & MPI_MAX,WHAM_COMM,IERROR)
189 call MPI_AllReduce(rgymin,rgymin_t,1,MPI_DOUBLE_PRECISION,
190 & MPI_MIN,WHAM_COMM,IERROR)
191 call MPI_AllReduce(rgymax,rgymax_t,1,MPI_DOUBLE_PRECISION,
192 & MPI_MAX,WHAM_COMM,IERROR)
198 rmsmin=deltrms*dint(rmsmin/deltrms)
199 rmsmax=deltrms*dint(rmsmax/deltrms)
200 rgymin=deltrms*dint(rgymin/deltrgy)
201 rgymax=deltrms*dint(rgymax/deltrgy)
202 nbin_rms=(rmsmax-rmsmin)/deltrms
203 nbin_rgy=(rgymax-rgymin)/deltrgy
204 write (iout,*) "rmsmin",rmsmin," rmsmax",rmsmax," rgymin",rgymin,
205 & " rgymax",rgymax," nbin_rms",nbin_rms," nbin_rgy",nbin_rgy
212 write (iout,*) "nFi",nFi
213 ! Compute the Boltzmann factor corresponing to restrain potentials in different
220 c write (9,'(3i5,f10.5)') i,(iparm,potE(i,iparm),iparm=1,nParmSet)
223 write (iout,'(2i5,21f8.2)') i,iparm,
224 & (enetb(k,i,iparm),k=1,21)
226 call restore_parm(iparm)
228 write (iout,*) wsc,wscp,welec,wvdwpp,wang,wtor,wscloc,
229 & wcorr,wcorr5,wcorr6,wturn4,wturn3,wturn6,wel_loc,
230 & wtor_d,wsccor,wbond
233 if (rescale_mode.eq.1) then
234 quot=1.0d0/(beta_h(ib,iparm)*1.987D-3*T0)
241 fT(l)=kfacl/(kfacl-1.0d0+quotl)
244 tt = 1.0d0/(beta_h(ib,iparm)*1.987D-3)
245 ft(6)=(320.0d0+80.0d0*dtanh((tt-320.0d0)/80.0d0))/320.0d0
247 ft(6)=1.0d0/(beta_h(ib,iparm)*1.987D-3*T0)
251 else if (rescale_mode.eq.2) then
252 quot=1.0d0/(T0*beta_h(ib,iparm)*1.987D-3)
256 fT(l)=1.12692801104297249644d0/
257 & dlog(dexp(quotl)+dexp(-quotl))
260 tt = 1.0d0/(beta_h(ib,iparm)*1.987D-3)
261 ft(6)=(320.0d0+80.0d0*dtanh((tt-320.0d0)/80.0d0))/320.0d0
263 ft(6)=1.0d0/(beta_h(ib,iparm)*1.987D-3*T0)
267 c write (iout,*) 1.0d0/(beta_h(ib,iparm)*1.987D-3),ft
268 else if (rescale_mode.eq.0) then
273 write (iout,*) "Error in WHAM_CALC: wrong RESCALE_MODE",
278 evdw=enetb(1,i,iparm)
279 evdw_t=enetb(21,i,iparm)
281 evdw2_14=enetb(17,i,iparm)
282 evdw2=enetb(2,i,iparm)+evdw2_14
284 evdw2=enetb(2,i,iparm)
289 evdw1=enetb(16,i,iparm)
294 ecorr=enetb(4,i,iparm)
295 ecorr5=enetb(5,i,iparm)
296 ecorr6=enetb(6,i,iparm)
297 eel_loc=enetb(7,i,iparm)
298 eello_turn3=enetb(8,i,iparm)
299 eello_turn4=enetb(9,i,iparm)
300 eturn6=enetb(10,i,iparm)
301 ebe=enetb(11,i,iparm)
302 escloc=enetb(12,i,iparm)
303 etors=enetb(13,i,iparm)
304 etors_d=enetb(14,i,iparm)
305 ehpb=enetb(15,i,iparm)
306 estr=enetb(18,i,iparm)
307 esccor=enetb(19,i,iparm)
308 edihcnstr=enetb(20,i,iparm)
310 write (iout,'(3i5,6f5.2,14f12.3)') i,ib,iparm,(ft(l),l=1,6),
311 & evdw+evdw_t,evdw2,ees,evdw1,ecorr,eel_loc,estr,ebe,escloc,
312 & etors,etors_d,eello_turn3,eello_turn4,esccor
316 etot=wsc*(evdw+ft(6)*evdw_t)+wscp*evdw2+ft(1)*welec*ees
318 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
319 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
320 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
321 & +ft(2)*wturn3*eello_turn3
322 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc
323 & +edihcnstr+ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
326 etot=wsc*(evdw+ft(6)*evdw_t)+wscp*evdw2
327 & +ft(1)*welec*(ees+evdw1)
328 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
329 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
330 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
331 & +ft(2)*wturn3*eello_turn3
332 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc+edihcnstr
333 & +ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
337 write (iout,*) i,iparm,1.0d0/(beta_h(ib,iparm)*1.987D-3),
341 if (iparm.eq.1 .and. ib.eq.1) then
342 write (iout,*)"Conformation",i
345 energia(k)=enetb(k,i,iparm)
347 call enerprint(energia(0),fT)
354 Econstr=Econstr+Kh(j,kk,ib,iparm)
355 & *(dd-q0(j,kk,ib,iparm))**2
358 & -beta_h(ib,iparm)*(etot+Econstr)
360 write (iout,'(4i5,4e15.5)') i,kk,ib,iparm,
361 & etot,v(i,kk,ib,iparm)
367 ! Simple iteration to calculate free energies corresponding to all simulation
371 ! Compute new free-energy values corresponding to the righ-hand side of the
372 ! equation and their derivatives.
373 write (iout,*) "------------------------fi"
384 vf=v(t,l,k,i)+f(l,k,i)
385 if (vf.gt.vmax) vmax=vf
393 aux=f(l,k,i)+v(t,l,k,i)-vmax
395 & denom=denom+snk(l,k,i,islice)*dexp(aux)
399 entfac(t)=-dlog(denom)-vmax
400 if (entfac(t).lt.entfac_min) entfac_min=entfac(t)
402 write (iout,*) t,"vmax",vmax," denom",denom,"entfac",entfac(t)
406 c write (iout,*) "entfac_min before AllReduce",entfac_min
407 c call MPI_AllReduce(entfac_min,entfac_min_t,1,
408 c & MPI_DOUBLE_PRECISION,MPI_MAX,WHAM_COMM,IERROR)
409 c entfac_min=entfac_min_t
410 c write (iout,*) "entfac_min after AllReduce",entfac_min
414 c entfac(t)=entfac(t)-entfac_min
417 c do t=1,ntot(islice)
418 c entfac(t)=entfac(t)-entfac_min
423 do ii=1,nR(iib,iparm)
425 fi_p_min(ii,iib,iparm)=-1.0d10
427 aux=v(t,ii,iib,iparm)+entfac(t)
428 if (aux.gt.fi_p_min(ii,iib,iparm))
429 & fi_p_min(ii,iib,iparm)=aux
433 aux=v(t,ii,iib,iparm)+entfac(t)
434 if (aux.gt.fi_min(ii,iib,iparm))
435 & fi_min(ii,iib,iparm)=aux
443 write (iout,*) "fi_min before AllReduce"
446 write (iout,*) (i,j,k,fi_p_min(k,j,i),k=1,nR(j,i))
450 call MPI_AllReduce(fi_p_min,fi_min,MaxR*MaxT_h*nParmSet,
451 & MPI_DOUBLE_PRECISION,MPI_MAX,WHAM_COMM,IERROR)
453 write (iout,*) "fi_min after AllReduce"
456 write (iout,*) (i,j,k,fi_min(k,j,i),k=1,nR(j,i))
463 do ii=1,nR(iib,iparm)
465 fi_p(ii,iib,iparm)=0.0d0
467 fi_p(ii,iib,iparm)=fi_p(ii,iib,iparm)
468 & +dexp(v(t,ii,iib,iparm)+entfac(t)-fi_min(ii,iib,iparm))
470 write (iout,'(4i5,4e15.5)') t,ii,iib,iparm,
471 & v(t,ii,iib,iparm),entfac(t),fi_min(ii,iib,iparm),
476 fi(ii,iib,iparm)=0.0d0
478 fi(ii,iib,iparm)=fi(ii,iib,iparm)
479 & +dexp(v(t,ii,iib,iparm)+entfac(t)-fi_min(ii,iib,iparm))
488 write (iout,*) "fi before MPI_Reduce me",me,' master',master
491 write (iout,*) "iparm",iparm," ib",ib
492 write (iout,*) "beta=",beta_h(ib,iparm)
493 write (iout,'(8e15.5)') (fi_p(i,ib,iparm),i=1,nR(ib,iparm))
498 write (iout,*) "REDUCE size",maxR,MaxT_h,nParmSet,
499 & maxR*MaxT_h*nParmSet
500 write (iout,*) "MPI_COMM_WORLD",MPI_COMM_WORLD,
501 & " WHAM_COMM",WHAM_COMM
503 call MPI_Reduce(fi_p(1,1,1),fi(1,1,1),maxR*MaxT_h*nParmSet,
504 & MPI_DOUBLE_PRECISION,
505 & MPI_SUM,Master,WHAM_COMM,IERROR)
507 write (iout,*) "fi after MPI_Reduce nparmset",nparmset
509 write (iout,*) "iparm",iparm
511 write (iout,*) "beta=",beta_h(ib,iparm)
512 write (iout,'(8e15.5)') (fi(i,ib,iparm),i=1,nR(ib,iparm))
516 if (me1.eq.Master) then
522 fi(i,ib,iparm)=-dlog(fi(i,ib,iparm))-fi_min(i,ib,iparm)
523 avefi=avefi+fi(i,ib,iparm)
529 write (iout,*) "Parameter set",iparm
531 write (iout,*) "beta=",beta_h(ib,iparm)
533 fi(i,ib,iparm)=fi(i,ib,iparm)-avefi
535 write (iout,'(8f10.5)') (fi(i,ib,iparm),i=1,nR(ib,iparm))
536 write (iout,'(8f10.5)') (f(i,ib,iparm),i=1,nR(ib,iparm))
540 ! Compute the norm of free-energy increments.
545 finorm=finorm+dabs(fi(i,ib,iparm)-f(i,ib,iparm))
546 f(i,ib,iparm)=fi(i,ib,iparm)
551 write (iout,*) 'Iteration',iter,' finorm',finorm
555 call MPI_Bcast(f(1,1,1),MaxR*MaxT_h*nParmSet,
556 & MPI_DOUBLE_PRECISION,Master,
558 call MPI_Bcast(finorm,1,MPI_DOUBLE_PRECISION,Master,
561 ! Exit, if the increment norm is smaller than pre-assigned tolerance.
562 if (finorm.lt.fimin) then
563 write (iout,*) 'Iteration converged'
570 ! Now, put together the histograms from all simulations, in order to get the
571 ! unbiased total histogram.
573 C Determine the minimum free energies
579 c write (9,'(3i5,f10.5)') i,(iparm,potE(i,iparm),iparm=1,nParmSet)
582 write (iout,'(2i5,21f8.2)') i,iparm,
583 & (enetb(k,i,iparm),k=1,21)
585 call restore_parm(iparm)
587 write (iout,*) wsc,wscp,welec,wvdwpp,wang,wtor,wscloc,
588 & wcorr,wcorr5,wcorr6,wturn4,wturn3,wturn6,wel_loc,
589 & wtor_d,wsccor,wbond
592 if (rescale_mode.eq.1) then
593 quot=1.0d0/(beta_h(ib,iparm)*1.987D-3*T0)
600 fT(l)=kfacl/(kfacl-1.0d0+quotl)
603 tt = 1.0d0/(beta_h(ib,iparm)*1.987D-3)
604 ft(6)=(320.0d0+80.0d0*dtanh((tt-320.0d0)/80.0d0))/320.0d0
606 ft(6)=1.0d0/(beta_h(ib,iparm)*1.987D-3*T0)
610 else if (rescale_mode.eq.2) then
611 quot=1.0d0/(T0*beta_h(ib,iparm)*1.987D-3)
615 fT(l)=1.12692801104297249644d0/
616 & dlog(dexp(quotl)+dexp(-quotl))
619 tt = 1.0d0/(beta_h(ib,iparm)*1.987D-3)
620 ft(6)=(320.0d0+80.0d0*dtanh((tt-320.0d0)/80.0d0))/320.0d0
622 ft(6)=1.0d0/(beta_h(ib,iparm)*1.987D-3*T0)
626 c write (iout,*) 1.0d0/(beta_h(ib,iparm)*1.987D-3),ft
627 else if (rescale_mode.eq.0) then
632 write (iout,*) "Error in WHAM_CALC: wrong RESCALE_MODE",
637 evdw=enetb(1,i,iparm)
638 evdw_t=enetb(21,i,iparm)
640 evdw2_14=enetb(17,i,iparm)
641 evdw2=enetb(2,i,iparm)+evdw2_14
643 evdw2=enetb(2,i,iparm)
648 evdw1=enetb(16,i,iparm)
653 ecorr=enetb(4,i,iparm)
654 ecorr5=enetb(5,i,iparm)
655 ecorr6=enetb(6,i,iparm)
656 eel_loc=enetb(7,i,iparm)
657 eello_turn3=enetb(8,i,iparm)
658 eello_turn4=enetb(9,i,iparm)
659 eturn6=enetb(10,i,iparm)
660 ebe=enetb(11,i,iparm)
661 escloc=enetb(12,i,iparm)
662 etors=enetb(13,i,iparm)
663 etors_d=enetb(14,i,iparm)
664 ehpb=enetb(15,i,iparm)
665 estr=enetb(18,i,iparm)
666 esccor=enetb(19,i,iparm)
667 edihcnstr=enetb(20,i,iparm)
669 write (iout,'(3i5,6f5.2,14f12.3)') i,ib,iparm,(ft(l),l=1,6),
670 & evdw+evdw_t,evdw2,ees,evdw1,ecorr,eel_loc,estr,ebe,escloc,
671 & etors,etors_d,eello_turn3,eello_turn4,esccor,edihcnstr
675 etot=wsc*(evdw+ft(6)*evdw_t)+wscp*evdw2+ft(1)*welec*ees
677 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
678 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
679 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
680 & +ft(2)*wturn3*eello_turn3
681 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc
682 & +edihcnstr+ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
685 etot=wsc*(evdw+ft(6)*evdw_t)+wscp*evdw2
686 & +ft(1)*welec*(ees+evdw1)
687 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
688 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
689 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
690 & +ft(2)*wturn3*eello_turn3
691 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc+edihcnstr
692 & +ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
695 c write (iout,*) "i",i," ib",ib,
696 c & " temp",1.0d0/(1.987d-3*beta_h(ib,iparm))," etot",etot,
697 c & " entfac",entfac(i)
698 etot=etot-entfac(i)/beta_h(ib,iparm)
699 if(etot.lt.potEmin_all(ib,iparm)) potEmin_all(ib,iparm)=etot
700 c write (iout,*) "efree",etot," potEmin",potEmin_all(ib,iparm)
705 write (iout,*) "The potEmin array before reduction"
707 write (iout,*) "Parameter set",i
709 write (iout,*) j,PotEmin_all(j,i)
712 write (iout,*) "potEmin_min",potEmin_min
715 C Determine the minimum energes for all parameter sets and temperatures
716 call MPI_AllReduce(potEmin_all(1,1),potEmin_t_all(1,1),
717 & maxT_h*nParmSet,MPI_DOUBLE_PRECISION,MPI_MIN,WHAM_COMM,IERROR)
720 potEmin_all(j,i)=potEmin_t_all(j,i)
724 potEmin_min=potEmin_all(1,1)
727 if (potEmin_all(j,i).lt.potEmin_min)
728 & potEmin_min=potEmin_all(j,i)
732 write (iout,*) "The potEmin array"
734 write (iout,*) "Parameter set",i
736 write (iout,*) j,PotEmin_all(j,i)
739 write (iout,*) "potEmin_min",potEmin_min
751 write (iout,*) "--------------hist"
755 sumW_p(i,iparm)=0.0d0
756 sumE_p(i,iparm)=0.0d0
757 sumEbis_p(i,iparm)=0.0d0
758 sumEsq_p(i,iparm)=0.0d0
760 sumQ_p(j,i,iparm)=0.0d0
761 sumQsq_p(j,i,iparm)=0.0d0
762 sumEQ_p(j,i,iparm)=0.0d0
772 sumEbis(i,iparm)=0.0d0
773 sumEsq(i,iparm)=0.0d0
775 sumQ(j,i,iparm)=0.0d0
776 sumQsq(j,i,iparm)=0.0d0
777 sumEQ(j,i,iparm)=0.0d0
783 c 8/26/05 entropy distribution
788 c ent=-dlog(entfac(t))
790 if (ent.lt.entmin_p) entmin_p=ent
791 if (ent.gt.entmax_p) entmax_p=ent
793 write (iout,*) "entmin",entmin_p," entmax",entmax_p
795 call MPI_Allreduce(entmin_p,entmin,1,MPI_DOUBLE_PRECISION,MPI_MIN,
797 call MPI_Allreduce(entmax_p,entmax,1,MPI_DOUBLE_PRECISION,MPI_MAX,
799 ientmax=entmax-entmin
800 if (ientmax.gt.2000) ientmax=2000
801 write (iout,*) "entmin",entmin," entmax",entmax," ientmax",ientmax
804 c ient=-dlog(entfac(t))-entmin
805 ient=entfac(t)-entmin
806 if (ient.le.2000) histent_p(ient)=histent_p(ient)+1
808 call MPI_Allreduce(histent_p(0),histent(0),ientmax+1,MPI_INTEGER,
809 & MPI_SUM,WHAM_COMM,IERROR)
810 if (me1.eq.Master) then
811 write (iout,*) "Entropy histogram"
813 write(iout,'(f15.4,i10)') entmin+i,histent(i)
821 if (ent.lt.entmin) entmin=ent
822 if (ent.gt.entmax) entmax=ent
824 ientmax=-dlog(entmax)-entmin
825 if (ientmax.gt.2000) ientmax=2000
827 ient=entfac(t)-entmin
828 if (ient.le.2000) histent(ient)=histent(ient)+1
830 write (iout,*) "Entropy histogram"
832 write(iout,'(2f15.4)') entmin+i,histent(i)
837 c write (iout,*) "me1",me1," scount",scount(me1)
863 hrmsrgy(j,i,ib)=0.0d0
865 hrmsrgy_p(j,i,ib)=0.0d0
877 hfin_ent_p(ind)=hfin_ent_p(ind)+dexp(entfac(t))
879 hfin_ent(ind)=hfin_ent(ind)+dexp(entfac(t))
881 call restore_parm(iparm)
882 evdw=enetb(21,t,iparm)
883 evdw_t=enetb(1,t,iparm)
885 evdw2_14=enetb(17,t,iparm)
886 evdw2=enetb(2,t,iparm)+evdw2_14
888 evdw2=enetb(2,t,iparm)
893 evdw1=enetb(16,t,iparm)
898 ecorr=enetb(4,t,iparm)
899 ecorr5=enetb(5,t,iparm)
900 ecorr6=enetb(6,t,iparm)
901 eel_loc=enetb(7,t,iparm)
902 eello_turn3=enetb(8,t,iparm)
903 eello_turn4=enetb(9,t,iparm)
904 eturn6=enetb(10,t,iparm)
905 ebe=enetb(11,t,iparm)
906 escloc=enetb(12,t,iparm)
907 etors=enetb(13,t,iparm)
908 etors_d=enetb(14,t,iparm)
909 ehpb=enetb(15,t,iparm)
910 estr=enetb(18,t,iparm)
911 esccor=enetb(19,t,iparm)
912 edihcnstr=enetb(20,t,iparm)
914 betaT=startGridT+k*delta_T
918 if (rescale_mode.eq.1) then
926 denom=kfacl-1.0d0+quotl
928 ftprim(l)=-l*ft(l)*quotl1/(T0*denom)
929 ftbis(l)=l*kfacl*quotl1*
930 & (2*l*quotl-(l-1)*denom)/(quot*t0*t0*denom**3)
933 ft(6)=(320.0d0+80.0d0*dtanh((betaT-320.0d0)/80.0d0))/
935 ftprim(6)=1.0d0/(320.0d0*dcosh((betaT-320.0d0)/80.0d0)**2)
936 ftbis(6)=-2.0d0*dtanh((betaT-320.0d0)/80.0d0)
937 & /(320.0d0*80.0d0*dcosh((betaT-320.0d0)/80.0d0)**3)
947 else if (rescale_mode.eq.2) then
955 logfac=1.0d0/dlog(eplus+eminus)
956 tanhT=(eplus-eminus)/(eplus+eminus)
957 fT(l)=1.12692801104297249644d0*logfac
958 ftprim(l)=-l*quotl1*ft(l)*tanhT*logfac/T0
959 ftbis(l)=(l-1)*ftprim(l)/(quot*T0)-
960 & 2*l*quotl1/T0*logfac*
961 & (2*l*quotl1*ft(l)/(T0*(eplus+eminus)**2)
965 ft(6)=(320.0d0+80.0d0*dtanh((betaT-320.0d0)/80.0d0))/
967 ftprim(6)=1.0d0/(320.0d0*dcosh((betaT-320.0d0)/80.0d0)**2)
968 ftbis(6)=-2.0d0*dtanh((betaT-320.0d0)/80.0d0)
969 & /(320.0d0*80.0d0*dcosh((betaT-320.0d0)/80.0d0)**3)
979 else if (rescale_mode.eq.0) then
985 write (iout,*) "Error in WHAM_CALC: wrong RESCALE_MODE",
990 c write (iout,*) "ftprim",ftprim
991 c write (iout,*) "ftbis",ftbis
992 betaT=1.0d0/(1.987D-3*betaT)
993 if (betaT.ge.beta_h(1,iparm)) then
994 potEmin=potEmin_all(1,iparm)
995 c write(iout,*) "first",temper,potEmin
996 else if (betaT.lt.beta_h(nT_h(iparm),iparm)) then
997 potEmin=potEmin_all(nT_h(iparm),iparm)
998 c write (iout,*) "last",temper,potEmin
1000 do l=1,nT_h(iparm)-1
1001 if (betaT.le.beta_h(l,iparm) .and.
1002 & betaT.gt.beta_h(l+1,iparm)) then
1003 potEmin=potEmin_all(l,iparm)
1004 c write (iout,*) "l",l,
1005 c & betaT,1.0d0/(1.987D-3*beta_h(l,iparm)),
1006 c & 1.0d0/(1.987D-3*beta_h(l+1,iparm)),temper,potEmin
1011 c write (iout,*) ib," PotEmin",potEmin
1013 etot=wsc*(evdw+ft(6)*evdw_t)+wscp*evdw2+ft(1)*welec*ees
1015 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
1016 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
1017 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
1018 & +ft(2)*wturn3*eello_turn3
1019 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc
1020 & +edihcnstr+ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
1022 eprim=ftprim(6)*evdw_t+ftprim(1)*welec*ees
1023 & +ftprim(1)*wtor*etors+
1024 & ftprim(3)*wcorr*ecorr+ftprim(4)*wcorr5*ecorr5+
1025 & ftprim(5)*wcorr6*ecorr6+ftprim(3)*wturn4*eello_turn4+
1026 & ftprim(2)*wturn3*eello_turn3+ftprim(5)*wturn6*eturn6+
1027 & ftprim(2)*wel_loc*eel_loc+ftprim(2)*wtor_d*etors_d+
1028 & ftprim(1)*wsccor*esccor
1029 ebis=ftbis(1)*welec*ees+ftbis(1)*wtor*etors+
1030 & ftbis(3)*wcorr*ecorr+ftbis(4)*wcorr5*ecorr5+
1031 & ftbis(5)*wcorr6*ecorr6+ftbis(3)*wturn4*eello_turn4+
1032 & ftbis(2)*wturn3*eello_turn3+ftbis(5)*wturn6*eturn6+
1033 & ftbis(2)*wel_loc*eel_loc+ftbis(2)*wtor_d*etors_d+
1034 & ftbis(1)*wsccor*esccor
1036 etot=wsc*(evdw+ft(6)*evdw_t)+wscp*evdw2
1037 & +ft(1)*welec*(ees+evdw1)
1038 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
1039 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
1040 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
1041 & +ft(2)*wturn3*eello_turn3
1042 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc+edihcnstr
1043 & +ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
1045 eprim=ftprim(6)*evdw_t+ftprim(1)*welec*(ees+evdw1)
1046 & +ftprim(1)*wtor*etors+
1047 & ftprim(3)*wcorr*ecorr+ftprim(4)*wcorr5*ecorr5+
1048 & ftprim(5)*wcorr6*ecorr6+ftprim(3)*wturn4*eello_turn4+
1049 & ftprim(2)*wturn3*eello_turn3+ftprim(5)*wturn6*eturn6+
1050 & ftprim(2)*wel_loc*eel_loc+ftprim(2)*wtor_d*etors_d+
1051 & ftprim(1)*wsccor*esccor
1052 ebis=ftbis(1)*welec*(ees+evdw1)+ftbis(1)*wtor*etors+
1053 & ftbis(3)*wcorr*ecorr+ftbis(4)*wcorr5*ecorr5+
1054 & ftbis(5)*wcorr6*ecorr6+ftbis(3)*wturn4*eello_turn4+
1055 & ftbis(2)*wturn3*eello_turn3+ftbis(5)*wturn6*eturn6+
1056 & ftbis(2)*wel_loc*eel_loc+ftbis(2)*wtor_d*etors_d+
1057 & ftprim(1)*wsccor*esccor
1059 weight=dexp(-betaT*(etot-potEmin)+entfac(t))
1061 write (iout,*) "iparm",iparm," t",t," temper",temper,
1062 & " etot",etot," entfac",entfac(t),
1063 & " efree",etot-entfac(t)/betaT," potEmin",potEmin,
1064 & " boltz",-betaT*(etot-potEmin)+entfac(t),
1065 & " weight",weight," ebis",ebis
1067 etot=etot-temper*eprim
1069 sumW_p(k,iparm)=sumW_p(k,iparm)+weight
1070 sumE_p(k,iparm)=sumE_p(k,iparm)+etot*weight
1071 sumEbis_p(k,iparm)=sumEbis_p(k,iparm)+ebis*weight
1072 sumEsq_p(k,iparm)=sumEsq_p(k,iparm)+etot**2*weight
1074 sumQ_p(j,k,iparm)=sumQ_p(j,k,iparm)+q(j,t)*weight
1075 sumQsq_p(j,k,iparm)=sumQsq_p(j,k,iparm)+q(j,t)**2*weight
1076 sumEQ_p(j,k,iparm)=sumEQ_p(j,k,iparm)
1077 & +etot*q(j,t)*weight
1080 sumW(k,iparm)=sumW(k,iparm)+weight
1081 sumE(k,iparm)=sumE(k,iparm)+etot*weight
1082 sumEbis(k,iparm)=sumEbis(k,iparm)+ebis*weight
1083 sumEsq(k,iparm)=sumEsq(k,iparm)+etot**2*weight
1085 sumQ(j,k,iparm)=sumQ(j,k,iparm)+q(j,t)*weight
1086 sumQsq(j,k,iparm)=sumQsq(j,k,iparm)+q(j,t)**2*weight
1087 sumEQ(j,k,iparm)=sumEQ(j,k,iparm)
1088 & +etot*q(j,t)*weight
1092 indE = aint(potE(t,iparm)-aint(potEmin))
1093 if (indE.ge.0 .and. indE.le.maxinde) then
1094 if (indE.gt.upindE_p) upindE_p=indE
1095 histE_p(indE)=histE_p(indE)+dexp(-entfac(t))
1099 potEmin=potEmin_all(ib,iparm)
1100 expfac=dexp(-beta_h(ib,iparm)*(etot-potEmin)+entfac(t))
1101 hfin_p(ind,ib)=hfin_p(ind,ib)+
1102 & dexp(-beta_h(ib,iparm)*(etot-potEmin)+entfac(t))
1104 indrgy=dint((q(nQ+2,t)-rgymin)/deltrgy)
1105 indrms=dint((q(nQ+1,t)-rmsmin)/deltrms)
1106 hrmsrgy_p(indrgy,indrms,ib)=
1107 & hrmsrgy_p(indrgy,indrms,ib)+expfac
1112 potEmin=potEmin_all(ib,iparm)
1113 expfac=dexp(-beta_h(ib,iparm)*(etot-potEmin)+entfac(t))
1114 hfin(ind,ib)=hfin(ind,ib)+
1115 & dexp(-beta_h(ib,iparm)*(etot-potEmin)+entfac(t))
1117 indrgy=dint((q(nQ+2,t)-rgymin)/deltrgy)
1118 indrms=dint((q(nQ+1,t)-rmsmin)/deltrms)
1119 hrmsrgy(indrgy,indrms,ib)=
1120 & hrmsrgy(indrgy,indrms,ib)+expfac
1126 if (histout) call MPI_Reduce(hfin_p(0,ib),hfin(0,ib),nbin,
1127 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
1129 call MPI_Reduce(hrmsrgy_p(0,0,ib),hrmsrgy(0,0,ib),
1130 & (MaxBinRgy+1)*(nbin_rms+1),MPI_DOUBLE_PRECISION,MPI_SUM,Master,
1134 call MPI_Reduce(upindE_p,upindE,1,
1135 & MPI_INTEGER,MPI_MAX,Master,WHAM_COMM,IERROR)
1136 call MPI_Reduce(histE_p(0),histE(0),maxindE,
1137 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
1139 if (me1.eq.master) then
1143 write (iout,'(6x,$)')
1144 write (iout,'(f20.2,$)') (1.0d0/(1.987D-3*beta_h(ib,iparm)),
1148 write (iout,'(/a)') 'Final histograms'
1150 if (nslice.eq.1) then
1151 if (separate_parset) then
1152 write(licz3,"(bz,i3.3)") myparm
1153 histname=prefix(:ilen(prefix))//'_par'//licz3//'.hist'
1155 histname=prefix(:ilen(prefix))//'.hist'
1158 if (separate_parset) then
1159 write(licz3,"(bz,i3.3)") myparm
1160 histname=prefix(:ilen(prefix))//'_par'//licz3//
1161 & '_slice_'//licz2//'.hist'
1163 histname=prefix(:ilen(prefix))//'_slice_'//licz2//'.hist'
1166 #if defined(AIX) || defined(PGI)
1167 open (ihist,file=histname,position='append')
1169 open (ihist,file=histname,access='append')
1177 sumH=sumH+hfin(t,ib)
1179 if (sumH.gt.0.0d0) then
1181 jj = mod(liczba,nbin1)
1183 write (iout,'(f6.3,$)') dmin+(jj+0.5d0)*delta
1185 & write (ihist,'(f6.3,$)') dmin+(jj+0.5d0)*delta
1188 write (iout,'(e20.10,$)') hfin(t,ib)
1189 if (histfile) write (ihist,'(e20.10,$)') hfin(t,ib)
1191 write (iout,'(i5)') iparm
1192 if (histfile) write (ihist,'(i5)') iparm
1199 if (nslice.eq.1) then
1200 if (separate_parset) then
1201 write(licz3,"(bz,i3.3)") myparm
1202 histname=prefix(:ilen(prefix))//"_par"//licz3//'.ent'
1204 histname=prefix(:ilen(prefix))//'.ent'
1207 if (separate_parset) then
1208 write(licz3,"(bz,i3.3)") myparm
1209 histname=prefix(:ilen(prefix))//'par_'//licz3//
1210 & '_slice_'//licz2//'.ent'
1212 histname=prefix(:ilen(prefix))//'_slice_'//licz2//'.ent'
1215 #if defined(AIX) || defined(PGI)
1216 open (ihist,file=histname,position='append')
1218 open (ihist,file=histname,access='append')
1220 write (ihist,'(a)') "# Microcanonical entropy"
1222 write (ihist,'(f8.0,$)') dint(potEmin)+i
1223 if (histE(i).gt.0.0e0) then
1224 write (ihist,'(f15.5,$)') dlog(histE(i))
1226 write (ihist,'(f15.5,$)') 0.0d0
1232 write (iout,*) "Microcanonical entropy"
1234 write (iout,'(f8.0,$)') dint(potEmin)+i
1235 if (histE(i).gt.0.0e0) then
1236 write (iout,'(f15.5,$)') dlog(histE(i))
1238 write (iout,'(f15.5,$)') 0.0d0
1243 if (nslice.eq.1) then
1244 if (separate_parset) then
1245 write(licz3,"(bz,i3.3)") myparm
1246 histname=prefix(:ilen(prefix))//'_par'//licz3//'.rmsrgy'
1248 histname=prefix(:ilen(prefix))//'.rmsrgy'
1251 if (separate_parset) then
1252 write(licz3,"(bz,i3.3)") myparm
1253 histname=prefix(:ilen(prefix))//'_par'//licz3//
1254 & '_slice_'//licz2//'.rmsrgy'
1256 histname=prefix(:ilen(prefix))//'_slice_'//licz2//'.rmsrgy'
1259 #if defined(AIX) || defined(PGI)
1260 open (ihist,file=histname,position='append')
1262 open (ihist,file=histname,access='append')
1266 write(ihist,'(2f8.2,$)')
1267 & rgymin+deltrgy*j,rmsmin+deltrms*i
1269 if (hrmsrgy(j,i,ib).gt.0.0d0) then
1270 write(ihist,'(e14.5,$)')
1271 & -dlog(hrmsrgy(j,i,ib))/beta_h(ib,iparm)
1274 write(ihist,'(e14.5,$)') 1.0d6
1277 write (ihist,'(i2)') iparm
1285 call MPI_Reduce(hfin_ent_p(0),hfin_ent(0),nbin,
1286 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
1287 call MPI_Reduce(sumW_p(0,1),sumW(0,1),(nGridT+1)*nParmSet,
1288 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
1289 call MPI_Reduce(sumE_p(0,1),sumE(0,1),(nGridT+1)*nParmSet,
1290 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
1291 call MPI_Reduce(sumEbis_p(0,1),sumEbis(0,1),(nGridT+1)*nParmSet,
1292 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
1293 call MPI_Reduce(sumEsq_p(0,1),sumEsq(0,1),(nGridT+1)*nParmSet,
1294 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
1295 call MPI_Reduce(sumQ_p(1,0,1),sumQ(1,0,1),
1296 & MaxQ1*(nGridT+1)*nParmSet,MPI_DOUBLE_PRECISION,MPI_SUM,Master,
1298 call MPI_Reduce(sumQsq_p(1,0,1),sumQsq(1,0,1),
1299 & MaxQ1*(nGridT+1)*nParmSet,MPI_DOUBLE_PRECISION,MPI_SUM,Master,
1301 call MPI_Reduce(sumEQ_p(1,0,1),sumEQ(1,0,1),
1302 & MaxQ1*(nGridT+1)*nParmSet,MPI_DOUBLE_PRECISION,MPI_SUM,Master,
1304 if (me.eq.master) then
1306 write (iout,'(/a)') 'Thermal characteristics of folding'
1307 if (nslice.eq.1) then
1310 nazwa=prefix(:ilen(prefix))//"_slice_"//licz2
1313 if (nparmset.eq.1 .and. .not.separate_parset) then
1314 nazwa=nazwa(:iln)//".thermal"
1315 else if (nparmset.eq.1 .and. separate_parset) then
1316 write(licz3,"(bz,i3.3)") myparm
1317 nazwa=nazwa(:iln)//"_par_"//licz3//".thermal"
1320 if (nparmset.gt.1) then
1321 write(licz3,"(bz,i3.3)") iparm
1322 nazwa=nazwa(:iln)//"_par_"//licz3//".thermal"
1325 if (separate_parset) then
1326 write (iout,'(a,i3)') "Parameter set",myparm
1328 write (iout,'(a,i3)') "Parameter set",iparm
1331 betaT=1.0d0/(1.987D-3*(startGridT+i*delta_T))
1332 if (betaT.ge.beta_h(1,iparm)) then
1333 potEmin=potEmin_all(1,iparm)
1334 else if (betaT.lt.beta_h(nT_h(iparm),iparm)) then
1335 potEmin=potEmin_all(nT_h(iparm),iparm)
1337 do l=1,nT_h(iparm)-1
1338 if (betaT.le.beta_h(l,iparm) .and.
1339 & betaT.gt.beta_h(l+1,iparm)) then
1340 potEmin=potEmin_all(l,iparm)
1345 sumE(i,iparm)=sumE(i,iparm)/sumW(i,iparm)
1346 sumEbis(i,iparm)=(startGridT+i*delta_T)*sumEbis(i,iparm)/
1348 sumEsq(i,iparm)=(sumEsq(i,iparm)/sumW(i,iparm)
1349 & -sumE(i,iparm)**2)/(1.987D-3*(startGridT+i*delta_T)**2)
1351 sumQ(j,i,iparm)=sumQ(j,i,iparm)/sumW(i,iparm)
1352 sumQsq(j,i,iparm)=sumQsq(j,i,iparm)/sumW(i,iparm)
1353 & -sumQ(j,i,iparm)**2
1354 sumEQ(j,i,iparm)=sumEQ(j,i,iparm)/sumW(i,iparm)
1355 & -sumQ(j,i,iparm)*sumE(i,iparm)
1357 sumW(i,iparm)=-dlog(sumW(i,iparm))*(1.987D-3*
1358 & (startGridT+i*delta_T))+potEmin
1359 write (iout,'(f7.1,2f15.5,$)') startGridT+i*delta_T,
1360 & sumW(i,iparm),sumE(i,iparm)
1361 write (iout,'(f10.5,$)') (sumQ(j,i,iparm),j=1,nQ+2)
1362 write (iout,'(e15.5,$)') sumEsq(i,iparm)-sumEbis(i,iparm),
1363 & (sumQsq(j,i,iparm),j=1,nQ+2),(sumEQ(j,i,iparm),j=1,nQ+2)
1365 write (34,'(f7.1,2f15.5,$)') startGridT+i*delta_T,
1366 & sumW(i,iparm),sumE(i,iparm)
1367 write (34,'(f10.5,$)') (sumQ(j,i,iparm),j=1,nQ+2)
1368 write (34,'(e15.5,$)') sumEsq(i,iparm)-sumEbis(i,iparm),
1369 & (sumQsq(j,i,iparm),j=1,nQ+2),(sumEQ(j,i,iparm),j=1,nQ+2)
1376 if (hfin_ent(t).gt.0.0d0) then
1378 jj = mod(liczba,nbin1)
1379 write (iout,'(f6.3,e20.10," ent")') dmin+(jj+0.5d0)*delta,
1381 if (histfile) write (ihist,'(f6.3,e20.10," ent")')
1382 & dmin+(jj+0.5d0)*delta,
1386 if (histfile) close(ihist)
1390 ! Write data for zscore
1391 if (nslice.eq.1) then
1392 zscname=prefix(:ilen(prefix))//".zsc"
1394 zscname=prefix(:ilen(prefix))//"_slice_"//licz2//".zsc"
1396 #if defined(AIX) || defined(PGI)
1397 open (izsc,file=prefix(:ilen(prefix))//'.zsc',position='append')
1399 open (izsc,file=prefix(:ilen(prefix))//'.zsc',access='append')
1401 write (izsc,'("NQ=",i1," NPARM=",i1)') nQ,nParmSet
1403 write (izsc,'("NT=",i1)') nT_h(iparm)
1405 write (izsc,'("TEMP=",f6.1," NR=",i2," SNK=",$)')
1406 & 1.0d0/(beta_h(ib,iparm)*1.987D-3),nR(ib,iparm)
1407 jj = min0(nR(ib,iparm),7)
1408 write (izsc,'(i8,$)') (snk(i,ib,iparm,islice),i=1,jj)
1409 write (izsc,'(a1,$)') (" ",i=22+8*jj+1,79)
1410 write (izsc,'("&")')
1411 if (nR(ib,iparm).gt.7) then
1412 do ii=8,nR(ib,iparm),9
1413 jj = min0(nR(ib,iparm),ii+8)
1414 write (izsc,'(i8,$)') (snk(i,ib,iparm,islice),i=ii,jj)
1415 write (izsc,'(a1,$') (" ",i=(jj-ii+1)*8+1,79)
1416 write (izsc,'("&")')
1419 write (izsc,'("FI=",$)')
1420 jj=min0(nR(ib,iparm),7)
1421 write (izsc,'(f10.5,$)') (fi(i,ib,iparm),i=1,jj)
1422 write (izsc,'(a1,$)') (" ",i=3+10*jj+1,79)
1423 write (izsc,'("&")')
1424 if (nR(ib,iparm).gt.7) then
1425 do ii=8,nR(ib,iparm),9
1426 jj = min0(nR(ib,iparm),ii+8)
1427 write (izsc,'(f10.5,$)') (fi(i,ib,iparm),i=ii,jj)
1428 if (jj.eq.nR(ib,iparm)) then
1431 write (izsc,'(a1,$)') (" ",i=10*(jj-ii+1)+1,79)
1432 write (izsc,'(t80,"&")')
1437 write (izsc,'("KH=",$)')
1438 write (izsc,'(f7.2,$)') (Kh(j,i,ib,iparm),j=1,nQ)
1439 write (izsc,'(" Q0=",$)')
1440 write (izsc,'(f7.5,$)') (q0(j,i,ib,iparm),j=1,nQ)