source: sans/Dev/trunk/NCNR_User_Procedures/Analysis/Models/NewModels_2008/Core_and_NShells_v40.ipf @ 510

Last change on this file since 510 was 510, checked in by srkline, 14 years ago

Simple change in all of the model function files to include the name of the parameter wave in the Keyword=list that is generated when a model is plotted. This is becoming an issue where the proper parameter wave can't be deduced from just the suffix, then there is nothing to put in the table.

I should have added this when I initially wrote the wrapper...

File size: 21.7 KB
Line 
1#pragma rtGlobals=1             // Use modern global access method.
2#pragma IgorVersion=6.0
3
4////////////////////////////////////////////////
5//
6// this function is for the form factor of a sphere with some
7// number of shells around a central core (currently 1-2-3)
8//
9// monodisperse and polydisperse (and smeared) versions are included
10// - for the polydisperse models, only polydispersity of the core is taken
11// into account, and done numerically. for a Schulz distribution, this
12// should be possible to do analytically, whith a great savings in computation
13// time.
14//
15// It may also be useful to think of scenarios where the layers as well are
16// polydisperse - to break up the very regular spacing of the layers, which
17// is not a very natural structure.
18//
19// 03 MAR 04 SRK
20////////////////////////////////////////////////
21// Four shell model added and polyCore versions written in separate ipf (polyCore_and_NShells.ipf)
22//
23// Four shell model added at request of H Wacklin to model lipid vesicles with varying deuteration
24//
25// DEC 2008 AJJ
26//////////////////////////////////////////////
27
28//this macro sets up all the necessary parameters and waves that are
29//needed to calculate the model function.
30//
31Proc PlotOneShell(num,qmin,qmax)
32        Variable num=200, qmin=0.001, qmax=0.7
33        Prompt num "Enter number of data points for model: "
34        Prompt qmin "Enter minimum q-value (^-1) for model: "
35        Prompt qmax "Enter maximum q-value (^-1) for model: "
36//
37        Make/O/D/n=(num) xwave_OneShell, ywave_OneShell
38        xwave_OneShell =  alog(log(qmin) + x*((log(qmax)-log(qmin))/num))
39        Make/O/D coef_OneShell = {1.,60,6.4e-6,10,1e-6,6.4e-6,0.001}
40        make/o/t parameters_OneShell = {"scale","core radius (A)","Core SLD (A-2)","Shell thickness (A)","Shell SLD (A-2)","Solvent SLD (A-2)","bkg (cm-1)"}
41        Edit parameters_OneShell, coef_OneShell
42       
43        Variable/G root:g_OneShell
44        g_OneShell := OneShell(coef_OneShell, ywave_OneShell, xwave_OneShell)
45        Display ywave_OneShell vs xwave_OneShell
46        ModifyGraph marker=29, msize=2, mode=4
47        ModifyGraph log=1,grid=1,mirror=2
48        Label bottom "q (\\S-1\\M) "
49        Label left "I(q) (cm\\S-1\\M)"
50        AutoPositionWindow/M=1/R=$(WinName(0,1)) $WinName(0,2)
51       
52        AddModelToStrings("OneShell","coef_OneShell","parameters_OneShell","OneShell")
53//
54End
55
56Proc PlotTwoShell(num,qmin,qmax)
57        Variable num=200, qmin=0.001, qmax=0.7
58        Prompt num "Enter number of data points for model: "
59        Prompt qmin "Enter minimum q-value (^-1) for model: "
60        Prompt qmax "Enter maximum q-value (^-1) for model: "
61//
62        Make/O/D/n=(num) xwave_TwoShell, ywave_TwoShell
63        xwave_TwoShell =  alog(log(qmin) + x*((log(qmax)-log(qmin))/num))
64        Make/O/D coef_TwoShell = {1.,60,6.4e-6,10,1e-6,10,2e-6,6.4e-6,0.001}
65        make/o/t parameters_TwoShell = {"scale","core radius (A)","Core SLD (A-2)","Shell 1 thickness","Shell 1 SLD (A-2)","Shell 2 thickness","Shell 2 SLD (A-2)","Solvent SLD (A-2)","bkg (cm-1)"}
66        Edit parameters_TwoShell, coef_TwoShell
67       
68        Variable/G root:g_TwoShell
69        g_TwoShell := TwoShell(coef_TwoShell, ywave_TwoShell, xwave_TwoShell)
70        Display ywave_TwoShell vs xwave_TwoShell
71        ModifyGraph marker=29, msize=2, mode=4
72        ModifyGraph log=1,grid=1,mirror=2
73        Label bottom "q (\\S-1\\M) "
74        Label left "I(q) (cm\\S-1\\M)"
75        AutoPositionWindow/M=1/R=$(WinName(0,1)) $WinName(0,2)
76       
77        AddModelToStrings("TwoShell","coef_TwoShell","parameters_TwoShell","TwoShell")
78//
79End
80
81Proc PlotThreeShell(num,qmin,qmax)
82        Variable num=200, qmin=0.001, qmax=0.7
83        Prompt num "Enter number of data points for model: "
84        Prompt qmin "Enter minimum q-value (^-1) for model: "
85        Prompt qmax "Enter maximum q-value (^-1) for model: "
86//
87        Make/O/D/n=(num) xwave_ThreeShell, ywave_ThreeShell
88        xwave_ThreeShell =  alog(log(qmin) + x*((log(qmax)-log(qmin))/num))
89        Make/O/D coef_ThreeShell ={1.,60,6.4e-6,10,1e-6,10,2e-6,10,3e-6,6.4e-6,0.001}
90        make/o/t parameters_ThreeShell = {"scale","core radius (A)","Core SLD (A-2)","Shell 1 thickness","Shell 1 SLD (A-2)","Shell 2 thickness","Shell 2 SLD (A-2)","Shell 3 thickness","Shell 3 SLD (A-2)","Solvent SLD (A-2)","bkg (cm-1)"}
91        Edit parameters_ThreeShell, coef_ThreeShell
92       
93        Variable/G root:g_ThreeShell
94        g_ThreeShell := ThreeShell(coef_ThreeShell, ywave_ThreeShell, xwave_ThreeShell)
95        Display ywave_ThreeShell vs xwave_ThreeShell
96        ModifyGraph marker=29, msize=2, mode=4
97        ModifyGraph log=1,grid=1,mirror=2
98        Label bottom "q (\\S-1\\M) "
99        Label left "I(q) (cm\\S-1\\M)"
100        AutoPositionWindow/M=1/R=$(WinName(0,1)) $WinName(0,2)
101       
102        AddModelToStrings("ThreeShell","coef_ThreeShell","parameters_ThreeShell","ThreeShell")
103//
104End
105
106Proc PlotFourShell(num,qmin,qmax)
107        Variable num=200, qmin=0.001, qmax=0.7
108        Prompt num "Enter number of data points for model: "
109        Prompt qmin "Enter minimum q-value (^-1) for model: "
110        Prompt qmax "Enter maximum q-value (^-1) for model: "
111//
112        Make/O/D/n=(num) xwave_FourShell, ywave_FourShell
113        xwave_FourShell =  alog(log(qmin) + x*((log(qmax)-log(qmin))/num))
114        Make/O/D coef_FourShell ={1.,60,6.4e-6,10,1e-6,10,2e-6,10,3e-6,10,4e-6,6.4e-6,0.001}
115        make/o/t parameters_FourShell = {"scale","core radius (A)","Core SLD (A-2)","Shell 1 thickness","Shell 1 SLD (A-2)","Shell 2 thickness","Shell 2 SLD (A-2)","Shell 3 thickness","Shell 3 SLD (A-2)","Shell 4 thickness","Shell 4 SLD (A-2)","Solvent SLD (A-2)","bkg (cm-1)"}
116        Edit parameters_FourShell, coef_FourShell
117       
118        Variable/G root:g_FourShell
119        g_FourShell := FourShell(coef_FourShell, ywave_FourShell, xwave_FourShell)
120        Display ywave_FourShell vs xwave_FourShell
121        ModifyGraph marker=29, msize=2, mode=4
122        ModifyGraph log=1,grid=1,mirror=2
123        Label bottom "q (\\S-1\\M) "
124        Label left "I(q) (cm\\S-1\\M)"
125        AutoPositionWindow/M=1/R=$(WinName(0,1)) $WinName(0,2)
126       
127        AddModelToStrings("FourShell","coef_FourShell","parameters_FourShell","FourShell")
128//
129End
130
131
132
133
134//
135//this macro sets up all the necessary parameters and waves that are
136//needed to calculate the  smeared model function.
137//
138//no input parameters are necessary, it MUST use the experimental q-values
139// from the experimental data read in from an AVE/QSIG data file
140////////////////////////////////////////////////////
141// - sets up a dependency to a wrapper, not the actual SmearedModelFunction
142Proc PlotSmearedOneShell(str)                                                           
143        String str
144        Prompt str,"Pick the data folder containing the resolution you want",popup,getAList(4)
145       
146        // if any of the resolution waves are missing => abort
147        if(ResolutionWavesMissingDF(str))               //updated to NOT use global strings (in GaussUtils)
148                Abort
149        endif
150       
151        SetDataFolder $("root:"+str)
152       
153        // Setup parameter table for model function
154        Make/O/D smear_coef_OneShell =  {1.,60,6.4e-6,10,1e-6,6.4e-6,0.001}
155        make/o/t smear_parameters_OneShell =  {"scale","core radius (A)","Core SLD (A-2)","Shell thickness (A)","Shell SLD (A-2)","Solvent SLD (A-2)","bkg (cm-1)"}
156        Edit smear_parameters_OneShell,smear_coef_OneShell                                      //display parameters in a table
157       
158        // output smeared intensity wave, dimensions are identical to experimental QSIG values
159        // make extra copy of experimental q-values for easy plotting
160        Duplicate/O $(str+"_q") smeared_OneShell,smeared_qvals
161        SetScale d,0,0,"1/cm",smeared_OneShell
162                                       
163        Variable/G gs_OneShell=0
164        gs_OneShell := fSmearedOneShell(smear_coef_OneShell,smeared_OneShell,smeared_qvals)     //this wrapper fills the STRUCT
165       
166        Display smeared_OneShell vs smeared_qvals
167        ModifyGraph log=1,marker=29,msize=2,mode=4
168        Label bottom "q (\\S-1\\M)"
169        Label left "I(q) (cm\\S-1\\M)"
170        AutoPositionWindow/M=1/R=$(WinName(0,1)) $WinName(0,2)
171       
172        SetDataFolder root:
173        AddModelToStrings("SmearedOneShell","smear_coef_OneShell","smear_parameters_OneShell","OneShell")
174End
175
176Proc PlotSmearedTwoShell(str)                                                           
177        String str
178        Prompt str,"Pick the data folder containing the resolution you want",popup,getAList(4)
179       
180        // if any of the resolution waves are missing => abort
181        if(ResolutionWavesMissingDF(str))               //updated to NOT use global strings (in GaussUtils)
182                Abort
183        endif
184       
185        SetDataFolder $("root:"+str)
186       
187        // Setup parameter table for model function
188        Make/O/D smear_coef_TwoShell =  {1.,60,6.4e-6,10,1e-6,10,2e-6,6.4e-6,0.001}
189        make/o/t smear_parameters_TwoShell =  {"scale","core radius (A)","Core SLD (A-2)","Shell 1 thickness","Shell 1 SLD (A-2)","Shell 2 thickness","Shell 2 SLD (A-2)","Solvent SLD (A-2)","bkg (cm-1)"}
190        Edit smear_parameters_TwoShell,smear_coef_TwoShell                                      //display parameters in a table
191       
192        // output smeared intensity wave, dimensions are identical to experimental QSIG values
193        // make extra copy of experimental q-values for easy plotting
194        Duplicate/O $(str+"_q") smeared_TwoShell,smeared_qvals
195        SetScale d,0,0,"1/cm",smeared_TwoShell
196                                       
197        Variable/G gs_TwoShell=0
198        gs_TwoShell := fSmearedTwoShell(smear_coef_TwoShell,smeared_TwoShell,smeared_qvals)     //this wrapper fills the STRUCT
199       
200        Display smeared_TwoShell vs smeared_qvals
201        ModifyGraph log=1,marker=29,msize=2,mode=4
202        Label bottom "q (\\S-1\\M)"
203        Label left "I(q) (cm\\S-1\\M)"
204        AutoPositionWindow/M=1/R=$(WinName(0,1)) $WinName(0,2)
205       
206        SetDataFolder root:
207        AddModelToStrings("SmearedTwoShell","smear_coef_TwoShell","smear_parameters_TwoShell","TwoShell")
208End
209
210
211Proc PlotSmearedThreeShell(str)                                                         
212        String str
213        Prompt str,"Pick the data folder containing the resolution you want",popup,getAList(4)
214       
215        // if any of the resolution waves are missing => abort
216        if(ResolutionWavesMissingDF(str))               //updated to NOT use global strings (in GaussUtils)
217                Abort
218        endif
219       
220        SetDataFolder $("root:"+str)
221       
222        // Setup parameter table for model function
223        Make/O/D smear_coef_ThreeShell =  {1.,60,6.4e-6,10,1e-6,10,2e-6,10,3e-6,6.4e-6,0.001}
224        make/o/t smear_parameters_ThreeShell =  {"scale","core radius (A)","Core SLD (A-2)","Shell 1 thickness","Shell 1 SLD (A-2)","Shell 2 thickness","Shell 2 SLD (A-2)","Shell 3 thickness","Shell 3 SLD (A-2)","Solvent SLD (A-2)","bkg (cm-1)"}
225        Edit smear_parameters_ThreeShell,smear_coef_ThreeShell                                  //display parameters in a table
226       
227        // output smeared intensity wave, dimensions are identical to experimental QSIG values
228        // make extra copy of experimental q-values for easy plotting
229        Duplicate/O $(str+"_q") smeared_ThreeShell,smeared_qvals
230        SetScale d,0,0,"1/cm",smeared_ThreeShell
231                                       
232        Variable/G gs_ThreeShell=0
233        gs_ThreeShell := fSmearedThreeShell(smear_coef_ThreeShell,smeared_ThreeShell,smeared_qvals)     //this wrapper fills the STRUCT
234       
235        Display smeared_ThreeShell vs smeared_qvals
236        ModifyGraph log=1,marker=29,msize=2,mode=4
237        Label bottom "q (\\S-1\\M)"
238        Label left "I(q) (cm\\S-1\\M)"
239        AutoPositionWindow/M=1/R=$(WinName(0,1)) $WinName(0,2)
240       
241        SetDataFolder root:
242        AddModelToStrings("SmearedThreeShell","smear_coef_ThreeShell","smear_parameters_ThreeShell","ThreeShell")
243End
244
245Proc PlotSmearedFourShell(str)                                                         
246        String str
247        Prompt str,"Pick the data folder containing the resolution you want",popup,getAList(4)
248       
249        // if any of the resolution waves are missing => abort
250        if(ResolutionWavesMissingDF(str))               //updated to NOT use global strings (in GaussUtils)
251                Abort
252        endif
253       
254        SetDataFolder $("root:"+str)
255       
256        // Setup parameter table for model function
257        Make/O/D smear_coef_FourShell = {1.,60,6.4e-6,10,1e-6,10,2e-6,10,3e-6,10,4e-6,6.4e-6,0.001}
258        make/o/t smear_parameters_FourShell =  {"scale","core radius (A)","Core SLD (A-2)","Shell 1 thickness","Shell 1 SLD (A-2)","Shell 2 thickness","Shell 2 SLD (A-2)","Shell 3 thickness","Shell 3 SLD (A-2)","Shell 4 thickness","Shell 4 SLD (A-2)","Solvent SLD (A-2)","bkg (cm-1)"}
259        Edit smear_parameters_FourShell,smear_coef_FourShell                                    //display parameters in a table
260       
261        // output smeared intensity wave, dimensions are identical to experimental QSIG values
262        // make extra copy of experimental q-values for easy plotting
263        Duplicate/O $(str+"_q") smeared_FourShell,smeared_qvals
264        SetScale d,0,0,"1/cm",smeared_FourShell
265                                       
266        Variable/G gs_FourShell=0
267        gs_FourShell := fSmearedFourShell(smear_coef_FourShell,smeared_FourShell,smeared_qvals) //this wrapper fills the STRUCT
268       
269        Display smeared_FourShell vs smeared_qvals
270        ModifyGraph log=1,marker=29,msize=2,mode=4
271        Label bottom "q (\\S-1\\M)"
272        Label left "I(q) (cm\\S-1\\M)"
273        AutoPositionWindow/M=1/R=$(WinName(0,1)) $WinName(0,2)
274       
275        SetDataFolder root:
276        AddModelToStrings("SmearedFourShell","smear_coef_FourShell","smear_parameters_FourShell","FourShell")
277End
278
279
280// nothing to change here
281//
282//AAO version, uses XOP if available
283// simply calls the original single point calculation with
284// a wave assignment (this will behave nicely if given point ranges)
285Function OneShell(cw,yw,xw) : FitFunc
286        Wave cw,yw,xw
287       
288#if exists("OneShellX")
289        yw = OneShellX(cw,xw)
290#else
291        yw = fOneShell(cw,xw)
292#endif
293        return(0)
294End
295
296Function TwoShell(cw,yw,xw) : FitFunc
297        Wave cw,yw,xw
298       
299#if exists("TwoShellX")
300        yw = TwoShellX(cw,xw)
301#else
302        yw = fTwoShell(cw,xw)
303#endif
304        return(0)
305End
306
307Function ThreeShell(cw,yw,xw) : FitFunc
308        Wave cw,yw,xw
309       
310#if exists("ThreeShellX")
311        yw =ThreeShellX(cw,xw)
312#else
313        yw = fThreeShell(cw,xw)
314#endif
315        return(0)
316End
317
318Function FourShell(cw,yw,xw) : FitFunc
319        Wave cw,yw,xw
320       
321#if exists("FourShellX")
322        yw = FourShellX(cw,xw)
323#else
324        yw = fFourShell(cw,xw)
325#endif
326        return(0)
327End
328
329
330//
331// unsmeared model calculation
332//
333Function fOneShell(w,x) : FitFunc
334        Wave w
335        Variable x
336       
337        // variables are:
338        //[0] scale factor
339        //[1] radius of core []
340        //[2] SLD of the core   [-2]
341        //[3] thickness of the shell    []
342        //[4] SLD of the shell
343        //[5] SLD of the solvent
344        //[6] background        [cm-1]
345       
346        // All inputs are in ANGSTROMS
347        //OUTPUT is normalized by the particle volume, and converted to [cm-1]
348       
349       
350        Variable scale,rcore,thick,rhocore,rhoshel,rhosolv,bkg
351        scale = w[0]
352        rcore = w[1]
353        rhocore = w[2]
354        thick = w[3]
355        rhoshel = w[4]
356        rhosolv = w[5]
357        bkg = w[6]
358       
359        // calculates scale *( f^2 + bkg)
360        Variable bes,f,vol,qr,contr,f2
361       
362        // core first, then add in shell
363        qr=x*rcore
364        contr = rhocore-rhoshel
365        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
366        vol = 4*pi/3*rcore^3
367        f = vol*bes*contr
368        //now the shell
369        qr=x*(rcore+thick)
370        contr = rhoshel-rhosolv
371        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
372        vol = 4*pi/3*(rcore+thick)^3
373        f += vol*bes*contr
374       
375        // normalize to particle volume and rescale from [-1] to [cm-1]
376        f2 = f*f/vol*1.0e8
377       
378        //scale if desired
379        f2 *= scale
380        // then add in the background
381        f2 += bkg
382       
383        return (f2)
384End
385
386
387Function fTwoShell(w,x) : FitFunc
388        Wave w
389        Variable x
390       
391        // variables are:
392        //[0] scale factor
393        //[1] radius of core []
394        //[2] SLD of the core   [-2]
395        //[3] thickness of shell 1 []
396        //[4] SLD of shell 1
397        //[5] thickness of shell 2 []
398        //[6] SLD of shell 2
399        //[7] SLD of the solvent
400        //[8] background        [cm-1]
401       
402        // All inputs are in ANGSTROMS
403        //OUTPUT is normalized by the particle volume, and converted to [cm-1]
404       
405       
406        Variable scale,rcore,thick1,thick2,rhocore,rhoshel1,rhoshel2,rhosolv,bkg
407        scale = w[0]
408        rcore = w[1]
409        rhocore = w[2]
410        thick1 = w[3]
411        rhoshel1 = w[4]
412        thick2 = w[5]
413        rhoshel2 = w[6]
414        rhosolv = w[7]
415        bkg = w[8]
416       
417        // calculates scale *( f^2 + bkg)
418        Variable bes,f,vol,qr,contr,f2
419       
420        // core first, then add in shells
421        qr=x*rcore
422        contr = rhocore-rhoshel1
423        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
424        vol = 4*pi/3*rcore^3
425        f = vol*bes*contr
426        //now the shell (1)
427        qr=x*(rcore+thick1)
428        contr = rhoshel1-rhoshel2
429        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
430        vol = 4*pi/3*(rcore+thick1)^3
431        f += vol*bes*contr
432        //now the shell (2)
433        qr=x*(rcore+thick1+thick2)
434        contr = rhoshel2-rhosolv
435        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
436        vol = 4*pi/3*(rcore+thick1+thick2)^3
437        f += vol*bes*contr
438       
439        // normalize to particle volume and rescale from [-1] to [cm-1]
440        f2 = f*f/vol*1.0e8
441       
442        //scale if desired
443        f2 *= scale
444        // then add in the background
445        f2 += bkg
446       
447        return (f2)
448       
449End
450
451Function fThreeShell(w,x) : FitFunc
452        Wave w
453        Variable x
454       
455        // variables are:
456        //[0] scale factor
457        //[1] radius of core []
458        //[2] SLD of the core   [-2]
459        //[3] thickness of shell 1 []
460        //[4] SLD of shell 1
461        //[5] thickness of shell 2 []
462        //[6] SLD of shell 2
463        //[7] SLD of the solvent
464        //[8] background        [cm-1]
465       
466        // All inputs are in ANGSTROMS
467        //OUTPUT is normalized by the particle volume, and converted to [cm-1]
468       
469       
470        Variable scale,rcore,thick1,thick2,thick3,rhoshel1,rhoshel2,rhoshel3
471        Variable rhocore,rhosolv,bkg
472        scale = w[0]
473        rcore = w[1]
474        rhocore = w[2]
475        thick1 = w[3]
476        rhoshel1 = w[4]
477        thick2 = w[5]
478        rhoshel2 = w[6]
479        thick3 = w[7]
480        rhoshel3 = w[8]
481        rhosolv = w[9]
482        bkg = w[10]
483       
484        // calculates scale *( f^2 + bkg)
485        Variable bes,f,vol,qr,contr,f2
486       
487        // core first, then add in shells
488        qr=x*rcore
489        contr = rhocore-rhoshel1
490        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
491        vol = 4*pi/3*rcore^3
492        f = vol*bes*contr
493        //now the shell (1)
494        qr=x*(rcore+thick1)
495        contr = rhoshel1-rhoshel2
496        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
497        vol = 4*pi/3*(rcore+thick1)^3
498        f += vol*bes*contr
499        //now the shell (2)
500        qr=x*(rcore+thick1+thick2)
501        contr = rhoshel2-rhoshel3
502        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
503        vol = 4*pi/3*(rcore+thick1+thick2)^3
504        f += vol*bes*contr
505        //now the shell (3)
506        qr=x*(rcore+thick1+thick2+thick3)
507        contr = rhoshel3-rhosolv
508        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
509        vol = 4*pi/3*(rcore+thick1+thick2+thick3)^3
510        f += vol*bes*contr
511       
512        // normalize to particle volume and rescale from [-1] to [cm-1]
513        f2 = f*f/vol*1.0e8
514       
515        //scale if desired
516        f2 *= scale
517        // then add in the background
518        f2 += bkg
519       
520        return (f2)
521End
522
523
524Function fFourShell(w,x) : FitFunc
525        Wave w
526        Variable x
527       
528        // variables are:
529        //[0] scale factor
530        //[1] radius of core []
531        //[2] SLD of the core   [-2]
532        //[3] thickness of shell 1 []
533        //[4] SLD of shell 1
534        //[5] thickness of shell 2 []
535        //[6] SLD of shell 2
536        //[7] SLD of the solvent
537        //[8] background        [cm-1]
538       
539        // All inputs are in ANGSTROMS
540        //OUTPUT is normalized by the particle volume, and converted to [cm-1]
541       
542       
543        Variable scale,rcore,thick1,thick2,thick3,thick4
544        Variable rhoshel1,rhoshel2,rhoshel3,rhoshel4
545        Variable rhocore,rhosolv,bkg
546        scale = w[0]
547        rcore = w[1]
548        rhocore = w[2]
549        thick1 = w[3]
550        rhoshel1 = w[4]
551        thick2 = w[5]
552        rhoshel2 = w[6]
553        thick3 = w[7]
554        rhoshel3 = w[8]
555        thick4 = w[9]
556        rhoshel4 = w[10]
557        rhosolv = w[11]
558        bkg = w[12]
559       
560        // calculates scale *( f^2 + bkg)
561        Variable bes,f,vol,qr,contr,f2
562       
563        // core first, then add in shells
564        qr=x*rcore
565        contr = rhocore-rhoshel1
566        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
567        vol = 4*pi/3*rcore^3
568        f = vol*bes*contr
569        //now the shell (1)
570        qr=x*(rcore+thick1)
571        contr = rhoshel1-rhoshel2
572        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
573        vol = 4*pi/3*(rcore+thick1)^3
574        f += vol*bes*contr
575        //now the shell (2)
576        qr=x*(rcore+thick1+thick2)
577        contr = rhoshel2-rhoshel3
578        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
579        vol = 4*pi/3*(rcore+thick1+thick2)^3
580        f += vol*bes*contr
581        //now the shell (3)
582        qr=x*(rcore+thick1+thick2+thick3)
583        contr = rhoshel3-rhoshel4
584        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
585        vol = 4*pi/3*(rcore+thick1+thick2+thick3)^3
586        f += vol*bes*contr
587        //now the shell (4)
588        qr=x*(rcore+thick1+thick2+thick3+thick4)
589        contr = rhoshel4-rhosolv
590        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
591        vol = 4*pi/3*(rcore+thick1+thick2+thick3+thick4)^3
592        f += vol*bes*contr
593       
594       
595        // normalize to particle volume and rescale from [-1] to [cm-1]
596        f2 = f*f/vol*1.0e8
597       
598        //scale if desired
599        f2 *= scale
600        // then add in the background
601        f2 += bkg
602       
603        return (f2)
604End
605
606
607
608
609///////////////////////////////////////////////////////////////
610// smeared model calculation
611//
612// you don't need to do anything with this function, as long as
613// your OneShell works correctly, you get the resolution-smeared
614// version for free.
615//
616// this is all there is to the smeared model calculation!
617Function SmearedOneShell(s) : FitFunc
618        Struct ResSmearAAOStruct &s
619
620//      the name of your unsmeared model (AAO) is the first argument
621        Smear_Model_20(OneShell,s.coefW,s.xW,s.yW,s.resW)
622
623        return(0)
624End
625
626Function SmearedTwoShell(s) : FitFunc
627        Struct ResSmearAAOStruct &s
628
629//      the name of your unsmeared model (AAO) is the first argument
630        Smear_Model_20(TwoShell,s.coefW,s.xW,s.yW,s.resW)
631
632        return(0)
633End
634
635Function SmearedThreeShell(s) : FitFunc
636        Struct ResSmearAAOStruct &s
637
638//      the name of your unsmeared model (AAO) is the first argument
639        Smear_Model_20(ThreeShell,s.coefW,s.xW,s.yW,s.resW)
640
641        return(0)
642End
643
644Function SmearedFourShell(s) : FitFunc
645        Struct ResSmearAAOStruct &s
646
647//      the name of your unsmeared model (AAO) is the first argument
648        Smear_Model_20(FourShell,s.coefW,s.xW,s.yW,s.resW)
649
650        return(0)
651End
652
653
654
655///////////////////////////////////////////////////////////////
656
657
658// nothing to change here
659//
660//wrapper to calculate the smeared model as an AAO-Struct
661// fills the struct and calls the ususal function with the STRUCT parameter
662//
663// used only for the dependency, not for fitting
664//
665Function fSmearedOneShell(coefW,yW,xW)
666        Wave coefW,yW,xW
667       
668        String str = getWavesDataFolder(yW,0)
669        String DF="root:"+str+":"
670       
671        WAVE resW = $(DF+str+"_res")
672       
673        STRUCT ResSmearAAOStruct fs
674        WAVE fs.coefW = coefW   
675        WAVE fs.yW = yW
676        WAVE fs.xW = xW
677        WAVE fs.resW = resW
678       
679        Variable err
680        err = SmearedOneShell(fs)
681       
682        return (0)
683End
684
685//wrapper to calculate the smeared model as an AAO-Struct
686// fills the struct and calls the ususal function with the STRUCT parameter
687//
688// used only for the dependency, not for fitting
689//
690Function fSmearedTwoShell(coefW,yW,xW)
691        Wave coefW,yW,xW
692       
693        String str = getWavesDataFolder(yW,0)
694        String DF="root:"+str+":"
695       
696        WAVE resW = $(DF+str+"_res")
697       
698        STRUCT ResSmearAAOStruct fs
699        WAVE fs.coefW = coefW   
700        WAVE fs.yW = yW
701        WAVE fs.xW = xW
702        WAVE fs.resW = resW
703       
704        Variable err
705        err = SmearedTwoShell(fs)
706       
707        return (0)
708End
709
710//wrapper to calculate the smeared model as an AAO-Struct
711// fills the struct and calls the ususal function with the STRUCT parameter
712//
713// used only for the dependency, not for fitting
714//
715Function fSmearedThreeShell(coefW,yW,xW)
716        Wave coefW,yW,xW
717       
718        String str = getWavesDataFolder(yW,0)
719        String DF="root:"+str+":"
720       
721        WAVE resW = $(DF+str+"_res")
722       
723        STRUCT ResSmearAAOStruct fs
724        WAVE fs.coefW = coefW   
725        WAVE fs.yW = yW
726        WAVE fs.xW = xW
727        WAVE fs.resW = resW
728       
729        Variable err
730        err = SmearedThreeShell(fs)
731       
732        return (0)
733End
734
735//wrapper to calculate the smeared model as an AAO-Struct
736// fills the struct and calls the ususal function with the STRUCT parameter
737//
738// used only for the dependency, not for fitting
739//
740Function fSmearedFourShell(coefW,yW,xW)
741        Wave coefW,yW,xW
742       
743        String str = getWavesDataFolder(yW,0)
744        String DF="root:"+str+":"
745       
746        WAVE resW = $(DF+str+"_res")
747       
748        STRUCT ResSmearAAOStruct fs
749        WAVE fs.coefW = coefW   
750        WAVE fs.yW = yW
751        WAVE fs.xW = xW
752        WAVE fs.resW = resW
753       
754        Variable err
755        err = SmearedFourShell(fs)
756       
757        return (0)
758End
Note: See TracBrowser for help on using the repository browser.