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"
18 parameter (NGridT=400)
19 integer MaxBinRms,MaxBinRgy
20 parameter (MaxBinRms=100,MaxBinRgy=100)
22 c parameter (MaxHdim=200)
24 parameter (maxinde=200)
28 integer ierror,errcode,status(MPI_STATUS_SIZE)
30 include "COMMON.CONTROL"
31 include "COMMON.IOUNITS"
33 include "COMMON.ENERGIES"
34 include "COMMON.FFIELD"
35 include "COMMON.SBRIDGE"
37 include "COMMON.ENEPS"
38 include "COMMON.SHIELD"
39 integer MaxPoint,MaxPointProc
40 parameter (MaxPoint=MaxStr,
41 & MaxPointProc=MaxStr_Proc)
42 double precision finorm_max,potfac,entmin,entmax,expfac,vf
43 parameter (finorm_max=1.0d0)
45 integer i,ii,j,jj,k,kk,l,m,ind,iter,t,tmax,ient,ientmax,iln
46 integer start,end,iharm,ib,iib,nbin1,nbin,nbin_rms,nbin_rgy,
47 & nbin_rmsrgy,liczba,iparm,nFi,indrgy,indrms
48 integer htot(0:MaxHdim),histent(0:2000)
49 double precision v(MaxPointProc,MaxR,MaxT_h,Max_Parm)
50 double precision energia(0:max_ene)
52 integer tmax_t,upindE_p
53 double precision fi_p(MaxR,MaxT_h,Max_Parm)
54 double precision sumW_p(0:nGridT,Max_Parm),
55 & sumE_p(0:nGridT,Max_Parm),sumEsq_p(0:nGridT,Max_Parm),
56 & sumQ_p(MaxQ1,0:nGridT,Max_Parm),
57 & sumQsq_p(MaxQ1,0:nGridT,Max_Parm),
58 & sumEQ_p(MaxQ1,0:nGridT,Max_Parm),
59 & sumEprim_p(MaxQ1,0:nGridT,Max_Parm),
60 & sumEbis_p(0:nGridT,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,entmin_p,entmax_p
66 double precision ePMF,ePMF_q
67 integer histent_p(0:2000)
68 logical lprint /.true./
70 double precision delta_T /1.0d0/
71 double precision rgymin,rmsmin,rgymax,rmsmax
72 double precision sumW(0:NGridT,Max_Parm),sumE(0:NGridT,Max_Parm),
73 & sumEsq(0:NGridT,Max_Parm),sumQ(MaxQ1,0:NGridT,Max_Parm),
74 & sumQsq(MaxQ1,0:NGridT,Max_Parm),sumEQ(MaxQ1,0:NGridT,Max_Parm),
75 & sumEprim(0:NGridT,Max_Parm),sumEbis(0:NGridT,Max_Parm),betaT,
77 double precision fi(MaxR,maxT_h,Max_Parm),
78 & dd,dd1,dd2,hh,dmin,denom,finorm,avefi,pom,
79 & hfin(0:MaxHdim,maxT_h),histE(0:maxindE),
80 & hrmsrgy(0:MaxBinRgy,0:MaxBinRms,maxT_h),
82 & hfin_ent(0:MaxHdim),vmax,aux
83 double precision fT(6),fTprim(6),fTbis(6),quot,quotl1,quotl,kfacl,
84 & eprim,ebis,temper,kfac/2.4d0/,T0/300.0d0/,startGridT/200.0d0/,
85 & eplus,eminus,logfac,tanhT,tt
86 double precision etot,evdw,evdw_t,evdw2,ees,evdw1,ebe,etors,
87 & escloc,ehpb,ecorr,ecorr5,ecorr6,eello_turn4,eello_turn3,
88 & eturn6,eel_loc,edihcnstr,etors_d,estr,evdw2_14,esccor,
91 integer ind_point(maxpoint),upindE,indE
100 write(licz2,'(bz,i2.2)') islice
102 write (iout,'(//80(1h-)/"Solving WHAM equations for slice",
103 & i2/80(1h-)//)') islice
104 write (iout,*) "delta",delta," nbin1",nbin1
105 write (iout,*) "MaxN",MaxN," MaxQ",MaxQ," MaHdim",MaxHdim
123 c write (iout,*) "i",i," potE",(potE(i,j),j=1,nParmset)
125 if (potE(i,j).le.potEmin) potEmin=potE(i,j)
127 if (q(nQ+1,i).lt.rmsmin) rmsmin=q(nQ+1,i)
128 if (q(nQ+1,i).gt.rmsmax) rmsmax=q(nQ+1,i)
129 if (q(nQ+2,i).lt.rgymin) rgymin=q(nQ+2,i)
130 if (q(nQ+2,i).gt.rgymax) rgymax=q(nQ+2,i)
133 ind=(q(j,i)-dmin+1.0d-8)/delta
135 ind_point(i)=ind_point(i)+ind
137 ind_point(i)=ind_point(i)+nbin1**(j-1)*ind
139 c write (iout,*) "i",i," j",j," q",q(j,i)," ind_point",
142 if (ind_point(i).lt.0 .or. ind_point(i).gt.MaxHdim) then
143 write (iout,*) "Error - index exceeds range for point",i,
144 & " q=",q(j,i)," ind",ind_point(i)
146 write (iout,*) "Processor",me1
148 call MPI_Abort(MPI_COMM_WORLD, Ierror, Errcode )
153 if (ind_point(i).gt.tmax) tmax=ind_point(i)
154 htot(ind_point(i))=htot(ind_point(i))+1
156 write (iout,*) "i",i,"q",(q(j,i),j=1,nQ)," ind",ind_point(i),
157 & " htot",htot(ind_point(i))
162 write (iout,*) "potEmin before reduce",potEmin
164 write (iout,'(a)') "Numbers of counts in Q bins"
166 if (htot(t).gt.0) then
167 write (iout,'(i15,$)') t
170 jj = mod(liczba,nbin1)
172 write (iout,'(i5,$)') jj
174 write (iout,'(i8)') htot(t)
178 write (iout,'(a,i3)') "Number of data points for parameter set",
180 write (iout,'(i7,$)') ((snk(m,ib,iparm,islice),m=1,nr(ib,iparm)),
182 write (iout,'(i8)') stot(islice)
188 call MPI_AllReduce(tmax,tmax_t,1,MPI_INTEGER,MPI_MAX,
191 call MPI_AllReduce(potEmin,potEmin_t,1,MPI_DOUBLE_PRECISION,
192 & MPI_MIN,WHAM_COMM,IERROR)
193 call MPI_AllReduce(rmsmin,rmsmin_t,1,MPI_DOUBLE_PRECISION,
194 & MPI_MIN,WHAM_COMM,IERROR)
195 call MPI_AllReduce(rmsmax,rmsmax_t,1,MPI_DOUBLE_PRECISION,
196 & MPI_MAX,WHAM_COMM,IERROR)
197 call MPI_AllReduce(rgymin,rgymin_t,1,MPI_DOUBLE_PRECISION,
198 & MPI_MIN,WHAM_COMM,IERROR)
199 call MPI_AllReduce(rgymax,rgymax_t,1,MPI_DOUBLE_PRECISION,
200 & MPI_MAX,WHAM_COMM,IERROR)
201 c potEmin=potEmin_t/2
207 write (iout,*) "potEmin",potEmin
209 rmsmin=deltrms*dint(rmsmin/deltrms)
210 rmsmax=deltrms*dint(rmsmax/deltrms)
211 rgymin=deltrms*dint(rgymin/deltrgy)
212 rgymax=deltrms*dint(rgymax/deltrgy)
213 nbin_rms=(rmsmax-rmsmin)/deltrms
214 nbin_rgy=(rgymax-rgymin)/deltrgy
215 write (iout,*) "rmsmin",rmsmin," rmsmax",rmsmax," rgymin",rgymin,
216 & " rgymax",rgymax," nbin_rms",nbin_rms," nbin_rgy",nbin_rgy
223 write (iout,*) "nFi",nFi
224 ! Compute the Boltzmann factor corresponing to restrain potentials in different
231 c write (9,'(3i5,f10.5)') i,(iparm,potE(i,iparm),iparm=1,nParmSet)
234 write (iout,'(2i5,21f8.2)') i,iparm,
235 & (enetb(k,i,iparm),k=1,22)
237 call restore_parm(iparm)
239 write (iout,*) wsc,wscp,welec,wvdwpp,wang,wtor,wscloc,
240 & wcorr,wcorr5,wcorr6,wturn4,wturn3,wturn6,wel_loc,
241 & wtor_d,wsccor,wbond,wsaxs
244 if (rescale_mode.eq.1) then
245 quot=1.0d0/(beta_h(ib,iparm)*1.987D-3*T0)
252 fT(l)=kfacl/(kfacl-1.0d0+quotl)
255 tt = 1.0d0/(beta_h(ib,iparm)*1.987D-3)
256 ft(6)=(320.0d0+80.0d0*dtanh((tt-320.0d0)/80.0d0))/320.0d0
258 ft(6)=1.0d0/(beta_h(ib,iparm)*1.987D-3*T0)
262 else if (rescale_mode.eq.2) then
263 quot=1.0d0/(T0*beta_h(ib,iparm)*1.987D-3)
267 fT(l)=1.12692801104297249644d0/
268 & dlog(dexp(quotl)+dexp(-quotl))
271 tt = 1.0d0/(beta_h(ib,iparm)*1.987D-3)
272 ft(6)=(320.0d0+80.0d0*dtanh((tt-320.0d0)/80.0d0))/320.0d0
274 ft(6)=1.0d0/(beta_h(ib,iparm)*1.987D-3*T0)
278 c write (iout,*) 1.0d0/(beta_h(ib,iparm)*1.987D-3),ft
279 else if (rescale_mode.eq.0) then
284 write (iout,*) "Error in WHAM_CALC: wrong RESCALE_MODE",
289 evdw=enetb(1,i,iparm)
290 evdw_t=enetb(21,i,iparm)
292 evdw2_14=enetb(17,i,iparm)
293 evdw2=enetb(2,i,iparm)+evdw2_14
295 evdw2=enetb(2,i,iparm)
300 evdw1=enetb(16,i,iparm)
305 ecorr=enetb(4,i,iparm)
306 ecorr5=enetb(5,i,iparm)
307 ecorr6=enetb(6,i,iparm)
308 eel_loc=enetb(7,i,iparm)
309 eello_turn3=enetb(8,i,iparm)
310 eello_turn4=enetb(9,i,iparm)
311 eturn6=enetb(10,i,iparm)
312 ebe=enetb(11,i,iparm)
313 escloc=enetb(12,i,iparm)
314 etors=enetb(13,i,iparm)
315 etors_d=enetb(14,i,iparm)
316 ehpb=enetb(15,i,iparm)
317 estr=enetb(18,i,iparm)
318 esccor=enetb(19,i,iparm)
319 edihcnstr=enetb(20,i,iparm)
320 eliptran=enetb(22,i,iparm)
321 esaxs=enetb(26,i,iparm)
324 write (iout,'(3i5,6f5.2,14f12.3)') i,ib,iparm,(ft(l),l=1,6),
325 & evdw+evdw_t,evdw2,ees,evdw1,ecorr,eel_loc,estr,ebe,escloc,
326 & etors,etors_d,eello_turn3,eello_turn4,esccor,esaxs
330 if (shield_mode.gt.0) then
331 etot=ft(1)*wsc*(evdw+ft(6)*evdw_t)+ft(1)*wscp*evdw2
333 & +ft(1)*wvdwpp*evdw1
334 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
335 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
336 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
337 & +ft(2)*wturn3*eello_turn3
338 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc
339 & +edihcnstr+ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
340 & +wbond*estr+wliptran*eliptran+wsaxs*esaxs
342 etot=wsc*(evdw+ft(6)*evdw_t)+wscp*evdw2+ft(1)*welec*ees
344 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
345 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
346 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
347 & +ft(2)*wturn3*eello_turn3
348 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc
349 & +edihcnstr+ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
350 & +wbond*estr+wliptran*eliptran+wsaxs*esaxs
353 if (shield_mode.gt.0) then
354 etot=ft(1)*wsc*(evdw+ft(6)*evdw_t)+ft(1)*wscp*evdw2
355 & +ft(1)*welec*(ees+evdw1)
356 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
357 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
358 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
359 & +ft(2)*wturn3*eello_turn3
360 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc+edihcnstr
361 & +ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
362 & +wbond*estr+wliptran*eliptran+wsaxs*esaxs
364 etot=wsc*(evdw+ft(6)*evdw_t)+wscp*evdw2
365 & +ft(1)*welec*(ees+evdw1)
366 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
367 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
368 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
369 & +ft(2)*wturn3*eello_turn3
370 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc+edihcnstr
371 & +ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
372 & +wbond*estr+wliptran*eliptran+wsaxs*esaxs
377 write (iout,*) i,iparm,1.0d0/(beta_h(ib,iparm)*1.987D-3),
381 if (iparm.eq.1 .and. ib.eq.1) then
382 write (iout,*)"Conformation",i
385 energia(k)=enetb(k,i,iparm)
387 call enerprint(energia(0),fT)
394 Econstr=Econstr+Kh(j,kk,ib,iparm)
395 & *(dd-q0(j,kk,ib,iparm))**2
397 c Adaptive potential contribution
399 call PMF_energy(q(1,i),ib,kk,iparm,ePMF,ePMF_q)
403 & -beta_h(ib,iparm)*(etot-potEmin+Econstr)
405 write (iout,'(4i5,4e15.5)') i,kk,ib,iparm,
406 & etot,potEmin,etot-potEmin,v(i,kk,ib,iparm)
412 ! Simple iteration to calculate free energies corresponding to all simulation
416 ! Compute new free-energy values corresponding to the righ-hand side of the
417 ! equation and their derivatives.
418 write (iout,*) "------------------------fi"
428 vf=v(t,l,k,i)+f(l,k,i)
429 if (vf.gt.vmax) vmax=vf
437 aux=f(l,k,i)+v(t,l,k,i)-vmax
439 & denom=denom+snk(l,k,i,islice)*dexp(aux)
443 entfac(t)=-dlog(denom)-vmax
445 write (iout,*) t,"vmax",vmax," denom",denom,"entfac",entfac(t)
450 do ii=1,nR(iib,iparm)
452 fi_p(ii,iib,iparm)=0.0d0
454 fi_p(ii,iib,iparm)=fi_p(ii,iib,iparm)
455 & +dexp(v(t,ii,iib,iparm)+entfac(t))
457 write (iout,'(4i5,3e15.5)') t,ii,iib,iparm,
458 & v(t,ii,iib,iparm),entfac(t),fi_p(ii,iib,iparm)
462 fi(ii,iib,iparm)=0.0d0
464 fi(ii,iib,iparm)=fi(ii,iib,iparm)
465 & +dexp(v(t,ii,iib,iparm)+entfac(t))
474 write (iout,*) "fi before MPI_Reduce me",me,' master',master
476 do ib=1,nT_h(nparmset)
477 write (iout,*) "iparm",iparm," ib",ib
478 write (iout,*) "beta=",beta_h(ib,iparm)
479 write (iout,'(8e15.5)') (fi_p(i,ib,iparm),i=1,nR(ib,iparm))
483 c write (iout,*) "REDUCE size",maxR,MaxT_h,nParmSet,
484 c & maxR*MaxT_h*nParmSet
485 c write (iout,*) "MPI_COMM_WORLD",MPI_COMM_WORLD,
486 c & " WHAM_COMM",WHAM_COMM
487 call MPI_Reduce(fi_p(1,1,1),fi(1,1,1),maxR*MaxT_h*nParmSet,
488 & MPI_DOUBLE_PRECISION,
489 & MPI_SUM,Master,WHAM_COMM,IERROR)
491 write (iout,*) "fi after MPI_Reduce nparmset",nparmset
493 write (iout,*) "iparm",iparm
495 write (iout,*) "beta=",beta_h(ib,iparm)
496 write (iout,'(8e15.5)') (fi(i,ib,iparm),i=1,nR(ib,iparm))
500 if (me1.eq.Master) then
506 fi(i,ib,iparm)=-dlog(fi(i,ib,iparm))
507 avefi=avefi+fi(i,ib,iparm)
513 write (iout,*) "Parameter set",iparm
515 write (iout,*) "beta=",beta_h(ib,iparm)
517 fi(i,ib,iparm)=fi(i,ib,iparm)-avefi
519 write (iout,'(8f10.5)') (fi(i,ib,iparm),i=1,nR(ib,iparm))
520 write (iout,'(8f10.5)') (f(i,ib,iparm),i=1,nR(ib,iparm))
524 ! Compute the norm of free-energy increments.
529 finorm=finorm+dabs(fi(i,ib,iparm)-f(i,ib,iparm))
530 f(i,ib,iparm)=fi(i,ib,iparm)
535 write (iout,*) 'Iteration',iter,' finorm',finorm
539 call MPI_Bcast(f(1,1,1),MaxR*MaxT_h*nParmSet,
540 & MPI_DOUBLE_PRECISION,Master,
542 call MPI_Bcast(finorm,1,MPI_DOUBLE_PRECISION,Master,
545 ! Exit, if the increment norm is smaller than pre-assigned tolerance.
546 if (finorm.lt.fimin) then
547 write (iout,*) 'Iteration converged'
554 ! Now, put together the histograms from all simulations, in order to get the
555 ! unbiased total histogram.
565 write (iout,*) "--------------hist"
569 sumW_p(i,iparm)=0.0d0
570 sumE_p(i,iparm)=0.0d0
571 sumEbis_p(i,iparm)=0.0d0
572 sumEsq_p(i,iparm)=0.0d0
574 sumQ_p(j,i,iparm)=0.0d0
575 sumQsq_p(j,i,iparm)=0.0d0
576 sumEQ_p(j,i,iparm)=0.0d0
586 sumEbis(i,iparm)=0.0d0
587 sumEsq(i,iparm)=0.0d0
589 sumQ(j,i,iparm)=0.0d0
590 sumQsq(j,i,iparm)=0.0d0
591 sumEQ(j,i,iparm)=0.0d0
597 c 8/26/05 entropy distribution
602 c ent=-dlog(entfac(t))
604 if (ent.lt.entmin_p) entmin_p=ent
605 if (ent.gt.entmax_p) entmax_p=ent
607 write (iout,*) "entmin",entmin_p," entmax",entmax_p
609 call MPI_Allreduce(entmin_p,entmin,1,MPI_DOUBLE_PRECISION,MPI_MIN,
611 call MPI_Allreduce(entmax_p,entmax,1,MPI_DOUBLE_PRECISION,MPI_MAX,
613 ientmax=entmax-entmin
614 if (ientmax.gt.2000) ientmax=2000
615 write (iout,*) "entmin",entmin," entmax",entmax," ientmax",ientmax
618 c ient=-dlog(entfac(t))-entmin
619 ient=entfac(t)-entmin
620 if (ient.le.2000) histent_p(ient)=histent_p(ient)+1
622 call MPI_Allreduce(histent_p(0),histent(0),ientmax+1,MPI_INTEGER,
623 & MPI_SUM,WHAM_COMM,IERROR)
624 if (me1.eq.Master) then
625 write (iout,*) "Entropy histogram"
627 write(iout,'(f15.4,i10)') entmin+i,histent(i)
635 if (ent.lt.entmin) entmin=ent
636 if (ent.gt.entmax) entmax=ent
638 ientmax=-dlog(entmax)-entmin
639 if (ientmax.gt.2000) ientmax=2000
641 ient=entfac(t)-entmin
642 if (ient.le.2000) histent(ient)=histent(ient)+1
644 write (iout,*) "Entropy histogram"
646 write(iout,'(2f15.4)') entmin+i,histent(i)
651 c write (iout,*) "me1",me1," scount",scount(me1)
677 hrmsrgy(j,i,ib)=0.0d0
679 hrmsrgy_p(j,i,ib)=0.0d0
691 hfin_ent_p(ind)=hfin_ent_p(ind)+dexp(entfac(t))
693 hfin_ent(ind)=hfin_ent(ind)+dexp(entfac(t))
695 c write (iout,'(2i5,20f8.2)') t,t,(enetb(k,t,iparm),k=1,18)
696 call restore_parm(iparm)
697 evdw=enetb(21,t,iparm)
698 evdw_t=enetb(1,t,iparm)
700 evdw2_14=enetb(17,t,iparm)
701 evdw2=enetb(2,t,iparm)+evdw2_14
703 evdw2=enetb(2,t,iparm)
708 evdw1=enetb(16,t,iparm)
713 ecorr=enetb(4,t,iparm)
714 ecorr5=enetb(5,t,iparm)
715 ecorr6=enetb(6,t,iparm)
716 eel_loc=enetb(7,t,iparm)
717 eello_turn3=enetb(8,t,iparm)
718 eello_turn4=enetb(9,t,iparm)
719 eturn6=enetb(10,t,iparm)
720 ebe=enetb(11,t,iparm)
721 escloc=enetb(12,t,iparm)
722 etors=enetb(13,t,iparm)
723 etors_d=enetb(14,t,iparm)
724 ehpb=enetb(15,t,iparm)
725 estr=enetb(18,t,iparm)
726 esccor=enetb(19,t,iparm)
727 edihcnstr=enetb(20,t,iparm)
728 esaxs=enetb(26,i,iparm)
730 betaT=startGridT+k*delta_T
734 if (rescale_mode.eq.1) then
742 denom=kfacl-1.0d0+quotl
744 ftprim(l)=-l*ft(l)*quotl1/(T0*denom)
745 ftbis(l)=l*kfacl*quotl1*
746 & (2*l*quotl-(l-1)*denom)/(quot*t0*t0*denom**3)
749 ft(6)=(320.0d0+80.0d0*dtanh((betaT-320.0d0)/80.0d0))/
751 ftprim(6)=1.0d0/(320.0d0*dcosh((betaT-320.0d0)/80.0d0)**2)
752 ftbis(6)=-2.0d0*dtanh((betaT-320.0d0)/80.0d0)
753 & /(320.0d0*80.0d0*dcosh((betaT-320.0d0)/80.0d0)**3)
763 else if (rescale_mode.eq.2) then
771 logfac=1.0d0/dlog(eplus+eminus)
772 tanhT=(eplus-eminus)/(eplus+eminus)
773 fT(l)=1.12692801104297249644d0*logfac
774 ftprim(l)=-l*quotl1*ft(l)*tanhT*logfac/T0
775 ftbis(l)=(l-1)*ftprim(l)/(quot*T0)-
776 & 2*l*quotl1/T0*logfac*
777 & (2*l*quotl1*ft(l)/(T0*(eplus+eminus)**2)
781 ft(6)=(320.0d0+80.0d0*dtanh((betaT-320.0d0)/80.0d0))/
783 ftprim(6)=1.0d0/(320.0d0*dcosh((betaT-320.0d0)/80.0d0)**2)
784 ftbis(6)=-2.0d0*dtanh((betaT-320.0d0)/80.0d0)
785 & /(320.0d0*80.0d0*dcosh((betaT-320.0d0)/80.0d0)**3)
795 else if (rescale_mode.eq.0) then
801 write (iout,*) "Error in WHAM_CALC: wrong RESCALE_MODE",
806 c write (iout,*) "ftprim",ftprim
807 c write (iout,*) "ftbis",ftbis
808 betaT=1.0d0/(1.987D-3*betaT)
810 if (shield_mode.gt.0) then
811 etot=ft(1)*wsc*(evdw+ft(6)*evdw_t)+ft(1)*wscp*evdw2
813 & +ft(1)*wvdwpp*evdw1
814 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
815 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
816 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
817 & +ft(2)*wturn3*eello_turn3
818 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc
819 & +edihcnstr+ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
820 & +wbond*estr+wliptran*eliptran+wsaxs*esaxs
821 eprim=ftprim(1)*(ft(6)*evdw_t+evdw)
822 C & +ftprim(6)*evdw_t
823 & +ftprim(1)*wscp*evdw2
824 & +ftprim(1)*welec*ees
825 & +ftprim(1)*wvdwpp*evdw1
826 & +ftprim(1)*wtor*etors+
827 & ftprim(3)*wcorr*ecorr+ftprim(4)*wcorr5*ecorr5+
828 & ftprim(5)*wcorr6*ecorr6+ftprim(3)*wturn4*eello_turn4+
829 & ftprim(2)*wturn3*eello_turn3+ftprim(5)*wturn6*eturn6+
830 & ftprim(2)*wel_loc*eel_loc+ftprim(2)*wtor_d*etors_d+
831 & ftprim(1)*wsccor*esccor
832 ebis=ftbis(1)*wsc*(evdw+ft(6)*evdw_t)
833 & +ftbis(1)*wscp*evdw2+
835 & +ftbis(1)*wvdwpp*evdw
836 & +ftbis(1)*wtor*etors+
837 & ftbis(3)*wcorr*ecorr+ftbis(4)*wcorr5*ecorr5+
838 & ftbis(5)*wcorr6*ecorr6+ftbis(3)*wturn4*eello_turn4+
839 & ftbis(2)*wturn3*eello_turn3+ftbis(5)*wturn6*eturn6+
840 & ftbis(2)*wel_loc*eel_loc+ftbis(2)*wtor_d*etors_d+
841 & ftbis(1)*wsccor*esccor
843 etot=wsc*(evdw+ft(6)*evdw_t)+wscp*evdw2+ft(1)*welec*ees
845 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
846 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
847 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
848 & +ft(2)*wturn3*eello_turn3
849 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc
850 & +edihcnstr+ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
851 & +wbond*estr+wliptran*eliptran+wsaxs*esaxs
852 eprim=ftprim(6)*evdw_t+ftprim(1)*welec*ees
853 & +ftprim(1)*wtor*etors+
854 & ftprim(3)*wcorr*ecorr+ftprim(4)*wcorr5*ecorr5+
855 & ftprim(5)*wcorr6*ecorr6+ftprim(3)*wturn4*eello_turn4+
856 & ftprim(2)*wturn3*eello_turn3+ftprim(5)*wturn6*eturn6+
857 & ftprim(2)*wel_loc*eel_loc+ftprim(2)*wtor_d*etors_d+
858 & ftprim(1)*wsccor*esccor
859 ebis=ftbis(1)*welec*ees+ftbis(1)*wtor*etors+
860 & ftbis(3)*wcorr*ecorr+ftbis(4)*wcorr5*ecorr5+
861 & ftbis(5)*wcorr6*ecorr6+ftbis(3)*wturn4*eello_turn4+
862 & ftbis(2)*wturn3*eello_turn3+ftbis(5)*wturn6*eturn6+
863 & ftbis(2)*wel_loc*eel_loc+ftbis(2)*wtor_d*etors_d+
864 & ftbis(1)*wsccor*esccor
867 if (shield_mode.gt.0) then
868 etot=ft(1)*wsc*(evdw+ft(6)*evdw_t)+ft(1)*wscp*evdw2
869 & +ft(1)*welec*(ees+evdw1)
870 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
871 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
872 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
873 & +ft(2)*wturn3*eello_turn3
874 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc+edihcnstr
875 & +ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
876 & +wbond*estr+wliptran*eliptran+wsaxs*esaxs
877 eprim=ftprim(1)*(evdw+ft(6)*evdw_t)
878 & +ftprim(1)*welec*(ees+evdw1)
879 & +ftprim(1)*wtor*etors+
880 & ftprim(1)*wscp*evdw2+
881 & ftprim(3)*wcorr*ecorr+ftprim(4)*wcorr5*ecorr5+
882 & ftprim(5)*wcorr6*ecorr6+ftprim(3)*wturn4*eello_turn4+
883 & ftprim(2)*wturn3*eello_turn3+ftprim(5)*wturn6*eturn6+
884 & ftprim(2)*wel_loc*eel_loc+ftprim(2)*wtor_d*etors_d+
885 & ftprim(1)*wsccor*esccor
886 ebis= ftbis(1)*(evdw+ft(6)*evdw_t)
887 & +ftbis(1)*wscp*evdw2
888 & +ftbis(1)*welec*(ees+evdw1)+ftbis(1)*wtor*etors+
889 & ftbis(3)*wcorr*ecorr+ftbis(4)*wcorr5*ecorr5+
890 & ftbis(5)*wcorr6*ecorr6+ftbis(3)*wturn4*eello_turn4+
891 & ftbis(2)*wturn3*eello_turn3+ftbis(5)*wturn6*eturn6+
892 & ftbis(2)*wel_loc*eel_loc+ftbis(2)*wtor_d*etors_d+
893 & ftprim(1)*wsccor*esccor
895 etot=wsc*(evdw+ft(6)*evdw_t)+wscp*evdw2
896 & +ft(1)*welec*(ees+evdw1)
897 & +wang*ebe+ft(1)*wtor*etors+wscloc*escloc
898 & +wstrain*ehpb+nss*ebr+ft(3)*wcorr*ecorr+ft(4)*wcorr5*ecorr5
899 & +ft(5)*wcorr6*ecorr6+ft(3)*wturn4*eello_turn4
900 & +ft(2)*wturn3*eello_turn3
901 & +ft(5)*wturn6*eturn6+ft(2)*wel_loc*eel_loc+edihcnstr
902 & +ft(2)*wtor_d*etors_d+ft(1)*wsccor*esccor
903 & +wbond*estr+wliptran*eliptran+wsaxs*esaxs
904 eprim=ftprim(6)*evdw_t+ftprim(1)*welec*(ees+evdw1)
905 & +ftprim(1)*wtor*etors+
906 & ftprim(3)*wcorr*ecorr+ftprim(4)*wcorr5*ecorr5+
907 & ftprim(5)*wcorr6*ecorr6+ftprim(3)*wturn4*eello_turn4+
908 & ftprim(2)*wturn3*eello_turn3+ftprim(5)*wturn6*eturn6+
909 & ftprim(2)*wel_loc*eel_loc+ftprim(2)*wtor_d*etors_d+
910 & ftprim(1)*wsccor*esccor
911 ebis=ftbis(1)*welec*(ees+evdw1)+ftbis(1)*wtor*etors+
912 & ftbis(3)*wcorr*ecorr+ftbis(4)*wcorr5*ecorr5+
913 & ftbis(5)*wcorr6*ecorr6+ftbis(3)*wturn4*eello_turn4+
914 & ftbis(2)*wturn3*eello_turn3+ftbis(5)*wturn6*eturn6+
915 & ftbis(2)*wel_loc*eel_loc+ftbis(2)*wtor_d*etors_d+
916 & ftprim(1)*wsccor*esccor
921 weight=dexp(-betaT*(etot-potEmin)+entfac(t))
923 write (iout,*) "iparm",iparm," t",t," betaT",betaT,
924 & " etot",etot," entfac",entfac(t),
925 & " weight",weight," ebis",ebis
927 etot=etot-temper*eprim
929 sumW_p(k,iparm)=sumW_p(k,iparm)+weight
930 sumE_p(k,iparm)=sumE_p(k,iparm)+etot*weight
931 sumEbis_p(k,iparm)=sumEbis_p(k,iparm)+ebis*weight
932 sumEsq_p(k,iparm)=sumEsq_p(k,iparm)+etot**2*weight
934 sumQ_p(j,k,iparm)=sumQ_p(j,k,iparm)+q(j,t)*weight
935 sumQsq_p(j,k,iparm)=sumQsq_p(j,k,iparm)+q(j,t)**2*weight
936 sumEQ_p(j,k,iparm)=sumEQ_p(j,k,iparm)
937 & +etot*q(j,t)*weight
940 sumW(k,iparm)=sumW(k,iparm)+weight
941 sumE(k,iparm)=sumE(k,iparm)+etot*weight
942 sumEbis(k,iparm)=sumEbis(k,iparm)+ebis*weight
943 sumEsq(k,iparm)=sumEsq(k,iparm)+etot**2*weight
945 sumQ(j,k,iparm)=sumQ(j,k,iparm)+q(j,t)*weight
946 sumQsq(j,k,iparm)=sumQsq(j,k,iparm)+q(j,t)**2*weight
947 sumEQ(j,k,iparm)=sumEQ(j,k,iparm)
948 & +etot*q(j,t)*weight
952 indE = aint(potE(t,iparm)-aint(potEmin))
953 if (indE.ge.0 .and. indE.le.maxinde) then
954 if (indE.gt.upindE_p) upindE_p=indE
955 histE_p(indE)=histE_p(indE)+dexp(-entfac(t))
959 expfac=dexp(-beta_h(ib,iparm)*(etot-potEmin)+entfac(t))
960 hfin_p(ind,ib)=hfin_p(ind,ib)+
961 & dexp(-beta_h(ib,iparm)*(etot-potEmin)+entfac(t))
963 indrgy=dint((q(nQ+2,t)-rgymin)/deltrgy)
964 indrms=dint((q(nQ+1,t)-rmsmin)/deltrms)
965 hrmsrgy_p(indrgy,indrms,ib)=
966 & hrmsrgy_p(indrgy,indrms,ib)+expfac
971 expfac=dexp(-beta_h(ib,iparm)*(etot-potEmin)+entfac(t))
972 hfin(ind,ib)=hfin(ind,ib)+
973 & dexp(-beta_h(ib,iparm)*(etot-potEmin)+entfac(t))
975 indrgy=dint((q(nQ+2,t)-rgymin)/deltrgy)
976 indrms=dint((q(nQ+1,t)-rmsmin)/deltrms)
977 hrmsrgy(indrgy,indrms,ib)=
978 & hrmsrgy(indrgy,indrms,ib)+expfac
984 if (histout) call MPI_Reduce(hfin_p(0,ib),hfin(0,ib),nbin,
985 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
987 call MPI_Reduce(hrmsrgy_p(0,0,ib),hrmsrgy(0,0,ib),
988 & (MaxBinRgy+1)*(nbin_rms+1),MPI_DOUBLE_PRECISION,MPI_SUM,Master,
992 call MPI_Reduce(upindE_p,upindE,1,
993 & MPI_INTEGER,MPI_MAX,Master,WHAM_COMM,IERROR)
994 call MPI_Reduce(histE_p(0),histE(0),maxindE,
995 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
997 if (me1.eq.master) then
1001 write (iout,'(6x,$)')
1002 write (iout,'(f20.2,$)') (1.0d0/(1.987D-3*beta_h(ib,iparm)),
1006 write (iout,'(/a)') 'Final histograms'
1008 if (nslice.eq.1) then
1009 if (separate_parset) then
1010 write(licz3,"(bz,i3.3)") myparm
1011 histname=prefix(:ilen(prefix))//'_par'//licz3//'.hist'
1013 histname=prefix(:ilen(prefix))//'.hist'
1016 if (separate_parset) then
1017 write(licz3,"(bz,i3.3)") myparm
1018 histname=prefix(:ilen(prefix))//'_par'//licz3//
1019 & '_slice_'//licz2//'.hist'
1021 histname=prefix(:ilen(prefix))//'_slice_'//licz2//'.hist'
1024 #if defined(AIX) || defined(PGI)
1025 open (ihist,file=histname,position='append')
1027 open (ihist,file=histname,access='append')
1035 sumH=sumH+hfin(t,ib)
1037 if (sumH.gt.0.0d0) then
1039 jj = mod(liczba,nbin1)
1041 write (iout,'(f6.3,$)') dmin+(jj+0.5d0)*delta
1043 & write (ihist,'(f6.3,$)') dmin+(jj+0.5d0)*delta
1046 write (iout,'(e20.10,$)') hfin(t,ib)
1047 if (histfile) write (ihist,'(e20.10,$)') hfin(t,ib)
1049 write (iout,'(i5)') iparm
1050 if (histfile) write (ihist,'(i5)') iparm
1057 if (nslice.eq.1) then
1058 if (separate_parset) then
1059 write(licz3,"(bz,i3.3)") myparm
1060 histname=prefix(:ilen(prefix))//"_par"//licz3//'.ent'
1062 histname=prefix(:ilen(prefix))//'.ent'
1065 if (separate_parset) then
1066 write(licz3,"(bz,i3.3)") myparm
1067 histname=prefix(:ilen(prefix))//'par_'//licz3//
1068 & '_slice_'//licz2//'.ent'
1070 histname=prefix(:ilen(prefix))//'_slice_'//licz2//'.ent'
1073 #if defined(AIX) || defined(PGI)
1074 open (ihist,file=histname,position='append')
1076 open (ihist,file=histname,access='append')
1078 write (ihist,'(a)') "# Microcanonical entropy"
1080 write (ihist,'(f8.0,$)') dint(potEmin)+i
1081 if (histE(i).gt.0.0e0) then
1082 write (ihist,'(f15.5,$)') dlog(histE(i))
1084 write (ihist,'(f15.5,$)') 0.0d0
1090 write (iout,*) "Microcanonical entropy"
1092 write (iout,'(f8.0,$)') dint(potEmin)+i
1093 if (histE(i).gt.0.0e0) then
1094 write (iout,'(f15.5,$)') dlog(histE(i))
1096 write (iout,'(f15.5,$)') 0.0d0
1101 if (nslice.eq.1) then
1102 if (separate_parset) then
1103 write(licz3,"(bz,i3.3)") myparm
1104 histname=prefix(:ilen(prefix))//'_par'//licz3//'.rmsrgy'
1106 histname=prefix(:ilen(prefix))//'.rmsrgy'
1109 if (separate_parset) then
1110 write(licz3,"(bz,i3.3)") myparm
1111 histname=prefix(:ilen(prefix))//'_par'//licz3//
1112 & '_slice_'//licz2//'.rmsrgy'
1114 histname=prefix(:ilen(prefix))//'_slice_'//licz2//'.rmsrgy'
1117 #if defined(AIX) || defined(PGI)
1118 open (ihist,file=histname,position='append')
1120 open (ihist,file=histname,access='append')
1124 write(ihist,'(2f8.2,$)')
1125 & rgymin+deltrgy*j,rmsmin+deltrms*i
1127 if (hrmsrgy(j,i,ib).gt.0.0d0) then
1128 write(ihist,'(e14.5,$)')
1129 & -dlog(hrmsrgy(j,i,ib))/beta_h(ib,iparm)
1132 write(ihist,'(e14.5,$)') 1.0d6
1135 write (ihist,'(i2)') iparm
1143 call MPI_Reduce(hfin_ent_p(0),hfin_ent(0),nbin,
1144 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
1145 call MPI_Reduce(sumW_p(0,1),sumW(0,1),(nGridT+1)*nParmSet,
1146 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
1147 call MPI_Reduce(sumE_p(0,1),sumE(0,1),(nGridT+1)*nParmSet,
1148 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
1149 call MPI_Reduce(sumEbis_p(0,1),sumEbis(0,1),(nGridT+1)*nParmSet,
1150 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
1151 call MPI_Reduce(sumEsq_p(0,1),sumEsq(0,1),(nGridT+1)*nParmSet,
1152 & MPI_DOUBLE_PRECISION,MPI_SUM,Master,WHAM_COMM,IERROR)
1153 call MPI_Reduce(sumQ_p(1,0,1),sumQ(1,0,1),
1154 & MaxQ1*(nGridT+1)*nParmSet,MPI_DOUBLE_PRECISION,MPI_SUM,Master,
1156 call MPI_Reduce(sumQsq_p(1,0,1),sumQsq(1,0,1),
1157 & MaxQ1*(nGridT+1)*nParmSet,MPI_DOUBLE_PRECISION,MPI_SUM,Master,
1159 call MPI_Reduce(sumEQ_p(1,0,1),sumEQ(1,0,1),
1160 & MaxQ1*(nGridT+1)*nParmSet,MPI_DOUBLE_PRECISION,MPI_SUM,Master,
1162 if (me.eq.master) then
1164 write (iout,'(/a)') 'Thermal characteristics of folding'
1165 if (nslice.eq.1) then
1168 nazwa=prefix(:ilen(prefix))//"_slice_"//licz2
1171 if (nparmset.eq.1 .and. .not.separate_parset) then
1172 nazwa=nazwa(:iln)//".thermal"
1173 else if (nparmset.eq.1 .and. separate_parset) then
1174 write(licz3,"(bz,i3.3)") myparm
1175 nazwa=nazwa(:iln)//"_par_"//licz3//".thermal"
1178 if (nparmset.gt.1) then
1179 write(licz3,"(bz,i3.3)") iparm
1180 nazwa=nazwa(:iln)//"_par_"//licz3//".thermal"
1183 if (separate_parset) then
1184 write (iout,'(a,i3)') "Parameter set",myparm
1186 write (iout,'(a,i3)') "Parameter set",iparm
1189 sumE(i,iparm)=sumE(i,iparm)/sumW(i,iparm)
1190 sumEbis(i,iparm)=(startGridT+i*delta_T)*sumEbis(i,iparm)/
1192 sumEsq(i,iparm)=(sumEsq(i,iparm)/sumW(i,iparm)
1193 & -sumE(i,iparm)**2)/(1.987D-3*(startGridT+i*delta_T)**2)
1195 sumQ(j,i,iparm)=sumQ(j,i,iparm)/sumW(i,iparm)
1196 sumQsq(j,i,iparm)=sumQsq(j,i,iparm)/sumW(i,iparm)
1197 & -sumQ(j,i,iparm)**2
1198 sumEQ(j,i,iparm)=sumEQ(j,i,iparm)/sumW(i,iparm)
1199 & -sumQ(j,i,iparm)*sumE(i,iparm)
1201 sumW(i,iparm)=-dlog(sumW(i,iparm))*(1.987D-3*
1202 & (startGridT+i*delta_T))+potEmin
1203 write (iout,'(f7.1,2f15.5,$)') startGridT+i*delta_T,
1204 & sumW(i,iparm),sumE(i,iparm)
1205 write (iout,'(f10.5,$)') (sumQ(j,i,iparm),j=1,nQ+2)
1206 write (iout,'(e15.5,$)') sumEsq(i,iparm)-sumEbis(i,iparm),
1207 & (sumQsq(j,i,iparm),j=1,nQ+2),(sumEQ(j,i,iparm),j=1,nQ+2)
1209 write (34,'(f7.1,2f15.5,$)') startGridT+i*delta_T,
1210 & sumW(i,iparm),sumE(i,iparm)
1211 write (34,'(f10.5,$)') (sumQ(j,i,iparm),j=1,nQ+2)
1212 write (34,'(e15.5,$)') sumEsq(i,iparm)-sumEbis(i,iparm),
1213 & (sumQsq(j,i,iparm),j=1,nQ+2),(sumEQ(j,i,iparm),j=1,nQ+2)
1221 if (hfin_ent(t).gt.0.0d0) then
1223 jj = mod(liczba,nbin1)
1224 write (iout,'(f6.3,e20.10," ent")') dmin+(jj+0.5d0)*delta,
1226 if (histfile) write (ihist,'(f6.3,e20.10," ent")')
1227 & dmin+(jj+0.5d0)*delta,
1231 if (histfile) close(ihist)
1235 ! Write data for zscore
1236 if (nslice.eq.1) then
1237 zscname=prefix(:ilen(prefix))//".zsc"
1239 zscname=prefix(:ilen(prefix))//"_slice_"//licz2//".zsc"
1241 #if defined(AIX) || defined(PGI)
1242 open (izsc,file=prefix(:ilen(prefix))//'.zsc',position='append')
1244 open (izsc,file=prefix(:ilen(prefix))//'.zsc',access='append')
1246 write (izsc,'("NQ=",i1," NPARM=",i1)') nQ,nParmSet
1248 write (izsc,'("NT=",i1)') nT_h(iparm)
1250 write (izsc,'("TEMP=",f6.1," NR=",i2," SNK=",$)')
1251 & 1.0d0/(beta_h(ib,iparm)*1.987D-3),nR(ib,iparm)
1252 jj = min0(nR(ib,iparm),7)
1253 write (izsc,'(i8,$)') (snk(i,ib,iparm,islice),i=1,jj)
1254 write (izsc,'(a1,$)') (" ",i=22+8*jj+1,79)
1255 write (izsc,'("&")')
1256 if (nR(ib,iparm).gt.7) then
1257 do ii=8,nR(ib,iparm),9
1258 jj = min0(nR(ib,iparm),ii+8)
1259 write (izsc,'(i8,$)') (snk(i,ib,iparm,islice),i=ii,jj)
1260 write (izsc,'(a1,$') (" ",i=(jj-ii+1)*8+1,79)
1261 write (izsc,'("&")')
1264 write (izsc,'("FI=",$)')
1265 jj=min0(nR(ib,iparm),7)
1266 write (izsc,'(f10.5,$)') (fi(i,ib,iparm),i=1,jj)
1267 write (izsc,'(a1,$)') (" ",i=3+10*jj+1,79)
1268 write (izsc,'("&")')
1269 if (nR(ib,iparm).gt.7) then
1270 do ii=8,nR(ib,iparm),9
1271 jj = min0(nR(ib,iparm),ii+8)
1272 write (izsc,'(f10.5,$)') (fi(i,ib,iparm),i=ii,jj)
1273 if (jj.eq.nR(ib,iparm)) then
1276 write (izsc,'(a1,$)') (" ",i=10*(jj-ii+1)+1,79)
1277 write (izsc,'(t80,"&")')
1282 write (izsc,'("KH=",$)')
1283 write (izsc,'(f7.2,$)') (Kh(j,i,ib,iparm),j=1,nQ)
1284 write (izsc,'(" Q0=",$)')
1285 write (izsc,'(f7.5,$)') (q0(j,i,ib,iparm),j=1,nQ)