source: sasmodels/doc/guide/plugin.rst @ a6d3f46

Last change on this file since a6d3f46 was 7e6bc45e, checked in by Paul Kienzle <pkienzle@…>, 7 years ago

update orientation docs

  • Property mode set to 100644
File size: 45.5 KB
RevLine 
[990d8df]1.. _Writing_a_Plugin:
2
3Writing a Plugin Model
4======================
5
6Overview
7^^^^^^^^
8
9In addition to the models provided with the sasmodels package, you are free to
10create your own models.
11
12Models can be of three types:
13
14- A pure python model : Example -
15  `broadpeak.py <https://github.com/SasView/sasmodels/blob/master/sasmodels/models/broad_peak.py>`_
16
17- A python model with embedded C : Example -
18  `sphere.py <https://github.com/SasView/sasmodels/blob/master/sasmodels/models/sphere.py>`_
19
20- A python wrapper with separate C code : Example -
21  `cylinder.py <https://github.com/SasView/sasmodels/blob/master/sasmodels/models/cylinder.py>`_,
22  `cylinder.c <https://github.com/SasView/sasmodels/blob/master/sasmodels/models/cylinder.c>`_
23
24When using SasView, plugin models should be saved to the SasView
25*plugin_models* folder *C:\\Users\\{username}\\.sasview\\plugin_models*
26(on Windows) or */Users/{username}/.sasview\\plugin_models* (on Mac).
27The next time SasView is started it will compile the plugin and add
28it to the list of *Plugin Models* in a FitPage.  Scripts can load
29the models from anywhere.
30
31The built-in modules are available in the *models* subdirectory
32of the sasmodels package.  For SasView on Windows, these will
33be found in *C:\\Program Files (x86)\\SasView\\sasmodels-data\\models*.
34On Mac OSX, these will be within the application bundle as
35*/Applications/SasView 4.0.app/Contents/Resources/sasmodels-data/models*.
36
37Other models are available for download from the
38`Model Marketplace <http://marketplace.sasview.org/>`_. You can contribute your
39own models to the Marketplace as well.
40
41Create New Model Files
42^^^^^^^^^^^^^^^^^^^^^^
43
44Copy the appropriate files to your plugin models directory (we recommend
45using the examples above as templates) as mymodel.py (and mymodel.c, etc)
46as required, where "mymodel" is the name for the model you are creating.
47
48*Please follow these naming rules:*
49
50- No capitalization and thus no CamelCase
51- If necessary use underscore to separate words (i.e. barbell not BarBell or
52  broad_peak not BroadPeak)
53- Do not include "model" in the name (i.e. barbell not BarBellModel)
54
55
56Edit New Model Files
57^^^^^^^^^^^^^^^^^^^^
58
59Model Contents
60..............
61
62The model interface definition is in the .py file.  This file contains:
63
64- a **model name**:
65   - this is the **name** string in the *.py* file
66   - titles should be:
67
68    - all in *lower* case
69    - without spaces (use underscores to separate words instead)
70    - without any capitalization or CamelCase
71    - without incorporating the word "model"
72    - examples: *barbell* **not** *BarBell*; *broad_peak* **not** *BroadPeak*;
73      *barbell* **not** *BarBellModel*
74
75- a **model title**:
76   - this is the **title** string in the *.py* file
77   - this is a one or two line description of the model, which will appear
78     at the start of the model documentation and as a tooltip in the SasView GUI
79
[3048ec6]80- a **short description**:
[990d8df]81   - this is the **description** string in the *.py* file
82   - this is a medium length description which appears when you click
83     *Description* on the model FitPage
84
85- a **parameter table**:
86   - this will be auto-generated from the *parameters* in the *.py* file
87
88- a **long description**:
89   - this is ReStructuredText enclosed between the r""" and """ delimiters
90     at the top of the *.py* file
91   - what you write here is abstracted into the SasView help documentation
92   - this is what other users will refer to when they want to know what
93     your model does; so please be helpful!
94
95- a **definition** of the model:
96   - as part of the **long description**
97
98- a **formula** defining the function the model calculates:
99   - as part of the **long description**
100
101- an **explanation of the parameters**:
102   - as part of the **long description**
103   - explaining how the symbols in the formula map to the model parameters
104
105- a **plot** of the function, with a **figure caption**:
106   - this is automatically generated from your default parameters
107
108- at least one **reference**:
109   - as part of the **long description**
110   - specifying where the reader can obtain more information about the model
111
112- the **name of the author**
113   - as part of the **long description**
114   - the *.py* file should also contain a comment identifying *who*
115     converted/created the model file
116
117Models that do not conform to these requirements will *never* be incorporated
118into the built-in library.
119
120
121Model Documentation
122...................
123
124The *.py* file starts with an r (for raw) and three sets of quotes
125to start the doc string and ends with a second set of three quotes.
126For example::
127
128    r"""
129    Definition
130    ----------
131
132    The 1D scattering intensity of the sphere is calculated in the following
133    way (Guinier, 1955)
134
135    .. math::
136
137        I(q) = \frac{\text{scale}}{V} \cdot \left[
138            3V(\Delta\rho) \cdot \frac{\sin(qr) - qr\cos(qr))}{(qr)^3}
139            \right]^2 + \text{background}
140
141    where *scale* is a volume fraction, $V$ is the volume of the scatterer,
142    $r$ is the radius of the sphere and *background* is the background level.
143    *sld* and *sld_solvent* are the scattering length densities (SLDs) of the
144    scatterer and the solvent respectively, whose difference is $\Delta\rho$.
145
146    You can included figures in your documentation, as in the following
147    figure for the cylinder model.
148
149    .. figure:: img/cylinder_angle_definition.jpg
150
151        Definition of the angles for oriented cylinders.
152
153    References
154    ----------
155
156    A Guinier, G Fournet, *Small-Angle Scattering of X-Rays*,
157    John Wiley and Sons, New York, (1955)
158    """
159
160This is where the FULL documentation for the model goes (to be picked up by
161the automatic documentation system).  Although it feels odd, you
162should start the documentation immediately with the **definition**---the model
163name, a brief description and the parameter table are automatically inserted
164above the definition, and the a plot of the model is automatically inserted
165before the **reference**.
166
167Figures can be included using the *figure* command, with the name
168of the *.png* file containing the figure and a caption to appear below the
169figure.  Figure numbers will be added automatically.
170
171See this `Sphinx cheat sheet <http://matplotlib.org/sampledoc/cheatsheet.html>`_
172for a quick guide to the documentation layout commands, or the
173`Sphinx Documentation <http://www.sphinx-doc.org/en/stable/>`_ for
174complete details.
175
176The model should include a **formula** written using LaTeX markup.
177The example above uses the *math* command to make a displayed equation.  You
178can also use *\$formula\$* for an inline formula. This is handy for defining
179the relationship between the model parameters and formula variables, such
180as the phrase "\$r\$ is the radius" used above.  The live demo MathJax
181page `<http://www.mathjax.org/>`_ is handy for checking that the equations
182will look like you intend.
183
184Math layout uses the `amsmath <http://www.ams.org/publications/authors/tex/amslatex>`_
185package for aligning equations (see amsldoc.pdf on that page for complete
186documentation). You will automatically be in an aligned environment, with
187blank lines separating the lines of the equation.  Place an ampersand before
188the operator on which to align.  For example::
189
190    .. math::
191
192      x + y &= 1 \\
193      y &= x - 1
194
195produces
196
197.. math::
198
199      x + y &= 1 \\
200      y &= x - 1
201
202If you need more control, use::
203
204    .. math::
205        :nowrap:
206
207
208Model Definition
209................
210
211Following the documentation string, there are a series of definitions::
212
213    name = "sphere"  # optional: defaults to the filename without .py
214
215    title = "Spheres with uniform scattering length density"
216
217    description = """\
218    P(q)=(scale/V)*[3V(sld-sld_solvent)*(sin(qr)-qr cos(qr))
219                    /(qr)^3]^2 + background
220        r: radius of sphere
221        V: The volume of the scatter
222        sld: the SLD of the sphere
223        sld_solvent: the SLD of the solvent
224    """
225
226    category = "shape:sphere"
227
228    single = True   # optional: defaults to True
229
230    opencl = False  # optional: defaults to False
231
232    structure_factor = False  # optional: defaults to False
233
234**name = "mymodel"** defines the name of the model that is shown to the user.
[3048ec6]235If it is not provided it will use the name of the model file. The name must
236be a valid variable name, starting with a letter and contains only letters,
237numbers or underscore.  Spaces, dashes, and other symbols are not permitted.
[990d8df]238
239**title = "short description"** is short description of the model which
240is included after the model name in the automatically generated documentation.
241The title can also be used for a tooltip.
242
243**description = """doc string"""** is a longer description of the model. It
244shows up when you press the "Description" button of the SasView FitPage.
245It should give a brief description of the equation and the parameters
246without the need to read the entire model documentation. The triple quotes
247allow you to write the description over multiple lines. Keep the lines
248short since the GUI will wrap each one separately if they are too long.
249**Make sure the parameter names in the description match the model definition!**
250
251**category = "shape:sphere"** defines where the model will appear in the
252model documentation.  In this example, the model will appear alphabetically
253in the list of spheroid models in the *Shape* category.
254
255**single = True** indicates that the model can be run using single
256precision floating point values.  Set it to False if the numerical
257calculation for the model is unstable, which is the case for about 20 of
258the built in models.  It is worthwhile modifying the calculation to support
259single precision, allowing models to run up to 10 times faster.  The
260section `Test_Your_New_Model`_  describes how to compare model values for
261single vs. double precision so you can decide if you need to set
262single to False.
263
264**opencl = False** indicates that the model should not be run using OpenCL.
265This may be because the model definition includes code that cannot be
266compiled for the GPU (for example, goto statements).  It can also be used
267for large models which can't run on most GPUs.  This flag has not been
268used on any of the built in models; models which were failing were
269streamlined so this flag was not necessary.
270
271**structure_factor = True** indicates that the model can be used as a
272structure factor to account for interactions between particles.  See
273`Form_Factors`_ for more details.
274
275Model Parameters
276................
277
278Next comes the parameter table.  For example::
279
280    # pylint: disable=bad-whitespace, line-too-long
281    #   ["name",        "units", default, [min, max], "type",    "description"],
282    parameters = [
283        ["sld",         "1e-6/Ang^2",  1, [-inf, inf], "sld",    "Layer scattering length density"],
284        ["sld_solvent", "1e-6/Ang^2",  6, [-inf, inf], "sld",    "Solvent scattering length density"],
285        ["radius",      "Ang",        50, [0, inf],    "volume", "Sphere radius"],
286    ]
287    # pylint: enable=bad-whitespace, line-too-long
288
289**parameters = [["name", "units", default, [min,max], "type", "tooltip"],...]**
290defines the parameters that form the model.
291
292**Note: The order of the parameters in the definition will be the order of the
293parameters in the user interface and the order of the parameters in Iq(),
[108e70e]294Iqac(), Iqabc() and form_volume(). And** *scale* **and** *background*
295**parameters are implicit to all models, so they do not need to be included
296in the parameter table.**
[990d8df]297
298- **"name"** is the name of the parameter shown on the FitPage.
299
[3048ec6]300  - the name must be a valid variable name, starting with a letter and
301    containing only letters, numbers and underscore.
302
[990d8df]303  - parameter names should follow the mathematical convention; e.g.,
304    *radius_core* not *core_radius*, or *sld_solvent* not *solvent_sld*.
305
306  - model parameter names should be consistent between different models,
307    so *sld_solvent*, for example, should have exactly the same name
308    in every model.
309
310  - to see all the parameter names currently in use, type the following in the
311    python shell/editor under the Tools menu::
312
313       import sasmodels.list_pars
314       sasmodels.list_pars.list_pars()
315
316    *re-use* as many as possible!!!
317
318  - use "name[n]" for multiplicity parameters, where *n* is the name of
319    the parameter defining the number of shells/layers/segments, etc.
320
321- **"units"** are displayed along with the parameter name
322
323  - every parameter should have units; use "None" if there are no units.
324
325  - **sld's should be given in units of 1e-6/Ang^2, and not simply
326    1/Ang^2 to be consistent with the builtin models.  Adjust your formulas
327    appropriately.**
328
329  - fancy units markup is available for some units, including::
330
331        Ang, 1/Ang, 1/Ang^2, 1e-6/Ang^2, degrees, 1/cm, Ang/cm, g/cm^3, mg/m^2
332
333  - the list of units is defined in the variable *RST_UNITS* within
334    `sasmodels/generate.py <https://github.com/SasView/sasmodels/tree/master/sasmodels/generate.py>`_
335
336    - new units can be added using the macros defined in *doc/rst_prolog*
337      in the sasmodels source.
338    - units should be properly formatted using sub-/super-scripts
339      and using negative exponents instead of the / operator, though
340      the unit name should use the / operator for consistency.
341    - please post a message to the SasView developers mailing list with your changes.
342
343- **default** is the initial value for the parameter.
344
345  - **the parameter default values are used to auto-generate a plot of
346    the model function in the documentation.**
347
348- **[min, max]** are the lower and upper limits on the parameter.
349
350  - lower and upper limits can be any number, or *-inf* or *inf*.
351
352  - the limits will show up as the default limits for the fit making it easy,
353    for example, to force the radius to always be greater than zero.
354
355  - these are hard limits defining the valid range of parameter values;
356    polydisperity distributions will be truncated at the limits.
357
358- **"type"** can be one of: "", "sld", "volume", or "orientation".
359
360  - "sld" parameters can have magnetic moments when fitting magnetic models;
361    depending on the spin polarization of the beam and the $q$ value being
362    examined, the effective sld for that material will be used to compute the
363    scattered intensity.
364
[108e70e]365  - "volume" parameters are passed to Iq(), Iqac(), Iqabc() and form_volume(),
366    and have polydispersity loops generated automatically.
[990d8df]367
[108e70e]368  - "orientation" parameters are not passed, but instead are combined with
369    orientation dispersity to translate *qx* and *qy* to *qa*, *qb* and *qc*.
370    These parameters should appear at the end of the table with the specific
371    names *theta*, *phi* and for asymmetric shapes *psi*, in that order.
[990d8df]372
[9844c3a]373Some models will have integer parameters, such as number of pearls in the
374pearl necklace model, or number of shells in the multi-layer vesicle model.
375The optimizers in BUMPS treat all parameters as floating point numbers which
376can take arbitrary values, even for integer parameters, so your model should
377round the incoming parameter value to the nearest integer inside your model
378you should round to the nearest integer.  In C code, you can do this using::
379
380    static double
381    Iq(double q, ..., double fp_n, ...)
382    {
383        int n = (int)(fp_n + 0.5);
384        ...
385    }
386
387in python::
388
389    def Iq(q, ..., n, ...):
390        n = int(n+0.5)
391        ...
392
[3048ec6]393Derivative based optimizers such as Levenberg-Marquardt will not work
[9844c3a]394for integer parameters since the partial derivative is always zero, but
395the remaining optimizers (DREAM, differential evolution, Nelder-Mead simplex)
396will still function.
397
[990d8df]398Model Computation
399.................
400
401Models can be defined as pure python models, or they can be a mixture of
402python and C models.  C models are run on the GPU if it is available,
403otherwise they are compiled and run on the CPU.
404
405Models are defined by the scattering kernel, which takes a set of parameter
406values defining the shape, orientation and material, and returns the
407expected scattering. Polydispersity and angular dispersion are defined
408by the computational infrastructure.  Any parameters defined as "volume"
409parameters are polydisperse, with polydispersity defined in proportion
410to their value.  "orientation" parameters use angular dispersion defined
411in degrees, and are not relative to the current angle.
412
413Based on a weighting function $G(x)$ and a number of points $n$, the
414computed value is
415
416.. math::
417
418     \hat I(q)
419     = \frac{\int G(x) I(q, x)\,dx}{\int G(x) V(x)\,dx}
420     \approx \frac{\sum_{i=1}^n G(x_i) I(q,x_i)}{\sum_{i=1}^n G(x_i) V(x_i)}
421
[3048ec6]422That is, the individual models do not need to include polydispersity
[990d8df]423calculations, but instead rely on numerical integration to compute the
[108e70e]424appropriately smeared pattern.
[990d8df]425
426Python Models
427.............
428
429For pure python models, define the *Iq* function::
430
431      import numpy as np
432      from numpy import cos, sin, ...
433
434      def Iq(q, par1, par2, ...):
435          return I(q, par1, par2, ...)
436      Iq.vectorized = True
437
438The parameters *par1, par2, ...* are the list of non-orientation parameters
439to the model in the order that they appear in the parameter table.
[3048ec6]440**Note that the auto-generated model file uses** *x* **rather than** *q*.
[990d8df]441
442The *.py* file should import trigonometric and exponential functions from
443numpy rather than from math.  This lets us evaluate the model for the whole
444range of $q$ values at once rather than looping over each $q$ separately in
445python.  With $q$ as a vector, you cannot use if statements, but must instead
446do tricks like
447
448::
449
450     a = x*q*(q>0) + y*q*(q<=0)
451
452or
453
454::
455
456     a = np.empty_like(q)
457     index = q>0
458     a[index] = x*q[index]
459     a[~index] = y*q[~index]
460
461which sets $a$ to $q \cdot x$ if $q$ is positive or $q \cdot y$ if $q$
462is zero or negative. If you have not converted your function to use $q$
463vectors, you can set the following and it will only receive one $q$
464value at a time::
465
466    Iq.vectorized = False
467
468Return np.NaN if the parameters are not valid (e.g., cap_radius < radius in
469barbell).  If I(q; pars) is NaN for any $q$, then those parameters will be
470ignored, and not included in the calculation of the weighted polydispersity.
471
472Models should define *form_volume(par1, par2, ...)* where the parameter
473list includes the *volume* parameters in order.  This is used for a weighted
474volume normalization so that scattering is on an absolute scale.  If
475*form_volume* is not defined, then the default *form_volume = 1.0* will be
476used.
477
478Embedded C Models
479.................
480
481Like pure python models, inline C models need to define an *Iq* function::
482
483    Iq = """
484        return I(q, par1, par2, ...);
485    """
486
487This expands into the equivalent C code::
488
489    #include <math.h>
490    double Iq(double q, double par1, double par2, ...);
491    double Iq(double q, double par1, double par2, ...)
492    {
493        return I(q, par1, par2, ...);
494    }
495
496*form_volume* defines the volume of the shape. As in python models, it
497includes only the volume parameters.
498
499**source=['fn.c', ...]** includes the listed C source files in the
[108e70e]500program before *Iq* and *form_volume* are defined. This allows you to
[ef85a09]501extend the library of C functions available to your model.
502
503*c_code* includes arbitrary C code into your kernel, which can be
504handy for defining helper functions for *Iq* and *form_volume*. Note that
[108e70e]505you can put the full function definition for *Iq* and *form_volume*
[ef85a09]506(include function declaration) into *c_code* as well, or put them into an
507external C file and add that file to the list of sources.
[990d8df]508
509Models are defined using double precision declarations for the
510parameters and return values.  When a model is run using single
511precision or long double precision, each variable is converted
512to the target type, depending on the precision requested.
513
514**Floating point constants must include the decimal point.**  This allows us
515to convert values such as 1.0 (double precision) to 1.0f (single precision)
516so that expressions that use these values are not promoted to double precision
517expressions.  Some graphics card drivers are confused when functions
518that expect floating point values are passed integers, such as 4*atan(1); it
519is safest to not use integers in floating point expressions.  Even better,
520use the builtin constant M_PI rather than 4*atan(1); it is faster and smaller!
521
522The C model operates on a single $q$ value at a time.  The code will be
523run in parallel across different $q$ values, either on the graphics card
524or the processor.
525
526Rather than returning NAN from Iq, you must define the *INVALID(v)*.  The
527*v* parameter lets you access all the parameters in the model using
528*v.par1*, *v.par2*, etc. For example::
529
530    #define INVALID(v) (v.bell_radius < v.radius)
531
[ef85a09]532The INVALID define can go into *Iq*, or *c_code*, or an external C file
533listed in *source*.
534
[108e70e]535Oriented Shapes
536...............
537
538If the scattering is dependent on the orientation of the shape, then you
539will need to include *orientation* parameters *theta*, *phi* and *psi*
[7e6bc45e]540at the end of the parameter table.  As described in the section
541:ref:`orientation`, the individual $(q_x, q_y)$ points on the detector will
542be rotated into $(q_a, q_b, q_c)$ points relative to the sample in its
543canonical orientation with $a$-$b$-$c$ aligned with $x$-$y$-$z$ in the
544laboratory frame and beam travelling along $-z$.
545
546The oriented C model is called using *Iqabc(qa, qb, qc, par1, par2, ...)* where
[108e70e]547*par1*, etc. are the parameters to the model.  If the shape is rotationally
548symmetric about *c* then *psi* is not needed, and the model is called
549as *Iqac(qab, qc, par1, par2, ...)*.  In either case, the orientation
550parameters are not included in the function call.
551
552For 1D oriented shapes, an integral over all angles is usually needed for
[b85227d]553the *Iq* function. Given symmetry and the substitution $u = \cos(\alpha)$,
[108e70e]554$du = -\sin(\alpha)\,d\alpha$ this becomes
555
556.. math::
557
[b85227d]558    I(q) &= \frac{1}{4\pi} \int_{-\pi/2}^{pi/2} \int_{-pi}^{pi}
559            F(q_a, q_b, q_c)^2 \sin(\alpha)\,d\beta\,d\alpha \\
560        &= \frac{8}{4\pi} \int_{0}^{pi/2} \int_{0}^{\pi/2}
561            F^2 \sin(\alpha)\,d\beta\,d\alpha \\
562        &= \frac{8}{4\pi} \int_1^0 \int_{0}^{\pi/2} - F^2 \,d\beta\,du \\
563        &= \frac{8}{4\pi} \int_0^1 \int_{0}^{\pi/2} F^2 \,d\beta\,du
564
565for
566
567.. math::
568
569    q_a &= q \sin(\alpha)\sin(\beta) = q \sqrt{1-u^2} \sin(\beta) \\
570    q_b &= q \sin(\alpha)\cos(\beta) = q \sqrt{1-u^2} \cos(\beta) \\
571    q_c &= q \cos(\alpha) = q u
[108e70e]572
573Using the $z, w$ values for Gauss-Legendre integration in "lib/gauss76.c", the
574numerical integration is then::
575
576    double outer_sum = 0.0;
577    for (int i = 0; i < GAUSS_N; i++) {
578        const double cos_alpha = 0.5*GAUSS_Z[i] + 0.5;
579        const double sin_alpha = sqrt(1.0 - cos_alpha*cos_alpha);
580        const double qc = cos_alpha * q;
581        double inner_sum = 0.0;
582        for (int j = 0; j < GAUSS_N; j++) {
583            const double beta = M_PI_4 * GAUSS_Z[j] + M_PI_4;
584            double sin_beta, cos_beta;
585            SINCOS(beta, sin_beta, cos_beta);
586            const double qa = sin_alpha * sin_beta * q;
[b85227d]587            const double qb = sin_alpha * cos_beta * q;
588            const double form = Fq(qa, qb, qc, ...);
589            inner_sum += GAUSS_W[j] * form * form;
[108e70e]590        }
591        outer_sum += GAUSS_W[i] * inner_sum;
592    }
593    outer_sum *= 0.25; // = 8/(4 pi) * outer_sum * (pi/2) / 4
594
595The *z* values for the Gauss-Legendre integration extends from -1 to 1, so
596the double sum of *w[i]w[j]* explains the factor of 4.  Correcting for the
597average *dz[i]dz[j]* gives $(1-0) \cdot (\pi/2-0) = \pi/2$.  The $8/(4 \pi)$
598factor comes from the integral over the quadrant.  With less symmetry (eg.,
599in the bcc and fcc paracrystal models), then an integral over the entire
600sphere may be necessary.
601
602For simpler models which are rotationally symmetric a single integral
603suffices:
604
605.. math::
606
[b85227d]607    I(q) &= \frac{1}{\pi}\int_{-\pi/2}^{\pi/2}
608            F(q_{ab}, q_c)^2 \sin(\alpha)\,d\alpha/\pi \\
609        &= \frac{2}{\pi} \int_0^1 F^2\,du
610
611for
612
613.. math::
614
615    q_{ab} &= q \sin(\alpha) = q \sqrt{1 - u^2} \\
616    q_c &= q \cos(\alpha) = q u
617
[108e70e]618
619with integration loop::
620
621    double sum = 0.0;
622    for (int i = 0; i < GAUSS_N; i++) {
623        const double cos_alpha = 0.5*GAUSS_Z[i] + 0.5;
624        const double sin_alpha = sqrt(1.0 - cos_alpha*cos_alpha);
625        const double qab = sin_alpha * q;
[b85227d]626        const double qc = cos_alpha * q;
627        const double form = Fq(qab, qc, ...);
628        sum += GAUSS_W[j] * form * form;
[108e70e]629    }
630    sum *= 0.5; // = 2/pi * sum * (pi/2) / 2
631
632Magnetism
633.........
634
635Magnetism is supported automatically for all shapes by modifying the
636effective SLD of particle according to the Halpern-Johnson vector
[c654160]637describing the interaction between neutron spin and magnetic field.  All
[108e70e]638parameters marked as type *sld* in the parameter table are treated as
639possibly magnetic particles with magnitude *M0* and direction
640*mtheta* and *mphi*.  Polarization parameters are also provided
641automatically for magnetic models to set the spin state of the measurement.
642
643For more complicated systems where magnetism is not uniform throughout
644the individual particles, you will need to write your own models.
645You should not mark the nuclear sld as type *sld*, but instead leave
646them unmarked and provide your own magnetism and polarization parameters.
647For 2D measurements you will need $(q_x, q_y)$ values for the measurement
648to compute the proper magnetism and orientation, which you can implement
649using *Iqxy(qx, qy, par1, par2, ...)*.
650
[990d8df]651Special Functions
652.................
653
654The C code follows the C99 standard, with the usual math functions,
655as defined in
656`OpenCL <https://www.khronos.org/registry/cl/sdk/1.1/docs/man/xhtml/mathFunctions.html>`_.
657This includes the following:
658
659    M_PI, M_PI_2, M_PI_4, M_SQRT1_2, M_E:
660        $\pi$, $\pi/2$, $\pi/4$, $1/\sqrt{2}$ and Euler's constant $e$
[d0dc9a3]661    exp, log, pow(x,y), expm1, log1p, sqrt, cbrt:
662        Power functions $e^x$, $\ln x$, $x^y$, $e^x - 1$, $\ln 1 + x$,
663        $\sqrt{x}$, $\sqrt[3]{x}$. The functions expm1(x) and log1p(x)
664        are accurate across all $x$, including $x$ very close to zero.
[990d8df]665    sin, cos, tan, asin, acos, atan:
666        Trigonometry functions and inverses, operating on radians.
667    sinh, cosh, tanh, asinh, acosh, atanh:
668        Hyperbolic trigonometry functions.
669    atan2(y,x):
670        Angle from the $x$\ -axis to the point $(x,y)$, which is equal to
671        $\tan^{-1}(y/x)$ corrected for quadrant.  That is, if $x$ and $y$ are
672        both negative, then atan2(y,x) returns a value in quadrant III where
673        atan(y/x) would return a value in quadrant I. Similarly for
674        quadrants II and IV when $x$ and $y$ have opposite sign.
[d0dc9a3]675    fabs(x), fmin(x,y), fmax(x,y), trunc, rint:
[990d8df]676        Floating point functions.  rint(x) returns the nearest integer.
677    NAN:
678        NaN, Not a Number, $0/0$.  Use isnan(x) to test for NaN.  Note that
679        you cannot use :code:`x == NAN` to test for NaN values since that
[d0dc9a3]680        will always return false.  NAN does not equal NAN!  The alternative,
681        :code:`x != x` may fail if the compiler optimizes the test away.
[990d8df]682    INFINITY:
683        $\infty, 1/0$.  Use isinf(x) to test for infinity, or isfinite(x)
684        to test for finite and not NaN.
685    erf, erfc, tgamma, lgamma:  **do not use**
686        Special functions that should be part of the standard, but are missing
687        or inaccurate on some platforms. Use sas_erf, sas_erfc and sas_gamma
688        instead (see below). Note: lgamma(x) has not yet been tested.
689
690Some non-standard constants and functions are also provided:
691
692    M_PI_180, M_4PI_3:
693        $\frac{\pi}{180}$, $\frac{4\pi}{3}$
694    SINCOS(x, s, c):
695        Macro which sets s=sin(x) and c=cos(x). The variables *c* and *s*
696        must be declared first.
697    square(x):
698        $x^2$
699    cube(x):
700        $x^3$
701    sas_sinx_x(x):
702        $\sin(x)/x$, with limit $\sin(0)/0 = 1$.
703    powr(x, y):
704        $x^y$ for $x \ge 0$; this is faster than general $x^y$ on some GPUs.
705    pown(x, n):
706        $x^n$ for $n$ integer; this is faster than general $x^n$ on some GPUs.
707    FLOAT_SIZE:
708        The number of bytes in a floating point value.  Even though all
709        variables are declared double, they may be converted to single
710        precision float before running. If your algorithm depends on
711        precision (which is not uncommon for numerical algorithms), use
712        the following::
713
714            #if FLOAT_SIZE>4
715            ... code for double precision ...
716            #else
717            ... code for single precision ...
718            #endif
719    SAS_DOUBLE:
720        A replacement for :code:`double` so that the declared variable will
721        stay double precision; this should generally not be used since some
722        graphics cards do not support double precision.  There is no provision
723        for forcing a constant to stay double precision.
724
725The following special functions and scattering calculations are defined in
726`sasmodels/models/lib <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib>`_.
727These functions have been tuned to be fast and numerically stable down
728to $q=0$ even in single precision.  In some cases they work around bugs
729which appear on some platforms but not others, so use them where needed.
730Add the files listed in :code:`source = ["lib/file.c", ...]` to your *model.py*
731file in the order given, otherwise these functions will not be available.
732
733    polevl(x, c, n):
734        Polynomial evaluation $p(x) = \sum_{i=0}^n c_i x^i$ using Horner's
735        method so it is faster and more accurate.
736
737        $c = \{c_n, c_{n-1}, \ldots, c_0 \}$ is the table of coefficients,
738        sorted from highest to lowest.
739
740        :code:`source = ["lib/polevl.c", ...]` (`link to code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/polevl.c>`_)
741
742    p1evl(x, c, n):
743        Evaluation of normalized polynomial $p(x) = x^n + \sum_{i=0}^{n-1} c_i x^i$
744        using Horner's method so it is faster and more accurate.
745
746        $c = \{c_{n-1}, c_{n-2} \ldots, c_0 \}$ is the table of coefficients,
747        sorted from highest to lowest.
748
749        :code:`source = ["lib/polevl.c", ...]`
[870a2f4]750        (`polevl.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/polevl.c>`_)
[990d8df]751
752    sas_gamma(x):
[30b60d2]753        Gamma function sas_gamma\ $(x) = \Gamma(x)$.
[990d8df]754
755        The standard math function, tgamma(x) is unstable for $x < 1$
756        on some platforms.
757
[870a2f4]758        :code:`source = ["lib/sas_gamma.c", ...]`
759        (`sas_gamma.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_gamma.c>`_)
[990d8df]760
761    sas_erf(x), sas_erfc(x):
762        Error function
[30b60d2]763        sas_erf\ $(x) = \frac{2}{\sqrt\pi}\int_0^x e^{-t^2}\,dt$
[990d8df]764        and complementary error function
[30b60d2]765        sas_erfc\ $(x) = \frac{2}{\sqrt\pi}\int_x^{\infty} e^{-t^2}\,dt$.
[990d8df]766
767        The standard math functions erf(x) and erfc(x) are slower and broken
768        on some platforms.
769
770        :code:`source = ["lib/polevl.c", "lib/sas_erf.c", ...]`
[870a2f4]771        (`sas_erf.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_erf.c>`_)
[990d8df]772
773    sas_J0(x):
[30b60d2]774        Bessel function of the first kind sas_J0\ $(x)=J_0(x)$ where
[990d8df]775        $J_0(x) = \frac{1}{\pi}\int_0^\pi \cos(x\sin(\tau))\,d\tau$.
776
777        The standard math function j0(x) is not available on all platforms.
778
779        :code:`source = ["lib/polevl.c", "lib/sas_J0.c", ...]`
[870a2f4]780        (`sas_J0.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_J0.c>`_)
[990d8df]781
782    sas_J1(x):
[30b60d2]783        Bessel function of the first kind  sas_J1\ $(x)=J_1(x)$ where
[990d8df]784        $J_1(x) = \frac{1}{\pi}\int_0^\pi \cos(\tau - x\sin(\tau))\,d\tau$.
785
786        The standard math function j1(x) is not available on all platforms.
787
788        :code:`source = ["lib/polevl.c", "lib/sas_J1.c", ...]`
[870a2f4]789        (`sas_J1.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_J1.c>`_)
[990d8df]790
791    sas_JN(n, x):
[30b60d2]792        Bessel function of the first kind and integer order $n$,
793        sas_JN\ $(n, x) =J_n(x)$ where
[990d8df]794        $J_n(x) = \frac{1}{\pi}\int_0^\pi \cos(n\tau - x\sin(\tau))\,d\tau$.
[30b60d2]795        If $n$ = 0 or 1, it uses sas_J0($x$) or sas_J1($x$), respectively.
[990d8df]796
797        The standard math function jn(n, x) is not available on all platforms.
798
799        :code:`source = ["lib/polevl.c", "lib/sas_J0.c", "lib/sas_J1.c", "lib/sas_JN.c", ...]`
[870a2f4]800        (`sas_JN.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_JN.c>`_)
[990d8df]801
802    sas_Si(x):
[30b60d2]803        Sine integral Si\ $(x) = \int_0^x \tfrac{\sin t}{t}\,dt$.
[990d8df]804
805        This function uses Taylor series for small and large arguments:
806
807        For large arguments,
808
809        .. math::
810
811             \text{Si}(x) \sim \frac{\pi}{2}
812             - \frac{\cos(x)}{x}\left(1 - \frac{2!}{x^2} + \frac{4!}{x^4} - \frac{6!}{x^6} \right)
813             - \frac{\sin(x)}{x}\left(\frac{1}{x} - \frac{3!}{x^3} + \frac{5!}{x^5} - \frac{7!}{x^7}\right)
814
815        For small arguments,
816
817        .. math::
818
819           \text{Si}(x) \sim x
820           - \frac{x^3}{3\times 3!} + \frac{x^5}{5 \times 5!} - \frac{x^7}{7 \times 7!}
821           + \frac{x^9}{9\times 9!} - \frac{x^{11}}{11\times 11!}
822
823        :code:`source = ["lib/Si.c", ...]`
[870a2f4]824        (`Si.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/Si.c>`_)
[990d8df]825
826    sas_3j1x_x(x):
827        Spherical Bessel form
[30b60d2]828        sph_j1c\ $(x) = 3 j_1(x)/x = 3 (\sin(x) - x \cos(x))/x^3$,
[990d8df]829        with a limiting value of 1 at $x=0$, where $j_1(x)$ is the spherical
830        Bessel function of the first kind and first order.
831
832        This function uses a Taylor series for small $x$ for numerical accuracy.
833
834        :code:`source = ["lib/sas_3j1x_x.c", ...]`
[870a2f4]835        (`sas_3j1x_x.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_3j1x_x.c>`_)
[990d8df]836
837
838    sas_2J1x_x(x):
[30b60d2]839        Bessel form sas_J1c\ $(x) = 2 J_1(x)/x$, with a limiting value
[990d8df]840        of 1 at $x=0$, where $J_1(x)$ is the Bessel function of first kind
841        and first order.
842
843        :code:`source = ["lib/polevl.c", "lib/sas_J1.c", ...]`
[870a2f4]844        (`sas_J1.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_J1.c>`_)
[990d8df]845
846
847    Gauss76Z[i], Gauss76Wt[i]:
848        Points $z_i$ and weights $w_i$ for 76-point Gaussian quadrature, respectively,
849        computing $\int_{-1}^1 f(z)\,dz \approx \sum_{i=1}^{76} w_i\,f(z_i)$.
850
851        Similar arrays are available in :code:`gauss20.c` for 20-point
852        quadrature and in :code:`gauss150.c` for 150-point quadrature.
[d0dc9a3]853        The macros :code:`GAUSS_N`, :code:`GAUSS_Z` and :code:`GAUSS_W` are
854        defined so that you can change the order of the integration by
855        selecting an different source without touching the C code.
[990d8df]856
857        :code:`source = ["lib/gauss76.c", ...]`
[870a2f4]858        (`gauss76.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/gauss76.c>`_)
[990d8df]859
860
861
862Problems with C models
863......................
864
865The graphics processor (GPU) in your computer is a specialized computer tuned
866for certain kinds of problems.  This leads to strange restrictions that you
867need to be aware of.  Your code may work fine on some platforms or for some
868models, but then return bad values on other platforms.  Some examples of
869particular problems:
870
871  **(1) Code is too complex, or uses too much memory.** GPU devices only
872  have a limited amount of memory available for each processor. If you run
873  programs which take too much memory, then rather than running multiple
874  values in parallel as it usually does, the GPU may only run a single
875  version of the code at a time, making it slower than running on the CPU.
876  It may fail to run on some platforms, or worse, cause the screen to go
877  blank or the system to reboot.
878
879  **(2) Code takes too long.** Because GPU devices are used for the computer
880  display, the OpenCL drivers are very careful about the amount of time they
881  will allow any code to run. For example, on OS X, the model will stop
882  running after 5 seconds regardless of whether the computation is complete.
883  You may end up with only some of your 2D array defined, with the rest
884  containing random data. Or it may cause the screen to go blank or the
885  system to reboot.
886
887  **(3) Memory is not aligned**. The GPU hardware is specialized to operate
888  on multiple values simultaneously. To keep the GPU simple the values in
889  memory must be aligned with the different GPU compute engines. Not
890  following these rules can lead to unexpected values being loaded into
891  memory, and wrong answers computed. The conclusion from a very long and
892  strange debugging session was that any arrays that you declare in your
893  model should be a multiple of four. For example::
894
895      double Iq(q, p1, p2, ...)
896      {
897          double vector[8];  // Only going to use seven slots, but declare 8
898          ...
899      }
900
901The first step when your model is behaving strangely is to set
902**single=False**. This automatically restricts the model to only run on the
903CPU, or on high-end GPU cards. There can still be problems even on high-end
904cards, so you can force the model off the GPU by setting **opencl=False**.
905This runs the model as a normal C program without any GPU restrictions so
906you know that strange results are probably from your code rather than the
907environment. Once the code is debugged, you can compare your output to the
908output on the GPU.
909
910Although it can be difficult to get your model to work on the GPU, the reward
911can be a model that runs 1000x faster on a good card.  Even your laptop may
912show a 50x improvement or more over the equivalent pure python model.
913
914
915.. _Form_Factors:
916
917Form Factors
918............
919
920Away from the dilute limit you can estimate scattering including
921particle-particle interactions using $I(q) = P(q)*S(q)$ where $P(q)$
922is the form factor and $S(q)$ is the structure factor.  The simplest
923structure factor is the *hardsphere* interaction, which
924uses the effective radius of the form factor as an input to the structure
925factor model.  The effective radius is the average radius of the
926form averaged over all the polydispersity values.
927
928::
929
930    def ER(radius, thickness):
931        """Effective radius of a core-shell sphere."""
932        return radius + thickness
933
934Now consider the *core_shell_sphere*, which has a simple effective radius
935equal to the radius of the core plus the thickness of the shell, as
936shown above. Given polydispersity over *(r1, r2, ..., rm)* in radius and
937*(t1, t2, ..., tn)* in thickness, *ER* is called with a mesh
938grid covering all possible combinations of radius and thickness.
939That is, *radius* is *(r1, r2, ..., rm, r1, r2, ..., rm, ...)*
940and *thickness* is *(t1, t1, ... t1, t2, t2, ..., t2, ...)*.
941The *ER* function returns one effective radius for each combination.
942The effective radius calculator weights each of these according to
943the polydispersity distributions and calls the structure factor
944with the average *ER*.
945
946::
947
948    def VR(radius, thickness):
949        """Sphere and shell volumes for a core-shell sphere."""
950        whole = 4.0/3.0 * pi * (radius + thickness)**3
951        core = 4.0/3.0 * pi * radius**3
952        return whole, whole - core
953
954Core-shell type models have an additional volume ratio which scales
955the structure factor.  The *VR* function returns the volume of
956the whole sphere and the volume of the shell. Like *ER*, there is
957one return value for each point in the mesh grid.
958
959*NOTE: we may be removing or modifying this feature soon. As of the
960time of writing, core-shell sphere returns (1., 1.) for VR, giving a volume
961ratio of 1.0.*
962
963Unit Tests
964..........
965
966THESE ARE VERY IMPORTANT. Include at least one test for each model and
967PLEASE make sure that the answer value is correct (i.e. not a random number).
968
969::
970
971    tests = [
972        [{}, 0.2, 0.726362],
973        [{"scale": 1., "background": 0., "sld": 6., "sld_solvent": 1.,
974          "radius": 120., "radius_pd": 0.2, "radius_pd_n":45},
975         0.2, 0.228843],
976        [{"radius": 120., "radius_pd": 0.2, "radius_pd_n":45}, "ER", 120.],
977        [{"radius": 120., "radius_pd": 0.2, "radius_pd_n":45}, "VR", 1.],
978    ]
979
980
981**tests=[[{parameters}, q, result], ...]** is a list of lists.
982Each list is one test and contains, in order:
983
984- a dictionary of parameter values. This can be *{}* using the default
985  parameters, or filled with some parameters that will be different from the
986  default, such as *{"radius":10.0, "sld":4}*. Unlisted parameters will
987  be given the default values.
988- the input $q$ value or tuple of $(q_x, q_y)$ values.
989- the output $I(q)$ or $I(q_x,q_y)$ expected of the model for the parameters
990  and input value given.
991- input and output values can themselves be lists if you have several
992  $q$ values to test for the same model parameters.
993- for testing *ER* and *VR*, give the inputs as "ER" and "VR" respectively;
994  the output for *VR* should be the sphere/shell ratio, not the individual
995  sphere and shell values.
996
997.. _Test_Your_New_Model:
998
999Test Your New Model
1000^^^^^^^^^^^^^^^^^^^
1001
1002Minimal Testing
1003...............
1004
1005From SasView either open the Python shell (*Tools* > *Python Shell/Editor*)
1006or the plugin editor (*Fitting* > *Plugin Model Operations* > *Advanced
1007Plugin Editor*), load your model, and then select *Run > Check Model* from
1008the menu bar. An *Info* box will appear with the results of the compilation
1009and a check that the model runs.
1010
1011If you are not using sasmodels from SasView, skip this step.
1012
1013Recommended Testing
1014...................
1015
1016If the model compiles and runs, you can next run the unit tests that
1017you have added using the **test =** values.
1018
1019From SasView, switch to the *Shell* tab and type the following::
1020
1021    from sasmodels.model_test import run_one
1022    run_one("~/.sasview/plugin_models/model.py")
1023
1024This should print::
1025
1026    test_model_python (sasmodels.model_test.ModelTestCase) ... ok
1027
1028To check whether single precision is good enough, type the following::
1029
1030    from sasmodels.compare import main as compare
1031    compare("~/.sasview/plugin_models/model.py")
1032
1033This will pop up a plot showing the difference between single precision
1034and double precision on a range of $q$ values.
1035
1036::
1037
1038  demo = dict(scale=1, background=0,
1039              sld=6, sld_solvent=1,
1040              radius=120,
1041              radius_pd=.2, radius_pd_n=45)
1042
1043**demo={'par': value, ...}** in the model file sets the default values for
1044the comparison. You can include polydispersity parameters such as
1045*radius_pd=0.2, radius_pd_n=45* which would otherwise be zero.
1046
1047These commands can also be run directly in the python interpreter:
1048
1049    $ python -m sasmodels.model_test -v ~/.sasview/plugin_models/model.py
1050    $ python -m sasmodels.compare ~/.sasview/plugin_models/model.py
1051
1052The options to compare are quite extensive; type the following for help::
1053
1054    compare()
1055
1056Options will need to be passed as separate strings.
1057For example to run your model with a random set of parameters::
1058
1059    compare("-random", "-pars", "~/.sasview/plugin_models/model.py")
1060
1061For the random models,
1062
1063- *sld* will be in the range (-0.5,10.5),
1064- angles (*theta, phi, psi*) will be in the range (-180,180),
1065- angular dispersion will be in the range (0,45),
1066- polydispersity will be in the range (0,1)
1067- other values will be in the range (0, 2\ *v*), where *v* is the value
1068  of the parameter in demo.
1069
1070Dispersion parameters *n*\, *sigma* and *type* will be unchanged from
1071demo so that run times are more predictable (polydispersity calculated
1072across multiple parameters can be very slow).
1073
[3048ec6]1074If your model has 2D orientation calculation, then you should also
[990d8df]1075test with::
1076
1077    compare("-2d", "~/.sasview/plugin_models/model.py")
1078
1079Check The Docs
1080^^^^^^^^^^^^^^
1081
1082You can get a rough idea of how the documentation will look using the
1083following::
1084
1085    compare("-help", "~/.sasview/plugin_models/model.py")
1086
1087This does not use the same styling as the rest of the docs, but it will
1088allow you to check that your ReStructuredText and LaTeX formatting.
1089Here are some tools to help with the inevitable syntax errors:
1090
1091- `Sphinx cheat sheet <http://matplotlib.org/sampledoc/cheatsheet.html>`_
1092- `Sphinx Documentation <http://www.sphinx-doc.org/en/stable/>`_
1093- `MathJax <http://www.mathjax.org/>`_
1094- `amsmath <http://www.ams.org/publications/authors/tex/amslatex>`_
1095
1096There is also a neat online WYSIWYG ReStructuredText editor at
1097http://rst.ninjs.org\ .
1098
1099
1100Clean Lint - (Developer Version Only)
1101^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1102
1103**NB: For now we are not providing pylint with the installer version
1104of SasView; so unless you have a SasView build environment available,
1105you can ignore this section!**
1106
1107Run the lint check with::
1108
1109    python -m pylint --rcfile=extra/pylint.rc ~/.sasview/plugin_models/model.py
1110
1111We are not aiming for zero lint just yet, only keeping it to a minimum.
1112For now, don't worry too much about *invalid-name*. If you really want a
1113variable name *Rg* for example because $R_g$ is the right name for the model
1114parameter then ignore the lint errors.  Also, ignore *missing-docstring*
[108e70e]1115for standard model functions *Iq*, *Iqac*, etc.
[990d8df]1116
1117We will have delinting sessions at the SasView Code Camps, where we can
1118decide on standards for model files, parameter names, etc.
1119
1120For now, you can tell pylint to ignore things.  For example, to align your
1121parameters in blocks::
1122
1123    # pylint: disable=bad-whitespace,line-too-long
1124    #   ["name",                  "units", default, [lower, upper], "type", "description"],
1125    parameters = [
1126        ["contrast_factor",       "barns",    10.0,  [-inf, inf], "", "Contrast factor of the polymer"],
1127        ["bjerrum_length",        "Ang",       7.1,  [0, inf],    "", "Bjerrum length"],
1128        ["virial_param",          "1/Ang^2",  12.0,  [-inf, inf], "", "Virial parameter"],
1129        ["monomer_length",        "Ang",      10.0,  [0, inf],    "", "Monomer length"],
1130        ["salt_concentration",    "mol/L",     0.0,  [-inf, inf], "", "Concentration of monovalent salt"],
1131        ["ionization_degree",     "",          0.05, [0, inf],    "", "Degree of ionization"],
1132        ["polymer_concentration", "mol/L",     0.7,  [0, inf],    "", "Polymer molar concentration"],
1133        ]
1134    # pylint: enable=bad-whitespace,line-too-long
1135
1136Don't put in too many pylint statements, though, since they make the code ugly.
1137
1138Share Your Model!
1139^^^^^^^^^^^^^^^^^
1140
1141Once compare and the unit test(s) pass properly and everything is done,
1142consider adding your model to the
1143`Model Marketplace <http://marketplace.sasview.org/>`_ so that others may use it!
1144
1145.. ZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZ
1146
1147*Document History*
1148
1149| 2016-10-25 Steve King
[c654160]1150| 2017-05-07 Paul Kienzle - Moved from sasview to sasmodels docs
Note: See TracBrowser for help on using the repository browser.