source: sasmodels/sasmodels/models/stacked_disks.py @ 98ce141

core_shell_microgelscostrafo411magnetic_modelticket-1257-vesicle-productticket_1156ticket_1265_superballticket_822_more_unit_tests
Last change on this file since 98ce141 was 98ce141, checked in by butler, 7 years ago

correct normalization problem in stacked_disk fixes #789. Also heavy
cleanup and correction of equations to match code as well as
standardization of the model documentation addressing #646.

  • Property mode set to 100644
File size: 8.3 KB
Line 
1r"""
2Definition
3----------
4
5This model provides the form factor, $P(q)$, for stacked discs (tactoids)
6with a core/layer structure which is constructed itself as $P(q) S(Q)$
7multiplying a $P(q)$ for individual core/layer disks by a structure factor
8$S(q)$ proposed by Kratky and Porod in 1949\ [#CIT1949]_ assuming the next
9neighbor distance (d-spacing) in the stack of parallel discs obeys a Gaussian
10distribution. As such the normalization of this "composite" form factor is
11relative to the individual disk volume, not the volume of the stack of disks.
12This model is appropriate for example for non non exfoliated clay particles such
13as Laponite.
14
15.. figure:: img/stacked_disks_geometry.png
16
17   Geometry of a single core/layer disk
18
19The scattered intensity $I(q)$ is calculated as
20
21.. math::
22
23    I(q) = N\int_{0}^{\pi /2}\left[ \Delta \rho_t
24    \left( V_t f_t(q,\alpha) - V_c f_c(q,\alpha)\right) + \Delta
25    \rho_c V_c f_c(q,\alpha)\right]^2 S(q,\alpha)\sin{\alpha}\ d\alpha
26    + \text{background}
27
28where the contrast
29
30.. math::
31
32    \Delta \rho_i = \rho_i - \rho_\text{solvent}
33
34and $N$ is the number of individual (single) discs per unit volume, $\alpha$ is
35the angle between the axis of the disc and $q$, and $V_t$ and $V_c$ are the
36total volume and the core volume of a single disc, respectively, and
37
38.. math::
39
40    f_t(q,\alpha) =
41    \left(\frac{\sin(q(d+h)\cos{\alpha})}{q(d+h)\cos{\alpha}}\right)
42    \left(\frac{2J_1(qR\sin{\alpha})}{qR\sin{\alpha}} \right)
43
44    f_c(q,\alpha) =
45    \left(\frac{\sin(qh)\cos{\alpha})}{qh\cos{\alpha}}\right)
46    \left(\frac{2J_1(qR\sin{\alpha})}{qR\sin{\alpha}}\right)
47
48where $d$ = thickness of the layer (*thick_layer*),
49$2h$ = core thickness (*thick_core*), and $R$ = radius of the disc (*radius*).
50
51.. math::
52
53    S(q,\alpha) = 1 + \frac{1}{2}\sum_{k=1}^n(n-k)\cos{(kDq\cos{\alpha})}
54    \exp\left[ -k(q)^2(D\cos{\alpha}~\sigma_d)^2/2\right]
55
56where $n$ is the total number of the disc stacked (*n_stacking*),
57$D = 2(d+h)$ is the next neighbor center-to-center distance (d-spacing),
58and $\sigma_d$ = the Gaussian standard deviation of the d-spacing (*sigma_d*).
59Note that $D\cos(\alpha)$ is the component of $D$ parallel to $q$ and the last
60term in the equation above is effectively a Debye-Waller factor term.
61
62.. note::
63
64    1. Each assembly in the stack is layer/core/layer, so the spacing of the
65    cores is core plus two layers. The 2nd virial coefficient of the cylinder
66    is calculated based on the *radius* and *length*
67    = *n_stacking* * (*thick_core* + 2 * *thick_layer*)
68    values, and used as the effective radius for $S(Q)$ when $P(Q) * S(Q)$
69    is applied.
70
71    2. the resolution smearing calculation uses 76 Gaussian quadrature points
72    to properly smear the model since the function is HIGHLY oscillatory,
73    especially around the q-values that correspond to the repeat distance of
74    the layers.
75
76To provide easy access to the orientation of the stacked disks, we define
77the axis of the cylinder using two angles $\theta$ and $\varphi$.
78
79.. figure:: img/cylinder_angle_definition.jpg
80
81    Examples of the angles against the detector plane.
82
83
84Our model is derived from the form factor calculations implemented in a
85c-library provided by the NIST Center for Neutron Research
86(Kline, 2006)\ [#CIT_Kline]_
87
88References
89----------
90
91.. [#CIT1949] O Kratky and G Porod, *J. Colloid Science*, 4, (1949) 35
92.. [#CIT_Kline] S R Kline, *J Appl. Cryst.*, 39 (2006) 895
93.. [#] J S Higgins and H C Benoit, *Polymers and Neutron Scattering*,
94   Clarendon, Oxford, 1994
95.. [#] A Guinier and G Fournet, *Small-Angle Scattering of X-Rays*,
96   John Wiley and Sons, New York, 1955
97
98Authorship and Verification
99----------------------------
100
101* **Author:** NIST IGOR/DANSE **Date:** pre 2010
102* **Last Modified by:** Paul Butler and Paul Kienzle **on:** November 26, 2016
103* **Last Reviewed by:** Paul Butler and Paul Kienzle **on:** November 26, 2016
104"""
105
106from numpy import inf
107
108name = "stacked_disks"
109title = "Form factor for a stacked set of non exfoliated core/shell disks"
110description = """\
111    One layer of disk consists of a core, a top layer, and a bottom layer.
112    radius =  the radius of the disk
113    thick_core = thickness of the core
114    thick_layer = thickness of a layer
115    sld_core = the SLD of the core
116    sld_layer = the SLD of the layers
117    n_stacking = the number of the disks
118    sigma_d =  Gaussian STD of d-spacing
119    sld_solvent = the SLD of the solvent
120    """
121category = "shape:cylinder"
122
123# pylint: disable=bad-whitespace, line-too-long
124#   ["name", "units", default, [lower, upper], "type","description"],
125parameters = [
126    ["thick_core",  "Ang",        10.0, [0, inf],    "volume",      "Thickness of the core disk"],
127    ["thick_layer", "Ang",        10.0, [0, inf],    "volume",      "Thickness of layer each side of core"],
128    ["radius",      "Ang",        15.0, [0, inf],    "volume",      "Radius of the stacked disk"],
129    ["n_stacking",  "",            1.0, [0, inf],    "volume",      "Number of stacked layer/core/layer disks"],
130    ["sigma_d",     "Ang",         0,   [0, inf],    "",            "Sigma of nearest neighbor spacing"],
131    ["sld_core",    "1e-6/Ang^2",  4,   [-inf, inf], "sld",         "Core scattering length density"],
132    ["sld_layer",   "1e-6/Ang^2",  0.0, [-inf, inf], "sld",         "Layer scattering length density"],
133    ["sld_solvent", "1e-6/Ang^2",  5.0, [-inf, inf], "sld",         "Solvent scattering length density"],
134    ["theta",       "degrees",     0,   [-inf, inf], "orientation", "Orientation of the stacked disk axis w/respect incoming beam"],
135    ["phi",         "degrees",     0,   [-inf, inf], "orientation", "Orientation of the stacked disk in the plane of the detector"],
136    ]
137# pylint: enable=bad-whitespace, line-too-long
138
139source = ["lib/polevl.c", "lib/sas_J1.c", "lib/gauss76.c", "stacked_disks.c"]
140
141demo = dict(background=0.001,
142            scale=0.01,
143            thick_core=10.0,
144            thick_layer=10.0,
145            radius=15.0,
146            n_stacking=1,
147            sigma_d=0,
148            sld_core=4,
149            sld_layer=0.0,
150            sld_solvent=5.0,
151            theta=0,
152            phi=0)
153#After redefinition to spherical coordinates find new reasonable test values
154#tests = [
155#    # Accuracy tests based on content in test/utest_extra_models.py.
156#    # Added 2 tests with n_stacked = 5 using SasView 3.1.2 - PDB
157#    [{'thick_core': 10.0,
158#      'thick_layer': 15.0,
159#      'radius': 3000.0,
160#      'n_stacking': 1.0,
161#      'sigma_d': 0.0,
162#      'sld_core': 4.0,
163#      'sld_layer': -0.4,
164#      'solvent_sd': 5.0,
165#      'theta': 0.0,
166#      'phi': 0.0,
167#      'scale': 0.01,
168#      'background': 0.001,
169#     }, 0.001, 5075.12],
170
171#    [{'thick_core': 10.0,
172#      'thick_layer': 15.0,
173#      'radius': 3000.0,
174#      'n_stacking': 5.0,
175#      'sigma_d': 0.0,
176#      'sld_core': 4.0,
177#      'sld_layer': -0.4,
178#      'solvent_sd': 5.0,
179#      'theta': 0.0,
180#      'phi': 0.0,
181#      'scale': 0.01,
182#      'background': 0.001,
183#     }, 0.001, 5065.12793824],
184
185#    [{'thick_core': 10.0,
186#      'thick_layer': 15.0,
187#      'radius': 3000.0,
188#      'n_stacking': 5.0,
189#      'sigma_d': 0.0,
190#      'sld_core': 4.0,
191#      'sld_layer': -0.4,
192#      'solvent_sd': 5.0,
193#      'theta': 0.0,
194#      'phi': 0.0,
195#      'scale': 0.01,
196#      'background': 0.001,
197#     }, 0.164, 0.0127673597265],
198
199#    [{'thick_core': 10.0,
200#      'thick_layer': 15.0,
201#      'radius': 3000.0,
202#      'n_stacking': 1.0,
203#      'sigma_d': 0.0,
204#      'sld_core': 4.0,
205#      'sld_layer': -0.4,
206#      'solvent_sd': 5.0,
207#      'theta': 0.0,
208#      'phi': 0.0,
209#      'scale': 0.01,
210#      'background': 0.001,
211#     }, [0.001, 90.0], [5075.12, 0.001]],
212
213#    [{'thick_core': 10.0,
214#      'thick_layer': 15.0,
215#      'radius': 3000.0,
216#      'n_stacking': 1.0,
217#      'sigma_d': 0.0,
218#      'sld_core': 4.0,
219#      'sld_layer': -0.4,
220#      'solvent_sd': 5.0,
221#      'theta': 0.0,
222#      'phi': 0.0,
223#      'scale': 0.01,
224#      'background': 0.001,
225#     }, ([0.4, 0.5]), [0.00105074, 0.00121761]],
226
227#    [{'thick_core': 10.0,
228#      'thick_layer': 15.0,
229#      'radius': 3000.0,
230#      'n_stacking': 1.0,
231#      'sigma_d': 0.0,
232#      'sld_core': 4.0,
233#      'sld_layer': -0.4,
234#     'solvent_sd': 5.0,
235#      'theta': 0.0,
236#      'phi': 0.0,
237#      'scale': 0.01,
238#      'background': 0.001,
239#     }, ([1.3, 1.57]), [0.0010039, 0.0010038]],
240#    ]
241# 21Mar2016   RKH notes that unit tests all have n_stacking=1, ought to test other values
242
Note: See TracBrowser for help on using the repository browser.