source: sasmodels/explore/realspace.py @ 3db96b0

core_shell_microgelsmagnetic_modelticket-1257-vesicle-productticket_1156ticket_1265_superballticket_822_more_unit_tests
Last change on this file since 3db96b0 was 3db96b0, checked in by Paul Kienzle <pkienzle@…>, 6 years ago

improve plotting for 2d monte carlo calculator

  • Property mode set to 100644
File size: 27.1 KB
Line 
1from __future__ import division, print_function
2
3import time
4from copy import copy
5import os
6import argparse
7from collections import OrderedDict
8
9import numpy as np
10from numpy import pi, radians, sin, cos, sqrt
11from numpy.random import poisson, uniform, randn, rand
12from numpy.polynomial.legendre import leggauss
13from scipy.integrate import simps
14from scipy.special import j1 as J1
15
16try:
17    import numba
18    USE_NUMBA = True
19except ImportError:
20    USE_NUMBA = False
21
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
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
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
93    def dims(self):
94        # type: () -> float, float, float
95        raise NotImplementedError()
96
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])
131        self.volume = sum(shape.volume for shape in self.shapes)
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)
146        self.dims = a, b, c
147        self.volume = a*b*c
148
149    def sample(self, density):
150        num_points = poisson(density*self.volume)
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)
164        self.dims = 2*ra, 2*rb, length
165        self.volume = pi*ra*rb*length
166
167    def sample(self, density):
168        # randomly sample from a box of side length 2*r, excluding anything
169        # not in the cylinder
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)
186        self.dims = 2*ra, 2*rb, 2*rc
187        self.volume = 4*pi/3 * ra * rb * rc
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)
205        helix_length = helix_pitch * tube_length/sqrt(helix_radius**2 + helix_pitch**2)
206        total_radius = self.helix_radius + self.tube_radius
207        self.helix_radius, self.helix_pitch = helix_radius, helix_pitch
208        self.tube_radius, self.tube_length = tube_radius, tube_length
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
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
214        self.volume = pi*self.tube_radius**2*self.tube_length
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
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))
263    shape.dims = 2*da+a, 2*db+b, 2*dc+c
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
274    from numba import njit
275    @njit("f8[:](f8[:],f8[:],f8[:],f8[:],f8[:],f8[:],f8[:])")
276    def _Iqxy(values, x, y, z, qa, qb, qc):
277        Iq = np.zeros_like(qa)
278        for j in range(len(Iq)):
279            total = 0. + 0j
280            for k in range(len(values)):
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)):
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
290    x, y, z = points.T
291    values, x, y, z, qa, qb, qc = [np.asarray(v, 'd')
292                                   for v in (values, x, y, z, qa, qb, qc)]
293
294    # I(q) = |sum V(r) rho(r) e^(1j q.r)|^2 / sum V(r)
295    Iq = _Iqxy(values, x, y, z, qa.flatten(), qb.flatten(), qc.flatten())
296    return np.asarray(Iq).reshape(qx.shape) / np.sum(volume)
297
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]
317    return Pr
318
319def _calc_Pr_uniform(r, rho, points):
320    # Make Pr a little be bigger than necessary so that only distances
321    # min < d < max end up in Pr
322    dr, n_max = r[0], len(r)
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:]
329        index = np.minimum(np.asarray(distances/dr, 'i'), n_max)
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
341    return Pr
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
364if USE_NUMBA:
365    # Override simple numpy solution with numba if available
366    @njit("f8[:](f8[:], f8[:], f8[:,:])")
367    def _calc_Pr_uniform(r, rho, points):
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
389    r, rho, points = [np.asarray(v, 'd') for v in (r, rho, points)]
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:
393        Pr = _calc_Pr_uniform(r, rho, points)
394    return Pr / Pr.max()
395
396
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
424# -------------- plotters ----------------
425def plot_calc(r, Pr, q, Iq, theory=None, title=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    if title is not None:
432        plt.title(title)
433    plt.subplot(212)
434    plt.loglog(q, Iq, '-', label='from Pr')
435    plt.xlabel('q (1/A')
436    plt.ylabel('Iq')
437    if theory is not None:
438        plt.loglog(theory[0], theory[1]/theory[1][0], '-', label='analytic')
439        plt.legend()
440
441def plot_calc_2d(qx, qy, Iqxy, theory=None, title=None):
442    import matplotlib.pyplot as plt
443    qx, qy = bin_edges(qx), bin_edges(qy)
444    #qx, qy = np.meshgrid(qx, qy)
445    if theory is not None:
446        plt.subplot(121)
447    #plt.pcolor(qx, qy, np.log10(Iqxy))
448    extent = [qx[0], qx[-1], qy[0], qy[-1]]
449    plt.imshow(np.log10(Iqxy), extent=extent, interpolation="nearest",
450               origin='lower')
451    plt.xlabel('qx (1/A)')
452    plt.ylabel('qy (1/A)')
453    plt.axis('equal')
454    plt.axis(extent)
455    #plt.grid(True)
456    if title is not None:
457        plt.title(title)
458    if theory is not None:
459        plt.subplot(122)
460        plt.imshow(np.log10(theory), extent=extent, interpolation="nearest",
461                   origin='lower')
462        plt.axis('equal')
463        plt.axis(extent)
464        plt.xlabel('qx (1/A)')
465
466def plot_points(rho, points):
467    import mpl_toolkits.mplot3d
468    import matplotlib.pyplot as plt
469
470    ax = plt.axes(projection='3d')
471    try:
472        ax.axis('square')
473    except Exception:
474        pass
475    n = len(points)
476    #print("len points", n)
477    index = np.random.choice(n, size=500) if n > 500 else slice(None, None)
478    ax.scatter(points[index, 0], points[index, 1], points[index, 2], c=rho[index])
479    #low, high = points.min(axis=0), points.max(axis=0)
480    #ax.axis([low[0], high[0], low[1], high[1], low[2], high[2]])
481    ax.autoscale(True)
482
483# ----------- Analytic models --------------
484def sas_sinx_x(x):
485    with np.errstate(all='ignore'):
486        retvalue = sin(x)/x
487    retvalue[x == 0.] = 1.
488    return retvalue
489
490def sas_2J1x_x(x):
491    with np.errstate(all='ignore'):
492        retvalue = 2*J1(x)/x
493    retvalue[x == 0] = 1.
494    return retvalue
495
496def sas_3j1x_x(x):
497    """return 3*j1(x)/x"""
498    with np.errstate(all='ignore'):
499        retvalue = 3*(sin(x) - x*cos(x))/x**3
500    retvalue[x == 0.] = 1.
501    return retvalue
502
503def cylinder_Iq(q, radius, length):
504    z, w = leggauss(76)
505    cos_alpha = (z+1)/2
506    sin_alpha = sqrt(1.0 - cos_alpha**2)
507    Iq = np.empty_like(q)
508    for k, qk in enumerate(q):
509        qab, qc = qk*sin_alpha, qk*cos_alpha
510        Fq = sas_2J1x_x(qab*radius) * sas_sinx_x(qc*length/2)
511        Iq[k] = np.sum(w*Fq**2)
512    Iq = Iq
513    return Iq
514
515def cylinder_Iqxy(qx, qy, radius, length, view=(0, 0, 0)):
516    qa, qb, qc = invert_view(qx, qy, view)
517    qab = sqrt(qa**2 + qb**2)
518    Fq = sas_2J1x_x(qab*radius) * sas_sinx_x(qc*length/2)
519    Iq = Fq**2
520    return Iq.reshape(qx.shape)
521
522def sphere_Iq(q, radius):
523    Iq = sas_3j1x_x(q*radius)**2
524    return Iq
525
526def box_Iq(q, a, b, c):
527    z, w = leggauss(76)
528    outer_sum = np.zeros_like(q)
529    for cos_alpha, outer_w in zip((z+1)/2, w):
530        sin_alpha = sqrt(1.0-cos_alpha*cos_alpha)
531        qc = q*cos_alpha
532        siC = c*sas_sinx_x(c*qc/2)
533        inner_sum = np.zeros_like(q)
534        for beta, inner_w in zip((z + 1)*pi/4, w):
535            qa, qb = q*sin_alpha*sin(beta), q*sin_alpha*cos(beta)
536            siA = a*sas_sinx_x(a*qa/2)
537            siB = b*sas_sinx_x(b*qb/2)
538            Fq = siA*siB*siC
539            inner_sum += inner_w * Fq**2
540        outer_sum += outer_w * inner_sum
541    Iq = outer_sum / 4  # = outer*um*zm*8.0/(4.0*M_PI)
542    return Iq
543
544def box_Iqxy(qx, qy, a, b, c, view=(0, 0, 0)):
545    qa, qb, qc = invert_view(qx, qy, view)
546    sia = sas_sinx_x(qa*a/2)
547    sib = sas_sinx_x(qb*b/2)
548    sic = sas_sinx_x(qc*c/2)
549    Fq = sia*sib*sic
550    Iq = Fq**2
551    return Iq.reshape(qx.shape)
552
553def csbox_Iq(q, a, b, c, da, db, dc, slda, sldb, sldc, sld_core):
554    z, w = leggauss(76)
555
556    sld_solvent = 0
557    overlapping = False
558    dr0 = sld_core - sld_solvent
559    drA, drB, drC = slda-sld_solvent, sldb-sld_solvent, sldc-sld_solvent
560    tA, tB, tC = a + 2*da, b + 2*db, c + 2*dc
561
562    outer_sum = np.zeros_like(q)
563    for cos_alpha, outer_w in zip((z+1)/2, w):
564        sin_alpha = sqrt(1.0-cos_alpha*cos_alpha)
565        qc = q*cos_alpha
566        siC = c*sas_sinx_x(c*qc/2)
567        siCt = tC*sas_sinx_x(tC*qc/2)
568        inner_sum = np.zeros_like(q)
569        for beta, inner_w in zip((z + 1)*pi/4, w):
570            qa, qb = q*sin_alpha*sin(beta), q*sin_alpha*cos(beta)
571            siA = a*sas_sinx_x(a*qa/2)
572            siB = b*sas_sinx_x(b*qb/2)
573            siAt = tA*sas_sinx_x(tA*qa/2)
574            siBt = tB*sas_sinx_x(tB*qb/2)
575            if overlapping:
576                Fq = (dr0*siA*siB*siC
577                      + drA*(siAt-siA)*siB*siC
578                      + drB*siAt*(siBt-siB)*siC
579                      + drC*siAt*siBt*(siCt-siC))
580            else:
581                Fq = (dr0*siA*siB*siC
582                      + drA*(siAt-siA)*siB*siC
583                      + drB*siA*(siBt-siB)*siC
584                      + drC*siA*siB*(siCt-siC))
585            inner_sum += inner_w * Fq**2
586        outer_sum += outer_w * inner_sum
587    Iq = outer_sum / 4  # = outer*um*zm*8.0/(4.0*M_PI)
588    return Iq/Iq[0]
589
590def csbox_Iqxy(qx, qy, a, b, c, da, db, dc, slda, sldb, sldc, sld_core, view=(0,0,0)):
591    qa, qb, qc = invert_view(qx, qy, view)
592
593    sld_solvent = 0
594    overlapping = False
595    dr0 = sld_core - sld_solvent
596    drA, drB, drC = slda-sld_solvent, sldb-sld_solvent, sldc-sld_solvent
597    tA, tB, tC = a + 2*da, b + 2*db, c + 2*dc
598    siA = a*sas_sinx_x(a*qa/2)
599    siB = b*sas_sinx_x(b*qb/2)
600    siC = c*sas_sinx_x(c*qc/2)
601    siAt = tA*sas_sinx_x(tA*qa/2)
602    siBt = tB*sas_sinx_x(tB*qb/2)
603    siCt = tC*sas_sinx_x(tC*qc/2)
604    Fq = (dr0*siA*siB*siC
605          + drA*(siAt-siA)*siB*siC
606          + drB*siA*(siBt-siB)*siC
607          + drC*siA*siB*(siCt-siC))
608    Iq = Fq**2
609    return Iq.reshape(qx.shape)
610
611# --------- Test cases -----------
612
613def build_cylinder(radius=25, length=125, rho=2.):
614    shape = EllipticalCylinder(radius, radius, length, rho)
615    fn = lambda q: cylinder_Iq(q, radius, length)*rho**2
616    fn_xy = lambda qx, qy, view: cylinder_Iqxy(qx, qy, radius, length, view=view)*rho**2
617    return shape, fn, fn_xy
618
619def build_sphere(radius=125, rho=2):
620    shape = TriaxialEllipsoid(radius, radius, radius, rho)
621    fn = lambda q: sphere_Iq(q, radius)*rho**2
622    fn_xy = lambda qx, qy, view: sphere_Iq(np.sqrt(qx**2+qy**2), radius)*rho**2
623    return shape, fn, fn_xy
624
625def build_box(a=10, b=20, c=30, rho=2.):
626    shape = Box(a, b, c, rho)
627    fn = lambda q: box_Iq(q, a, b, c)*rho**2
628    fn_xy = lambda qx, qy, view: box_Iqxy(qx, qy, a, b, c, view=view)*rho**2
629    return shape, fn, fn_xy
630
631def build_csbox(a=10, b=20, c=30, da=1, db=2, dc=3, slda=1, sldb=2, sldc=3, sld_core=4):
632    shape = csbox(a, b, c, da, db, dc, slda, sldb, sldc, sld_core)
633    fn = lambda q: csbox_Iq(q, a, b, c, da, db, dc, slda, sldb, sldc, sld_core)
634    fn_xy = lambda qx, qy, view: csbox_Iqxy(qx, qy, a, b, c, da, db, dc,
635                                            slda, sldb, sldc, sld_core, view=view)
636    return shape, fn, fn_xy
637
638def build_cubic_lattice(shape, nx=1, ny=1, nz=1, dx=2, dy=2, dz=2,
639                  shuffle=0, rotate=0):
640    a, b, c = shape.dims
641    shapes = [copy(shape)
642              .shift((ix+(randn() if shuffle < 0.3 else rand())*shuffle)*dx*a,
643                     (iy+(randn() if shuffle < 0.3 else rand())*shuffle)*dy*b,
644                     (iz+(randn() if shuffle < 0.3 else rand())*shuffle)*dz*c)
645              .rotate(*((randn(3) if rotate < 30 else rand(3))*rotate))
646              for ix in range(nx)
647              for iy in range(ny)
648              for iz in range(nz)]
649    lattice = Composite(shapes)
650    return lattice
651
652
653SHAPE_FUNCTIONS = OrderedDict([
654    ("cylinder", build_cylinder),
655    ("sphere", build_sphere),
656    ("box", build_box),
657    ("csbox", build_csbox),
658])
659SHAPES = list(SHAPE_FUNCTIONS.keys())
660
661def check_shape(title, shape, fn=None, show_points=False,
662                mesh=100, qmax=1.0, r_step=0.01, samples=5000):
663    rho_solvent = 0
664    qmin = qmax/100.
665    q = np.logspace(np.log10(qmin), np.log10(qmax), mesh)
666    r = shape.r_bins(q, r_step=r_step)
667    sampling_density = samples / shape.volume
668    rho, points = shape.sample(sampling_density)
669    t0 = time.time()
670    Pr = calc_Pr(r, rho-rho_solvent, points)
671    print("calc Pr time", time.time() - t0)
672    Iq = calc_Iq(q, r, Pr)
673    theory = (q, fn(q)) if fn is not None else None
674
675    import pylab
676    if show_points:
677         plot_points(rho, points); pylab.figure()
678    plot_calc(r, Pr, q, Iq, theory=theory, title=title)
679    pylab.gcf().canvas.set_window_title(title)
680    pylab.show()
681
682def check_shape_2d(title, shape, fn=None, view=(0, 0, 0), show_points=False,
683                   mesh=100, qmax=1.0, samples=5000):
684    rho_solvent = 0
685    qx = np.linspace(0.0, qmax, mesh)
686    qy = np.linspace(0.0, qmax, mesh)
687    Qx, Qy = np.meshgrid(qx, qy)
688    sampling_density = samples / shape.volume
689    t0 = time.time()
690    rho, points = shape.sample(sampling_density)
691    print("point generation time", time.time() - t0)
692    t0 = time.time()
693    Iqxy = calc_Iqxy(Qx, Qy, rho, points, view=view)
694    print("calc Iqxy time", time.time() - t0)
695    theory = fn(Qx, Qy, view) if fn is not None else None
696    Iqxy += 0.001 * Iqxy.max()
697    if theory is not None:
698        theory += 0.001 * theory.max()
699
700    import pylab
701    if show_points:
702        plot_points(rho, points); pylab.figure()
703    plot_calc_2d(qx, qy, Iqxy, theory=theory, title=title)
704    pylab.gcf().canvas.set_window_title(title)
705    pylab.show()
706
707def main():
708    parser = argparse.ArgumentParser(
709        description="Compute scattering from realspace sampling",
710        formatter_class=argparse.ArgumentDefaultsHelpFormatter,
711        )
712    parser.add_argument('-d', '--dim', type=int, default=1,
713                        help='dimension 1 or 2')
714    parser.add_argument('-m', '--mesh', type=int, default=100,
715                        help='number of mesh points')
716    parser.add_argument('-s', '--samples', type=int, default=5000,
717                        help="number of sample points")
718    parser.add_argument('-q', '--qmax', type=float, default=0.5,
719                        help='max q')
720    parser.add_argument('-v', '--view', type=str, default='0,0,0',
721                        help='theta,phi,psi angles')
722    parser.add_argument('-n', '--lattice', type=str, default='1,1,1',
723                        help='lattice size')
724    parser.add_argument('-z', '--spacing', type=str, default='2,2,2',
725                        help='lattice spacing')
726    parser.add_argument('-r', '--rotate', type=float, default=0.,
727                        help="rotation relative to lattice, gaussian < 30 degrees, uniform otherwise")
728    parser.add_argument('-w', '--shuffle', type=float, default=0.,
729                        help="position relative to lattice, gaussian < 0.3, uniform otherwise")
730    parser.add_argument('-p', '--plot', action='store_true',
731                        help='plot points')
732    parser.add_argument('shape', choices=SHAPES, nargs='?', default=SHAPES[0],
733                        help='oriented shape')
734    parser.add_argument('pars', type=str, nargs='*', help='shape parameters')
735    opts = parser.parse_args()
736    pars = {key: float(value) for p in opts.pars for key, value in [p.split('=')]}
737    nx, ny, nz = [int(v) for v in opts.lattice.split(',')]
738    dx, dy, dz = [float(v) for v in opts.spacing.split(',')]
739    shuffle, rotate = opts.shuffle, opts.rotate
740    shape, fn, fn_xy = SHAPE_FUNCTIONS[opts.shape](**pars)
741    if nx > 1 or ny > 1 or nz > 1:
742        shape = build_cubic_lattice(shape, nx, ny, nz, dx, dy, dz, shuffle, rotate)
743    title = "%s(%s)" % (opts.shape, " ".join(opts.pars))
744    if opts.dim == 1:
745        check_shape(title, shape, fn, show_points=opts.plot,
746                    mesh=opts.mesh, qmax=opts.qmax, samples=opts.samples)
747    else:
748        view = tuple(float(v) for v in opts.view.split(','))
749        check_shape_2d(title, shape, fn_xy, view=view, show_points=opts.plot,
750                       mesh=opts.mesh, qmax=opts.qmax, samples=opts.samples)
751
752
753if __name__ == "__main__":
754    main()
Note: See TracBrowser for help on using the repository browser.