source: sasmodels/explore/realspace.py @ 97be877

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

update authorship/verification for core-shell parallelepiped

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