source: sans/Release/trunk/NCNR_User_Procedures/SANS/Analysis/Models/NewModels_2006/Vesicle_UL_v40.ipf @ 345

Last change on this file since 345 was 345, checked in by srkline, 15 years ago

Merging NewModels_2006 back into where it belongs

File size: 7.1 KB
Line 
1#pragma rtGlobals=1             // Use modern global access method.
2#pragma IgorVersion = 6.0
3
4////////////////////////////////////////////////
5// GaussUtils.proc and PlotUtils.proc MUST be included for the smearing calculation to compile
6// Adopting these into the experiment will insure that they are always present
7////////////////////////////////////////////////
8//
9// this function is for the form factor of a unilamellar vesicle
10//
11// the "scale" or "volume fraction" factor is the "material" volume fraction
12// - i.e. the volume fraction of surfactant added. NOT the excluded volume
13// of the vesicles, which can be much larger. See the Vesicle_Volume_N_Rg macro
14//
15// this excluded volume is accounted for in the structure factor calculations.
16//
17// a macro is also provided to calculate the number density, I(q=0)
18// the Rg, and all of the volumes of the particle.
19//
20// 13 JUL 04 SRK
21////////////////////////////////////////////////
22
23Proc PlotVesicleForm(num,qmin,qmax)
24        Variable num=128,qmin=0.001,qmax=0.7
25        Prompt num "Enter number of data points for model: "
26        Prompt qmin "Enter minimum q-value (^-1) for model: "
27        Prompt qmax "Enter maximum q-value (^-1) for model: "
28       
29        make/o/d/n=(num) xwave_vesicle,ywave_vesicle
30        xwave_vesicle =alog(log(qmin) + x*((log(qmax)-log(qmin))/num))
31        make/o/d coef_vesicle = {1.,100,30,6.36e-6,0.5e-6,0}
32        make/o/t parameters_vesicle = {"scale","core radius (A)","shell thickness (A)","Core and Solvent SLD (A-2)","Shell SLD (A-2)","bkg (cm-1)"}
33        Edit parameters_vesicle,coef_vesicle
34       
35        Variable/G root:g_vesicle
36        g_vesicle := VesicleForm(coef_vesicle,ywave_vesicle,xwave_vesicle)
37        Display ywave_vesicle vs xwave_vesicle
38        ModifyGraph log=1,marker=29,msize=2,mode=4
39        Label bottom "q (\\S-1\\M)"
40        Label left "Intensity (cm\\S-1\\M)"
41        AutoPositionWindow/M=1/R=$(WinName(0,1)) $WinName(0,2)
42       
43        AddModelToStrings("VesicleForm","coef_vesicle","vesicle")
44End
45
46///////////////////////////////////////////////////////////
47
48// - sets up a dependency to a wrapper, not the actual SmearedModelFunction
49Proc PlotSmearedVesicleForm(str)                                                               
50        String str
51        Prompt str,"Pick the data folder containing the resolution you want",popup,getAList(4)
52       
53        // if any of the resolution waves are missing => abort
54        if(ResolutionWavesMissingDF(str))               //updated to NOT use global strings (in GaussUtils)
55                Abort
56        endif
57       
58        SetDataFolder $("root:"+str)
59       
60        // Setup parameter table for model function
61        make/o/d smear_coef_vesicle = {1.,100,30,6.36e-6,0.5e-6,0}
62        make/o/t smear_parameters_vesicle = {"scale","core radius (A)","shell thickness (A)","Core and Solvent SLD (A-2)","Shell SLD (A-2)","bkg (cm-1)"}
63        Edit smear_parameters_vesicle,smear_coef_vesicle
64       
65        // output smeared intensity wave, dimensions are identical to experimental QSIG values
66        // make extra copy of experimental q-values for easy plotting
67       
68        Duplicate/O $(str+"_q") smeared_vesicle,smeared_qvals                           
69        SetScale d,0,0,"1/cm",smeared_vesicle                                                   
70                                       
71        Variable/G gs_vesicle=0
72        gs_vesicle := fSmearedVesicleForm(smear_coef_vesicle,smeared_vesicle,smeared_qvals)     //this wrapper fills the STRUCT
73       
74        Display smeared_vesicle vs smeared_qvals                                                                       
75        ModifyGraph log=1,marker=29,msize=2,mode=4
76        Label bottom "q (\\S-1\\M)"
77        Label left "Intensity (cm\\S-1\\M)"
78        AutoPositionWindow/M=1/R=$(WinName(0,1)) $WinName(0,2)
79       
80        SetDataFolder root:
81        AddModelToStrings("SmearedVesicleForm","smear_coef_vesicle","vesicle")
82End
83
84
85//AAO version, uses XOP if available
86// simply calls the original single point calculation with
87// a wave assignment (this will behave nicely if given point ranges)
88Function VesicleForm(cw,yw,xw) : FitFunc
89        Wave cw,yw,xw
90       
91#if exists("VesicleFormX")
92        yw = VesicleFormX(cw,xw)
93#else
94        yw = fVesicleForm(cw,xw)
95#endif
96        return(0)
97End
98
99///////////////////////////////////////////////////////////////
100// unsmeared model calculation
101///////////////////////////
102Function fVesicleForm(w,x) : FitFunc
103        Wave w
104        Variable x
105       
106        // variables are:
107        //[0] scale factor
108        //[1] radius of core []
109        //[2] thickness of the shell    []
110        //[3] SLD of the core and solvent[-2]
111        //[4] SLD of the shell
112        //[5] background        [cm-1]
113       
114        // All inputs are in ANGSTROMS
115        //OUTPUT is normalized by the particle volume, and converted to [cm-1]
116       
117       
118        Variable scale,rcore,thick,rhocore,rhoshel,rhosolv,bkg
119        scale = w[0]
120        rcore = w[1]
121        thick = w[2]
122        rhocore = w[3]
123        rhosolv = rhocore
124        rhoshel = w[4]
125        bkg = w[5]
126       
127        // calculates scale *( f^2 + bkg)
128        Variable bes,f,vol,qr,contr,f2
129       
130        // core first, then add in shell
131        qr=x*rcore
132        contr = rhocore-rhoshel
133        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
134        vol = 4*pi/3*rcore^3
135        f = vol*bes*contr
136        //now the shell
137        qr=x*(rcore+thick)
138        contr = rhoshel-rhosolv
139        bes = 3*(sin(qr)-qr*cos(qr))/qr^3
140        vol = 4*pi/3*(rcore+thick)^3
141        f += vol*bes*contr
142       
143        // normalize to the particle volume and rescale from [-1] to [cm-1]
144        //note that for the vesicle model, the volume is ONLY the shell volume
145        vol = 4*pi/3*((rcore+thick)^3-rcore^3)
146        f2 = f*f/vol*1.0e8
147       
148        //scale if desired
149        f2 *= scale
150        // then add in the background
151        f2 += bkg
152       
153        return (f2)
154End
155
156//wrapper to calculate the smeared model as an AAO-Struct
157// fills the struct and calls the ususal function with the STRUCT parameter
158//
159// used only for the dependency, not for fitting
160//
161Function fSmearedVesicleForm(coefW,yW,xW)
162        Wave coefW,yW,xW
163       
164        String str = getWavesDataFolder(yW,0)
165        String DF="root:"+str+":"
166       
167        WAVE resW = $(DF+str+"_res")
168       
169        STRUCT ResSmearAAOStruct fs
170        WAVE fs.coefW = coefW   
171        WAVE fs.yW = yW
172        WAVE fs.xW = xW
173        WAVE fs.resW = resW
174       
175        Variable err
176        err = SmearedVesicleForm(fs)
177       
178        return (0)
179End
180
181// this is all there is to the smeared calculation!
182Function SmearedVesicleForm(s) :FitFunc
183        Struct ResSmearAAOStruct &s
184
185//      the name of your unsmeared model (AAO) is the first argument
186        Smear_Model_20(VesicleForm,s.coefW,s.xW,s.yW,s.resW)
187
188        return(0)
189End
190
191Macro Vesicle_Volume_N_Rg()
192        Variable totVol,core,shell,i0,nden,rhoCore,rhoShell,rhoSolvent
193        Variable phi
194       
195        if(Exists("coef_vesicle")!=1)
196                abort "You need to plot the unsmeared vesicle model first to create the coefficient table"
197        Endif
198        totvol=4*pi/3*(coef_vesicle[1]+coef_vesicle[2])^3
199        core=4*pi/3*(coef_vesicle[1])^3
200        shell = totVol-core
201       
202//      nden = phi/(shell volume) or phi/Vtotal
203        nden = coef_vesicle[0]/shell
204        rhoCore = coef_vesicle[3]
205        rhoShell = coef_vesicle[4]
206        rhoSolvent = rhoCore
207       
208        i0 = nden*shell*shell*(rhoCore-rhoShell)^2*1e8
209        Print "Total Volume [A^3] = ",totVol
210        Print "Core Volume [A^3] = ",core
211        Print "Shell Volume [A^3] = ",shell
212        Print "Material volume fraction = ",coef_vesicle[0]
213        Print "Excluded volume fraction = ",nden*totvol
214//      Print "I(q=0) = ",i0
215        Print "I(Q=0) = n Vshell^2(DR)^2 [1/cm] = ",i0
216        Print "Number Density [1/A^3]= ",nden
217//      Print "model I(0) = ",ywave_vesicle[0]
218//      Print "model/limit = ",ywave_vesicle[0]/i0
219       
220        CalcRg_Vesicle(coef_vesicle)
221End
222
223
224Function CalcRg_Vesicle(coef_vesicle)
225        Wave coef_vesicle
226
227        Variable Rc,Rsh,r1,r2,rs,ans
228       
229        Rc = coef_vesicle[1]
230        Rsh = Rc + coef_vesicle[2]
231        r1 = coef_vesicle[3]
232        r2 = coef_vesicle[4]
233        rs = coef_vesicle[3]
234       
235//      ans = 0
236//      ans = ( (r1-r2)/(r2-rs) )*Rc^5/Rsh^5 - 1
237//      ans /= ( (r1-r2)/(r2-rs) )*Rc^3/Rsh^3 - 1
238//      ans *= 3/5*Rsh^2
239//      Print "Rg of vesicle [A] = ",sqrt(ans)
240       
241        ans = 0
242        ans = Rc^5/Rsh^5 + 1
243        ans /= Rc^3/Rsh^3 + 1
244        ans *= 3/5*Rsh^2
245       
246        Print "Rg of vesicle [A] = ",sqrt(ans)
247End
Note: See TracBrowser for help on using the repository browser.