source: sasmodels/explore/realspace.py @ e309e23

core_shell_microgelsmagnetic_modelticket-1257-vesicle-productticket_1156ticket_1265_superballticket_822_more_unit_tests
Last change on this file since e309e23 was 032646d, checked in by Paul Kienzle <pkienzle@…>, 7 years ago

explore/realspace.py: separate shape from lattice spec

  • Property mode set to 100644
File size: 26.6 KB
RevLine 
[cfa28d3]1from __future__ import division, print_function
2
3import time
4from copy import copy
[7d0afa3]5import os
6import argparse
7from collections import OrderedDict
[cfa28d3]8
9import numpy as np
10from numpy import pi, radians, sin, cos, sqrt
[032646d]11from numpy.random import poisson, uniform, randn, rand
[a1c32c2]12from numpy.polynomial.legendre import leggauss
[cfa28d3]13from scipy.integrate import simps
14from scipy.special import j1 as J1
15
[7d0afa3]16try:
17    import numba
18    USE_NUMBA = True
19except ImportError:
20    USE_NUMBA = False
21
[cfa28d3]22# Definition of rotation matrices comes from wikipedia:
23#    https://en.wikipedia.org/wiki/Rotation_matrix#Basic_rotations
24def Rx(angle):
25    """Construct a matrix to rotate points about *x* by *angle* degrees."""
26    a = radians(angle)
27    R = [[1, 0, 0],
28         [0, +cos(a), -sin(a)],
29         [0, +sin(a), +cos(a)]]
30    return np.matrix(R)
31
32def Ry(angle):
33    """Construct a matrix to rotate points about *y* by *angle* degrees."""
34    a = radians(angle)
35    R = [[+cos(a), 0, +sin(a)],
36         [0, 1, 0],
37         [-sin(a), 0, +cos(a)]]
38    return np.matrix(R)
39
40def Rz(angle):
41    """Construct a matrix to rotate points about *z* by *angle* degrees."""
42    a = radians(angle)
43    R = [[+cos(a), -sin(a), 0],
44         [+sin(a), +cos(a), 0],
45         [0, 0, 1]]
46    return np.matrix(R)
47
48def rotation(theta, phi, psi):
49    """
50    Apply the jitter transform to a set of points.
51
52    Points are stored in a 3 x n numpy matrix, not a numpy array or tuple.
53    """
54    return Rx(phi)*Ry(theta)*Rz(psi)
55
[226473d]56def apply_view(points, view):
57    """
58    Apply the view transform (theta, phi, psi) to a set of points.
59
60    Points are stored in a 3 x n numpy array.
61
62    View angles are in degrees.
63    """
64    theta, phi, psi = view
65    return np.asarray((Rz(phi)*Ry(theta)*Rz(psi))*np.matrix(points.T)).T
66
67
68def invert_view(qx, qy, view):
69    """
70    Return (qa, qb, qc) for the (theta, phi, psi) view angle at detector
71    pixel (qx, qy).
72
73    View angles are in degrees.
74    """
75    theta, phi, psi = view
76    q = np.vstack((qx.flatten(), qy.flatten(), 0*qx.flatten()))
77    return np.asarray((Rz(-psi)*Ry(-theta)*Rz(-phi))*np.matrix(q))
78
79
[cfa28d3]80class Shape:
81    rotation = np.matrix([[1., 0, 0], [0, 1, 0], [0, 0, 1]])
82    center = np.array([0., 0., 0.])[:, None]
83    r_max = None
84
85    def volume(self):
86        # type: () -> float
87        raise NotImplementedError()
88
89    def sample(self, density):
90        # type: (float) -> np.ndarray[N], np.ndarray[N, 3]
91        raise NotImplementedError()
92
[032646d]93    def dims(self):
94        # type: () -> float, float, float
95        raise NotImplementedError()
96
[cfa28d3]97    def rotate(self, theta, phi, psi):
98        self.rotation = rotation(theta, phi, psi) * self.rotation
99        return self
100
101    def shift(self, x, y, z):
102        self.center = self.center + np.array([x, y, z])[:, None]
103        return self
104
105    def _adjust(self, points):
106        points = np.asarray(self.rotation * np.matrix(points.T)) + self.center
107        return points.T
108
109    def r_bins(self, q, over_sampling=1, r_step=0.):
110        r_max = min(2 * pi / q[0], self.r_max)
111        if r_step == 0.:
112            r_step = 2 * pi / q[-1] / over_sampling
113        #r_step = 0.01
114        return np.arange(r_step, r_max, r_step)
115
116class Composite(Shape):
117    def __init__(self, shapes, center=(0, 0, 0), orientation=(0, 0, 0)):
118        self.shapes = shapes
119        self.rotate(*orientation)
120        self.shift(*center)
121
122        # Find the worst case distance between any two points amongst a set
123        # of shapes independent of orientation.  This could easily be a
124        # factor of two worse than necessary, e.g., a pair of thin rods
125        # end-to-end vs the same pair side-by-side.
126        distances = [((s1.r_max + s2.r_max)/2
127                      + sqrt(np.sum((s1.center - s2.center)**2)))
128                     for s1 in shapes
129                     for s2 in shapes]
130        self.r_max = max(distances + [s.r_max for s in shapes])
[032646d]131        self.volume = sum(shape.volume for shape in self.shapes)
[cfa28d3]132
133    def sample(self, density):
134        values, points = zip(*(shape.sample(density) for shape in self.shapes))
135        return np.hstack(values), self._adjust(np.vstack(points))
136
137class Box(Shape):
138    def __init__(self, a, b, c,
139                 value, center=(0, 0, 0), orientation=(0, 0, 0)):
140        self.value = np.asarray(value)
141        self.rotate(*orientation)
142        self.shift(*center)
143        self.a, self.b, self.c = a, b, c
144        self._scale = np.array([a/2, b/2, c/2])[None, :]
145        self.r_max = sqrt(a**2 + b**2 + c**2)
[032646d]146        self.dims = a, b, c
147        self.volume = a*b*c
[cfa28d3]148
149    def sample(self, density):
[032646d]150        num_points = poisson(density*self.volume)
[cfa28d3]151        points = self._scale*uniform(-1, 1, size=(num_points, 3))
152        values = self.value.repeat(points.shape[0])
153        return values, self._adjust(points)
154
155class EllipticalCylinder(Shape):
156    def __init__(self, ra, rb, length,
157                 value, center=(0, 0, 0), orientation=(0, 0, 0)):
158        self.value = np.asarray(value)
159        self.rotate(*orientation)
160        self.shift(*center)
161        self.ra, self.rb, self.length = ra, rb, length
162        self._scale = np.array([ra, rb, length/2])[None, :]
163        self.r_max = sqrt(4*max(ra, rb)**2 + length**2)
[032646d]164        self.dims = 2*ra, 2*rb, length
165        self.volume = pi*ra*rb*length
[cfa28d3]166
167    def sample(self, density):
[032646d]168        # randomly sample from a box of side length 2*r, excluding anything
169        # not in the cylinder
[cfa28d3]170        num_points = poisson(density*4*self.ra*self.rb*self.length)
171        points = uniform(-1, 1, size=(num_points, 3))
172        radius = points[:, 0]**2 + points[:, 1]**2
173        points = self._scale*points[radius <= 1]
174        values = self.value.repeat(points.shape[0])
175        return values, self._adjust(points)
176
177class TriaxialEllipsoid(Shape):
178    def __init__(self, ra, rb, rc,
179                 value, center=(0, 0, 0), orientation=(0, 0, 0)):
180        self.value = np.asarray(value)
181        self.rotate(*orientation)
182        self.shift(*center)
183        self.ra, self.rb, self.rc = ra, rb, rc
184        self._scale = np.array([ra, rb, rc])[None, :]
185        self.r_max = 2*max(ra, rb, rc)
[032646d]186        self.dims = 2*ra, 2*rb, 2*rc
187        self.volume = 4*pi/3 * ra * rb * rc
[cfa28d3]188
189    def sample(self, density):
190        # randomly sample from a box of side length 2*r, excluding anything
191        # not in the ellipsoid
192        num_points = poisson(density*8*self.ra*self.rb*self.rc)
193        points = uniform(-1, 1, size=(num_points, 3))
194        radius = np.sum(points**2, axis=1)
195        points = self._scale*points[radius <= 1]
196        values = self.value.repeat(points.shape[0])
197        return values, self._adjust(points)
198
199class Helix(Shape):
200    def __init__(self, helix_radius, helix_pitch, tube_radius, tube_length,
201                 value, center=(0, 0, 0), orientation=(0, 0, 0)):
202        self.value = np.asarray(value)
203        self.rotate(*orientation)
204        self.shift(*center)
[032646d]205        helix_length = helix_pitch * tube_length/sqrt(helix_radius**2 + helix_pitch**2)
206        total_radius = self.helix_radius + self.tube_radius
[cfa28d3]207        self.helix_radius, self.helix_pitch = helix_radius, helix_pitch
208        self.tube_radius, self.tube_length = tube_radius, tube_length
[032646d]209        self.r_max = sqrt(4*total_radius + (helix_length + 2*tube_radius)**2)
210        self.dims = 2*total_radius, 2*total_radius, helix_length
[cfa28d3]211        # small tube radius approximation; for larger tubes need to account
212        # for the fact that the inner length is much shorter than the outer
213        # length
[032646d]214        self.volume = pi*self.tube_radius**2*self.tube_length
[cfa28d3]215
216    def points(self, density):
217        num_points = poisson(density*4*self.tube_radius**2*self.tube_length)
218        points = uniform(-1, 1, size=(num_points, 3))
219        radius = points[:, 0]**2 + points[:, 1]**2
220        points = points[radius <= 1]
221
222        # Based on math stackexchange answer by Jyrki Lahtonen
223        #     https://math.stackexchange.com/a/461637
224        # with helix along z rather than x [so tuples in answer are (z, x, y)]
225        # and with random points in the cross section (p1, p2) rather than
226        # uniform points on the surface (cos u, sin u).
227        a, R = self.tube_radius, self.helix_radius
228        h = self.helix_pitch
229        scale = 1/sqrt(R**2 + h**2)
230        t = points[:, 3] * (self.tube_length * scale/2)
231        cos_t, sin_t = cos(t), sin(t)
232
233        # rx = R*cos_t
234        # ry = R*sin_t
235        # rz = h*t
236        # nx = -a * cos_t * points[:, 1]
237        # ny = -a * sin_t * points[:, 1]
238        # nz = 0
239        # bx = (a * h/scale) * sin_t * points[:, 2]
240        # by = (-a * h/scale) * cos_t * points[:, 2]
241        # bz = a*R/scale
242        # x = rx + nx + bx
243        # y = ry + ny + by
244        # z = rz + nz + bz
245        u, v = (R - a*points[:, 1]), (a * h/scale)*points[:, 2]
246        x = u * cos_t + v * sin_t
247        y = u * sin_t - v * cos_t
248        z = a*R/scale + h * t
249
250        points = np.hstack((x, y, z))
251        values = self.value.repeat(points.shape[0])
252        return values, self._adjust(points)
253
[7d0afa3]254def csbox(a=10, b=20, c=30, da=1, db=2, dc=3, slda=1, sldb=2, sldc=3, sld_core=4):
255    core = Box(a, b, c, sld_core)
256    side_a = Box(da, b, c, slda, center=((a+da)/2, 0, 0))
257    side_b = Box(a, db, c, sldb, center=(0, (b+db)/2, 0))
258    side_c = Box(a, b, dc, sldc, center=(0, 0, (c+dc)/2))
259    side_a2 = copy(side_a).shift(-a-da, 0, 0)
260    side_b2 = copy(side_b).shift(0, -b-db, 0)
261    side_c2 = copy(side_c).shift(0, 0, -c-dc)
262    shape = Composite((core, side_a, side_b, side_c, side_a2, side_b2, side_c2))
[032646d]263    shape.dims = 2*da+a, 2*db+b, 2*dc+c
[7d0afa3]264    return shape
265
266def _Iqxy(values, x, y, z, qa, qb, qc):
267    """I(q) = |sum V(r) rho(r) e^(1j q.r)|^2 / sum V(r)"""
268    Iq = [abs(np.sum(values*np.exp(1j*(qa_k*x + qb_k*y + qc_k*z))))**2
269            for qa_k, qb_k, qc_k in zip(qa.flat, qb.flat, qc.flat)]
270    return Iq
271
272if USE_NUMBA:
273    # Override simple numpy solution with numba if available
[8fb2a94]274    from numba import njit
275    @njit("f8[:](f8[:],f8[:],f8[:],f8[:],f8[:],f8[:],f8[:])")
[8bd379a]276    def _Iqxy(values, x, y, z, qa, qb, qc):
[8fb2a94]277        Iq = np.zeros_like(qa)
278        for j in range(len(Iq)):
279            total = 0. + 0j
[8bd379a]280            for k in range(len(values)):
[8fb2a94]281                total += values[k]*np.exp(1j*(qa[j]*x[k] + qb[j]*y[k] + qc[j]*z[k]))
282            Iq[j] = abs(total)**2
283        return Iq
284
285def calc_Iqxy(qx, qy, rho, points, volume=1.0, view=(0, 0, 0)):
[226473d]286    qx, qy = np.broadcast_arrays(qx, qy)
287    qa, qb, qc = invert_view(qx, qy, view)
288    rho, volume = np.broadcast_arrays(rho, volume)
289    values = rho*volume
[8fb2a94]290    x, y, z = points.T
[8bd379a]291    values, x, y, z, qa, qb, qc = [np.asarray(v, 'd')
292                                   for v in (values, x, y, z, qa, qb, qc)]
[226473d]293
294    # I(q) = |sum V(r) rho(r) e^(1j q.r)|^2 / sum V(r)
[7d0afa3]295    Iq = _Iqxy(values, x, y, z, qa.flatten(), qb.flatten(), qc.flatten())
[8fb2a94]296    return np.asarray(Iq).reshape(qx.shape) / np.sum(volume)
[226473d]297
[cfa28d3]298def _calc_Pr_nonuniform(r, rho, points):
299    # Make Pr a little be bigger than necessary so that only distances
300    # min < d < max end up in Pr
301    n_max = len(r)+1
302    extended_Pr = np.zeros(n_max+1, 'd')
303    # r refers to bin centers; find corresponding bin edges
304    bins = bin_edges(r)
305    t_next = time.clock() + 3
306    for k, rho_k in enumerate(rho[:-1]):
307        distance = sqrt(np.sum((points[k] - points[k+1:])**2, axis=1))
308        weights = rho_k * rho[k+1:]
309        index = np.searchsorted(bins, distance)
310        # Note: indices may be duplicated, so "Pr[index] += w" will not work!!
311        extended_Pr += np.bincount(index, weights, n_max+1)
312        t = time.clock()
313        if t > t_next:
314            t_next = t + 3
315            print("processing %d of %d"%(k, len(rho)-1))
316    Pr = extended_Pr[1:-1]
[8fb2a94]317    return Pr
[cfa28d3]318
[8fb2a94]319def _calc_Pr_uniform(r, rho, points):
[cfa28d3]320    # Make Pr a little be bigger than necessary so that only distances
321    # min < d < max end up in Pr
[8fb2a94]322    dr, n_max = r[0], len(r)
[cfa28d3]323    extended_Pr = np.zeros(n_max+1, 'd')
324    t0 = time.clock()
325    t_next = t0 + 3
326    for k, rho_k in enumerate(rho[:-1]):
327        distances = sqrt(np.sum((points[k] - points[k+1:])**2, axis=1))
328        weights = rho_k * rho[k+1:]
[8fb2a94]329        index = np.minimum(np.asarray(distances/dr, 'i'), n_max)
[cfa28d3]330        # Note: indices may be duplicated, so "Pr[index] += w" will not work!!
331        extended_Pr += np.bincount(index, weights, n_max+1)
332        t = time.clock()
333        if t > t_next:
334            t_next = t + 3
335            print("processing %d of %d"%(k, len(rho)-1))
336    #print("time py:", time.clock() - t0)
337    Pr = extended_Pr[:-1]
338    # Make Pr independent of sampling density.  The factor of 2 comes because
339    # we are only accumulating the upper triangular distances.
340    #Pr = Pr * 2 / len(rho)**2
[8fb2a94]341    return Pr
[cfa28d3]342
343    # Can get an additional 2x by going to C.  Cuda/OpenCL will allow even
344    # more speedup, though still bounded by the n^2 cose.
345    """
346void pdfcalc(int n, const double *pts, const double *rho,
347             int nPr, double *Pr, double rstep)
348{
349  int i,j;
350
351  for (i=0; i<n-2; i++) {
352    for (j=i+1; j<=n-1; j++) {
353      const double dxx=pts[3*i]-pts[3*j];
354      const double dyy=pts[3*i+1]-pts[3*j+1];
355      const double dzz=pts[3*i+2]-pts[3*j+2];
356      const double d=sqrt(dxx*dxx+dyy*dyy+dzz*dzz);
357      const int k=rint(d/rstep);
358      if (k < nPr) Pr[k]+=rho[i]*rho[j];
359    }
360  }
361}
362"""
363
[7d0afa3]364if USE_NUMBA:
365    # Override simple numpy solution with numba if available
[8fb2a94]366    @njit("f8[:](f8[:], f8[:], f8[:,:])")
[7d0afa3]367    def _calc_Pr_uniform(r, rho, points):
[8fb2a94]368        dr = r[0]
369        n_max = len(r)
370        Pr = np.zeros_like(r)
371        for j in range(len(rho) - 1):
372            x, y, z = points[j, 0], points[j, 1], points[j, 2]
373            for k in range(j+1, len(rho)):
374                distance = sqrt((x - points[k, 0])**2
375                                + (y - points[k, 1])**2
376                                + (z - points[k, 2])**2)
377                index = int(distance/dr)
378                if index < n_max:
379                    Pr[index] += rho[j] * rho[k]
380        # Make Pr independent of sampling density.  The factor of 2 comes because
381        # we are only accumulating the upper triangular distances.
382        #Pr = Pr * 2 / len(rho)**2
383        return Pr
384
385
386def calc_Pr(r, rho, points):
387    # P(r) with uniform steps in r is 3x faster; check if we are uniform
388    # before continuing
[8bd379a]389    r, rho, points = [np.asarray(v, 'd') for v in (r, rho, points)]
[8fb2a94]390    if np.max(np.abs(np.diff(r) - r[0])) > r[0]*0.01:
391        Pr = _calc_Pr_nonuniform(r, rho, points)
392    else:
[7d0afa3]393        Pr = _calc_Pr_uniform(r, rho, points)
[8fb2a94]394    return Pr / Pr.max()
395
396
[cfa28d3]397def j0(x):
398    return np.sinc(x/np.pi)
399
400def calc_Iq(q, r, Pr):
401    Iq = np.array([simps(Pr * j0(qk*r), r) for qk in q])
402    #Iq = np.array([np.trapz(Pr * j0(qk*r), r) for qk in q])
403    Iq /= Iq[0]
404    return Iq
405
406# NOTE: copied from sasmodels/resolution.py
407def bin_edges(x):
408    """
409    Determine bin edges from bin centers, assuming that edges are centered
410    between the bins.
411
412    Note: this uses the arithmetic mean, which may not be appropriate for
413    log-scaled data.
414    """
415    if len(x) < 2 or (np.diff(x) < 0).any():
416        raise ValueError("Expected bins to be an increasing set")
417    edges = np.hstack([
418        x[0]  - 0.5*(x[1]  - x[0]),  # first point minus half first interval
419        0.5*(x[1:] + x[:-1]),        # mid points of all central intervals
420        x[-1] + 0.5*(x[-1] - x[-2]), # last point plus half last interval
421        ])
422    return edges
423
[7d0afa3]424# -------------- plotters ----------------
[cfa28d3]425def plot_calc(r, Pr, q, Iq, theory=None):
426    import matplotlib.pyplot as plt
427    plt.subplot(211)
428    plt.plot(r, Pr, '-', label="Pr")
429    plt.xlabel('r (A)')
430    plt.ylabel('Pr (1/A^2)')
431    plt.subplot(212)
432    plt.loglog(q, Iq, '-', label='from Pr')
433    plt.xlabel('q (1/A')
434    plt.ylabel('Iq')
435    if theory is not None:
[7d0afa3]436        plt.loglog(theory[0], theory[1]/theory[1][0], '-', label='analytic')
[cfa28d3]437        plt.legend()
438
[226473d]439def plot_calc_2d(qx, qy, Iqxy, theory=None):
440    import matplotlib.pyplot as plt
441    qx, qy = bin_edges(qx), bin_edges(qy)
442    #qx, qy = np.meshgrid(qx, qy)
443    if theory is not None:
444        plt.subplot(121)
445    plt.pcolormesh(qx, qy, np.log10(Iqxy))
446    plt.xlabel('qx (1/A)')
447    plt.ylabel('qy (1/A)')
448    if theory is not None:
449        plt.subplot(122)
450        plt.pcolormesh(qx, qy, np.log10(theory))
451        plt.xlabel('qx (1/A)')
452
[cfa28d3]453def plot_points(rho, points):
454    import mpl_toolkits.mplot3d
455    import matplotlib.pyplot as plt
456
457    ax = plt.axes(projection='3d')
458    try:
459        ax.axis('square')
460    except Exception:
461        pass
462    n = len(points)
[2ab331f]463    #print("len points", n)
464    index = np.random.choice(n, size=500) if n > 500 else slice(None, None)
[cfa28d3]465    ax.scatter(points[index, 0], points[index, 1], points[index, 2], c=rho[index])
466    #low, high = points.min(axis=0), points.max(axis=0)
467    #ax.axis([low[0], high[0], low[1], high[1], low[2], high[2]])
468    ax.autoscale(True)
469
[7d0afa3]470# ----------- Analytic models --------------
[97be877]471def sas_sinx_x(x):
472    with np.errstate(all='ignore'):
473        retvalue = sin(x)/x
474    retvalue[x == 0.] = 1.
475    return retvalue
476
[cfa28d3]477def sas_2J1x_x(x):
478    with np.errstate(all='ignore'):
479        retvalue = 2*J1(x)/x
480    retvalue[x == 0] = 1.
481    return retvalue
482
483def sas_3j1x_x(x):
484    """return 3*j1(x)/x"""
485    with np.errstate(all='ignore'):
486        retvalue = 3*(sin(x) - x*cos(x))/x**3
487    retvalue[x == 0.] = 1.
488    return retvalue
489
490def cylinder_Iq(q, radius, length):
[a1c32c2]491    z, w = leggauss(76)
492    cos_alpha = (z+1)/2
493    sin_alpha = sqrt(1.0 - cos_alpha**2)
[cfa28d3]494    Iq = np.empty_like(q)
495    for k, qk in enumerate(q):
[a1c32c2]496        qab, qc = qk*sin_alpha, qk*cos_alpha
[7d0afa3]497        Fq = sas_2J1x_x(qab*radius) * sas_sinx_x(qc*length/2)
[a1c32c2]498        Iq[k] = np.sum(w*Fq**2)
[032646d]499    Iq = Iq
[cfa28d3]500    return Iq
501
[226473d]502def cylinder_Iqxy(qx, qy, radius, length, view=(0, 0, 0)):
503    qa, qb, qc = invert_view(qx, qy, view)
[7d0afa3]504    qab = sqrt(qa**2 + qb**2)
505    Fq = sas_2J1x_x(qab*radius) * sas_sinx_x(qc*length/2)
[226473d]506    Iq = Fq**2
507    return Iq.reshape(qx.shape)
508
[cfa28d3]509def sphere_Iq(q, radius):
510    Iq = sas_3j1x_x(q*radius)**2
[032646d]511    return Iq
[cfa28d3]512
[7d0afa3]513def box_Iq(q, a, b, c):
514    z, w = leggauss(76)
515    outer_sum = np.zeros_like(q)
516    for cos_alpha, outer_w in zip((z+1)/2, w):
517        sin_alpha = sqrt(1.0-cos_alpha*cos_alpha)
518        qc = q*cos_alpha
519        siC = c*sas_sinx_x(c*qc/2)
520        inner_sum = np.zeros_like(q)
521        for beta, inner_w in zip((z + 1)*pi/4, w):
522            qa, qb = q*sin_alpha*sin(beta), q*sin_alpha*cos(beta)
523            siA = a*sas_sinx_x(a*qa/2)
524            siB = b*sas_sinx_x(b*qb/2)
525            Fq = siA*siB*siC
526            inner_sum += inner_w * Fq**2
527        outer_sum += outer_w * inner_sum
528    Iq = outer_sum / 4  # = outer*um*zm*8.0/(4.0*M_PI)
[032646d]529    return Iq
[7d0afa3]530
531def box_Iqxy(qx, qy, a, b, c, view=(0, 0, 0)):
532    qa, qb, qc = invert_view(qx, qy, view)
533    sia = sas_sinx_x(qa*a/2)
534    sib = sas_sinx_x(qb*b/2)
535    sic = sas_sinx_x(qc*c/2)
536    Fq = sia*sib*sic
537    Iq = Fq**2
538    return Iq.reshape(qx.shape)
539
[cfa28d3]540def csbox_Iq(q, a, b, c, da, db, dc, slda, sldb, sldc, sld_core):
[a1c32c2]541    z, w = leggauss(76)
542
[cfa28d3]543    sld_solvent = 0
544    overlapping = False
545    dr0 = sld_core - sld_solvent
546    drA, drB, drC = slda-sld_solvent, sldb-sld_solvent, sldc-sld_solvent
547    tA, tB, tC = a + 2*da, b + 2*db, c + 2*dc
548
[a1c32c2]549    outer_sum = np.zeros_like(q)
550    for cos_alpha, outer_w in zip((z+1)/2, w):
[cfa28d3]551        sin_alpha = sqrt(1.0-cos_alpha*cos_alpha)
552        qc = q*cos_alpha
[7d0afa3]553        siC = c*sas_sinx_x(c*qc/2)
554        siCt = tC*sas_sinx_x(tC*qc/2)
[a1c32c2]555        inner_sum = np.zeros_like(q)
556        for beta, inner_w in zip((z + 1)*pi/4, w):
[cfa28d3]557            qa, qb = q*sin_alpha*sin(beta), q*sin_alpha*cos(beta)
[7d0afa3]558            siA = a*sas_sinx_x(a*qa/2)
559            siB = b*sas_sinx_x(b*qb/2)
560            siAt = tA*sas_sinx_x(tA*qa/2)
561            siBt = tB*sas_sinx_x(tB*qb/2)
[cfa28d3]562            if overlapping:
[a1c32c2]563                Fq = (dr0*siA*siB*siC
564                      + drA*(siAt-siA)*siB*siC
565                      + drB*siAt*(siBt-siB)*siC
566                      + drC*siAt*siBt*(siCt-siC))
[cfa28d3]567            else:
[a1c32c2]568                Fq = (dr0*siA*siB*siC
569                      + drA*(siAt-siA)*siB*siC
570                      + drB*siA*(siBt-siB)*siC
571                      + drC*siA*siB*(siCt-siC))
572            inner_sum += inner_w * Fq**2
573        outer_sum += outer_w * inner_sum
574    Iq = outer_sum / 4  # = outer*um*zm*8.0/(4.0*M_PI)
[cfa28d3]575    return Iq/Iq[0]
576
[97be877]577def csbox_Iqxy(qx, qy, a, b, c, da, db, dc, slda, sldb, sldc, sld_core, view=(0,0,0)):
578    qa, qb, qc = invert_view(qx, qy, view)
[226473d]579
[97be877]580    sld_solvent = 0
581    overlapping = False
582    dr0 = sld_core - sld_solvent
583    drA, drB, drC = slda-sld_solvent, sldb-sld_solvent, sldc-sld_solvent
584    tA, tB, tC = a + 2*da, b + 2*db, c + 2*dc
585    siA = a*sas_sinx_x(a*qa/2)
586    siB = b*sas_sinx_x(b*qb/2)
587    siC = c*sas_sinx_x(c*qc/2)
588    siAt = tA*sas_sinx_x(tA*qa/2)
589    siBt = tB*sas_sinx_x(tB*qb/2)
590    siCt = tC*sas_sinx_x(tC*qc/2)
591    Fq = (dr0*siA*siB*siC
592          + drA*(siAt-siA)*siB*siC
593          + drB*siA*(siBt-siB)*siC
594          + drC*siA*siB*(siCt-siC))
595    Iq = Fq**2
596    return Iq.reshape(qx.shape)
[226473d]597
[7d0afa3]598# --------- Test cases -----------
[cfa28d3]599
[032646d]600def build_cylinder(radius=25, length=125, rho=2.):
[226473d]601    shape = EllipticalCylinder(radius, radius, length, rho)
[7d0afa3]602    fn = lambda q: cylinder_Iq(q, radius, length)*rho**2
603    fn_xy = lambda qx, qy, view: cylinder_Iqxy(qx, qy, radius, length, view=view)*rho**2
604    return shape, fn, fn_xy
[226473d]605
[032646d]606def build_sphere(radius=125, rho=2):
[cfa28d3]607    shape = TriaxialEllipsoid(radius, radius, radius, rho)
[032646d]608    fn = lambda q: sphere_Iq(q, radius)*rho**2
[7d0afa3]609    fn_xy = lambda qx, qy, view: sphere_Iq(np.sqrt(qx**2+qy**2), radius)*rho**2
610    return shape, fn, fn_xy
611
[032646d]612def build_box(a=10, b=20, c=30, rho=2.):
[7d0afa3]613    shape = Box(a, b, c, rho)
614    fn = lambda q: box_Iq(q, a, b, c)*rho**2
615    fn_xy = lambda qx, qy, view: box_Iqxy(qx, qy, a, b, c, view=view)*rho**2
616    return shape, fn, fn_xy
617
[032646d]618def build_csbox(a=10, b=20, c=30, da=1, db=2, dc=3, slda=1, sldb=2, sldc=3, sld_core=4):
[7d0afa3]619    shape = csbox(a, b, c, da, db, dc, slda, sldb, sldc, sld_core)
620    fn = lambda q: csbox_Iq(q, a, b, c, da, db, dc, slda, sldb, sldc, sld_core)
621    fn_xy = lambda qx, qy, view: csbox_Iqxy(qx, qy, a, b, c, da, db, dc,
622                                            slda, sldb, sldc, sld_core, view=view)
623    return shape, fn, fn_xy
624
[032646d]625def build_cubic_lattice(shape, nx=1, ny=1, nz=1, dx=2, dy=2, dz=2,
626                  shuffle=0, rotate=0):
627    a, b, c = shape.dims
628    shapes = [copy(shape)
629              .shift((ix+(randn() if shuffle < 0.3 else rand())*shuffle)*dx*a,
630                     (iy+(randn() if shuffle < 0.3 else rand())*shuffle)*dy*b,
631                     (iz+(randn() if shuffle < 0.3 else rand())*shuffle)*dz*c)
632              .rotate(*((randn(3) if rotate < 30 else rand(3))*rotate))
633              for ix in range(nx)
634              for iy in range(ny)
635              for iz in range(nz)]
636    lattice = Composite(shapes)
637    return lattice
638
[7d0afa3]639
640SHAPE_FUNCTIONS = OrderedDict([
[032646d]641    ("cylinder", build_cylinder),
642    ("sphere", build_sphere),
643    ("box", build_box),
644    ("csbox", build_csbox),
[7d0afa3]645])
646SHAPES = list(SHAPE_FUNCTIONS.keys())
647
648def check_shape(title, shape, fn=None, show_points=False,
649                mesh=100, qmax=1.0, r_step=0.01, samples=5000):
650    rho_solvent = 0
[032646d]651    qmin = qmax/100.
[7d0afa3]652    q = np.logspace(np.log10(qmin), np.log10(qmax), mesh)
653    r = shape.r_bins(q, r_step=r_step)
[032646d]654    sampling_density = samples / shape.volume
[7d0afa3]655    rho, points = shape.sample(sampling_density)
656    t0 = time.time()
657    Pr = calc_Pr(r, rho-rho_solvent, points)
658    print("calc Pr time", time.time() - t0)
659    Iq = calc_Iq(q, r, Pr)
660    theory = (q, fn(q)) if fn is not None else None
661
662    import pylab
663    if show_points:
664         plot_points(rho, points); pylab.figure()
665    plot_calc(r, Pr, q, Iq, theory=theory)
666    pylab.gcf().canvas.set_window_title(title)
667    pylab.show()
668
669def check_shape_2d(title, shape, fn=None, view=(0, 0, 0), show_points=False,
670                   mesh=100, qmax=1.0, samples=5000):
671    rho_solvent = 0
672    qx = np.linspace(0.0, qmax, mesh)
673    qy = np.linspace(0.0, qmax, mesh)
674    Qx, Qy = np.meshgrid(qx, qy)
[032646d]675    sampling_density = samples / shape.volume
676    t0 = time.time()
[7d0afa3]677    rho, points = shape.sample(sampling_density)
[032646d]678    print("point generation time", time.time() - t0)
[7d0afa3]679    t0 = time.time()
680    Iqxy = calc_Iqxy(Qx, Qy, rho, points, view=view)
681    print("calc Iqxy time", time.time() - t0)
682    theory = fn(Qx, Qy, view) if fn is not None else None
683    Iqxy += 0.001 * Iqxy.max()
684    if theory is not None:
685        theory += 0.001 * theory.max()
686
687    import pylab
688    if show_points:
689        plot_points(rho, points); pylab.figure()
690    plot_calc_2d(qx, qy, Iqxy, theory=theory)
691    pylab.gcf().canvas.set_window_title(title)
692    pylab.show()
693
694def main():
695    parser = argparse.ArgumentParser(
696        description="Compute scattering from realspace sampling",
697        formatter_class=argparse.ArgumentDefaultsHelpFormatter,
698        )
[032646d]699    parser.add_argument('-d', '--dim', type=int, default=1,
700                        help='dimension 1 or 2')
701    parser.add_argument('-m', '--mesh', type=int, default=100,
702                        help='number of mesh points')
703    parser.add_argument('-s', '--samples', type=int, default=5000,
704                        help="number of sample points")
705    parser.add_argument('-q', '--qmax', type=float, default=0.5,
706                        help='max q')
707    parser.add_argument('-v', '--view', type=str, default='0,0,0',
708                        help='theta,phi,psi angles')
709    parser.add_argument('-n', '--lattice', type=str, default='1,1,1',
710                        help='lattice size')
711    parser.add_argument('-z', '--spacing', type=str, default='2,2,2',
712                        help='lattice spacing')
713    parser.add_argument('-r', '--rotate', type=float, default=0.,
714                        help="rotation relative to lattice, gaussian < 30 degrees, uniform otherwise")
715    parser.add_argument('-w', '--shuffle', type=float, default=0.,
716                        help="position relative to lattice, gaussian < 0.3, uniform otherwise")
717    parser.add_argument('-p', '--plot', action='store_true',
718                        help='plot points')
719    parser.add_argument('shape', choices=SHAPES, nargs='?', default=SHAPES[0],
720                        help='oriented shape')
[7d0afa3]721    parser.add_argument('pars', type=str, nargs='*', help='shape parameters')
722    opts = parser.parse_args()
723    pars = {key: float(value) for p in opts.pars for key, value in [p.split('=')]}
[032646d]724    nx, ny, nz = [int(v) for v in opts.lattice.split(',')]
725    dx, dy, dz = [float(v) for v in opts.spacing.split(',')]
726    shuffle, rotate = opts.shuffle, opts.rotate
[7d0afa3]727    shape, fn, fn_xy = SHAPE_FUNCTIONS[opts.shape](**pars)
[032646d]728    if nx > 1 or ny > 1 or nz > 1:
729        shape = build_cubic_lattice(shape, nx, ny, nz, dx, dy, dz, shuffle, rotate)
[7d0afa3]730    title = "%s(%s)" % (opts.shape, " ".join(opts.pars))
731    if opts.dim == 1:
732        check_shape(title, shape, fn, show_points=opts.plot,
733                    mesh=opts.mesh, qmax=opts.qmax, samples=opts.samples)
734    else:
735        view = tuple(float(v) for v in opts.view.split(','))
736        check_shape_2d(title, shape, fn_xy, view=view, show_points=opts.plot,
737                       mesh=opts.mesh, qmax=opts.qmax, samples=opts.samples)
[cfa28d3]738
[97be877]739
[cfa28d3]740if __name__ == "__main__":
[7d0afa3]741    main()
Note: See TracBrowser for help on using the repository browser.