6. formfactor (ff)¶
Particle form factors
The scattering intensity I(q) of a single particle with real scattering length densities is calculated. If the scattering length density is not defined as e.g. for beaucage model the normalized particle form factor F(q) is calculated.
The scattering per particle is
I(q)= I_0 F(q)
with particle form factor
F(q)=<F_a(q)F^*_a(q)>=<|F_a(q)|^2>
and forward scattering per particle
I_0=V_p^2(\rho-\rho_{solvent})^2
Here V_p is particle volume, \rho the average scattering length density and <> indicates the ensemble average.
The particle scattering amplitude is
F_a(q)= \int_V b(r) e^{iqr} \mathrm{d}r / \int_V b(r) \mathrm{d}r = \sum_N b_i e^{iqr} / \sum_N b_i
In this module units for I(q) and I_0 are nm^2=10^{-14} cm^2 per particle.
The scattering of particles with concentration c in mol/liter in units of \frac{1}{cm} is I_{[1/cm]}(q)=N_A \frac{c}{1000} 10^{-14} I_{[nm^2]}(q).
The scattering of arbitrary shaped particles can be calculated by cloudScattering()
as a cloud of points representing the desired shape.
In the same way distributions of particles as e.g. clusters of particles or nanocrystals can be calculated.
Oriented scattering of e.g. oriented nanoclusters can be calculated by
orientedCloudScattering()
.
Methods to build clouds of scatterers e.g. a cube decorated with spheres at the corners can be found in Lattice with examples. The advantage here is that there is no double counted overlap.
- Some scattering length densities as guide to choose realistic values for SLD and solventSLD :
- neutron scattering unit nm-2:
- protonated polyethylene glycol = 0.71e-6 A-2 = 0.71e-4 nm-2
- protonated polyethylene =-0.315e-6 A-2 =-0.315e-4 nm-2
- SiO2 = 4.185e-6 A-2 = 4.185e-4 nm-2
- D2O = 6.335e-6 A-2 = 6.335e-4 nm-2
- H2O =-0.560e-6 A-2 =-0.560e-4 nm-2
- protein ≈ 2.0e-6 A-2 ≈ 2.0e-4 nm-2
- gold = 4.500e-6 A-2 = 4.500e-4 nm-2
- Xray scattering unit nm^-2:
- polyethylene glycol = 1.11e-3 nm-2 = 396 e/nm3
- polyethylene = 0.85e-3 nm-2 = 302 e/nm3
- SiO2 = 2.25e-3 nm-2 = 796 e/nm3
- D2O = 0.94e-3 nm-2 = 332 e/nm3
- H2O = 0.94e-3 nm-2 = 333 e/nm3
- protein ≈ 1.20e-3 nm-2 ≈ 430 e/nm3
- gold = 13.1e-3 nm-2 =4662 e/nm3
Density SiO2 = 2.65 g/ml quartz; ≈ 2.2 g/ml quartz glass
Return values are dataArrays were useful. To get only Y values use .Y
6.1. Form Factors¶
beaucage (q[, Rg, G, d]) |
Beaucage introduced a model based on the polymer fractal model. |
genGuinier (q[, Rg, A, alpha]) |
Generalized Guinier approximation for low wavevector q scattering q*Rg< 1-1.3 |
guinier (q[, Rg, A]) |
Classical Guinier |
guinierPorod (q, Rg, s, G, dd) |
Generalized Guinier-Porod Model with high Q power law. |
guinierPorod3d (q, Rg1, s1, Rg2, s2, G2, dd) |
Generalized Guinier-Porod Model with high Q power law with 3 length scales. |
gaussianChain (q, Rg[, nu]) |
General formfactor of a gaussian polymer chain with excluded volume parameter. |
ringPolymer (q, Rg) |
General formfactor of a polymer ring in theta solvent. |
sphere (q, radius[, contrast]) |
Scattering of a single homogeneous sphere. |
sphereCoreShell (q, Rc, Rs, bc, bs[, solventSLD]) |
Scattering of a spherical core shell particle. |
ellipsoid (q, Ra, Rb[, SLD, solventSLD, …]) |
Form factor for a simple ellipsoid (ellipsoid of revolution). |
disc (q, R, D, SLD[, solventSLD, alpha]) |
Disc form factor . |
multilayer (q, layerd, SLD[, gaussw, solventSLD]) |
Form factor of a multilayer with rectangular/Gaussian density profiles perpendicular to the layer. |
multiShellSphere (q, shellthickness, shellSLD) |
Scattering of spherical multi shell particle including linear contrast variation in subshells. |
multiShellEllipsoid (q, poleshells, …[, …]) |
Scattering of multi shell ellipsoidal particle with varying shell thickness at pole and equator. |
multiShellDisc (q, radialthickness, …[, …]) |
Multi shell disc in solvent averaged over axis orientations. |
multiShellCylinder (q, L, shellthickness, …) |
Multi shell cylinder with caps in solvent averaged over axis orientations. |
multilamellarVesicles (Q, R, N, phi[, …]) |
Scattering intensity of a multilamellar vesicle with random displacements of the inner vesicles [1]. |
cuboid (q, a[, b, c, SLD, solventSLD, N]) |
Formfactor of cuboid with different edge lengths. |
sphereFuzzySurface (q, R, sigmasurf, contrast) |
Scattering of a sphere with a fuzzy interface. |
sphereGaussianCorona (q, R, Rg, Ncoil, coilequR) |
Scattering of a sphere surrounded by gaussian coils as model for grafted polymers on particle e.g. |
sphereCoreShellGaussianCorona (q, Rc, Rs, Rg, …) |
Scattering of a core-shell particle surrounded by gaussian coils as model for grafted polymers on particle. |
pearlNecklace (Q, Rc, l, N[, A1, A2, A3, ms, mr]) |
Formfactor of a pearl necklace (freely jointed chain of pearls connected by rods) |
wormlikeChain (q, N, a[, R, SLD, solventSLD, …]) |
Scattering of a wormlike chain, which correctly reproduces the rigid-rod and random-coil limits. |
teubnerStrey (q, xi, d[, eta2]) |
Phenomenological model for the scattering intensity of a two-component system using the Teubner-Strey model [1]. |
ellipsoidFilledCylinder ([q, R, L, Ra, Rb, …]) |
Scattering of a single cylinder filled with ellipsoidal particles . |
superball (q, R, p[, SLD, solventSLD, nGrid, …]) |
A superball is a general mathematical shape that can be used to describe rounded cubes, sphere and octahedron’s. |
6.2. Cloud of scatterers¶
cloudScattering (q, cloud[, relError, V, …]) |
Scattering of a cloud of scatterers with variable scattering length. |
orientedCloudScattering (qxzw, cloud[, rms, …]) |
2D scattering of an oriented cloud of scatterers with equal or variable scattering length. |
6.3. Distribution¶
scatteringFromSizeDistribution (q, …[, …]) |
Scattering of objects with one multimodal parameter as e.g. |
Particle form factors
The scattering intensity I(q) of a single particle with real scattering length densities is calculated. If the scattering length density is not defined as e.g. for beaucage model the normalized particle form factor F(q) is calculated.
The scattering per particle is
I(q)= I_0 F(q)
with particle form factor
F(q)=<F_a(q)F^*_a(q)>=<|F_a(q)|^2>
and forward scattering per particle
I_0=V_p^2(\rho-\rho_{solvent})^2
Here V_p is particle volume, \rho the average scattering length density and <> indicates the ensemble average.
The particle scattering amplitude is
F_a(q)= \int_V b(r) e^{iqr} \mathrm{d}r / \int_V b(r) \mathrm{d}r = \sum_N b_i e^{iqr} / \sum_N b_i
In this module units for I(q) and I_0 are nm^2=10^{-14} cm^2 per particle.
The scattering of particles with concentration c in mol/liter in units of \frac{1}{cm} is I_{[1/cm]}(q)=N_A \frac{c}{1000} 10^{-14} I_{[nm^2]}(q).
The scattering of arbitrary shaped particles can be calculated by cloudScattering()
as a cloud of points representing the desired shape.
In the same way distributions of particles as e.g. clusters of particles or nanocrystals can be calculated.
Oriented scattering of e.g. oriented nanoclusters can be calculated by
orientedCloudScattering()
.
Methods to build clouds of scatterers e.g. a cube decorated with spheres at the corners can be found in Lattice with examples. The advantage here is that there is no double counted overlap.
- Some scattering length densities as guide to choose realistic values for SLD and solventSLD :
- neutron scattering unit nm-2:
- protonated polyethylene glycol = 0.71e-6 A-2 = 0.71e-4 nm-2
- protonated polyethylene =-0.315e-6 A-2 =-0.315e-4 nm-2
- SiO2 = 4.185e-6 A-2 = 4.185e-4 nm-2
- D2O = 6.335e-6 A-2 = 6.335e-4 nm-2
- H2O =-0.560e-6 A-2 =-0.560e-4 nm-2
- protein ≈ 2.0e-6 A-2 ≈ 2.0e-4 nm-2
- gold = 4.500e-6 A-2 = 4.500e-4 nm-2
- Xray scattering unit nm^-2:
- polyethylene glycol = 1.11e-3 nm-2 = 396 e/nm3
- polyethylene = 0.85e-3 nm-2 = 302 e/nm3
- SiO2 = 2.25e-3 nm-2 = 796 e/nm3
- D2O = 0.94e-3 nm-2 = 332 e/nm3
- H2O = 0.94e-3 nm-2 = 333 e/nm3
- protein ≈ 1.20e-3 nm-2 ≈ 430 e/nm3
- gold = 13.1e-3 nm-2 =4662 e/nm3
Density SiO2 = 2.65 g/ml quartz; ≈ 2.2 g/ml quartz glass
Return values are dataArrays were useful. To get only Y values use .Y
-
jscatter.formfactor.
beaucage
(q, Rg=1, G=1, d=3)[source]¶ Beaucage introduced a model based on the polymer fractal model.
Beaucage used the numerical integration form (Benoit, 1957) although the analytical integral form was available [1]. This is an artificial connection of Guinier and Porod Regime . Better use the polymer fractal model [1] used in gaussianChain.
Parameters: - q : array
Wavevector
- Rg : float
Radius of gyration in 1/q units
- G : float
Guinier scaling factor, transition between Guinier and Porod
- d : float
Porod exponent for large wavevectors
Returns: - dataArray [q,Fq]
Notes
Equation 9+10 in [1]
I(q) &= G e^{-q^2 R_g^2 / 3.} + C q^{-d} \left[erf(qR_g / 6^{0.5})\right]^{3d} C &= \frac{G d}{R_g^d} \left[\frac{6d^2}{(2+d)(2+2d)}\right]^{d / 2.} \Gamma(d/2)
with the Gamma function \Gamma(x) .
Polymer fractals:
d = 5/3 fully swollen chains,d = 2 ideal Gaussian chains andd = 3 collapsed chains. (volume scattering)d = 4 surface scattering at a sharp interface/surfaced = 6-dim rough surface area with a dimensionality dim between 2-3 (rough surface)d = 3 Volume scatteringd < r mass fractals (eg gaussian chain)The Beaucage model is used to analyze small-angle scattering (SAS) data from fractal and particulate systems. It models the Guinier and Porod regions with a smooth transition between them and yields a radius of gyration and a Porod exponent. This model is an approximate form of an earlier polymer fractal model that has been generalized to cover a wider scope. The practice of allowing both the Guinier and the Porod scale factors to vary independently during nonlinear least-squares fits introduces undesired artefact’s in the fitting of SAS data to this model.
[1] (1, 2, 3) Analysis of the Beaucage model Boualem Hammouda J. Appl. Cryst. (2010). 43, 1474–1478 http://dx.doi.org/10.1107/S0021889810033856
-
jscatter.formfactor.
cloudScattering
(q, cloud, relError=50, V=0, formfactor=None, rms=0, ffpolydispersity=0, ncpu=0)[source]¶ Scattering of a cloud of scatterers with variable scattering length. Uses multiprocessing.
Cloud can represent any object described by a cloud of scatterers with scattering amplitudes as constant, sphere scattering amplitude, Gaussian scattering amplitude or explicitly given one. The result is normalized by I_0=(\sum b_i)^2 to equal one for q=0 (except for polydispersity).
- .I0 represents the forward scattering if b_i=b_vV_{unit cell} with b_v as scattering length density in the unit cell.
- Remember that the atomic bond length are on the order 0.1-0.2 nm.
- Methods to build clouds of scatterers e.g. a cube decorated with spheres at the corners can be found in Lattice with examples.
- By default explicit spherical average is done. If rms and polydispersity are not needed the Debye-function can be used (for particle numbers<500 it is faster).
Parameters: - q : array, ndim= Nx1
wavevectors in 1/nm
- cloud : array Nx3 or Nx4
- Center of mass positions (in nm) of the N scatterers in the cloud.
- If given cloud[3] is the scattering length b_i at positions cloud[:3], otherwise b=1.
- To compare with material scattering length density b_v use b=b_vV_{unit cell} with b_v as scattering length density and V_{unit cell} as cloud unit cell volume.
- relError : float
- Determines calculation method.
- relError>1 Explicit calculation of spherical average with Fibonacci lattice on sphere
- of 2*relError+1 points. Already 150 gives good results (see Examples)
- 0<relError<1 Monte Carlo integration on sphere until changes in successive iterations
- become smaller than relError. (Monte carlo integration with pseudo random numbers, see sphereAverage). This might take long for too small error.
- relError=0 The Debye equation is used (no asymmetry factor beta, no rms, no ffpolydispersity).
- Computation is of order N^2 opposite to above which is order N. For about 1000 particles same computing time,for 500 Debye is 4 times faster than above. If beta, rms or polydispersity is needed use above.
- rms : float, default=0
Root mean square displacement =<u**2>**0.5 of the positions in cloud as random (Gaussian) displacements in nm. Displacement u is random for each orientation in sphere scattering. rms can be used to simulate a Debye-Waller factor.
- V : float, default=0
Volume of the scatterers for scattering amplitude (see formfactor).
- formfactor : None,’gauss’,’sphere’,’cube’
- Gridpoint scattering amplitudes F(q) are described by:
None : const scattering amplitude.
- ‘sphere’: Sphere scattering amplitude according to [3].
The sphere radius is R=(\frac{3V}{4\pi})^{1/3}
- ‘gauss’ : Gaussian function b_i(q)=b V exp(-\pi V^{2/3}q^2) according to [2].
The Gaussian shows no artificial minima compared to the sphere.
Explicit isotropic form factor ff as array with [q,ff] e.g. from multishell. The normalized scattering amplitude fa for each gridpoint is calculated as fa=ff**0.5/fa(0). Missing values are linear interpolated (np.interp), q values outside interval are mapped to qmin or qmax. If explicit formfactor is from an asymmetric object (e.g. cube) it is implicated that the explicit ff orientation is isotropic.
- ffpolydispersity : float
Polydispersity of the gridpoints in relative units for sphere, gauss, explicit. Assuming F(q*R) for each gridpoint F is scaled as F(q*f*R) with f as normal distribution around 1 and standard deviation ffpolydispersity. The scattering length b is scaled according to the respective volume change by f**3. (f<0 is set to zero . ) This results in a change of the forward scattering because of the stronger weight of larger objects.
- ncpu : int, default 0
- Number of cpus used in the pool for multiprocessing.
- not given or 0 : all cpus are used
- int>0 : min(ncpu, mp.cpu_count)
- int<0 : ncpu not to use
- 1 : single core usage for testing or comparing speed to Debye
Returns: - dataArray with columns [q, Pq, beta]
- .I0 : =I(q=0)=(\sum_N b_i)^2
- .sumblength : =\sum_N b_i
- .formfactoramplitude : formfactor amplitude of cloudpoints according to type for all q values.
- .formfactoramplitude_q : corresponding q values
Notes
We calculate the scattering amplitude F(q) for N particles in a volume V with scattering length density b(r)
F(q)= \int_V b(r) e^{i\mathbf{qr}} \mathrm{d}r / \int_V b(r) \mathrm{d}r = \sum_N b_i e^{i\mathbf{qr}} / \sum_N b_i
with the form factor P(Q) after explicit orientational average <>
P(Q)=< F(q) \cdot F^*(q) >=< |F(q)|^2 >
The scattering intensity of a single object represented by the cloud is
I(Q)=P(Q) \cdot (\int_V b(r) \mathrm{d}r)^2
beta is the asymmetry factor [1] beta =|< F(q) >|^2 / < |F(q)|^2 >
One has to expect a peak at q \approx 2\pi/d_{NN} with d_{NN} as the next neighbour distance between particles.
One might replace b_i = b_i(q) to include a formfactor amplitude of the particles as e.g. q dependent Xray scattering amplitude or the formfactor of the cloud particles.
Random displacements u_i lead to r_i=r_i+u_i and to the the Debye-Waller factor for Bragg peaks and diffusive scattering at higher q.
The explicit orientational average can be simplified using the Debye scattering equation [4]
I(Q)=\sum_i \sum_j b_i b_j \frac{\sin(qr_{ij})}{qr_{ij}} =\sum_i b_i^2 + 2\sum_i \sum_{j>i} b_i b_j \frac{\sin(qr_{ij})}{qr_{ij}}
Here no rms or ffpolydispersity are included. The calculation of beta requires an additional calculation.
The scattering of a cloud can represent the scattering of a cluster of particles with polydispersity and position distortion according to root mean square displacements (rms). Polydispersity and rms displacements are randomly changed within the orientational average to represent an ensemble average (opposite to the time average of a single cluster). See
latticeStructureFactor()
for nanocubes and example A nano cube build of different lattices .References
[1] (1, 2) - Kotlarchyk and S.-H. Chen, J. Chem. Phys. 79, 2461 (1983).1
[2] (1, 2) An improved method for calculating the contribution of solvent to the X-ray diffraction pattern of biological molecules Fraser R MacRae T Suzuki E IUCr Journal of Applied Crystallography 1978 vol: 11 (6) pp: 693-694 [3] (1, 2) X-ray diffuse scattering by proteins in solution. Consideration of solvent influence B. A. Fedorov, O. B. Ptitsyn and L. A. Voronin J. Appl. Cryst. (1974). 7, 181-186 doi: 10.1107/S0021889874009137 [4] (1, 2) Zerstreuung von Röntgenstrahlen Debye P. Annalen der Physik 1915 vol: 351 (6) pp: 809-823 DOI: 10.1002/andp.19153510606 Examples
The example compares to the analytic solution for an ellipsoid then for a cube. For other shapes the grid may be better rotated away from the object symmetry or a random grid should be used. The example shows a good approximation with NN=20. Because of the grid peak at q=2\pi/d_{NN} the grid scatterer distance d_{NN} should be d_{NN} < \frac{1}{3} 2\pi/q_{max} .
Inspecting A nano cube build of different lattices shows other possibilities building a grid. Also a pseudo random grid can be used
pseudoRandomLattice()
.# ellipsoid with grid build by mgrid import jscatter as js import numpy as np import matplotlib.pyplot as plt from mpl_toolkits.mplot3d import Axes3D # cubic grid points R=3;NN=20;relError=50 grid= np.mgrid[-R:R:1j*NN, -R:R:1j*NN,-2*R:2*R:2j*NN].reshape(3,-1).T # points inside of sphere with radius R p=1;p2=1*2 # p defines a superball with 1->sphere p=inf cuboid .... inside=lambda xyz,R1,R2,R3:(np.abs(xyz[:,0])/R1)**p2+(np.abs(xyz[:,1])/R2)**p2+(np.abs(xyz[:,2])/R3)**p2<=1 insidegrid=grid[inside(grid,R,R,2*R)] q=np.r_[0:5:0.1] p=js.grace() p.title('compare form factors of an ellipsoid') ffe=js.ff.cloudScattering(q,insidegrid,relError=relError) p.plot(ffe,legend='cloud ff explicit') ffa=js.ff.ellipsoid(q,2*R,R) p.plot(ffa.X,ffa.Y/ffa.I0,li=1,sy=0,legend='analytic formula') p.legend() # show only each 20th point js.mpl.scatter3d(insidegrid[::10,:])
# cube # grid points generated by cubic grid import jscatter as js import numpy as np q=np.r_[0.1:5:0.1] p=js.grace() R=3;N=10;relError=0.01 # random points on sphere grid= js.sf.scLattice(R/N,N) ffe=js.ff.cloudScattering(q,grid,relError=relError) p.plot(ffe,legend='cloud ff explicit 10') # each point has a cube around it including the border ffa=js.ff.cuboid(q,2*R+R/N) p.plot(ffa.X,ffa.Y/ffa.I0,li=1,sy=0,legend='analytic formula') p.yaxis(scale='l') p.title('compare form factors of an cube') p.legend(x=2,y=0.1)
An objekt with explicit given formfactor for each gridpoint.
# 5 coreshell particles in line with polydispersity rod0 = np.zeros([5, 3]) rod0[:, 1] = np.r_[0, 1, 2, 3, 4] * 4 q = js.loglist(0.01, 7, 100) cs = js.ff.sphereCoreShell(q=q, Rc=1, Rs=2, bc=0.1, bs=1, solventSLD=0) ffe = js.ff.cloudScattering(q, rod0, formfactor=cs,relError=100,ffpolydispersity=0.1) p=js.grace() p.plot(ffe)
Using cloudScattering as fit model.
We have to define a model that parametrizes the building of the cloud that we get fit parameters. As example we use two overlapping spheres. The model can be used to fit some data. The build of the model is important as it describes how the overlap is treated e.g. as average.
- We have to consider some points:
- It is important that the model is continuous in its parameters to avoid steps as any fit algorithm cannot handle this.
- We have to limit some parameters that make giant grids. Fit algorithm make first a small step then a large one to estimate a good step size for parameter changes. If in the dumbbell example the radii R1 or R2 is increased to >1000 then the grid size burst the RAM and we get a Memory Error. Use hard limits for the radii to a reasonable value as shown below (see setlimit).
- The argument “factor” limits the initial step size. Reduce it (default 100 -> [0.1..100]).
- In the below example the first fit is fast but bad as we find a local minimum. A global fit algorithm takes quite long but finds the correct solution.
import jscatter as js import numpy as np # #: test if distance from point on X axis isInside=lambda x,A,R:((x-np.r_[A,0,0])**2).sum(axis=1)**0.5<R #: model def dumbbell(q,A,R1,b1,bgr=0,dx=0.3,relError=50): # D sphere distance # R1, R2 radii # b1,b2 scattering length # bgr background # dx grid distance not a fit parameter!! R2=R1 b2=b1 mR=max(R1,R2) # xyz coordinates grid=np.mgrid[-A/2-mR:A/2+mR:dx,-mR:mR:dx,-mR:mR:dx].reshape(3,-1).T insidegrid=grid[isInside(grid,-A/2.,R1) | isInside(grid,A/2.,R2)] # add blength column insidegrid=np.c_[insidegrid,insidegrid[:,0]*0] # set the corresponding blength; the order is important as here b2 overwrites b1 insidegrid[isInside(insidegrid[:,:3],-A/2.,R1),3]=b1 insidegrid[isInside(insidegrid[:,:3],A/2.,R2),3]=b2 # and maybe a mix ; this depends on your model insidegrid[isInside(insidegrid[:,:3],-A/2.,R1) & isInside(insidegrid[:,:3],A/2.,R2),3]=(b2+b1)/2. # calc the scattering result=js.ff.cloudScattering(q,insidegrid,relError=relError) result.Y=result.Y*result.I0+bgr # add attributes for later usage result.A=A result.R1=R1 result.b1=b1 result.dx=dx result.insidegrid=insidegrid return result # # test it q=np.r_[0.01:5:0.02] data=dumbbell(q,3,2,1) # show result configuration js.mpl.scatter3d(data.insidegrid[:,0],data.insidegrid[:,1],data.insidegrid[:,2]) # # Fit your data like this. # It may be a good idea to use not the highest resolution in the beginning because of speed. # If you have a good set of starting parameters you can decrease dx. data2=data.prune(number=100) data2.makeErrPlot(yscale='l') data2=data.prune(number=100) data2.makeErrPlot(yscale='l') data2.setLimit(R1=[None,None,1,4],A=[None,None,1,10]) # this results in a fast but bad fit result # a local minima is found but the basics is working. data2.fit(model=dumbbell, freepar={'A':3,'R1':2.4,'b1':1}, fixpar={'dx':0.3,'bgr':0}, mapNames={'q':'X'},factor=1) # To get a good result we need to find the global minimum by a different algorithm ('differential evolution') # The limits are used as border to search in an limited area. # The fit takes about 3500 iterations (1000s on Ryzen 1600X 6 cores) data2.fit(model=dumbbell,method='differential_evolution', freepar={'A':3,'R1':2.4,'b1':1}, fixpar={'dx':0.3,'bgr':0}, mapNames={'q':'X'})
Fit a sphere formfactor.
The quality of the grid approximation (number of gridpoints) may improve the correct description of higher order minima.
import numpy as np import jscatter as js # a function to discriminate what is inside of the sphere # basically a superball p2=2 is a sphere inside=lambda xyz,R1,p2:(np.abs(xyz[:,0]))**p2+(np.abs(xyz[:,1]))**p2+(np.abs(xyz[:,2]))**p2<=R1**2 def test(q,R,b,p2=2,relError=20): # make cubic grid with right size (increase NN for better approximation) NN=20 grid= np.mgrid[-R:R:1j*NN, -R:R:1j*NN,-R:R:1j*NN].reshape(3,-1).T # cut the edges to get a sphere insidegrid=grid[inside(grid,R,p2)] # add scattering length for points # the average scattering length density is sum(b)/sphereVolume insidegrid=np.c_[insidegrid,insidegrid[:,0]*0] insidegrid[:,3]=b # calc formfactor (normalised) for a single sphere ffs=js.ff.cloudScattering(q,insidegrid,relError=relError) # the total scattering is sumblength**2 ffs.Y*=ffs.sumblength**2 # save radius and the grid for later ffs.R=R ffs.insidegrid=insidegrid return ffs ####main q=np.r_[0:3:0.01] sp=js.formfactor.sphere(q,3,1) sp.makeErrPlot(yscale='l') # show intermediate results sp.setlimit(R=[0.3,10]) # set some reasonable limits for R sp.fit(model=test, freepar={'b':6,'R':2.1}, fixpar={}, mapNames={'q':'X'}) # show the resulting sphere grid resultgrid=sp.lastfit.insidegrid js.mpl.scatter3d(resultgrid[:,0],resultgrid[:,1],resultgrid[:,2])
Here we compare explicit calculation with the Debye equation as the later gets quite slow for larger numbers.
import jscatter as js import numpy as np R=6;NN=20 q=np.r_[0:5:0.1] grid=js.formel.randomPointsInCube(10000)*R-R/2 ffe=js.ff.cloudScattering(q,grid,relError=150) # takes about 1.3 s on six core ffd=js.ff.cloudScattering(q,grid,relError=0) # takes about 11.4 s on six core grid=js.formel.randomPointsInCube(500)*R-R/2 ffe=js.ff.cloudScattering(q,grid,relError=150) # takes about 132 ms on six core ffd=js.ff.cloudScattering(q,grid,relError=0) # takes about 33 ms on six core p=js.grace() p.plot(ffe) p.plot(ffd) p.yaxis(scale='l')
-
jscatter.formfactor.
cuboid
(q, a, b=None, c=None, SLD=1, solventSLD=0, N=30)[source]¶ Formfactor of cuboid with different edge lengths.
Parameters: - q : array
Wavevector in 1/nm
- a,b,c : float, None
Edge length, for a=b=c its a cube, Units in nm. If b=None b=a. If c=None c=b.
- SLD : float, default =1
Scattering length density of cuboid.unit nm^-2 e.g. SiO2 = 4.186*1e-6 A^-2 = 4.186*1e-4 nm^-2 for neutrons
- solventSLD : float, default =0
Scattering length density of solvent. unit nm^-2 e.g. D2O = 6.335*1e-6 A^-2 = 6.335*1e-4 nm^-2 for neutrons
- N : int
Order for Gaussian integration over both phi and theta.
Returns: - dataArray [q,Iq]
- .I0 forward scattering
- .edges
- .contrast
Notes
I(q)=\rho^2V_{cube}^2 \int_{0}^{2\pi}\int_{0}^{\pi} \lvert sinc(q_xa/2 ) sinc(q_yb/2) sinc(q_zc/2)\rvert^2 \sin\theta d\theta d\phi
with q = (q_x,q_y,q_z) = (q\sin\theta\cos\phi,q\sin\theta\sin\phi,q\cos\theta) and contrast \rho [1].
In [1] the edge length is only half of it.
References
[1] (1, 2, 3) Analysis of small-angle scattering data from colloids and polymer solutions: modeling and least-squares fitting Pedersen, Jan Skov Advances in Colloid and Interface Science 70, 171 (1997) http://dx.doi.org/10.1016/S0001-8686(97)00312-6 Examples
import jscatter as js import numpy as np q=np.r_[0.1:5:0.01] p=js.grace() p.plot(js.ff.cuboid(q,60,4,6)) p.plot(js.ff.cuboid(q,10,4,60)) p.plot(js.ff.cuboid(q,11,11,11),li=1) p.yaxis(scale='l') p.xaxis(scale='l')
-
jscatter.formfactor.
cylinder
(q, L, radius, SLD, solventSLD=0, alpha=None)[source]¶ Cylinder form factor (open cap).
See multiShellCylinder
-
jscatter.formfactor.
disc
(q, R, D, SLD, solventSLD=0, alpha=None)[source]¶ Disc form factor .
Parameters: - q : array
Wavevectors, units 1/nm
- R : float
Radius in nm
- D : float
Thickness of shell
- SLD,solventSLD : float
Scattering length density in nm^-2.
- alpha : float, [float,float] , unit rad
Orientation, angle between the cylinder axis and the scattering vector q. 0 means parallel, pi/2 is perpendicular If alpha =[start,end] is integrated between start,end start > 0, end < pi/2
Notes
See multiShellCylinder
-
jscatter.formfactor.
ellipsoid
(q, Ra, Rb, SLD=1, solventSLD=0, alpha=None, tol=1e-06)[source]¶ Form factor for a simple ellipsoid (ellipsoid of revolution).
Parameters: - q : float
Scattering vector unit e.g. 1/A or 1/nm 1/Ra
- Ra : float
Radius rotation axis units in 1/unit(q)
- Rb : float
Radius rotated axis units in 1/unit(q)
- SLD : float, default =1
Scattering length density of unit nm^-2 e.g. SiO2 = 4.186*1e-6 A^-2 = 4.186*1e-4 nm^-2 for neutrons
- solventSLD : float, default =0
Scattering length density of solvent. unit nm^-2 e.g. D2O = 6.335*1e-6 A^-2 = 6.335*1e-4 nm^-2 for neutrons
- alpha : [float,float] , default [0,90]
Angle between rotation axis Ra and scattering vector q in unit grad Between these angles orientation is averaged alpha=0 axis and q are parallel, other orientation is averaged
- tol : float
relative tolerance for integration between alpha
Returns: - dataArray with columns [q; Iq; beta ]
- .RotationAxisRadius
- .RotatedAxisRadius
- .EllipsoidVolume
- .I0 forward scattering q=0
- beta is asymmetry factor according to [3]. \beta = |<F(Q)>|^2/<|F(Q)|^2> with scattering amplitude F(Q) and form factor P(Q)=<|F(Q)|²>
References
[1] Structure Analysis by Small-Angle X-Ray and Neutron Scattering Feigin, L. A, and D. I. Svergun, Plenum Press, New York, (1987). [2] http://www.ncnr.nist.gov/resources/sansmodels/Ellipsoid.html [3] (1, 2) - Kotlarchyk and S.-H. Chen, J. Chem. Phys. 79, 2461 (1983).
-
jscatter.formfactor.
ellipsoidFilledCylinder
(q=1, R=10, L=0, Ra=1, Rb=2, eta=0.1, SLDcylinder=0.1, SLDellipsoid=1, SLDmatrix=0, alpha=90, epsilon=None, fPY=1, dim=3)[source]¶ Scattering of a single cylinder filled with ellipsoidal particles .
A cylinder filled with ellipsoids of revolution with cylinder formfactor and ellipsoid scattering as described by Siefker [1]. Ellipsoids have a fluid like distribution and hard core interaction leading to Percus-Yevick structure factor between ellipsoids. Ellipsoids can be oriented along cylinder axis. If cylinders are in a lattice, the ellipsoid scattering (column 2) is observed in the diffusive scattering and the dominating cylinder contributes only to the bragg peaks as a form factor.
Parameters: - q : array
Wavevectors in units 1/nm
- R : float
Cylinder radius in nm
- L : float
Length of the cylinder in nm If zero infinite length is assumed, but absolute intensity is not valid, only relative intensity.
- Ra : float
Radius rotation axis units in nm
- Rb : float
Radius rotated axis units in nm
- eta : float
Volume fraction of ellipsoids in cylinder for use in Percus-Yevick structure factor. Radius in PY corresponds to sphere with same Volume as the ellipsoid.
- SLDcylinder : float,default 1
Scattering length density cylinder material in nm**-2
- SLDellipsoid : float,default 1
Scattering length density of ellipsoids in cylinder in nm**-2
- SLDmatrix : float
Scattering length density of the matrix outside the cylinder in nm**-2
- alpha : float, default 90
Orientation of the cylinder axis to wavevector in degrees
- epsilon : [float,float], default [0,90]
Orientation range of ellipsoids rotation axis relative to cylinder axis in degrees.
- fPY : float
Factor between radius of ellipsoids Rv (equivalent volume) and radius used in structure factor Rpy Rpy=fPY*(Ra*Rb*Rb)**(1/3)
- dim : 3,1, default 3
Dimensionality of the Percus-Yevick structure factor 1 is one dimensional stricture factor, anything else is 3 dimensional (normal PY)
Returns: - dataArray : [q,n*conv(ellipsoids,cylinder)*sf_b + cylinder,
n *conv(ellipsoids,cylinder)*sf_b, cylinder, n * ellipsoids, sf, beta_ellipsoids]
- Each contributing formfactor is given with its absolute contribution V^2contrast^2 (NOT normalized to 1)
- The observed structurefactor is sf\_b = S_{\beta}(q)=1+\beta (S(q)-1).
- beta_ellipsoids =\beta(q) is the asymmetry factor of Kotlarchyk and Chen [2].
- conv(ellipsoids,cylinder) -> ellipsoid formfactor convoluted with cylinder formfactor
- .ellipsoidNumberDensity -> n ellipsoid number density in cylinder volume
- .cylinderRadius
- .cylinderLength
- .cylinderVolume
- .ellipsoidRa
- .ellipsoidRb
- .ellipsoidRg
- .ellipsoidVolume
- .ellipsoidVolumefraction
- .ellipsoidNumberDensity unit 1/nm**3
- .alpha orientation range
- .ellipdoidAxisOrientation
References
[1] (1, 2) Confinement Facilitated Protein Stabilization As Investigated by Small-Angle Neutron Scattering. Siefker, J., Biehl, R., Kruteva, M., Feoktystov, A., & Coppens, M. O. (2018) Journal of the American Chemical Society, 140(40), 12720–12723. https://doi.org/10.1021/jacs.8b08454 [2] (1, 2) - Kotlarchyk and S.-H. Chen, J. Chem. Phys. 79, 2461 (1983).
Examples
import jscatter as js p=js.grace() q=js.loglist(0.01,5,800) ff=js.ff.ellipsoidFilledCylinder(q,L=100,R=5.4,Ra=1.63,Rb=1.63,eta=0.4,alpha=90,epsilon=[0,90]) p.plot(ff.X,ff[2],legend='convolution cylinder x ellipsoids') p.plot(ff.X,ff[3],legend='cylinder only') p.plot(ff.X,ff[4],legend='ellipsoid only') p.plot(ff.X,ff[5],legend='structure factor ellipsoids') p.plot(ff.X,ff.Y,legend='conv. ellipsoid + filled cylinder') p.legend() p.yaxis(scale='l',label='I(q)') p.xaxis(scale='n',label='q / nm\S-1') # an angular averaged formfactor def averageEFC(q,R,L,Ra,Rb,eta,alpha=[alpha0,alpha1],fPY=fPY): res=js.dL() alphas=np.deg2rad(np.r_[alpha0:alpha1:13j]) for alpha in alphas: ffe=js.ff.ellipsoidFilledCylinder(q,R=R,L=L,Ra=Ra,Rb=Rb,eta=ata,alpha=alpha,) res.append(ffe) result=res[0].copy() result.Y=scipy.integrate.simps(res.Y,alphas)/(alpha1-alpha0) return result
-
jscatter.formfactor.
gaussianChain
(q, Rg, nu=0.5)[source]¶ General formfactor of a gaussian polymer chain with excluded volume parameter.
For nu=0.5 this is the Debye model for Gaussian chain in theta solvent. nu>0.5 for good solvents,nu<0.5 for bad solvents.
Parameters: - q : array
Scattering vector, unit eg 1/A or 1/nm
- Rg : float
Radius of gyration, units in 1/unit(q)
- nu : float, default=0.5
ν is the excluded volume parameter, which is related to the Porod exponent d as ν = 1/d and [5/3 <= d <= 3].
Returns: - dataArray [q,Fq]
- .radiusOfGyration
- .nu excluded volume parameter
Notes
Rg^2=l^2 N^{2\nu} with monomer length l and monomer number N.
calcs
F(Q) = \frac{1}{\nu U^{\frac{1}{2\nu}}} \gamma_{inc}(\frac{1}{2\nu}, U) - \frac{1}{\nu U^{\frac{1}{\nu}}} \gamma_{inc}(\frac{1}{\nu}, U)
with U=(qR_g)^2 and \gamma_{inc} as lower incomplete gamma function.
The absolute scattering is proportional to b^2 N^2=b^2 (R_g/l)^{1/\nu} with monomer number N and monomer scattering length b.
From [1]: “Note that this model describing polymer chains with excluded volume applies only in the mass fractal range ([5/3 <= d <= 3]) and does not apply to surface fractals ([3 < d < 4]). It does not reproduce the rigid-rod limit (d = 1) because it assumes chain flexibility from the outset, nor does it describe semi-flexible chains ([1 < d < 5/3]). “
This model should be favoured compared to the Beaucage model as it is not an artificial connection between two regimes.
References
[1] (1, 2) Analysis of the Beaucage model Boualem Hammouda J. Appl. Cryst. (2010). 43, 1474–1478 http://dx.doi.org/10.1107/S0021889810033856 [2] SANS from homogeneous polymer mixtures: A unified overview. Hammouda, B. in Polymer Characteristics 87–133 (Springer-Verlag, 1993). doi:10.1007/BFb0025862 Examples
import jscatter as js import numpy as np q=js.loglist(0.1,8,100) p=js.grace() for nu in np.r_[0.3:0.61:0.05]: iq=js.ff.gaussianChain(q,2,nu) p.plot(iq,le='nu= $nu') p.yaxis(scale='l') p.xaxis(scale='l') p.legend(x=0.2,y=0.5)
-
jscatter.formfactor.
genGuinier
(q, Rg=1, A=1, alpha=0)[source]¶ Generalized Guinier approximation for low wavevector q scattering q*Rg< 1-1.3
Parameters: - q : array of float
Wavevector
- Rg : float
Radius of gyration in units=1/q
- alpha : float
Shape [α = 0] spheroid, [α = 1] rod-like [α = 2] plane
- A : float
Amplitudes
Returns: - dataArray [q,Fq]
Notes
- Quantitative analysis of particle size and shape starts with the Guinier approximations.
- For three-dimensional objects the Guinier approximation is given by I(q) = A e^{− Rg^2q^2/3}
- This approximation can be extended also to rod-like and plane objects by I(q) =(\alpha \pi q^{-\alpha}) A e^{ − Rg^2q^2/(3-\alpha) }
If the particle has one dimension of length L that is much larger than the others (i.e., elongated, rod-like, or worm-like), then there is a q range such that qR_c < 1 << qL, where α = 1.
References
[1] Form and structure of self-assembling particles in monoolein-bile salt mixtures Rex P. Hjelm, Claudio Schteingart, Alan F. Hofmann, and Devinderjit S. Sivia J. Phys. Chem., 99:16395–16406, 1995
-
jscatter.formfactor.
guinier
(q, Rg=1, A=1)[source]¶ Classical Guinier
see genGuinier with alpha=0
Parameters: - q :array
- A : float
- Rg : float
-
jscatter.formfactor.
guinierPorod
(q, Rg, s, G, dd)[source]¶ Generalized Guinier-Porod Model with high Q power law.
Parameters: - q : float
Wavevector in units of 1/nm
- Rg : float
Radii of gyration in units nm.
- s : float
Dimensionality parameter describing the low Q region.
- dd : float
Porod exponent describing the high Q slope.
- G : float
intensity
Returns: - dataArray [q,Iq]
- Iq scattering intensity
References
[1] A new Guinier/Porod Model B. Hammouda J. Appl. Cryst. (2010) 43, 716-719 Author M. Kruteva JCNS 2019
Examples
import jscatter as js q=js.loglist(0.01,5,300) I=js.ff.guinierPorod(q,s=0,Rg=5,G=1,dd=4) p=js.grace() p.plot(I) p.xaxis(scale='l',label='q / nm\S-1') p.yaxis(scale='l',label='I(q) / a.u.')
-
jscatter.formfactor.
guinierPorod3d
(q, Rg1, s1, Rg2, s2, G2, dd)[source]¶ Generalized Guinier-Porod Model with high Q power law with 3 length scales.
The model represents the most general case containing three Guinier regions [1].
Parameters: - q : float
Wavevector in units of 1/nm
- Rg1 : float
Radii of gyration for the short size of scattering object in units nm.
- Rg2 : float
Radii of gyration for the overall size of scattering object in units nm.
- s1 : float
Dimensionality parameter for the short size of scattering object (s1=1 for a cylinder)
- s2 : float
dimensionality parameter for the overall size of scattering object (s2=0 for a cylinder)
- G2 : float
Intensity for q=0.
- dd : float
Porod exponent
Returns: - dataArray [q,Iq]
- Iq scattering intensity
Notes
For a cylinder with length L and radius R (see [1]) R_{g2} = (L^2/12+R^2/2)^{\frac{1}{2}} and R_{g1}=R/\sqrt{2}
References
[1] (1, 2, 3) A new Guinier/Porod Model B. Hammouda J. Appl. Cryst. (2010) 43, 716-719 Author M. Kruteva JCNS 2019
Examples
import jscatter as js q=js.loglist(0.01,5,300) I=js.ff.guinierPorod3d(q,Rg1=1,s1=1,Rg2=10,s2=0,G2=1,dd=4) p=js.grace() p.plot(I) p.xaxis(scale='l',label='q / nm\S-1') p.yaxis(scale='l',label='I(q) / a.u.')
-
jscatter.formfactor.
multiShellCylinder
(q, L, shellthickness, shellSLD, solventSLD=0, alpha=None, h=None, nalpha=30, ncap=31)[source]¶ Multi shell cylinder with caps in solvent averaged over axis orientations.
Each shell has a constant SLD and may have a cap with same SLD sequence. Caps may be globular (barbell) or small (like lenses). For zero length L a lens shaped disc or a double sphere like shape is recovered.
Parameters: - q : array
Wavevectors, units 1/nm
- L : float
Length of cylinder, units nm L=0 infinite cylinder if h=None.
- shellthickness : list of float or float, all >0
Thickness of shells in sequence, units nm. radii r=cumulativeSum(shellthickness)
- shellSLD : list of float/list
Scattering length density of shells in nm^-2. A shell can be divided in sub shells if instead of a single float a list of floats is given. These list values are used as scattering length of equal thickness subshells. E.g. [1,2,[3,2,1]] results in the last shell with 3 subshell of equal thickness. The sum of subshell thickness is the thickness given in shellthickness. See second example. SiO2 = 4.186*1e-6 A^-2 = 4.186*1e-4 nm^-2
- solventSLD : float
Scattering length density of surrounding solvent in nm^-2. D2O = 6.335*1e-6 A^-2 = 6.335*1e-4 nm^-2
- h : float, default=None
Geometry of the caps with cap radii R=(r**2+h**2)**0.5 h is distance of cap center with radius R from the flat cylinder cap and r as radii of the cylinder shells.
- None No caps, flat ends as default.
- 0 cap radii equal cylinder radii (same shellthickness as cylinder shells)
- >0 cap radius larger cylinder radii as barbell
- <0 cap radius smaller cylinder radii as lens caps
- alpha : float, [float,float] , unit rad
Orientation, angle between the cylinder axis and the scattering vector q. 0 means parallel, pi/2 is perpendicular If alpha =[start,end] is integrated between start,end start > 0, end < pi/2
- nalpha : int, default 30
Number of points in Gauss integration along alpha.
- ncap : int, default=31
Number of points in Gauss integration for cap.
Returns: - dataArray [q ,Iq ]
- .outerCylinderVolume
- .Radius
- .cylinderLength
- .alpha
- .shellthickness
- .shellSLD
- .solventSLD
- .modelname
- .contrastprofile
- .capRadii
Notes
- Multishell of types:
- flat cap cylinder L>0, radii>0, h=None
- lens cap cylinder L>0, radii>0, h<0
- globular cap cylinder L>0, radii>0, h>0
- lens L=0, radii>0, h<0
- barbell no cylinder L=0, radii>0, h>0
- infinite flat disc L=0. h=None
References
Single cylinder
[1] Guinier, A. and G. Fournet, “Small-Angle Scattering of X-Rays”, John Wiley and Sons, New York, (1955) [2] http://www.ncnr.nist.gov/resources/sansmodels/Cylinder.html Double cylinder
[3] Use of viscous shear alignment to study anisotropic micellar structure by small-angle neutron scattering, J. B. Hayter and J. Penfold J. Phys. Chem., 88:4589–4593, 1984 [4] http://www.ncnr.nist.gov/resources/sansmodels/CoreShellCylinder.html Barbell, cylinder with small end-caps, circular lens
[5] Scattering from cylinders with globular end-caps Kaya (2004). J. Appl. Cryst. 37, 223-230] DOI: 10.1107/S0021889804000020 Scattering from capped cylinders. Addendum H. Kaya and Nicolas-Raphael de Souza J. Appl. Cryst. (2004). 37, 508-509 DOI: 10.1107/S0021889804005709 Examples
Alternating shells with different thickness 0.3 nm h2o and 0.2 nm d2o in vacuum:
import jscatter as js import numpy as np x=np.r_[0.0:10:0.01] ashell=js.ff.multiShellCylinder(x,20,[0.4,0.6]*5,[-0.56e-4,6.39e-4]*5) #plot it p=js.grace() p.multi(2,1) p[0].plot(ashell) p[1].plot(ashell.contrastprofile,li=1) # a contour of the SLDs
Double shell with exponential decreasing exterior shell to solvent scattering:
import jscatter as js import numpy as np x=np.r_[0.0:10:0.01] def doubleexpshells(q,L,d1,d2,e3,sd1,sd2,sol): # The third layer will have 9 subshells with combined thickness of e3. # The scattering length decays to e**(-3) in last subshell. return js.ff.multiShellCylinder(q,L,[d1,d2,e3],[sd1,sd2,((sd2-sol)*np.exp(-np.r_[0:3:9j])+sol)],solventSLD=sol) dde=doubleexpshells(x,10,0.5,0.5,3,1e-4,2e-4,0) #plot it p=js.grace() p.multi(2,1) p[0].plot(dde) p[1].plot(dde.contrastprofile,li=1) # a contour of the SLDs
Cylinder with cap:
x=np.r_[0.1:10:0.01] p=js.grace() p.title('Comparison of dumbbell cylinder with simple models') p.subtitle('thin lines correspond to simple models as sphere and dshell sphere') p.plot(js.ff.multiShellCylinder(x,0,[10],[1],h=0),sy=[1,0.5,2],le='simple sphere') p.plot(js.ff.sphere(x,10),sy=0,li=1) p.plot(js.ff.multiShellCylinder(x,0,[2,1],[1,2],h=0),sy=[1,0.5,3],le='double shell sphere') p.plot(js.ff.multiShellSphere(x,[2,1],[1,2]),sy=0,li=1) p.plot(js.ff.multiShellCylinder(x,10,[3],[20],h=-5),sy=[1,0.5,4],le='thin lens cap cylinder=flat cap cylinder') p.plot(js.ff.multiShellCylinder(x,10,[3],[20],h=None),sy=0,li=[1,2,1],le='flat cap cylinder') p.plot(js.ff.multiShellCylinder(x,10,[3],[20],h=-0.5),sy=0,li=[3,2,6],le='thick lens cap cylinder') p.yaxis(scale='l') p.xaxis(scale='l') p.legend(x=0.15,y=0.01)
-
jscatter.formfactor.
multiShellDisc
(q, radialthickness, shellthickness, shellSLD, solventSLD=0, alpha=None, nalpha=30)[source]¶ Multi shell disc in solvent averaged over axis orientations.
Parameters: - q : array
Wavevectors, units 1/nm
- radialthickness : float, all >0
Radial thickness of disc shells from inner to outer, units nm radii r=cumulativeSum(radialthickness)
- shellthickness : list of float or float, all >0
Thickness of shells from inner to outer, units nm. Innermost thickness is only once. total thickness = shellthickness[0]+2*cumulativeSum(shellthickness[1:])
- shellSLD : list of float/list
Scattering length density of shells in nm^-2. A shell can be divided in sub shells if instead of a single float a list of floats is given. These list values are used as scattering length of equal thickness subshells. E.g. [1,2,[3,2,1]] results in the last shell with 3 subshell of equal thickness. The sum of subshell thickness is the thickness given in shellthickness. See second example. SiO2 = 4.186*1e-6 A^-2 = 4.186*1e-4 nm^-2
- solventSLD : float
Scattering length density of surrounding solvent in nm^-2. D2O = 6.335*1e-6 A^-2 = 6.335*1e-4 nm^-2
- alpha : float, [float,float] , unit rad
Orientation, angle between the cylinder axis and the scattering vector q. 0 means parallel, pi/2 is perpendicular If alpha =[start,end] is integrated between start,end start > 0, end < pi/2
- nalpha : int, default 30
Number of points in Gauss integration along alpha.
Returns: - dataArray [q ,Iq ]
- .outerDiscVolume
- .radii
- .alpha
- .discthickness
- .shellSLD
- .solventSLD
- .modelname
References
[1] Guinier, A. and G. Fournet, “Small-Angle Scattering of X-Rays”, John Wiley and Sons, New York, (1955) Examples
Alternating shells with different thickness 0.3 nm h2o and 0.2 nm d2o in vacuum:
import jscatter as js import numpy as np x=np.r_[0.0:10:0.01] ashell=js.ff.multiShellDisc(x,[0.6,0.4]*2,[0.4,0.6]*2,[-0.56e-4,6.39e-4]*2) p=js.grace() p[0].plot(ashell) bshell=js.ff.multiShellDisc(x,2,2,6.39e-4) p[0].plot(bshell)
-
jscatter.formfactor.
multiShellEllipsoid
(q, poleshells, equatorshells, shellSLD, solventSLD=0, alpha=None, tol=1e-06)[source]¶ Scattering of multi shell ellipsoidal particle with varying shell thickness at pole and equator.
Shell thicknesses add up to form complex particles with any combination of axial ratios and shell thickness. A const axial ratio means different shell thickness at equator and pole.
Parameters: - q : array
Wavevectors, unit 1/nm
- equatorshells : list of float
Thickness of shells starting from inner most for rotated axis Re making the equator. unit nm. The absolute values are used.
- poleshells : list of float
Thickness of shells starting from inner most for rotating axis Rp pointing to pole. unit nm. The absolute values are used.
- shellSLD : list of float
List of scattering length densities of the shells in sequence corresponding to shellthickness. unit nm^-2.
- solventSLD : float, default=0
Scattering length density of the surrounding solvent. unit nm^-2
- alpha : [float,float], default [0,90]
Angular range of rotated axis to average over. Default is no preferred orientation.
- tol : float
Absolute tolerance for above adaptive integration of alpha.
Returns: - dataArray [q, Iq]
- Iq, scattering cross section in units nm**2
- .contrastprofile as radius and contrast values at edge points of equatorshells
- .equatorshellthicknes consecutive shell thickness
- .poleshellthickness
- .shellcontrast contrast of the shells to the solvent
- .equatorshellradii outer radius of the shells
- .poleshellradii
- .outerVolume Volume of complete sphere
- .I0 forward scattering for Q=0
References
[1] Structure Analysis by Small-Angle X-Ray and Neutron Scattering Feigin, L. A, and D. I. Svergun, Plenum Press, New York, (1987). [2] http://www.ncnr.nist.gov/resources/sansmodels/Ellipsoid.html [3] - Kotlarchyk and S.-H. Chen, J. Chem. Phys. 79, 2461 (1983).
Examples
Simple ellipsoid in vacuum:
x=np.r_[0.0:10:0.01] Rp=2. Re=1. ashell=js.ff.multiShellEllipsoid(x,Rp,Re,1) #plot it p=js.grace() p.multi(2,1) p[0].plot(ashell) p[1].plot(ashell.contrastprofile,li=1) # a contour of the SLDs
Alternating shells with thickness 0.3 nm h2o and 0.2 nm d2o in vacuum:
x=np.r_[0.0:10:0.01] shell=np.r_[[0.3,0.2]*3] sld=[-0.56e-4,6.39e-4]*3 # constant axial ratio for all shells but nonconstant shell thickness axialratio=2 ashell=js.ff.multiShellEllipsoid(x,axialratio*shell,shell,sld) # shell with constant shellthickness of one component and other const axialratio pshell=shell[:] pshell[0]=shell[0]*axialratio pshell[2]=shell[2]*axialratio pshell[4]=shell[4]*axialratio bshell=js.ff.multiShellEllipsoid(x,pshell,shell,sld) #plot it p=js.grace() p.multi(2,1) p[0].plot(ashell,le='const. axial ratio') p[1].plot(ashell.contrastprofile,li=2) # a contour of the SLDs p[0].plot(bshell,le='const shell thickness') p[1].plot(bshell.contrastprofile,li=2) # a contour of the SLDs p[0].legend()
double shell with exponential decreasing exterior shell to solvent scattering:
x=np.r_[0.0:10:0.01] def doubleexpshells(q,d1,ax,d2,e3,sd1,sd2,sol): shells =[d1 ,d2]+[e3]*9 shellsp=[d1*ax,d2]+[e3]*9 sld=[sd1,sd2]+list(((sd2-sol)*np.exp(-np.r_[0:3:9j]))) return js.ff.multiShellEllipsoid(q,shellsp,shells,sld,solventSLD=sol) dde=doubleexpshells(x,0.5,1,0.5,1,1e-4,2e-4,0) #plot it p=js.grace() p.multi(2,1) p[0].plot(dde) p[1].plot(dde.contrastprofile,li=1) # a countour of the SLDs
-
jscatter.formfactor.
multiShellSphere
(q, shellthickness, shellSLD, solventSLD=0)[source]¶ Scattering of spherical multi shell particle including linear contrast variation in subshells.
The results needs to be multiplied with the concentration to get the measured scattering. The resulting contrastprofile can be accessed as .contrastprofile
Parameters: - q : array
Wavevectors to calculate form factor, unit e.g. 1/nm.
- shellthickness : list of float
Thickness of shells starting from inner most, unit in nm. There is no limit for the number of shells.
- shellSLD : list of float or list
- List of scattering length densities of the shells in sequence corresponding to shellthickness. unit in nm**-2
- Innermost shell needs to be constant shell.
- If an element of the list is itself a list of SLD values it is interpreted as equal thick subshells with linear progress between SLD values in sum giving shellthickness. Here any shape can be approximated as sequence of linear pieces.
- If subshell list has only one float e.g. [1e.4] the second value is the SLD of the following shell.
- If empty list is given as [] the SLD of the previous and following shells are used as smooth transition.
- solventSLD : float, default=0
Scattering length density of the surrounding solvent. If equal to zero (default) then in profile the contrast is given. Unit in 1/nm**2
Returns: - dataArray [wavevector, Iq]
- Iq scattering cross section in units nm**2
- .contrastprofile as radius and contrast values at edge points
- .shellthickness consecutive shell thickness
- .shellcontrast contrast of the shells to the solvent
- .shellradii outer radius of the shells
- .slopes slope of linear increase of each shell
- .outerVolume Volume of complete sphere
- .I0 forward scattering for Q=0
Notes
The scattering intensity for a multishell particle with several subshells is
I(q) = \left( \sum_i f_a(q) \right)^2
The scattering amplitude of a subshell with inner and outer radius R_{i,o} is
f_a(q) = 4\pi\int_{R_i}^{R_o} \rho(r) \frac{sin(qr)}{qr}r^2dr
where we use always the scattering length density difference to the solvent (contrast) \rho(r) = \hat{\rho}(r) - \hat{\rho}_{solvent}.
For constant scattering length density \rho(r) = \rho we get
f_{a,const}(q) = \frac{4\pi}{3}r^3\rho \left. \frac{3(sin(qr)-qR cos(qr))}{(qr)^3}\right\rvert_{r=R_i}^{r=R_o}
with forward scattering contribution
f_{a,const}(q=0) = \frac{4\pi\rho}{3} (R_i^{3} - R_o^{3})
For a linear variation as \rho(r)=\Delta\rho(r-R_i)/d + \rho_i with \Delta\rho=\rho_o-\rho_i and thickness d=(R_o-R_i) we may sum a constant subshell as above with \rho(r)=\rho_i and contribution of the linear increase \rho(r)=\Delta\rho(r-R_i)/d resulting in
f_{a,lin}(q) =f_{a,const}(q) + \frac{4\pi\Delta\rho}{d} \left. \frac{(q(2r-R_i))sin(qr)-(q^2r(r-R_i)-2)cos(qr) }{q^4} \right\rvert_{r=R_i}^{r=R_o}
with the forward scattering contribution
f_{a,lin}(q=0)= f_{a,const}(q=0) + \frac{\pi \Delta\rho}{3 d} \left(R_{i} - R_{o}\right)^{2} \left(R_{i}^{2} + 2 R_{i} R_{o} + 3 R_{o}^{2}\right)
The solution is unstable (digital resolution) for really low QR values, which are set to the I0 scattering.
Examples
Alternating shells with 5 alternating thickness 0.4 nm and 0.6 nm with h2o, d2o scattering contrast in vacuum:
x=np.r_[0.0:10:0.01] ashell=js.ff.multiShellSphere(x,[0.4,0.6]*5,[-0.56e-4,6.39e-4]*5) #plot it p=js.grace() p.multi(2,1) p[0].plot(ashell) p[0].yaxis(scale='l') p[1].plot(ashell.contrastprofile,li=1) # a contour of the SLDs
Double shell with exponential decreasing exterior shell to solvent scattering:
x=np.r_[0.0:10:0.01] def doubleexpshells(q,d1,d2,e3,sd1,sd2,sol): return js.ff.multiShellSphere(q,[d1,d2,e3*3],[sd1,sd2,((sd2-sol)*np.exp(-np.r_[0:3:9j]))+sol],solventSLD=sol) dde=doubleexpshells(x,0.5,0.5,1,1e-4,2e-4,0) #plot it p=js.grace() p.multi(2,1) p[0].plot(dde) p[0].yaxis(scale='l') p[1].plot(dde.contrastprofile,li=1) # a contour of the SLDs
-
jscatter.formfactor.
multilamellarVesicles
(Q, R, N, phi, displace=0, dR=0, dN=0, layers=0, gaussw=None, ds=0, SLD=1, solventSLD=0, nGauss=100)[source]¶ Scattering intensity of a multilamellar vesicle with random displacements of the inner vesicles [1].
The result contains the full scattering, the structure factor of the lamella and a multilayer formfactor of the lamella layer structure. Other layer structures as mentioned in [2]. Multilayer formfactor is described in
multilayer()
.Parameters: - Q : float
Wavevector in 1/nm.
- R : float
Outer radius of the Vesicle in units nm.
- dR : float
Width of outer radius distribution in units nm.
- displace : float
Displacements of the vesicle centers in units nm. This describes the displacement steps in a random walk of the centers. displace=0 it is concentric, all have same center. displace< R/N.
- N : int
Number of lamella.
- dN : int, default=0
Width of distribution for number of lamella. (dN< 0.4 is single N) A zero truncated normal distribution is used with N>0 and N<R/displace. Check .Ndistribution and .Nweight = Nweight for the resulting distribution.
- phi : float
Volume fraction \phi of layers inside of vesicle.
- layers: list of float, default=0
Thickness of layer sequence in a lamella in units nm. Zero assumes infinite thin lamella with constant formfactor. List gives consecutive layer thickness from center to outside.
- gaussw : list of float, default None
Width of Gaussian in a layer. If float>0 a gaussian profile with this width is assumed. [1,None,1] is box in center and Gaussian width of 1 for outer layers. A gaussian needs a layers thickness of about the FWHM=2.35*sigma of the gaussian. See
multilayer()
.- SLD : list of float
Scattering length density of layers in lamella in nm^-2.
- ds : float
Thickness fluctuation of central layer in above lamella in nm. The variation is calculated as Gaussian distribution of the central layer (inner most one for even number of layers).
- solventSLD
Solvent scattering length density in nm^-2.
- nGauss : int, default 100
Number of Gaussian quadrature points in integration over dR distribution.
Returns: - dataArray with [q,I(q),S(q),F(q)]
- I(q)=S(q)F(q) scattering intensity
- S(q) multilamellar vesicle structure factor
- F(q) lamella formfactor
- .columnname=’q;Iq;Sq;Fq’
- .outerShellVolume
- .Ndistribution
- .Nweight
- .displace
- .phi
- .layerthickness
- .SLD
- .solventSLD
- .shellfluctuations=ds
- .preFactor=phi*Voutershell**2
Notes
The left shows a concentric lamellar structure. The right shows the random path of the consecutive centers of the spheres. See Multilamellar Vesicles for resulting scattering curves.
The function returns I(Q) as (see [1] equ. 17 )
I(Q)=\phi V_{outershell} S(Q) F(Q)
with the multishell structure factor S(Q) as described in [1]. For a single layer we have the formfactor F(Q)
F(Q)= ( \sum_i \rho_i d_i sinc( Q d_i) )^2
with \rho_i as scattering length density and thickness d_i. For a complex multilayer we find (see
multilayer()
)F(Q)= \sum_{i,j} \rho_i\rho_j d_i d_j sinc(qd_i)sinc(qd_j)cos(q(a_i-a_j))
with a_i as positions of the layers.
The amphiphile concentration phi is roughly given by phi = d/a, with d being the bilayer thickness and a being the spacing of the shells. The spacing of the shells is given by the scattering vector of the first correlation peak, i.e., a = 2pi/Q. Once the MLVs leave considerable space between each other then phi < d/a holds. This condition coincides with the assumption of dilution of the Guinier law. (from [1])
Structure factor part is normalized that S(0)=\sum_{j=1}^N (j/N)^2
To use a different shell form factor the structure factor is given explicitly.
Comparing a unilamellar vesicle (N=1) with multiShellSphere shows that R is located in the center of the shell:
Q=js.loglist(0.0001,5,1000)#np.r_[0.01:5:0.01] ffmV=js.ff.multilamellarVesicles p=js.grace() p.multi(1,2) # comparison single layer mV=ffmV(Q=Q, R=100., displace=0, dR=0,N=1,dN=0, phi=1,layers=6, SLD=1e-4) p[0].plot(mV) p[0].plot(js.ff.multiShellSphere(Q,[97,6],[0,1e-4]),li=[1,1,3],sy=0) # triple layer mV1=ffmV(Q=Q, R=100., displace=0, dR=0,N=1,dN=0, phi=1,layers=[1,4,1], SLD=[0.07e-3,0.6e-3,0.07e-3]) p[1].plot(mV1,sy=[1,0.5,2]) p[1].plot(js.ff.multiShellSphere(Q,[97,1,4,1],[0,0.07e-3,0.6e-3,0.07e-3]),li=[1,1,4],sy=0) p[1].yaxis(label='S(Q)',scale='l',min=1e-10,max=1e6,ticklabel=['power',0]) p[0].yaxis(label='S(Q)',scale='l',min=1e-10,max=1e6,ticklabel=['power',0]) p[1].xaxis(label='Q / nm\S-1',scale='l',min=1e-3,max=5,ticklabel=['power',0]) p[0].xaxis(label='Q / nm\S-1',scale='l',min=1e-3,max=5,ticklabel=['power',0])
References
[1] (1, 2, 3, 4, 5, 6) Small-angle scattering model for multilamellar vesicles H. Frielinghaus Physical Review E 76, 051603 (2007) [2] Small-Angle Scattering from Homogenous and Heterogeneous Lipid Bilayers N. Kučerka Advances in Planar Lipid Bilayers and Liposomes 12, 201-235 (2010) Examples
import jscatter as js import numpy as np ffmV=js.ff.multilamellarVesicles Q=js.loglist(0.01,5,500) dd=1.5 dR=5 nG=100 R=50 N=3 st=[1.5,3.5,1.5] p=js.grace(1,1) p.title('Lipid bilayer in SAXS/SANS') # SAXS sld=[0.07e-3,0.6e-3,0.07e-3] saxm=ffmV(Q=Q, R=R, displace=dd, dR=dR,N=N,dN=0, phi=0.2,layers=st, SLD=sld,solventSLD=0.94e-3,nGauss=nG) p.plot(saxm,sy=[1,0.3,1],le='SAXS multilamellar') saxu=ffmV(Q=Q, R=R, displace=0, dR=dR,N=1,dN=0, phi=0.2,layers=st, SLD=sld,solventSLD=0.94e-3,nGauss=100) p.plot(saxu,sy=0,li=[3,2,1],le='SAXS unilamellar') # SANS sld=[0.3e-4,1.5e-4,0.3e-4] sanm=ffmV(Q=Q, R=R, displace=dd, dR=dR,N=N,dN=0, phi=0.2,layers=st, SLD=sld,solventSLD=6.335e-4,nGauss=nG) p.plot( sanm,sy=[1,0.3,2],le='SANS multilamellar') sanu=ffmV(Q=Q, R=R, displace=0, dR=dR,N=1,dN=0, phi=0.2,layers=st, SLD=sld,solventSLD=6.335e-4,nGauss=100) p.plot(sanu,sy=0,li=[3,2,2],le='SANS unilamellar') # p.legend(x=0.015,y=1e-1) p.subtitle('R=50 nm, N=5, layerthickness=[1.5,3.5,1.5] nm, dR=5') p.yaxis(label='S(Q)',scale='l',min=1e-5,max=1e5,ticklabel=['power',0]) p.xaxis(label='Q / nm\S-1',scale='l',min=1e-2,max=5,ticklabel=['power',0])
-
jscatter.formfactor.
multilayer
(q, layerd, SLD, gaussw=None, solventSLD=0)[source]¶ Form factor of a multilayer with rectangular/Gaussian density profiles perpendicular to the layer.
To describe smeared interfaces or complex profiles use more layers.
Parameters: - q : array
Wavevectors in units 1/nm.
- layerd : list of float
Thickness of layers in units nm. List gives consecutive layer thickness from center to outside.
- gaussw : list of float, default None
Width of Gaussian in a layer. If float>0 a gaussian profile with this width is assumed. [1,None,1] is box in center and Gaussian width of 1 for outer layers.
- SLD : list of float
Scattering length density of layers in nm^-2.
- solventSLD : float, default=0
Solvent scattering length density in nm^-2.
Returns: - dataArray with [q,Fq]
- Fq is multilayer scattering per layer area.
- To get the scattering intensity of a volume the result needs to be multiplied with the layer area [2].
- .contrastprofile contrastprofile relative to layer center (contrast to solvent)
- .layerthickness
Notes
The scattering amplitude F_a is the Fourier transform of the density profile \rho(r)
F_a(q)=\int \rho(r)e^{iqr}dr
For a rectangular profile [1] of thickness d_i centered at a_i and layer scattering length density \rho_i we find
F_{a,box}(q)= \rho_i d_i sinc(qd_i/2)e^{iqa_i}
For a Gaussian profile [2] \rho(r)=\rho_ie^{-r^2/d_i^2/2} with same area as the rectangular profile we find
F_{a,gauss}(q)= \rho_i d_i e^{-(q^2d_i^2/2)}e^{iqa_i}
The scattering amplitude for a multi box/gauss profile is F_a(q)=\sum_i F_a(q,i)
The formfactor F(q) of this multi layer profile is in time average <>
F(q)=<\sum_{ij} F_a(q,i)F_a^*(q,j)>
resulting e.g. for a profile of rectangular boxes in
F_{box}(q)=\sum_{i,j} \rho_i\rho_j d_i d_j sinc(qd_i)sinc(qd_j)cos(q(a_i-a_j))
To get the 3D orientational average one has to add a Lorentz correction q^{⁻2} to describe the correct scattering in isotropic average (see [2]) or use multilamellarVesicles which includes a full structure factor and also size averaging.
- Restricting parameters for Fitting
- If the model is used during fits one has to consider dependencies between the parameters and restrict them somehow. For example the scattering length densities scale all with the particle concentration. This requires that maybe one is fixed. Symmetry in the layers may also be used to restrict the parameter space.
References
- Multi box profile
[1] Modelling X-ray or neutron scattering spectra of lyotropic lamellar phases : interplay between form and structure factors F. Nallet, R. Laversanne, D. Roux Journal de Physique II, EDP Sciences, 1993, 3 (4), pp.487-502 https://hal.archives-ouvertes.fr/jpa-00247849/document - Gaussian profile
[2] (1, 2) X-ray scattering from unilamellar lipid vesicles Brzustowicz and Brunger, J. Appl. Cryst. (2005). 38, 126–131 [3] Structural information from multilamellar liposomes at full hydration: Full q-range fitting with high quality X-ray data. Pabst, G., Rappolt, M., Amenitsch, H. & Laggner, P. Phys. Rev. E - Stat. Physics, Plasmas, Fluids, Relat. Interdiscip. Top. 62, 4000–4009 (2000).
Examples
Some symmetric membranes in comparison. Gaussian profiles correspond roughly to a box width of the Gaussian FWHM= 2\sqrt{2 log(2)}\sigma.
import jscatter as js import numpy as np q=np.r_[0.01:5:0.001] p=js.grace() p.multi(2,1) p[0].title('multilayer membranes') p[0].subtitle('fractions choosen to result in same forward scattering at q=0') profile=np.r_[2,1,2] # pf1=js.ff.multilayer(q,[1,5,1],profile) p[0].plot(pf1,sy=[1,0.3,1],le='box') p[1].plot(pf1.contrastprofile,li=[1,2,1],le='box') # f=2*(2*np.log(2))**0.5 # factor between sigma and FWHM pf2=js.ff.multilayer(q,np.r_[1,5,1],profile*f,gaussw=np.r_[1,5,1]/f) p[0].plot(pf2,sy=[2,0.3,2],le='gauss') p[1].plot(pf2.contrastprofile,li=[1,2,2],le='gauss') # pf3=js.ff.multilayer(q,np.r_[1,5,1],profile*np.r_[f,1,f],gaussw=np.r_[1/f,0,1/f]) p[0].plot(pf3,sy=[3,0.3,3],le='gauss-box-gauss') p[1].plot(pf3.contrastprofile,li=[1,2,3],le='gauss-box-gauss') # p[0].yaxis(scale='n',min=0.001,max=90,label='I(Q)',charsize=1)#,ticklabel=['power',0,1] p[0].xaxis(label='',charsize=1) p[1].yaxis(label='contrast profile ()') p[1].xaxis(label='Q / nm\S-1') p[0].legend(x=3,y=50)
How to use in a fit model. Due to the large number of possible models, smearing and so on one has to define what seems to be important. Also complex profiles are possible.
# define new model def membrane3(q,d1,d2,d3,s1,s2,s3,dd): # define your model #outer layers have half the scattering length density and thickness dd result=js.ff.multilayer(q,[dd,d1,d2+dd2,d3,dd],[s1/2,s1,s2,s3,s3/2])
-
jscatter.formfactor.
orientedCloudScattering
(qxzw, cloud, rms=0, coneangle=10, nCone=50, V=0, formfactor=None, ncpu=0)[source]¶ 2D scattering of an oriented cloud of scatterers with equal or variable scattering length. Using multiprocessing.
Cloud can represent an object described by a cloud of isotropic scatterers with orientation averaged in a cone. Scattering amplitudes may be constant, sphere scattering amplitude, Gaussian scattering amplitude or explicitly given form factor. Remember that the atomic bond length are on the order 0.1-0.2 nm and one expects Bragg peaks.
Parameters: - qxzw : array, ndim= Nx3
Wavevectors in 1/nm If 2D the 3rd dim is set to zero.
- cloud : array Nx3 or Nx4
- Center of mass positions (in nm) of the N scatterers in the cloud.
- If given cloud[3] is the scattering length b at positions cloud[:3], otherwise b=1.
- coneangle : float
Coneangle in units degrees.
- rms : float, default=0
Root mean square displacement =<u**2>**0.5 of the positions in cloud as random (Gaussian) displacements in nm. Displacement u is random for each orientation nCone. rms can be used to simulate a Debye-Waller factor. Larger nCone is advised to smooth data.
- nCone : int
Cone average as average over nCone Fibonacci lattice points in cone.
- V : float, default=0
Volume of the scatterers for formfactor ‘gauss’ and ‘sphere’.
- formfactor : ‘gauss’,’sphere’,array 2xN,default=None
- Gridpoint scattering amplitudes are described by:
- None : const scattering amplitude, point like particle.
- ‘sphere’: Sphere scattering amplitude according to [3].
- The sphere radius is R=(\frac{3V}{4\pi})^{1/3}
- ‘gauss’ : Gaussian function b_i(q)=b V exp(- \frac{V^{2/3.}}{4\pi}q^2) according to [2].
- explicit isotropic form factor f as array with [q,f] e.g. from multishell. The normalized scattering amplitude f_a for each gridpoint is calculated as f_a=f^{0.5}/f_a(0). Missing values are linear interpolated, q values outside interval are mapped to qmin or qmax.
- ncpu : int, default 0
- Number of cpus used in the pool for multiprocessing.
- not given or 0 : all cpus are used
- int>0 : min(ncpu, mp.cpu_count)
- int<0 : ncpu not to use
- 1 : single core usage for testing or comparing speed to Debye
Returns: - dataArray [qx, qz, qw, Sq]
- The forward scattering is Pq(q=0)= sumblength**2
- .sumblength : Sum of blength with sumblength**2
- .formfactoramplitude : formfactoramplitude of cloudpoints according to type for all q values.
- .formfactoramplitude_q :corresponding q values.
References
[1] - Kotlarchyk and S.-H. Chen, J. Chem. Phys. 79, 2461 (1983).1
[2] (1, 2) An improved method for calculating the contribution of solvent to the X-ray diffraction pattern of biological molecules Fraser R MacRae T Suzuki E IUCr Journal of Applied Crystallography 1978 vol: 11 (6) pp: 693-694 [3] (1, 2) X-ray diffuse scattering by proteins in solution. Consideration of solvent influence B. A. Fedorov, O. B. Ptitsyn and L. A. Voronin J. Appl. Cryst. (1974). 7, 181-186 doi: 10.1107/S0021889874009137 Examples
How to use orientedCloudScattering for fitting see last Example in cloudScattering.
Two points along y result in pattern independent of x but cos**2 for z with larger coneangle Ix becomes qx dependent.
import jscatter as js import numpy as np rod0=np.zeros([2,3]) rod0[:,1]=np.r_[0,np.pi] qxzw=np.mgrid[-6:6:50j, -6:6:50j].reshape(2,-1).T ffe=js.ff.orientedCloudScattering(qxzw,rod0,coneangle=5,nCone=10,rms=0) fig=js.mpl.surface(ffe.X,ffe.Z,ffe.Y) fig.axes[0].set_title(r'cos**2 for Z and slow decay for X due to 5 degree cone') fig.show() # noise in positions ffe=js.ff.orientedCloudScattering(qxzw,rod0,coneangle=5,nCone=100,rms=0.1) fig=js.mpl.surface(ffe.X,ffe.Z,ffe.Y) fig.axes[0].set_title('cos**2 for Y and slow decay for X with position noise') fig.show()
Two points along z result in symmetric pattern around zero symmetry reflects fibonacci lattice -> increase nCone
rod0=np.zeros([2,3]) rod0[:,2]=np.r_[0,np.pi] ffe=js.ff.orientedCloudScattering(qxzw,rod0,coneangle=45,nCone=10,rms=0.005) fig2=js.mpl.surface(ffe.X,ffe.Z,ffe.Y) fig2.axes[0].set_title('symmetric because of orientation along z; \n nCone needs to be larger for large cones') fig2.show()
5 spheres in line with small position distortion
rod0 = np.zeros([5, 3]) rod0[:, 1] = np.r_[0, 1, 2, 3, 4] * 3 qxzw = np.mgrid[-6:6:50j, -6:6:50j].reshape(2, -1).T ffe = js.ff.orientedCloudScattering(qxzw,rod0,formfactor='sphere',V=4/3.*np.pi*2**3,coneangle=20,nCone=30,rms=0.02) fig4 = js.mpl.surface(ffe.X, ffe.Z, np.log10(ffe.Y), colorMap='gnuplot') fig4.axes[0].set_title('5 spheres with R=2 along Z with noise (rms=0.02)') fig4.show()
5 core shell particles in line with small position distortion (Gaussian)
rod0 = np.zeros([5, 3]) rod0[:, 1] = np.r_[0, 1, 2, 3, 4] * 3 qxzw = np.mgrid[-6:6:50j, -6:6:50j].reshape(2, -1).T # only as demo : extract q from qxzw qxzw = np.c_[qxzw, np.zeros_like(qxzw[:, 0])] qrpt = js.formel.xyz2rphitheta(qxzw) q = np.unique(sorted(qrpt[:, 0])) # or use interpolation q = js.loglist(0.01, 7, 100) cs = js.ff.sphereCoreShell(q=q, Rc=1, Rs=2, bc=0.1, bs=1, solventSLD=0) ffe = js.ff.orientedCloudScattering(qxzw, rod0, formfactor=cs, coneangle=20, nCone=100, rms=0.05) fig4 = js.mpl.surface(ffe.X, ffe.Z, np.log10(ffe.Y), colorMap='gnuplot') fig4.axes[0].set_title('5 core shell particles with R=2 along Z with noise (rms=0.05)') fig4.show()
Make a slice for an angular region but with higher resolution to see the additional peaks due to alignment
# rod0 will be position of 5 points in a row rod0 = np.zeros([5, 3]) rod0[:, 1] = np.r_[0, 1, 2, 3, 4] * 3 qxzw = np.mgrid[-4:4:150j, -4:4:150j].reshape(2, -1).T # xz plane grid # only as demo : extract q from qxzw qxzw = np.c_[qxzw, np.zeros_like(qxzw[:, 0])] # add y=0 component qrpt = js.formel.xyz2rphitheta(qxzw) # spherical coordinates q = np.unique(sorted(qrpt[:, 0])) # or use interpolation; cs will be our formfactor q = js.loglist(0.01, 7, 100) cs = js.ff.sphereCoreShell(q=q, Rc=1, Rs=2, bc=0.1, bs=1, solventSLD=0) # calc scattering in plane qxzw ffe = js.ff.orientedCloudScattering(qxzw, rod0, formfactor=cs, coneangle=20, nCone=100, rms=0.05) # show it in surface plot fig4 = js.mpl.surface(ffe.X, ffe.Z, np.log10(ffe.Y), colorMap='gnuplot') fig4.axes[0].set_title('5 core shell particles with R=2 along Z with noise (rms=0.05)') fig4.show() # We do an explicit radial average # transform X,Z to spherical coordinates qphi=js.formel.xyz2rphitheta([ffe.X,ffe.Z,np.zeros_like(ffe.X)],transpose=True )[:,:2] # add qphi or use later rp[1] for selection ffb=ffe.addColumn(2,qphi.T) # select a portion of the phi angles phi=np.pi/2 dphi=0.2 ffn=ffb[:,(ffb[-1]<phi+dphi)&(ffb[-1]>phi-dphi)] ffn.isort(-2) # sort along radial q p=js.grace() p.plot(ffn[-2],ffn.Y,le='oriented spheres form factor') # compare to coreshell formfactor scaled p.plot(cs.X,cs.Y/cs.Y[0]*25,li=1,le='coreshell form factor') p.yaxis(label='F(Q,phi=90°+-11°)', scale='log') p.title('5 aligned core shell particle with additional interferences',size=1.) p.subtitle(' due to sphere alignment dependent on observation angle') # 2: direct way with 2D q in xz plane rod0 = np.zeros([5, 3]) rod0[:, 1] = np.r_[0, 1, 2, 3, 4] * 3 x=np.r_[0.0:6:0.05] qxzw = np.c_[x, x*0,x*0] for alpha in np.r_[0:91:30]: R=js.formel.rotationMatrix(np.r_[0,0,1],np.deg2rad(alpha)) # rotate around Z axis qa=np.dot(R,qxzw.T).T[:,:2] ffe = js.ff.orientedCloudScattering(qa, rod0, formfactor=cs, coneangle=20, nCone=100, rms=0.05) p.plot(x,ffe.Y,li=[1,2,-1],sy=0,le='alpha=%g' %alpha) p.xaxis(label=r'Q / nm\S-1') p.legend()
-
jscatter.formfactor.
pearlNecklace
(Q, Rc, l, N, A1=None, A2=None, A3=None, ms=None, mr=None)[source]¶ Formfactor of a pearl necklace (freely jointed chain of pearls connected by rods)
The formfactor is normalized that the pearls contribution equals 1.
Parameters: - Q : array
wavevector in nm
- Rc : float
pearl radius in nm
- N : float
number of pearls (homogeneous spheres)
- l : float
physical length of the rods
- A1, A2, A3 : float
Amplitudes of pearl-pearl, rod-rod and pearl-rod scattering. Can be calculated with the number of chemical monomers in a pearl ms and rod mr (see below for further information) If ms and mr are given A1,A2,A3 are calculated from these.
- ms : float, default None
number of chemical monomers in each pearl
- mr : float, default None
number of chemical monomers in rod like strings
Returns: - dataArray [q, Iq]
- .pearlRadius
- .A1
- .A2
- .A3
- .numberPearls
- .mr
- .ms
- .stringLength
Notes
- M : number of rod like strings (M=N-1)
- A1 = ms²/(M*mr+N*ms)²
- A2 = mr²/(M*mr+N*ms)²
- A3 = (mr*ms)/(M*mr+N*ms)²
References
[1] - Schweins, K. Huber, Macromol. Symp., 211, 25-42, 2004.
written by L. S. Fruhner, FZJ Juelich 2016
-
jscatter.formfactor.
ringPolymer
(q, Rg)[source]¶ General formfactor of a polymer ring in theta solvent.
Parameters: - q : array
Scattering vector, unit eg 1/A or 1/nm
- Rg : float
Radius of gyration, units in 1/unit(q)
Returns: - dataArray [q,Fq]
- .radiusOfGyration
References
[1] SANS from homogeneous polymer mixtures: A unified overview. Hammouda, B. in Polymer Characteristics 87–133 (Springer-Verlag, 1993). doi:10.1007/BFb0025862 Examples
import jscatter as js import numpy as np q=js.loglist(0.1,8,100) p=js.grace() p.multi(1,2) iq=js.ff.ringPolymer(q,5) p[0].plot(iq) p[0].yaxis(scale='l',label='I(q) / a.u.') p[0].xaxis(scale='l') p[1].plot(iq.X,iq.Y*iq.X**2) p[1].yaxis(scale='l',label=['I(q)q\S2\N / a.u.',1,'opposite'],ticklabel=['power',0,1,'opposite']) p[1].xaxis(scale='l',label='q / nm\S-1') p[1].legend(x=0.2,y=0.5)
-
jscatter.formfactor.
scatteringFromSizeDistribution
(q, sizedistribution, size=None, func=<function beaucage>, weight=None, **kwargs)[source]¶ Scattering of objects with one multimodal parameter as e.g. multimodal size distribution.
Distributions might be mixtures of small and large particles bi or multimodal. For predefined distributions see formel.parDistributedAverage with examples. The weighted average over given sizedistribution is calculated.
Parameters: - q : array of float;
Wavevectors to calculate scattering; unit = 1/unit(size distribution)
- sizedistribution : dataArray or array
Explicit given distribution of sizes as [ [list size],[list probability]]
- size : string
Name of the parameter describing the size (may be also something different than size).
- func : lambda or function
Function that describes the form factor with first arguments (q,size,…) and should return dataArray with .Y as result.
- kwargs :
Any additional keyword arguments passed to for func.
- weight : function
Weight function dependent on size. E.g. weight = lambda R:rho**2 * (4/3*np.pi*R**3)**2 with V= 4pi/3 R**3 for normalized form factors to account for forward scattering of volume objects of dimension 3.
Returns: - dataArray [q,I(q)]
Notes
We have to discriminate between formfactor normalized to 1 (e.g. beaucage) and form factors returning the absolute scattering (e.g. sphere) including the contrast. The later contains already \rho^2 V^2, the first not.
We need for normalized formfactors P(q) I(q) = n \rho^2 V^2 P(q) with n as number density \rho as difference in average scattering length (contrast), V as volume of particle (~r³ ~ mass) and use weight = \rho^2 V(R)^2
I(q)= \sum_{R_i} [ weight(R_i) * probability(R_i) * P(q, R_i , *kwargs).Y ]
For a gaussian chain with R_g^2=l^2 N^{2\nu} and monomer number N (nearly 2D object) we find N^2=(R_g/l)^{1/\nu} and the forward scattering as weight I_0=b^2 N^2=b^2 (R_g/l)^{1/\nu}
Examples
The contribution of different simple sizes to Beaucage
import jscatter as js q=js.loglist(0.01,6,100) p=js.grace() # bimodal with equal concentration bimodal=[[12,70],[1,1]] Iq=js.ff.scatteringFromSizeDistribution(q=q,sizedistribution=bimodal,d=3,weight=lambda r:(r/12)**6) p.plot(Iq,legend='with aggregates') Iq=js.ff.scatteringFromSizeDistribution(q=q,sizedistribution=bimodal,d=3) p.plot(Iq,legend='with aggregates') # 2:1 concentration bimodal=[[12,70],[1,2]] p.plot(js.ff.scatteringFromSizeDistribution(q=q,sizedistribution=bimodal,d=2.5),legend='no aggregates') p.yaxis(scale='l') p.xaxis(scale='l')
Three sphere sizes:
import jscatter as js q=js.loglist(0.001,6,1000) p=js.grace() # trimodal with equal concentration trimodal=[[10,50,500],[1,0.01,0.00001]] Iq=js.ff.scatteringFromSizeDistribution(q=q,sizedistribution=trimodal,size='radius',func=js.ff.sphere) p.plot(Iq,legend='with aggregates') p.yaxis(label='I(q)',scale='l',max=1e13,min=1) p.xaxis(scale='l',label='q / nm\S-1') p.text(r'minimum \nlargest',x=0.002,y=1e10) p.text(r'minimum \nmiddle',x=0.02,y=1e7) p.text(r'minimum \nsmallest',x=0.1,y=1e5) p.title('trimodal spheres') p.subtitle('first minima indicated')
-
jscatter.formfactor.
sphere
(q, radius, contrast=1)[source]¶ Scattering of a single homogeneous sphere.
Parameters: - q : float
Wavevector in units of 1/nm
- radius : float
Radius in units nm
- contrast : float, default=1
Difference in scattering length to the solvent = contrast
Returns: - dataArray [q,Iq]
- Iq scattering intensity
- .I0 forward scattering
Notes
I(q)= 4\pi\rho^2V^2\left[\frac{3(sin(qR) - qr cos(qR))}{(qR)^3}\right]^2
with contrast \rho and sphere volume V=\frac{4\pi}{3}R^3
The first minimum of the form factor is at qR=4.493
References
[1] Guinier, A. and G. Fournet, “Small-Angle Scattering of X-Rays”, John Wiley and Sons, New York, (1955).
-
jscatter.formfactor.
sphereCoreShell
(q, Rc, Rs, bc, bs, solventSLD=0)[source]¶ Scattering of a spherical core shell particle.
Parameters: - q : float
Wavevector in units of 1/(R units)
- Rc,Rs : float
Radius core and radius of shell Rs>Rc
- bc,bs : float
Contrast to solvent scattering length density of core and shell.
- solventSLD : float, default =0
Scattering length density of the surrounding solvent. If equal to zero (default) then in profile the contrast is given.
Returns: - dataArray [wavevector ,Iq ]
Notes
Calls multiShellSphere.
-
jscatter.formfactor.
sphereCoreShellGaussianCorona
(q, Rc, Rs, Rg, Ncoil, thicknessCoils, coilSLD=6.4e-05, coreSLD=0.0004186, shellSLD=0.0002, solventSLD=0.0006335, d=1)[source]¶ Scattering of a core-shell particle surrounded by gaussian coils as model for grafted polymers on particle.
The model is in analogy to the sphereGaussianCorona replacing the sphere by a core shell particle in [1]. The additional scattering from the coils is uniformly distributed at the surface, which might fail for lower aggregation numbers as 1, 2, 3. Instead of aggregation number in [1] we use volume of the gaussian coils collapsed to the surface.
Parameters: - q: array of float
Wavevectors in unit 1/nm.
- Rc,Rs : float
Radius of core and shell in unit nm.
- Rg : float
Radius of gyration of coils in unit nm.
- d : float, default 1
Coils centre located d*Rg away from the sphere surface This might be equivalent to Rg
- Ncoil : float
Number of coils at the surface (aggregation number)
- thicknessCoils : float
Thickness of a layer if all coils collapsed on the surface as additional shell in nm. Needed to calculate absolute scattering of the expanded coils. The densely packed coil shell volume is V_{coils}= 4/3\pi((R_{s}+thicknessCoils)^3-R_s^3) and the volume of a single polymer V_m =V_{coils} / Ncoils.
- coilSLD : float
Scattering length density of coil in bulk as if collapsed on surface unit nm^-2. default hPEG = 0.64*1e-6 A^-2 = 0.64*1e-4 nm^-2
- coreSLD,shellSLD : float
Scattering length density of core and shell in unit nm^-2. default :
- core SiO2 = 4.186*1e-4 nm^-2
- protein coating = 2*1e-4 nm^-2
- solventSLD : float
Scattering length density of solvent. unit nm^-2. default D2O = 6.335*1e-6 A^-2 = 6.335*1e-4 nm^-2
Returns: - dataArray [q,Iq]
- .coilRg
- .Radii
- .numberOfCoils
- .coildistancefactor
- .coilequVolume
- .coilSLD
- .coreshellSLD
- .solventSLD
Notes
- The defaults result in a silica sphere with protein coating and with hPEG grafted at the surface in D2O.
- Rg=N**0.5*b with N monomers of length b
- Vcoilsphere=N*monomerVolume=4/3.*np.pi*coilequR**3
- coilequR=(N*monomerVolume/(4/3.*np.pi))**(1/3.)
References
[1] (1, 2, 3) Form factors of block copolymer micelles with spherical, ellipsoidal and cylindrical cores Pedersen J Journal of Applied Crystallography 2000 vol: 33 (3) pp: 637-640 [2] Hammouda, B. (1992).J. Polymer Science B: Polymer Physics30 , 1387–1390 Examples
import jscatter as js q=js.loglist(0.01,5,300) p=js.grace() p.plot(js.ff.sphereCoreShellGaussianCorona(q,8,12,6,20,1.5),sy=0,li=1) p.yaxis(scale='l') p.xaxis(scale='l')
-
jscatter.formfactor.
sphereFuzzySurface
(q, R, sigmasurf, contrast)[source]¶ Scattering of a sphere with a fuzzy interface.
Parameters: - q : float
Wavevector in units of 1/(R units)
- R : float
The particle radius R represents the radius of the particle where the scattering length density profile decreased to 1/2 of the core density.
- sigmasurf : float
Sigmasurf is the width of the smeared particle surface.
- contrast : float
Difference in scattering length to the solvent = contrast
Returns: - dataArray [q, Iq]
- Iq scattering intensity related to sphere volume.
- .I0 forward scattering
Notes
The “fuzziness” of the interface is defined by the parameter sigmasurf. The particle radius R represents the radius of the particle where the scattering length density profile decreased to 1/2 of the core density. sigmasurf is the width of the smeared particle surface. The inner regions of the microgel that display a higher density are described by the radial box profile extending to a radius of approximately Rbox ~ R - 2(sigma). In dilute solution, the profile approaches zero as Rsans ~ R + 2(sigma).
References
[1] - Stieger, J. S. Pedersen, P. Lindner, W. Richtering, Langmuir 20 (2004) 7283-7292
-
jscatter.formfactor.
sphereGaussianCorona
(q, R, Rg, Ncoil, coilequR, coilSLD=6.4e-05, sphereSLD=0.0004186, solventSLD=0.0006335, d=1)[source]¶ Scattering of a sphere surrounded by gaussian coils as model for grafted polymers on particle e.g. a micelle.
The additional scattering is uniformly distributed at the surface, which might fail for lower aggregation numbers as 1, 2, 3. Instead of aggregation number in [1] we use sphere volume and a equivalent volume of the gaussian coils.
Parameters: - q: array of float
Wavevectors in unit 1/nm
- R : float
Sphere radius in unit nm
- Rg : float
Radius of gyration of coils in unit nm
- d : float, default 1
Coils centre located d*Rg away from the sphere surface
- Ncoil : float
Number of coils at the surface (aggregation number)
- coilequR : float
Equivalent radius to calc volume of one coil if densely packed as a sphere. Needed to calculate absolute scattering of the coil.
- coilSLD : float
Scattering length density of coil in bulk. unit nm^-2. default hPEG = 0.64*1e-6 A^-2 = 0.64*1e-4 nm^-2
- sphereSLD : float
Scattering length density of sphere.unit nm^-2. default SiO2 = 4.186*1e-6 A^-2 = 4.186*1e-4 nm^-2
- solventSLD : float
Scattering length density of solvent. unit nm^-2. default D2O = 6.335*1e-6 A^-2 = 6.335*1e-4 nm^-2
Returns: - dataArray [q,Iq]
- .coilRg
- .sphereRadius
- .numberOfCoils
- .coildistancefactor
- .coilequVolume
- .coilSLD
- .sphereSLD
- .solventSLD
Notes
- The defaults result in a silica sphere with hPEG grafted at the surface in D2O.
- Rg=N**0.5*b with N monomers of length b
- Vcoilsphere=N*monomerVolume=4/3.*np.pi*coilequR**3
- coilequR=(N*monomerVolume/(4/3.*np.pi))**(1/3.)
References
[1] (1, 2) Form factors of block copolymer micelles with spherical, ellipsoidal and cylindrical cores Pedersen J Journal of Applied Crystallography 2000 vol: 33 (3) pp: 637-640 [2] Hammouda, B. (1992).J. Polymer Science B: Polymer Physics30 , 1387–1390 Examples
import jscatter as js q=js.loglist(0.1,5,100) p=js.grace() p.plot(js.ff.sphereGaussianCorona(q,4.4,2,30,2))
-
jscatter.formfactor.
superball
(q, R, p, SLD=1, solventSLD=0, nGrid=15, returngrid=False)[source]¶ A superball is a general mathematical shape that can be used to describe rounded cubes, sphere and octahedron’s.
The numerical integration is done by a pseudorandom grid of scatterers. The integration is valid at low Q. Validity can be increased if nGrid is increased.
Parameters: - q : array
Wavevector in 1/nm
- R : float, None
2R = edge length
- p : float, 0<p<100
Parameter that describes shape | p=0 empty space | p<0.5 concave octahedron’s | p=0.5 octahedron | 0.5<p<1 convex octahedron’s | p=1 spheres | p>1 rounded cubes | p->inf cubes
- SLD : float, default =1
Scattering length density of cuboid.unit nm^-2
- solventSLD : float, default =0
Scattering length density of solvent. unit nm^-2
- nGrid : int
Number of gridpoints is nGrid**3. relError=nGrid*4 is used for Fibonacci lattice with 2*relError+1 orientations in spherical average.
- returngrid : bool
Return grid.
Returns: - dataArray [q,Iq, beta]
References
[1] Periodic lattices of arbitrary nano-objects: modeling and applications for self-assembled systems Yager, K.G.; Zhang, Y.; Lu, F.; Gang, O. Journal of Applied Crystallography 2014, 47, 118–129. doi: 10.1107/S160057671302832X [2] http://gisaxs.com/index.php/Form_Factor:Superball Examples
Compare to extreme cases of sphere (p=1) and cube (p->inf , use here 100)
import jscatter as js import numpy as np # q=np.r_[0:3.5:0.02] R=6; p=js.grace() p.multi(2,1) p[0].yaxis(scale='l') ss=js.ff.superball(q,R,p=1,returngrid=1) p[0].plot(ss,legend='superball p=1') ss=js.ff.superball(q,R,p=1,returngrid=1,nGrid=20) p[0].plot(ss,legend='superball p=1 nGrid=20') p[0].plot(js.ff.sphere(q,R),li=1,sy=0,legend='sphere ff') p[0].legend(x=2,y=2e5) # p[1].yaxis(scale='l') cc=js.ff.superball(q,R,p=100,nGrid=15) p[1].plot(cc,sy=[1,0.3,4],legend='superball p=100 nGrid=15') cc=js.ff.superball(q,R,p=100,nGrid=20) p[1].plot(cc,sy=[1,0.3,5],legend='superball p=100 nGrid=20') p[1].plot(js.ff.cuboid(q,2*R),li=4,sy=0,legend='cuboid') p[1].legend(x=2,y=2e5) # # visualisation import matplotlib.pyplot as plt from mpl_toolkits.mplot3d import Axes3D # cubic grid points fig = plt.figure() ax = fig.add_subplot(111, projection='3d') q=np.r_[0:5:0.1] R=3 pxyz=js.ff.superball(q,R,p=200,nGrid=10,returngrid=True).grid pxyz=pxyz[pxyz[:,0]>0] ax.scatter(pxyz[:,0],pxyz[:,1],pxyz[:,2],color="k",s=20) ax.set_xlim([-3,3]) ax.set_ylim([-3,3]) ax.set_zlim([-3,3]) ax.set_aspect("equal") plt.tight_layout() plt.show(block=False)
-
jscatter.formfactor.
teubnerStrey
(q, xi, d, eta2=1)[source]¶ Phenomenological model for the scattering intensity of a two-component system using the Teubner-Strey model [1].
Often used for bi-continuous micro-emulsions.
Parameters: - q : array
wavevectors
- xi : float
correlation length
- d : float
characteristic domain size, periodicity
- eta2 : float, default=1
squared scattering length density contrast
Returns: - dataArray [q, Iq]
Notes
- q_{max}=((2\pi/d)^2-\xi^{-2})^{1/2}
References
[1] (1, 2, 3) M. Teubner and R. Strey, Origin of the scattering peak in microemulsions, J. Chem. Phys., 87:3195, 1987 [2] K. V. Schubert, R. Strey, S. R. Kline, and E. W. Kaler, Small angle neutron scattering near lifshitz lines: Transition from weakly structured mixtures to microemulsions, J. Chem. Phys., 101:5343, 1994 Examples
Fit Teubner-Strey with background and a power law for low Q
#import jscatter as js #import numpy as np def tbpower(q,B,xi,dd,A,beta,bgr): # Model Teubner Strey + power law and background tb=js.ff.teubnerStrey(q=q,xi=xi,d=dd) # add power law and background tb.Y=B*tb.Y+A*q**beta+bgr tb.A=A tb.bgr=bgr tb.beta=beta return tb # simulate some data q=js.loglist(0.01,5,600) data=tbpower(q,1,10,20,0.002,-3,0.1) # or read them # data=js.dA('filename.chi') # plot data p=js.grace() p.plot(data,legend='simulated data') p.xaxis(scale='l',label=r'Q / nm\S-1') p.yaxis(scale='l',label='I(Q) / a.u.') p.title('TeubnerStrey model with power and background')
-
jscatter.formfactor.
wormlikeChain
(q, N, a, R=None, SLD=1, solventSLD=0, rtol=0.02)[source]¶ Scattering of a wormlike chain, which correctly reproduces the rigid-rod and random-coil limits.
The forward scattering is :math:´I0=V^2(SLD-solventSLD)^2´ volume :math:´V=piR^2N´.
Parameters: - q : array
wavevectors in 1/nm
- N : float
Chain length, units of 1/q
- a : float
Persistence length with l=2a l=Kuhn length (segment length), units of nm.
- R : float
Radius in units of nm.
- SLD : float
Scattering length density segments.
- solventSLD :
Solvent scattering length density.
- rtol : float
Maximum relative tolerance in integration.
Returns: - dataArray [q, Iq]
- .chainRadius
- .chainLength
- .persistenceLength
- .Rg
- .volume
- .contrast
Notes
- From [1] :
- The Kratky plot (Figure 4 ) is not the most convenient way to determine a as was pointed out in ref 20. Figure 5 provides an alternative way of measuring a by plotting the experimentally measurable combination Nk2S(k) versus a for fixed wavelength k. As Figure 5 indicates, this plot is rather insensitive to the chain length N and therefore is universal. The numerical analysis of eq 17 shows that this remains true for as long as k is not too small. Taking into account that the excluded-volume effects leave S(k) practically unchanged (e.g., see Figures 2 and 4 of ref 231, the plot of Figure 5 can serve as a useful alternative to the Kratky plot which, in addition, does not suffer from the polydispersity effects
References
[1] (1, 2) Analytical calculation of the scattering function for polymers of arbitrary flexibility using the dirac propagator A. L. Kholodenko, Macromolecules, 26:4179–4183, 1993 [2] (1, 2) The structure factor of a wormlike chain and the random-phase-approximation solution for the spinodal line of a diblock copolymer melt Zhang X et. al. Soft Matter 10, 5405 (2014) [3] (1, 2) Models of Polymer Chains Teraoka I. in Polymer Solutions: An Introduction to Physical Properties pp: 1-67, New York, John Wiley & Sons, Inc. Examples
import jscatter as js p=js.grace() p.multi(2,1) p.title('figure 3 (2 scaled) of ref Kholodenko Macromolecules 26, 4179 (1993)',size=1) q=js.loglist(0.01,10,100) for a in [1,2.5,5,20,50,1000]: ff=js.ff.wormlikeChain(q,200,a) p[0].plot(ff.X,200*ff.Y*ff.X**2,legend='a=%.4g' %ff.persistenceLength) p[1].plot(ff.X,ff.Y,legend='a=%.4g' %ff.persistenceLength) p[0].legend() p[0].yaxis(label='Nk\S2\NS(k)') p[1].xaxis(label='k',scale='l') p[1].yaxis(label='S(k)',scale='l') # p=js.grace() p.multi(2,1) p.title('figure 4 of ref Kholodenko Macromolecules 26, 4179 (1993)',size=1) # fig 4 seems to be wrong scale in [Re57b872e77e7-1]_ as for large N with a=1 fig 2 and 4 should have same plateau. a=1 q=js.loglist(0.01,4./a,100) for NN in [1,20,50,150,500]: ff=js.ff.wormlikeChain(q,NN,a) p[0].plot(ff.X*a,NN*a*ff.Y*ff.X**2,legend='N=%.4g' %ff.chainLength) p[1].plot(ff.X,ff.Y,legend='a=%.4g' %ff.persistenceLength) p[0].legend() p[0].yaxis(label='(N/a)(ka)\S2\NS(k)') p[0].xaxis(label='ka') p[1].xaxis(label='k',scale='l') p[1].yaxis(label='S(k)',scale='l')