Changeset 98ce141 in sasmodels
- Timestamp:
- Nov 27, 2016 4:25:24 PM (8 years ago)
- Branches:
- master, core_shell_microgels, costrafo411, magnetic_model, ticket-1257-vesicle-product, ticket_1156, ticket_1265_superball, ticket_822_more_unit_tests
- Children:
- 07300ea
- Parents:
- 797a8e3
- Location:
- sasmodels/models
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
sasmodels/models/stacked_disks.c
r3ac4e1b r98ce141 71 71 // loop for the structure factor S(q) 72 72 double qd_cos_alpha = q*d*cos_alpha; 73 //d*cos_alpha is the projection of d onto q (in other words the component 74 //of d that is parallel to q. 73 75 double debye_arg = -0.5*square(qd_cos_alpha*sigma_dnn); 74 76 double sq=0.0; … … 79 81 sq = 1.0 + 2.0*sq/n_stacking; 80 82 81 return pq * sq; 83 return pq * sq * n_stacking; 84 // volume normalization should be per disk not per stack but form_volume 85 // is per stack so correct here for now. Could change form_volume but 86 // if one ever wants to use P*S we need the ER based on the total volume 82 87 } 83 88 -
sasmodels/models/stacked_disks.py
ref5a314 r98ce141 1 1 r""" 2 This model provides the form factor, $P(q)$, for stacked discs (tactoids)3 with a core/layer structure where the form factor is normalized by the volume4 of the cylinder. Assuming the next neighbor distance (d-spacing) in a stack5 of parallel discs obeys a Gaussian distribution, a structure factor $S(q)$6 proposed by Kratky and Porod in 1949 is used in this function.7 8 Note that the resolution smearing calculation uses 76 Gauss quadrature points9 to properly smear the model since the function is HIGHLY oscillatory,10 especially around the q-values that correspond to the repeat distance of11 the layers.12 13 The 2D scattering intensity is the same as 1D, regardless of the orientation14 of the q vector which is defined as15 16 .. math:: q = \sqrt{q_x^2 + q_y^2}17 18 2 Definition 19 3 ---------- 20 4 5 This model provides the form factor, $P(q)$, for stacked discs (tactoids) 6 with a core/layer structure which is constructed itself as $P(q) S(Q)$ 7 multiplying 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 9 neighbor distance (d-spacing) in the stack of parallel discs obeys a Gaussian 10 distribution. As such the normalization of this "composite" form factor is 11 relative to the individual disk volume, not the volume of the stack of disks. 12 This model is appropriate for example for non non exfoliated clay particles such 13 as Laponite. 14 21 15 .. figure:: img/stacked_disks_geometry.png 22 16 17 Geometry of a single core/layer disk 18 23 19 The scattered intensity $I(q)$ is calculated as 24 20 … … 26 22 27 23 I(q) = N\int_{0}^{\pi /2}\left[ \Delta \rho_t 28 \left( V_t f_t(q) - V_c f_c(q)\right) + \Delta \rho_c V_c f_c(q) 29 \right]^2 S(q)\sin{\alpha}\ d\alpha + \text{background} 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} 30 27 31 28 where the contrast … … 35 32 \Delta \rho_i = \rho_i - \rho_\text{solvent} 36 33 37 and $N$ is the number of discs per unit volume, 38 $\alpha$ is the angle between the axis of the disc and $q$, 39 and $V_t$ and $V_c$ are the total volume and the core volume of 40 a single disc, respectively. 41 42 .. math:: 43 44 \left\langle f_{t}^2(q)\right\rangle_{\alpha} = 45 \int_{0}^{\pi/2}\left[ 34 and $N$ is the number of individual (single) discs per unit volume, $\alpha$ is 35 the angle between the axis of the disc and $q$, and $V_t$ and $V_c$ are the 36 total volume and the core volume of a single disc, respectively, and 37 38 .. math:: 39 40 f_t(q,\alpha) = 46 41 \left(\frac{\sin(q(d+h)\cos{\alpha})}{q(d+h)\cos{\alpha}}\right) 47 42 \left(\frac{2J_1(qR\sin{\alpha})}{qR\sin{\alpha}} \right) 48 \right]^2 \sin{\alpha}\ d\alpha 49 50 \left\langle f_{c}^2(q)\right\rangle_{\alpha} = 51 \int_{0}^{\pi/2}\left[ 43 44 f_c(q,\alpha) = 52 45 \left(\frac{\sin(qh)\cos{\alpha})}{qh\cos{\alpha}}\right) 53 46 \left(\frac{2J_1(qR\sin{\alpha})}{qR\sin{\alpha}}\right) 54 \right]^2 \sin{\alpha}\ d\alpha55 47 56 48 where $d$ = thickness of the layer (*thick_layer*), … … 59 51 .. math:: 60 52 61 S(q ) = 1 + \frac{1}{2}\sum_{k=1}^n(n-k)\cos{(kDq\cos{\alpha})}62 \exp\left[ -k(q \cos{\alpha})^2\sigma_d/2\right]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] 63 55 64 56 where $n$ is the total number of the disc stacked (*n_stacking*), 65 57 $D = 2(d+h)$ is the next neighbor center-to-center distance (d-spacing), 66 58 and $\sigma_d$ = the Gaussian standard deviation of the d-spacing (*sigma_d*). 59 Note that $D\cos(\alpha)$ is the component of $D$ parallel to $q$ and the last 60 term in the equation above is effectively a Debye-Waller factor term. 67 61 68 62 .. note:: 69 Each assembly in the stack is layer/core/layer, so the spacing of the 63 64 1. Each assembly in the stack is layer/core/layer, so the spacing of the 70 65 cores is core plus two layers. The 2nd virial coefficient of the cylinder 71 66 is calculated based on the *radius* and *length* … … 74 69 is applied. 75 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 76 76 To provide easy access to the orientation of the stacked disks, we define 77 77 the axis of the cylinder using two angles $\theta$ and $\varphi$. … … 79 79 .. figure:: img/cylinder_angle_definition.jpg 80 80 81 Examples of the angles against 82 the detector plane. 83 84 85 Our model uses the form factor calculations implemented in a c-library provided 86 by the NIST Center for Neutron Research (Kline, 2006) 81 Examples of the angles against the detector plane. 82 83 84 Our model is derived from the form factor calculations implemented in a 85 c-library provided by the NIST Center for Neutron Research 86 (Kline, 2006)\ [#CIT_Kline]_ 87 87 88 88 References 89 89 ---------- 90 90 91 A Guinier and G Fournet, *Small-Angle Scattering of X-Rays*, 92 John Wiley and Sons, New York, 195593 94 O Kratky and G Porod, *J. Colloid Science*, 4, (1949) 35 95 96 J S Higgins and H C Benoit, *Polymers and Neutron Scattering*, 97 Clarendon, Oxford, 1994 98 99 **Author:** NIST IGOR/DANSE **on:** pre 2010 100 101 * *Last Modified by:** Piotr Rozyczko **on:** February 18, 2016102 103 * *Last Reviewed by:** Richard Heenan **on:** March 22, 201691 .. [#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 98 Authorship 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 104 """ 105 105 … … 107 107 108 108 name = "stacked_disks" 109 title = " "109 title = "Form factor for a stacked set of non exfoliated core/shell disks" 110 110 description = """\ 111 111 One layer of disk consists of a core, a top layer, and a bottom layer.
Note: See TracChangeset
for help on using the changeset viewer.