def convolution(simulated, resolution, expdata, convolved, dak=None, norm2one=False): """Convolve a simulated S(Q,E) with a resolution file Arguments: simulated: Nexus file containing S(Q,E) from a simulation resolution: Nexus file containing the resolution. This will be used to produce a elastic line. convolved: Output Nexus file containing the convolution of the simulated S(Q,E) with the model beamline. expdata: Optional, experimental nexus file. Convolved will be binned as expdata. Returns: workspace for the convolution """ from mantid.simpleapi import (LoadNexus, Rebin, ConvertToHistogram, NormaliseToUnity, SaveNexus, SaveAscii, AddSampleLog) wss=LoadNexus(Filename=simulated,OutputWorkspace='simulated') width=wss.readX(0)[1]-wss.readX(0)[0] # rebin resolution as simulated wsr=LoadNexus(Filename=resolution,OutputWorkspace='resolution') #symmetrize the domain of the resolution function. Otherwise the #convolution results in a function with its peak shifted from the origin min=wsr.readX(0)[0] max=wsr.readX(0)[-1] delta=min+max if delta<0: wsr=Rebin(wsr, Params=(-max,width,max)) elif delta>0: wsr=Rebin(wsr, Params=(min,width,-min)) else: wsr=Rebin(wsr, Params=(min,width,max)) # convolve now, overwriting simulateds for i in range(wss.getNumberHistograms()): v=wsr.readY(i) w=wss.readY(i) x=camm_convolve(w,v,mode='same') wss.setY(i,x) wse=LoadNexus(Filename=expdata,OutputWorkspace='expdata') width=wse.readX(0)[1]-wse.readX(0)[0] # rebin simulated as expdata Rebin(InputWorkspace='simulated', Params=(wse.readX(0)[0],width,wse.readX(0)[-1]), OutputWorkspace='convolved') ConvertToHistogram(InputWorkspace='convolved', OutputWorkspace='convolved') # remember that output from sassena are point-data if norm2one: wsc=NormaliseToUnity(InputWorkspace='convolved', OutputWorkspace='convolved') AddSampleLog(Workspace='convolved',LogName='NormaliseToUnity',LogText=str(float(norm2one)),LogType='Number') if dak: dakota_vals = getParams(dak) # read in Dakota params file AddSampleLog(Workspace='convolved',LogName='FF1',LogText=str(dakota_vals["FF1"]),LogType='Number') from mantid.simpleapi import mtd SaveNexus(InputWorkspace='convolved', Filename=convolved) return
def modelB_EC_C(model, resolution, convolved, qvalues, assembled, expdata=None, costfile=None): """Assemble the Background, Elastic line and Convolution of the resolution with the simulated S(Q,E) This is a hard-coded model consisting of a linear background, and elastic line, and a convolution: b0+b1*E + +EC(Q)*e0*exp(-e1*Q^2)*Elastic(E) + c0*Resolution(E)xSimulated(Q,E) We load Resolution(E)xSimulated(Q,E) as Convolved(Q,E) EC(Q) is a fit to the Q-dependence of the integrated intensity of the empty can EC(Q) = 2.174495971 - 2.065826056*Q + 0.845367259*Q^2 Arguments: model: beamline model file is a single line, e.g, b0=1.3211; b1=0.00 e0=0.99; e1=1.9; c0=2.3 resolution: Nexus file containing the resolution. This will be used to produce a elastic line. convolved: Nexus file containing the convolution of the simulated S(Q,E) with the resolution. qvalues: single-column file containing list of Q-values assembled: output Nexus file containing the assembled S(Q,E) of the beamline model and the simulated S(Q,E) expdata: Optional, experimental nexus file. If passed, output convolved will be binned as expdata. costfile: Optional, file to store cost. If passed, the cost of comparing convolved and expdata will be saved. Returns: workspace containing the assembled S(Q,E) """ import numpy from mantid.simpleapi import (LoadNexus, ScaleX, ConvertToPointData, SaveNexus, DakotaChiSquared) EC = lambda Q: 2.174495971 - 2.065826056*Q + 0.845367259*Q*Q Q=[float(q) for q in open(qvalues,'r').read().split('\n')] p={} for pair in open(model,'r').readline().split(';'): key,val=pair.split('=') p[key.strip()]=float(val.strip()) wsr=LoadNexus(Filename=resolution,OutputWorkspace='resolution') wsr=ConvertToPointData(wsr) E=wsr.readX(0) wse=ScaleX(InputWorkspace=wsr, OutputWorkspace='elastic',factor=-1) # elastic line wsc=LoadNexus(Filename=convolved,OutputWorkspace='convolved') for i in range(wsc.getNumberHistograms()): elastic=wse.readY(i) # elastic spectrum at a given Q convolved=wsc.readY(i) # convolved spectrum at a given Q wsc.setY(i, (p['b0']+p['b1']*E) + (EC(Q[i])*p['e0']*numpy.exp(-p['e1']*Q[i])*elastic) + (p['c0']*convolved) ) # overwrite spectrum SaveNexus(InputWorkspace=wsc, Filename=assembled) if expdata and costfile: DakotaChiSquared(DataFile=assembled,CalculatedFile=expdata,OutputFile=costfile) return wsc
def modelB_freeE_C(modeltpl, elastic, convolved, expdata, initparfile=None): """Estimate initial beamline parameters for the modelB_freeE_C This is a hard-coded model consisting of a linear background, and elastic line, and a convolution: b0+b1*E + +e0(Q)*Elastic(E) + c0*Resolution(E)xSimulated(Q,E) We load Resolution(E)xSimulated(Q,E) as Convolved(Q,E) e0(Q) are a set of fitting parameters, one for each Q Initial values are estimated as follows: b0:0.0 b1:0.0 Evaluation of the model at E=0: e0*elastic(Q,0) + c0*convolved(Q,0) ~ experiment(Q,0) {Eq.1}, with 'convolved' the convolution of the experimental resolution and the Fourier transform of the simulated intermediate structure factor For the lowest Q, we assume contributions fromt the elastic line and simuation are equal. Thus: c0*convolved(Qmin,0) ~ 1/2*experiment(Qmin,0) ---> provides estimation for c0 e0(Qmin)*elastic(Qmin,0) ~ 1/2*experiment(Qmin,0) ---> provides estimation for e0 For the remaining Q, we use {Eq.1} substituting the c0 found above. Finally, eshift:0.0 Arguments: model: beamline template model file (xml format) elastic: Nexus file containing the elastic line convolved: Nexus file containing convolution of the resolution and simulated structure factor expdata: Nexus file containing the experimental data [initparfile]: Output the initial parameters as a string in file with name initparfile Returns: initparms: string with initial values for the parameters """ from simulation.src.molmec.ffupdate.ff_update import loadFFtpl,updateTemplate from mantid.simpleapi import LoadNexus wse=LoadNexus(Filename=elastic,OutputWorkspace='elastic') wsc=LoadNexus(Filename=convolved,OutputWorkspace='convolved') wsx=LoadNexus(Filename=expdata,OutputWorkspace='experiment') parl,template=loadFFtpl(modeltpl) pard={} for par in parl: pard[par._name]=par nhist=wsx.getNumberHistograms() le=len(wsx.readX(0)) for ws in wse,wsc: if ws.getNumberHistograms()!=nhist or len(ws.readX(0))!=le: error_message='%s %d histograms of length %d do not conform to those of experiment'%(ws.getName(),ws.getNumberHistograms(),len(ws.readX(0))) ws.getName()+' histograms do not conform to those of experiment' g_log.error(error_message) raise StandardError(error_message) pard['b0'].setValue(1e-10) # needs to be positive pard['b1'].setValue(0.) ezero=le/2 # assume E=0 in the middle of the histogram span pard['c0'].setValue(0.5*wsx.readY(0)[ezero]/wsc.readY(0)[ezero]) pard['e0.0'].setValue(0.5*wsx.readY(0)[ezero]/wse.readY(0)[ezero]) trace() for ihist in range(1,nhist): pard['e0.'+str(ihist)].setValue((wsx.readY(ihist)[ezero] - pard['c0']._value*wsc.readY(ihist)[ezero]) / wse.readY(ihist)[ezero]) pard['eshift'].setValue(0.) template=updateTemplate(template,parl) if initparfile: open(initparfile,'w').write(template) return template
def modelB_freeE_C(model, resolution, convolved, assembled, expdata=None, costfile=None, derivdata=None, derivexclude=[], doshift=None): """Assemble the Background, Elastic line and Convolution of the resolution with the simulated S(Q,E) This is a hard-coded model consisting of a linear background, and elastic line, and a convolution: b0+b1*E + +e0(Q)*Elastic(E) + c0*Resolution(E)xSimulated(Q,E) We load Resolution(E)xSimulated(Q,E) as Convolved(Q,E) e0(Q) are a set of fitting parameters, one for each Q Arguments: model: beamline model file is a single line, e.g, b0=1.3211; b1=0.00; e0.0=0.99; e0.1=0.99; e0.2=0.99;...e0.N=0.99; e1=1.9; c0=2.3 resolution: Nexus file containing the resolution. This will be used to produce a elastic line. convolved: Nexus file containing the convolution of the simulated S(Q,E) with the resolution. assembled: output Nexus file containing the assembled S(Q,E) of the beamline model and the simulated S(Q,E) expdata: Optional, experimental nexus file. If passed, output convolved will be binned as expdata. costfile: Optional, file to store cost. If passed, residuals and (optionally) partial derivatives will be stored derivdata: Optional, perform analytic derivatives (store in costfile if provided) derivexclude: list of fitting parameters for which partial derivatives will not be computed doshift: Optional, perform the shift of the model function Returns: wsm: workspace containing the assembled S(Q,E) gradients: dictionary of partial derivatives with respect to model parameters """ import numpy from copy import copy,deepcopy from mantid.simpleapi import (LoadNexus, ScaleX, ConvertToPointData, SaveNexus, DakotaChiSquared, AddSampleLog) from math import sqrt def shiftalongX(*kargs,**kwargs): """ Function to do the shift along the E-axis. By default, does nothing """ pass import interpX if doshift: # replace the dummy function with the real thing if doshift in dir(interpX): shiftalongX=getattr(__import__('interpX'), doshift) else: shiftalongX = getattr(__import__('interpX'), 'itp_simple') def computemodel(p,wse,wsc): """Assemble the model Arguments p: dictionary with parameter values wse: Mantid workspace holding the elastic line wsc: Mantid workspace holding the convolution of the resolution and the simulation Returns: wsm: Mantid workspace holding the resulting model """ from mantid.simpleapi import CloneWorkspace wsm=CloneWorkspace(wsc) E=wse.readX(0) # energy values, bins boundary values Eshifted=(E[1:]+E[:-1])/2 # energy values, center bin values for i in range(wsc.getNumberHistograms()): elastic=wse.readY(i) # elastic spectrum at a given Q convolved=wsc.readY(i) # convolved spectrum at a given Q wsm.setY(i, p['b0']+p['b1']*Eshifted + p['e0.'+str(i)]*elastic + p['c0']*convolved) # overwrite spectrum return wsm # init list of parameters names for which analytical derivative exists, same order as in the input model file derivparnames=[] # filled only if derivdata different than None p={} for pair in open(model,'r').readline().split(';'): key,val=[x.strip() for x in pair.split('=')] if derivdata and key not in derivexclude: derivparnames.append(key) p[key]=float(val) # read various inputs wsr=LoadNexus(Filename=resolution,OutputWorkspace='resolution') wse=ScaleX(InputWorkspace=wsr, OutputWorkspace='elastic',factor=-1) # elastic line wsc=LoadNexus(Filename=convolved,OutputWorkspace='convolved') E=wsr.readX(0) # energy values, bins boundary values de=E[1]-E[0] # assume all bins have same bin width Eshifted=(E[1:]+E[:-1])/2 # energy values, center bin values nhist=wsc.getNumberHistograms() nrsl=len(Eshifted)*nhist # number of residuals # calculate partial numerical derivative with respect to eshift gradients={} if 'eshift' in derivparnames: wsm=computemodel(p,wse,wsc) eshiftderiv=(shiftalongX(wsm,p['eshift']+0.5*de,newWorkspace='wsplus') - shiftalongX(wsm,p['eshift']-0.5*de,newWorkspace='wsminus')) / de # forward-backward difference with a 0.5*de step gradients['eshift']=numpy.zeros(0) for i in range(eshiftderiv.getNumberHistograms()): gradients['eshift'] = numpy.concatenate([gradients['eshift'], eshiftderiv.readY(i)]) # do eshift of component workspaces if doshift: Eshifted-=p['eshift'] wse=shiftalongX(wse,p['eshift']) # shift the spectrum, does nothing if shiftalongX is the dummy function wsc=shiftalongX(wsc,p['eshift']) # shift the spectrum, does nothing if shiftalongX is the dummy function # find difference in convolutions with FF1 changed if 'FF1' not in derivexclude: derivparnames.append('FF1') # difference in FF1 workspaces convolvedf=convolved.replace('.nxs','_1.nxs') wscf=LoadNexus(Filename=convolvedf,OutputWorkspace='convolvedf') convolvedb=convolved.replace('.nxs','_0.nxs') wscb=LoadNexus(Filename=convolvedb,OutputWorkspace='convolvedb') wksp_diff=wscf-wscb # calculate analytic partial derivatives with respect to the fit parameters if derivparnames: gradients['b0']=numpy.ones(nrsl) gradients['b1']=numpy.zeros(0) gradients['c0']=numpy.zeros(0) gradients['FF1']=numpy.zeros(0) for i in range(nhist): gradients['b1'] = numpy.concatenate([gradients['b1'], Eshifted]) gradients['c0'] = numpy.concatenate([gradients['c0'], wsc.readY(i)]) if 'FF1' not in derivexclude: gradients['FF1'] = numpy.concatenate([gradients['FF1'], wksp_diff.readY(i)]) gradients['e0.'+str(i)]=numpy.zeros(0) for j in range(nhist): if i==j: gradients['e0.'+str(i)] = numpy.concatenate([gradients['e0.'+str(i)], wse.readY(i)]) else: gradients['e0.'+str(i)] = numpy.concatenate([gradients['e0.'+str(i)], numpy.zeros(len(Eshifted))]) if 'FF1' not in derivexclude: FF1_f=wscf.getRun().getLogData('FF1').value FF1_b=wscb.getRun().getLogData('FF1').value gradients['FF1'] *= p['c0']/(FF1_f-FF1_b) # save model to file wsm=computemodel(p,wse,wsc) # add all parameters to assembled file for pair in open(model,'r').readline().split(';'): key,val=[x.strip() for x in pair.split('=')] AddSampleLog(Workspace=wsm,LogName=key,LogText=str(p[key]),LogType='Number') print key, "=", p[key] print "FF1 =", wsm.getRun().getLogData('FF1').value SaveNexus(InputWorkspace=wsm, Filename=assembled) # save residuals and partial derivatives buf='' if expdata and costfile: wex=LoadNexus(Filename=expdata,OutputWorkspace='experiment') chisq,wR=DakotaChiSquared(DataFile=expdata,CalculatedFile=assembled,OutputFile=costfile,ResidualsWorkspace='wR') Xe=numpy.zeros(0) # list of errors for each residual for i in range(nhist): Xe = numpy.concatenate([Xe, wex.readE(i)]) Ry=wR.readY(i) for j in range(len(Ry)): buf+=str(Ry[j])+" least_sq_term_"+str(i*len(Ry)+j+1)+"\n" for parname in derivparnames: gradients[parname]/=numpy.where(Xe>0,Xe,1) # divide by experimental error (with non-positive elements replaced by one) if derivparnames: for i in range(nrsl): buf+="[" for parname in derivparnames: buf+=" %.10e"%(-gradients[parname][i]) buf+=" ]\n" open(costfile,'w').write(buf) AddSampleLog(Workspace=wsm,LogName="chisq",LogText=str(chisq),LogType='Number') norm_chisq=chisq/(len(Ry)-len(derivparnames)) print costfile, " R = ", sqrt(norm_chisq) AddSampleLog(Workspace=wsm,LogName="norm_chisq",LogText=str(norm_chisq),LogType='Number') AddSampleLog(Workspace=wsm,LogName="norm_chi",LogText=str(sqrt(norm_chisq)),LogType='Number') SaveNexus(InputWorkspace=wsm, Filename=assembled) return {'model':wsm, 'gradients':gradients}