1 | from __future__ import division, print_function |
---|
2 | |
---|
3 | import time |
---|
4 | from copy import copy |
---|
5 | |
---|
6 | import numpy as np |
---|
7 | from numpy import pi, radians, sin, cos, sqrt |
---|
8 | from numpy.random import poisson, uniform |
---|
9 | from numpy.polynomial.legendre import leggauss |
---|
10 | from scipy.integrate import simps |
---|
11 | from 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 |
---|
15 | def 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 | |
---|
23 | def 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 | |
---|
31 | def 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 | |
---|
39 | def 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 | |
---|
47 | def 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 | |
---|
59 | def 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 | |
---|
71 | class 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 | |
---|
103 | class 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 | |
---|
126 | class 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 | |
---|
145 | class 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 | |
---|
167 | class 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 | |
---|
190 | class 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 | |
---|
246 | def 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 | |
---|
258 | def _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 | |
---|
279 | def 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 | """ |
---|
312 | void 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 | |
---|
330 | def j0(x): |
---|
331 | return np.sinc(x/np.pi) |
---|
332 | |
---|
333 | def 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 |
---|
340 | def 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 | |
---|
357 | def 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 | |
---|
371 | def 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 | |
---|
385 | def 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 | |
---|
402 | def sas_2J1x_x(x): |
---|
403 | with np.errstate(all='ignore'): |
---|
404 | retvalue = 2*J1(x)/x |
---|
405 | retvalue[x == 0] = 1. |
---|
406 | return retvalue |
---|
407 | |
---|
408 | def sas_3j1x_x(x): |
---|
409 | """return 3*j1(x)/x""" |
---|
410 | with np.errstate(all='ignore'): |
---|
411 | retvalue = 3*(sin(x) - x*cos(x))/x**3 |
---|
412 | retvalue[x == 0.] = 1. |
---|
413 | return retvalue |
---|
414 | |
---|
415 | def cylinder_Iq(q, radius, length): |
---|
416 | z, w = leggauss(76) |
---|
417 | cos_alpha = (z+1)/2 |
---|
418 | sin_alpha = sqrt(1.0 - cos_alpha**2) |
---|
419 | Iq = np.empty_like(q) |
---|
420 | for k, qk in enumerate(q): |
---|
421 | qab, qc = qk*sin_alpha, qk*cos_alpha |
---|
422 | Fq = sas_2J1x_x(qab*radius) * j0(qc*length/2) |
---|
423 | Iq[k] = np.sum(w*Fq**2) |
---|
424 | Iq = Iq/Iq[0] |
---|
425 | return Iq |
---|
426 | |
---|
427 | def cylinder_Iqxy(qx, qy, radius, length, view=(0, 0, 0)): |
---|
428 | qa, qb, qc = invert_view(qx, qy, view) |
---|
429 | qab = np.sqrt(qa**2 + qb**2) |
---|
430 | Fq = sas_2J1x_x(qab*radius) * j0(qc*length/2) |
---|
431 | Iq = Fq**2 |
---|
432 | return Iq.reshape(qx.shape) |
---|
433 | |
---|
434 | def sphere_Iq(q, radius): |
---|
435 | Iq = sas_3j1x_x(q*radius)**2 |
---|
436 | return Iq/Iq[0] |
---|
437 | |
---|
438 | def csbox_Iq(q, a, b, c, da, db, dc, slda, sldb, sldc, sld_core): |
---|
439 | z, w = leggauss(76) |
---|
440 | |
---|
441 | sld_solvent = 0 |
---|
442 | overlapping = False |
---|
443 | dr0 = sld_core - sld_solvent |
---|
444 | drA, drB, drC = slda-sld_solvent, sldb-sld_solvent, sldc-sld_solvent |
---|
445 | tA, tB, tC = a + 2*da, b + 2*db, c + 2*dc |
---|
446 | |
---|
447 | outer_sum = np.zeros_like(q) |
---|
448 | for cos_alpha, outer_w in zip((z+1)/2, w): |
---|
449 | sin_alpha = sqrt(1.0-cos_alpha*cos_alpha) |
---|
450 | qc = q*cos_alpha |
---|
451 | siC = c*j0(c*qc/2) |
---|
452 | siCt = tC*j0(tC*qc/2) |
---|
453 | inner_sum = np.zeros_like(q) |
---|
454 | for beta, inner_w in zip((z + 1)*pi/4, w): |
---|
455 | qa, qb = q*sin_alpha*sin(beta), q*sin_alpha*cos(beta) |
---|
456 | siA = a*j0(a*qa/2) |
---|
457 | siB = b*j0(b*qb/2) |
---|
458 | siAt = tA*j0(tA*qa/2) |
---|
459 | siBt = tB*j0(tB*qb/2) |
---|
460 | if overlapping: |
---|
461 | Fq = (dr0*siA*siB*siC |
---|
462 | + drA*(siAt-siA)*siB*siC |
---|
463 | + drB*siAt*(siBt-siB)*siC |
---|
464 | + drC*siAt*siBt*(siCt-siC)) |
---|
465 | else: |
---|
466 | Fq = (dr0*siA*siB*siC |
---|
467 | + drA*(siAt-siA)*siB*siC |
---|
468 | + drB*siA*(siBt-siB)*siC |
---|
469 | + drC*siA*siB*(siCt-siC)) |
---|
470 | inner_sum += inner_w * Fq**2 |
---|
471 | outer_sum += outer_w * inner_sum |
---|
472 | Iq = outer_sum / 4 # = outer*um*zm*8.0/(4.0*M_PI) |
---|
473 | return Iq/Iq[0] |
---|
474 | |
---|
475 | def check_shape(shape, fn=None): |
---|
476 | rho_solvent = 0 |
---|
477 | q = np.logspace(-3, 0, 200) |
---|
478 | r = shape.r_bins(q, r_step=0.01) |
---|
479 | sampling_density = 50000 / shape.volume() |
---|
480 | rho, points = shape.sample(sampling_density) |
---|
481 | Pr = calc_Pr(r, rho-rho_solvent, points) |
---|
482 | Iq = calc_Iq(q, r, Pr) |
---|
483 | theory = (q, fn(q)) if fn is not None else None |
---|
484 | |
---|
485 | import pylab |
---|
486 | #plot_points(rho, points); pylab.figure() |
---|
487 | plot_calc(r, Pr, q, Iq, theory=theory) |
---|
488 | pylab.show() |
---|
489 | |
---|
490 | def check_shape_2d(shape, fn=None, view=(0, 0, 0)): |
---|
491 | rho_solvent = 0 |
---|
492 | nq, qmax = 100, 0.1 |
---|
493 | nq, qmax = 100, 0.03 |
---|
494 | qx = np.linspace(0.0, qmax, nq) |
---|
495 | qy = np.linspace(0.0, qmax, nq) |
---|
496 | Qx, Qy = np.meshgrid(qx, qy) |
---|
497 | sampling_density = 15000 / shape.volume() |
---|
498 | #t0 = time.time() |
---|
499 | rho, points = shape.sample(sampling_density) |
---|
500 | #print("sample time", time.time() - t0) |
---|
501 | t0 = time.time() |
---|
502 | Iqxy = calc_Iqxy(Qx, Qy, rho, points, view=view) |
---|
503 | print("calc time", time.time() - t0) |
---|
504 | theory = fn(Qx, Qy) if fn is not None else None |
---|
505 | Iqxy += 0.001 * Iqxy.max() |
---|
506 | if theory is not None: |
---|
507 | theory += 0.001 * theory.max() |
---|
508 | |
---|
509 | import pylab |
---|
510 | #plot_points(rho, points); pylab.figure() |
---|
511 | plot_calc_2d(qx, qy, Iqxy, theory=theory) |
---|
512 | pylab.show() |
---|
513 | |
---|
514 | def check_cylinder(radius=25, length=125, rho=2.): |
---|
515 | shape = EllipticalCylinder(radius, radius, length, rho) |
---|
516 | fn = lambda q: cylinder_Iq(q, radius, length) |
---|
517 | check_shape(shape, fn) |
---|
518 | |
---|
519 | def check_cylinder_2d(radius=25, length=125, rho=2., view=(0, 0, 0)): |
---|
520 | shape = EllipticalCylinder(radius, radius, length, rho) |
---|
521 | fn = lambda qx, qy, view=view: cylinder_Iqxy(qx, qy, radius, length, view=view) |
---|
522 | check_shape_2d(shape, fn, view=view) |
---|
523 | |
---|
524 | def check_cylinder_2d_lattice(radius=25, length=125, rho=2., |
---|
525 | view=(0, 0, 0)): |
---|
526 | nx, dx = 1, 2*radius |
---|
527 | ny, dy = 30, 2*radius |
---|
528 | nz, dz = 30, length |
---|
529 | dx, dy, dz = 2*dx, 2*dy, 2*dz |
---|
530 | def center(*args): |
---|
531 | sigma = 0.333 |
---|
532 | space = 2 |
---|
533 | return [(space*n+np.random.randn()*sigma)*x for n, x in args] |
---|
534 | shapes = [EllipticalCylinder(radius, radius, length, rho, |
---|
535 | #center=(ix*dx, iy*dy, iz*dz) |
---|
536 | orientation=np.random.randn(3)*0, |
---|
537 | center=center((ix, dx), (iy, dy), (iz, dz)) |
---|
538 | ) |
---|
539 | for ix in range(nx) |
---|
540 | for iy in range(ny) |
---|
541 | for iz in range(nz)] |
---|
542 | shape = Composite(shapes) |
---|
543 | fn = lambda qx, qy, view=view: cylinder_Iqxy(qx, qy, radius, length, view=view) |
---|
544 | check_shape_2d(shape, fn, view=view) |
---|
545 | |
---|
546 | def check_sphere(radius=125, rho=2): |
---|
547 | shape = TriaxialEllipsoid(radius, radius, radius, rho) |
---|
548 | fn = lambda q: sphere_Iq(q, radius) |
---|
549 | check_shape(shape, fn) |
---|
550 | |
---|
551 | def check_csbox(a=10, b=20, c=30, da=1, db=2, dc=3, slda=1, sldb=2, sldc=3, sld_core=4): |
---|
552 | core = Box(a, b, c, sld_core) |
---|
553 | side_a = Box(da, b, c, slda, center=((a+da)/2, 0, 0)) |
---|
554 | side_b = Box(a, db, c, sldb, center=(0, (b+db)/2, 0)) |
---|
555 | side_c = Box(a, b, dc, sldc, center=(0, 0, (c+dc)/2)) |
---|
556 | side_a2 = copy(side_a).shift(-a-da, 0, 0) |
---|
557 | side_b2 = copy(side_b).shift(0, -b-db, 0) |
---|
558 | side_c2 = copy(side_c).shift(0, 0, -c-dc) |
---|
559 | shape = Composite((core, side_a, side_b, side_c, side_a2, side_b2, side_c2)) |
---|
560 | fn = lambda q: csbox_Iq(q, a, b, c, da, db, dc, slda, sldb, sldc, sld_core) |
---|
561 | check_shape(shape, fn) |
---|
562 | |
---|
563 | if __name__ == "__main__": |
---|
564 | check_cylinder(radius=10, length=40) |
---|
565 | #check_cylinder_2d(radius=10, length=40, view=(90,30,0)) |
---|
566 | #check_cylinder_2d_lattice(radius=10, length=50, view=(90,30,0)) |
---|
567 | #check_sphere() |
---|
568 | #check_csbox() |
---|
569 | #check_csbox(da=100, db=200, dc=300) |
---|