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

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

Merge branch 'master' into beta_approx

  • Property mode set to 100644
File size: 50.0 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
[9d8a027]275**model_info = ...** lets you define a model directly, for example, by
276loading and modifying existing models.  This is done implicitly by
277:func:`sasmodels.core.load_model_info`, which can create a mixture model
278from a pair of existing models.  For example::
279
280    from sasmodels.core import load_model_info
281    model_info = load_model_info('sphere+cylinder')
282
283See :class:`sasmodels.modelinfo.ModelInfo` for details about the model
284attributes that are defined.
285
[990d8df]286Model Parameters
287................
288
289Next comes the parameter table.  For example::
290
291    # pylint: disable=bad-whitespace, line-too-long
292    #   ["name",        "units", default, [min, max], "type",    "description"],
293    parameters = [
294        ["sld",         "1e-6/Ang^2",  1, [-inf, inf], "sld",    "Layer scattering length density"],
295        ["sld_solvent", "1e-6/Ang^2",  6, [-inf, inf], "sld",    "Solvent scattering length density"],
296        ["radius",      "Ang",        50, [0, inf],    "volume", "Sphere radius"],
297    ]
298    # pylint: enable=bad-whitespace, line-too-long
299
300**parameters = [["name", "units", default, [min,max], "type", "tooltip"],...]**
301defines the parameters that form the model.
302
303**Note: The order of the parameters in the definition will be the order of the
[31fc4ad]304parameters in the user interface and the order of the parameters in Fq(), Iq(),
305Iqac(), Iqabc(), form_volume() and shell_volume().
306And** *scale* **and** *background* **parameters are implicit to all models,
307so they do not need to be included in the parameter table.**
[990d8df]308
309- **"name"** is the name of the parameter shown on the FitPage.
310
[3048ec6]311  - the name must be a valid variable name, starting with a letter and
312    containing only letters, numbers and underscore.
313
[990d8df]314  - parameter names should follow the mathematical convention; e.g.,
315    *radius_core* not *core_radius*, or *sld_solvent* not *solvent_sld*.
316
317  - model parameter names should be consistent between different models,
318    so *sld_solvent*, for example, should have exactly the same name
319    in every model.
320
321  - to see all the parameter names currently in use, type the following in the
322    python shell/editor under the Tools menu::
323
324       import sasmodels.list_pars
325       sasmodels.list_pars.list_pars()
326
327    *re-use* as many as possible!!!
328
329  - use "name[n]" for multiplicity parameters, where *n* is the name of
330    the parameter defining the number of shells/layers/segments, etc.
331
332- **"units"** are displayed along with the parameter name
333
334  - every parameter should have units; use "None" if there are no units.
335
336  - **sld's should be given in units of 1e-6/Ang^2, and not simply
337    1/Ang^2 to be consistent with the builtin models.  Adjust your formulas
338    appropriately.**
339
340  - fancy units markup is available for some units, including::
341
342        Ang, 1/Ang, 1/Ang^2, 1e-6/Ang^2, degrees, 1/cm, Ang/cm, g/cm^3, mg/m^2
343
344  - the list of units is defined in the variable *RST_UNITS* within
345    `sasmodels/generate.py <https://github.com/SasView/sasmodels/tree/master/sasmodels/generate.py>`_
346
347    - new units can be added using the macros defined in *doc/rst_prolog*
348      in the sasmodels source.
349    - units should be properly formatted using sub-/super-scripts
350      and using negative exponents instead of the / operator, though
351      the unit name should use the / operator for consistency.
352    - please post a message to the SasView developers mailing list with your changes.
353
354- **default** is the initial value for the parameter.
355
356  - **the parameter default values are used to auto-generate a plot of
357    the model function in the documentation.**
358
359- **[min, max]** are the lower and upper limits on the parameter.
360
361  - lower and upper limits can be any number, or *-inf* or *inf*.
362
363  - the limits will show up as the default limits for the fit making it easy,
364    for example, to force the radius to always be greater than zero.
365
366  - these are hard limits defining the valid range of parameter values;
367    polydisperity distributions will be truncated at the limits.
368
369- **"type"** can be one of: "", "sld", "volume", or "orientation".
370
371  - "sld" parameters can have magnetic moments when fitting magnetic models;
372    depending on the spin polarization of the beam and the $q$ value being
373    examined, the effective sld for that material will be used to compute the
374    scattered intensity.
375
[31fc4ad]376  - "volume" parameters are passed to Fq(), Iq(), Iqac(), Iqabc(), form_volume()
377    and shell_volume(), and have polydispersity loops generated automatically.
[990d8df]378
[108e70e]379  - "orientation" parameters are not passed, but instead are combined with
380    orientation dispersity to translate *qx* and *qy* to *qa*, *qb* and *qc*.
381    These parameters should appear at the end of the table with the specific
382    names *theta*, *phi* and for asymmetric shapes *psi*, in that order.
[990d8df]383
[9844c3a]384Some models will have integer parameters, such as number of pearls in the
385pearl necklace model, or number of shells in the multi-layer vesicle model.
386The optimizers in BUMPS treat all parameters as floating point numbers which
387can take arbitrary values, even for integer parameters, so your model should
388round the incoming parameter value to the nearest integer inside your model
389you should round to the nearest integer.  In C code, you can do this using::
390
391    static double
392    Iq(double q, ..., double fp_n, ...)
393    {
394        int n = (int)(fp_n + 0.5);
395        ...
396    }
397
398in python::
399
400    def Iq(q, ..., n, ...):
401        n = int(n+0.5)
402        ...
403
[3048ec6]404Derivative based optimizers such as Levenberg-Marquardt will not work
[9844c3a]405for integer parameters since the partial derivative is always zero, but
406the remaining optimizers (DREAM, differential evolution, Nelder-Mead simplex)
407will still function.
408
[990d8df]409Model Computation
410.................
411
412Models can be defined as pure python models, or they can be a mixture of
413python and C models.  C models are run on the GPU if it is available,
414otherwise they are compiled and run on the CPU.
415
416Models are defined by the scattering kernel, which takes a set of parameter
417values defining the shape, orientation and material, and returns the
418expected scattering. Polydispersity and angular dispersion are defined
419by the computational infrastructure.  Any parameters defined as "volume"
420parameters are polydisperse, with polydispersity defined in proportion
421to their value.  "orientation" parameters use angular dispersion defined
422in degrees, and are not relative to the current angle.
423
424Based on a weighting function $G(x)$ and a number of points $n$, the
425computed value is
426
427.. math::
428
429     \hat I(q)
430     = \frac{\int G(x) I(q, x)\,dx}{\int G(x) V(x)\,dx}
431     \approx \frac{\sum_{i=1}^n G(x_i) I(q,x_i)}{\sum_{i=1}^n G(x_i) V(x_i)}
432
[3048ec6]433That is, the individual models do not need to include polydispersity
[990d8df]434calculations, but instead rely on numerical integration to compute the
[108e70e]435appropriately smeared pattern.
[990d8df]436
[2015f02]437Each .py file also contains a function::
438
439        def random():
440        ...
[31fc4ad]441
442This function provides a model-specific random parameter set which shows model
443features in the USANS to SANS range.  For example, core-shell sphere sets the
444outer radius of the sphere logarithmically in `[20, 20,000]`, which sets the Q
445value for the transition from flat to falling.  It then uses a beta distribution
446to set the percentage of the shape which is shell, giving a preference for very
447thin or very thick shells (but never 0% or 100%).  Using `-sets=10` in sascomp
448should show a reasonable variety of curves over the default sascomp q range.
449The parameter set is returned as a dictionary of `{parameter: value, ...}`.
450Any model parameters not included in the dictionary will default according to
[2015f02]451the code in the `_randomize_one()` function from sasmodels/compare.py.
452
[990d8df]453Python Models
454.............
455
456For pure python models, define the *Iq* function::
457
458      import numpy as np
459      from numpy import cos, sin, ...
460
461      def Iq(q, par1, par2, ...):
462          return I(q, par1, par2, ...)
463      Iq.vectorized = True
464
465The parameters *par1, par2, ...* are the list of non-orientation parameters
466to the model in the order that they appear in the parameter table.
[3048ec6]467**Note that the auto-generated model file uses** *x* **rather than** *q*.
[990d8df]468
469The *.py* file should import trigonometric and exponential functions from
470numpy rather than from math.  This lets us evaluate the model for the whole
471range of $q$ values at once rather than looping over each $q$ separately in
472python.  With $q$ as a vector, you cannot use if statements, but must instead
473do tricks like
474
475::
476
477     a = x*q*(q>0) + y*q*(q<=0)
478
479or
480
481::
482
483     a = np.empty_like(q)
484     index = q>0
485     a[index] = x*q[index]
486     a[~index] = y*q[~index]
487
488which sets $a$ to $q \cdot x$ if $q$ is positive or $q \cdot y$ if $q$
489is zero or negative. If you have not converted your function to use $q$
490vectors, you can set the following and it will only receive one $q$
491value at a time::
492
493    Iq.vectorized = False
494
495Return np.NaN if the parameters are not valid (e.g., cap_radius < radius in
496barbell).  If I(q; pars) is NaN for any $q$, then those parameters will be
497ignored, and not included in the calculation of the weighted polydispersity.
498
499Models should define *form_volume(par1, par2, ...)* where the parameter
500list includes the *volume* parameters in order.  This is used for a weighted
501volume normalization so that scattering is on an absolute scale.  If
502*form_volume* is not defined, then the default *form_volume = 1.0* will be
503used.
504
[31fc4ad]505Hollow shapes, where the volume fraction of particle corresponds to the
506material in the shell rather than the volume enclosed by the shape, must
507also define a *shell_volume(par1, par2, ...)* function.  The parameters
508are the same as for *form_volume*.  The *I(q)* calculation should use
509*shell_volume* squared as its scale factor for the volume normalization.
510The structure factor calculation needs *form_volume* in order to properly
511scale the volume fraction parameter, so both functions are required for
512hollow shapes.
513
514Note: Pure python models do not yet support direct computation of the
515average of $F(q)$ and $F^2(q)$.
516
[990d8df]517Embedded C Models
518.................
519
520Like pure python models, inline C models need to define an *Iq* function::
521
522    Iq = """
523        return I(q, par1, par2, ...);
524    """
525
526This expands into the equivalent C code::
527
528    double Iq(double q, double par1, double par2, ...);
529    double Iq(double q, double par1, double par2, ...)
530    {
531        return I(q, par1, par2, ...);
532    }
533
534*form_volume* defines the volume of the shape. As in python models, it
535includes only the volume parameters.
536
[31fc4ad]537*form_volume* defines the volume of the shell for hollow shapes. As in
538python models, it includes only the volume parameters.
539
[990d8df]540**source=['fn.c', ...]** includes the listed C source files in the
[108e70e]541program before *Iq* and *form_volume* are defined. This allows you to
[ef85a09]542extend the library of C functions available to your model.
543
544*c_code* includes arbitrary C code into your kernel, which can be
545handy for defining helper functions for *Iq* and *form_volume*. Note that
[108e70e]546you can put the full function definition for *Iq* and *form_volume*
[ef85a09]547(include function declaration) into *c_code* as well, or put them into an
548external C file and add that file to the list of sources.
[990d8df]549
550Models are defined using double precision declarations for the
551parameters and return values.  When a model is run using single
552precision or long double precision, each variable is converted
553to the target type, depending on the precision requested.
554
555**Floating point constants must include the decimal point.**  This allows us
556to convert values such as 1.0 (double precision) to 1.0f (single precision)
557so that expressions that use these values are not promoted to double precision
558expressions.  Some graphics card drivers are confused when functions
559that expect floating point values are passed integers, such as 4*atan(1); it
560is safest to not use integers in floating point expressions.  Even better,
561use the builtin constant M_PI rather than 4*atan(1); it is faster and smaller!
562
563The C model operates on a single $q$ value at a time.  The code will be
564run in parallel across different $q$ values, either on the graphics card
565or the processor.
566
567Rather than returning NAN from Iq, you must define the *INVALID(v)*.  The
568*v* parameter lets you access all the parameters in the model using
569*v.par1*, *v.par2*, etc. For example::
570
571    #define INVALID(v) (v.bell_radius < v.radius)
572
[ef85a09]573The INVALID define can go into *Iq*, or *c_code*, or an external C file
574listed in *source*.
575
[31fc4ad]576Structure Factors
577.................
578
579Structure factor calculations may need the underlying $<F(q)>$ and $<F^2(q)>$
580rather than $I(q)$.  This is used to compute $\beta = <F(q)>^2/<F^2(q)>$ in
581the decoupling approximation to the structure factor.
582
583Instead of defining the *Iq* function, models can define *Fq* as
584something like::
585
586    double Fq(double q, double *F1, double *F2, double par1, double par2, ...);
587    double Fq(double q, double *F1, double *F2, double par1, double par2, ...)
588    {
589        // Polar integration loop over all orientations.
590        ...
591        *F1 = 1e-2 * total_F1 * contrast * volume;
592        *F2 = 1e-4 * total_F2 * square(contrast * volume);
593        return I(q, par1, par2, ...);
594    }
595
596If the volume fraction scale factor is built into the model (as occurs for
597the vesicle model, for example), then scale *F1* by $\surd V_f$ so that
598$\beta$ is computed correctly.
599
600Structure factor calculations are not yet supported for oriented shapes.
601
602Note: only available as a separate C file listed in *source*, or within
603a *c_code* block within the python model definition file.
604
[108e70e]605Oriented Shapes
606...............
607
608If the scattering is dependent on the orientation of the shape, then you
609will need to include *orientation* parameters *theta*, *phi* and *psi*
[7e6bc45e]610at the end of the parameter table.  As described in the section
611:ref:`orientation`, the individual $(q_x, q_y)$ points on the detector will
612be rotated into $(q_a, q_b, q_c)$ points relative to the sample in its
613canonical orientation with $a$-$b$-$c$ aligned with $x$-$y$-$z$ in the
614laboratory frame and beam travelling along $-z$.
615
616The oriented C model is called using *Iqabc(qa, qb, qc, par1, par2, ...)* where
[108e70e]617*par1*, etc. are the parameters to the model.  If the shape is rotationally
618symmetric about *c* then *psi* is not needed, and the model is called
619as *Iqac(qab, qc, par1, par2, ...)*.  In either case, the orientation
620parameters are not included in the function call.
621
622For 1D oriented shapes, an integral over all angles is usually needed for
[b85227d]623the *Iq* function. Given symmetry and the substitution $u = \cos(\alpha)$,
[108e70e]624$du = -\sin(\alpha)\,d\alpha$ this becomes
625
626.. math::
627
[b85227d]628    I(q) &= \frac{1}{4\pi} \int_{-\pi/2}^{pi/2} \int_{-pi}^{pi}
629            F(q_a, q_b, q_c)^2 \sin(\alpha)\,d\beta\,d\alpha \\
630        &= \frac{8}{4\pi} \int_{0}^{pi/2} \int_{0}^{\pi/2}
631            F^2 \sin(\alpha)\,d\beta\,d\alpha \\
632        &= \frac{8}{4\pi} \int_1^0 \int_{0}^{\pi/2} - F^2 \,d\beta\,du \\
633        &= \frac{8}{4\pi} \int_0^1 \int_{0}^{\pi/2} F^2 \,d\beta\,du
634
635for
636
637.. math::
638
639    q_a &= q \sin(\alpha)\sin(\beta) = q \sqrt{1-u^2} \sin(\beta) \\
640    q_b &= q \sin(\alpha)\cos(\beta) = q \sqrt{1-u^2} \cos(\beta) \\
641    q_c &= q \cos(\alpha) = q u
[108e70e]642
643Using the $z, w$ values for Gauss-Legendre integration in "lib/gauss76.c", the
644numerical integration is then::
645
646    double outer_sum = 0.0;
647    for (int i = 0; i < GAUSS_N; i++) {
648        const double cos_alpha = 0.5*GAUSS_Z[i] + 0.5;
649        const double sin_alpha = sqrt(1.0 - cos_alpha*cos_alpha);
650        const double qc = cos_alpha * q;
651        double inner_sum = 0.0;
652        for (int j = 0; j < GAUSS_N; j++) {
653            const double beta = M_PI_4 * GAUSS_Z[j] + M_PI_4;
654            double sin_beta, cos_beta;
655            SINCOS(beta, sin_beta, cos_beta);
656            const double qa = sin_alpha * sin_beta * q;
[b85227d]657            const double qb = sin_alpha * cos_beta * q;
658            const double form = Fq(qa, qb, qc, ...);
659            inner_sum += GAUSS_W[j] * form * form;
[108e70e]660        }
661        outer_sum += GAUSS_W[i] * inner_sum;
662    }
663    outer_sum *= 0.25; // = 8/(4 pi) * outer_sum * (pi/2) / 4
664
665The *z* values for the Gauss-Legendre integration extends from -1 to 1, so
666the double sum of *w[i]w[j]* explains the factor of 4.  Correcting for the
667average *dz[i]dz[j]* gives $(1-0) \cdot (\pi/2-0) = \pi/2$.  The $8/(4 \pi)$
668factor comes from the integral over the quadrant.  With less symmetry (eg.,
669in the bcc and fcc paracrystal models), then an integral over the entire
670sphere may be necessary.
671
672For simpler models which are rotationally symmetric a single integral
673suffices:
674
675.. math::
676
[b85227d]677    I(q) &= \frac{1}{\pi}\int_{-\pi/2}^{\pi/2}
678            F(q_{ab}, q_c)^2 \sin(\alpha)\,d\alpha/\pi \\
679        &= \frac{2}{\pi} \int_0^1 F^2\,du
680
681for
682
683.. math::
684
685    q_{ab} &= q \sin(\alpha) = q \sqrt{1 - u^2} \\
686    q_c &= q \cos(\alpha) = q u
687
[108e70e]688
689with integration loop::
690
691    double sum = 0.0;
692    for (int i = 0; i < GAUSS_N; i++) {
693        const double cos_alpha = 0.5*GAUSS_Z[i] + 0.5;
694        const double sin_alpha = sqrt(1.0 - cos_alpha*cos_alpha);
695        const double qab = sin_alpha * q;
[b85227d]696        const double qc = cos_alpha * q;
697        const double form = Fq(qab, qc, ...);
698        sum += GAUSS_W[j] * form * form;
[108e70e]699    }
700    sum *= 0.5; // = 2/pi * sum * (pi/2) / 2
701
702Magnetism
703.........
704
705Magnetism is supported automatically for all shapes by modifying the
706effective SLD of particle according to the Halpern-Johnson vector
[c654160]707describing the interaction between neutron spin and magnetic field.  All
[108e70e]708parameters marked as type *sld* in the parameter table are treated as
709possibly magnetic particles with magnitude *M0* and direction
710*mtheta* and *mphi*.  Polarization parameters are also provided
711automatically for magnetic models to set the spin state of the measurement.
712
713For more complicated systems where magnetism is not uniform throughout
714the individual particles, you will need to write your own models.
715You should not mark the nuclear sld as type *sld*, but instead leave
716them unmarked and provide your own magnetism and polarization parameters.
717For 2D measurements you will need $(q_x, q_y)$ values for the measurement
718to compute the proper magnetism and orientation, which you can implement
719using *Iqxy(qx, qy, par1, par2, ...)*.
720
[990d8df]721Special Functions
722.................
723
724The C code follows the C99 standard, with the usual math functions,
725as defined in
726`OpenCL <https://www.khronos.org/registry/cl/sdk/1.1/docs/man/xhtml/mathFunctions.html>`_.
727This includes the following:
728
729    M_PI, M_PI_2, M_PI_4, M_SQRT1_2, M_E:
730        $\pi$, $\pi/2$, $\pi/4$, $1/\sqrt{2}$ and Euler's constant $e$
[d0dc9a3]731    exp, log, pow(x,y), expm1, log1p, sqrt, cbrt:
732        Power functions $e^x$, $\ln x$, $x^y$, $e^x - 1$, $\ln 1 + x$,
733        $\sqrt{x}$, $\sqrt[3]{x}$. The functions expm1(x) and log1p(x)
734        are accurate across all $x$, including $x$ very close to zero.
[990d8df]735    sin, cos, tan, asin, acos, atan:
736        Trigonometry functions and inverses, operating on radians.
737    sinh, cosh, tanh, asinh, acosh, atanh:
738        Hyperbolic trigonometry functions.
739    atan2(y,x):
740        Angle from the $x$\ -axis to the point $(x,y)$, which is equal to
741        $\tan^{-1}(y/x)$ corrected for quadrant.  That is, if $x$ and $y$ are
742        both negative, then atan2(y,x) returns a value in quadrant III where
743        atan(y/x) would return a value in quadrant I. Similarly for
744        quadrants II and IV when $x$ and $y$ have opposite sign.
[d0dc9a3]745    fabs(x), fmin(x,y), fmax(x,y), trunc, rint:
[990d8df]746        Floating point functions.  rint(x) returns the nearest integer.
747    NAN:
748        NaN, Not a Number, $0/0$.  Use isnan(x) to test for NaN.  Note that
749        you cannot use :code:`x == NAN` to test for NaN values since that
[d0dc9a3]750        will always return false.  NAN does not equal NAN!  The alternative,
751        :code:`x != x` may fail if the compiler optimizes the test away.
[990d8df]752    INFINITY:
753        $\infty, 1/0$.  Use isinf(x) to test for infinity, or isfinite(x)
754        to test for finite and not NaN.
755    erf, erfc, tgamma, lgamma:  **do not use**
756        Special functions that should be part of the standard, but are missing
[fba9ca0]757        or inaccurate on some platforms. Use sas_erf, sas_erfc, sas_gamma
758        and sas_lgamma instead (see below).
[990d8df]759
760Some non-standard constants and functions are also provided:
761
762    M_PI_180, M_4PI_3:
763        $\frac{\pi}{180}$, $\frac{4\pi}{3}$
764    SINCOS(x, s, c):
765        Macro which sets s=sin(x) and c=cos(x). The variables *c* and *s*
766        must be declared first.
767    square(x):
768        $x^2$
769    cube(x):
770        $x^3$
771    sas_sinx_x(x):
772        $\sin(x)/x$, with limit $\sin(0)/0 = 1$.
773    powr(x, y):
774        $x^y$ for $x \ge 0$; this is faster than general $x^y$ on some GPUs.
775    pown(x, n):
776        $x^n$ for $n$ integer; this is faster than general $x^n$ on some GPUs.
777    FLOAT_SIZE:
778        The number of bytes in a floating point value.  Even though all
779        variables are declared double, they may be converted to single
780        precision float before running. If your algorithm depends on
781        precision (which is not uncommon for numerical algorithms), use
782        the following::
783
784            #if FLOAT_SIZE>4
785            ... code for double precision ...
786            #else
787            ... code for single precision ...
788            #endif
789    SAS_DOUBLE:
790        A replacement for :code:`double` so that the declared variable will
791        stay double precision; this should generally not be used since some
792        graphics cards do not support double precision.  There is no provision
793        for forcing a constant to stay double precision.
794
795The following special functions and scattering calculations are defined in
796`sasmodels/models/lib <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib>`_.
797These functions have been tuned to be fast and numerically stable down
798to $q=0$ even in single precision.  In some cases they work around bugs
799which appear on some platforms but not others, so use them where needed.
800Add the files listed in :code:`source = ["lib/file.c", ...]` to your *model.py*
801file in the order given, otherwise these functions will not be available.
802
803    polevl(x, c, n):
804        Polynomial evaluation $p(x) = \sum_{i=0}^n c_i x^i$ using Horner's
805        method so it is faster and more accurate.
806
807        $c = \{c_n, c_{n-1}, \ldots, c_0 \}$ is the table of coefficients,
808        sorted from highest to lowest.
809
810        :code:`source = ["lib/polevl.c", ...]` (`link to code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/polevl.c>`_)
811
812    p1evl(x, c, n):
813        Evaluation of normalized polynomial $p(x) = x^n + \sum_{i=0}^{n-1} c_i x^i$
814        using Horner's method so it is faster and more accurate.
815
816        $c = \{c_{n-1}, c_{n-2} \ldots, c_0 \}$ is the table of coefficients,
817        sorted from highest to lowest.
818
819        :code:`source = ["lib/polevl.c", ...]`
[870a2f4]820        (`polevl.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/polevl.c>`_)
[990d8df]821
822    sas_gamma(x):
[30b60d2]823        Gamma function sas_gamma\ $(x) = \Gamma(x)$.
[990d8df]824
[fba9ca0]825        The standard math function, tgamma(x), is unstable for $x < 1$
[990d8df]826        on some platforms.
827
[870a2f4]828        :code:`source = ["lib/sas_gamma.c", ...]`
829        (`sas_gamma.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_gamma.c>`_)
[990d8df]830
[fba9ca0]831    sas_gammaln(x):
832        log gamma function sas_gammaln\ $(x) = \log \Gamma(|x|)$.
833
834        The standard math function, lgamma(x), is incorrect for single
835        precision on some platforms.
836
837        :code:`source = ["lib/sas_gammainc.c", ...]`
838        (`sas_gammainc.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_gammainc.c>`_)
839
840    sas_gammainc(a, x), sas_gammaincc(a, x):
841        Incomplete gamma function
842        sas_gammainc\ $(a, x) = \int_0^x t^{a-1}e^{-t}\,dt / \Gamma(a)$
843        and complementary incomplete gamma function
844        sas_gammaincc\ $(a, x) = \int_x^\infty t^{a-1}e^{-t}\,dt / \Gamma(a)$
845
846        :code:`source = ["lib/sas_gammainc.c", ...]`
847        (`sas_gammainc.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_gammainc.c>`_)
848
[990d8df]849    sas_erf(x), sas_erfc(x):
850        Error function
[30b60d2]851        sas_erf\ $(x) = \frac{2}{\sqrt\pi}\int_0^x e^{-t^2}\,dt$
[990d8df]852        and complementary error function
[30b60d2]853        sas_erfc\ $(x) = \frac{2}{\sqrt\pi}\int_x^{\infty} e^{-t^2}\,dt$.
[990d8df]854
855        The standard math functions erf(x) and erfc(x) are slower and broken
856        on some platforms.
857
858        :code:`source = ["lib/polevl.c", "lib/sas_erf.c", ...]`
[870a2f4]859        (`sas_erf.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_erf.c>`_)
[990d8df]860
861    sas_J0(x):
[30b60d2]862        Bessel function of the first kind sas_J0\ $(x)=J_0(x)$ where
[990d8df]863        $J_0(x) = \frac{1}{\pi}\int_0^\pi \cos(x\sin(\tau))\,d\tau$.
864
865        The standard math function j0(x) is not available on all platforms.
866
867        :code:`source = ["lib/polevl.c", "lib/sas_J0.c", ...]`
[870a2f4]868        (`sas_J0.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_J0.c>`_)
[990d8df]869
870    sas_J1(x):
[30b60d2]871        Bessel function of the first kind  sas_J1\ $(x)=J_1(x)$ where
[990d8df]872        $J_1(x) = \frac{1}{\pi}\int_0^\pi \cos(\tau - x\sin(\tau))\,d\tau$.
873
874        The standard math function j1(x) is not available on all platforms.
875
876        :code:`source = ["lib/polevl.c", "lib/sas_J1.c", ...]`
[870a2f4]877        (`sas_J1.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_J1.c>`_)
[990d8df]878
879    sas_JN(n, x):
[30b60d2]880        Bessel function of the first kind and integer order $n$,
881        sas_JN\ $(n, x) =J_n(x)$ where
[990d8df]882        $J_n(x) = \frac{1}{\pi}\int_0^\pi \cos(n\tau - x\sin(\tau))\,d\tau$.
[30b60d2]883        If $n$ = 0 or 1, it uses sas_J0($x$) or sas_J1($x$), respectively.
[990d8df]884
[57c609b]885        Warning: JN(n,x) can be very inaccurate (0.1%) for x not in [0.1, 100].
886
[990d8df]887        The standard math function jn(n, x) is not available on all platforms.
888
889        :code:`source = ["lib/polevl.c", "lib/sas_J0.c", "lib/sas_J1.c", "lib/sas_JN.c", ...]`
[870a2f4]890        (`sas_JN.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_JN.c>`_)
[990d8df]891
892    sas_Si(x):
[30b60d2]893        Sine integral Si\ $(x) = \int_0^x \tfrac{\sin t}{t}\,dt$.
[990d8df]894
[57c609b]895        Warning: Si(x) can be very inaccurate (0.1%) for x in [0.1, 100].
896
[990d8df]897        This function uses Taylor series for small and large arguments:
898
[57c609b]899        For large arguments use the following Taylor series,
[990d8df]900
901        .. math::
902
903             \text{Si}(x) \sim \frac{\pi}{2}
904             - \frac{\cos(x)}{x}\left(1 - \frac{2!}{x^2} + \frac{4!}{x^4} - \frac{6!}{x^6} \right)
905             - \frac{\sin(x)}{x}\left(\frac{1}{x} - \frac{3!}{x^3} + \frac{5!}{x^5} - \frac{7!}{x^7}\right)
906
[94bfa42]907        For small arguments,
[990d8df]908
909        .. math::
910
911           \text{Si}(x) \sim x
912           - \frac{x^3}{3\times 3!} + \frac{x^5}{5 \times 5!} - \frac{x^7}{7 \times 7!}
913           + \frac{x^9}{9\times 9!} - \frac{x^{11}}{11\times 11!}
914
915        :code:`source = ["lib/Si.c", ...]`
[f796469]916        (`Si.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_Si.c>`_)
[990d8df]917
918    sas_3j1x_x(x):
919        Spherical Bessel form
[30b60d2]920        sph_j1c\ $(x) = 3 j_1(x)/x = 3 (\sin(x) - x \cos(x))/x^3$,
[990d8df]921        with a limiting value of 1 at $x=0$, where $j_1(x)$ is the spherical
922        Bessel function of the first kind and first order.
923
924        This function uses a Taylor series for small $x$ for numerical accuracy.
925
926        :code:`source = ["lib/sas_3j1x_x.c", ...]`
[870a2f4]927        (`sas_3j1x_x.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_3j1x_x.c>`_)
[990d8df]928
929
930    sas_2J1x_x(x):
[30b60d2]931        Bessel form sas_J1c\ $(x) = 2 J_1(x)/x$, with a limiting value
[990d8df]932        of 1 at $x=0$, where $J_1(x)$ is the Bessel function of first kind
933        and first order.
934
935        :code:`source = ["lib/polevl.c", "lib/sas_J1.c", ...]`
[870a2f4]936        (`sas_J1.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_J1.c>`_)
[990d8df]937
938
939    Gauss76Z[i], Gauss76Wt[i]:
940        Points $z_i$ and weights $w_i$ for 76-point Gaussian quadrature, respectively,
941        computing $\int_{-1}^1 f(z)\,dz \approx \sum_{i=1}^{76} w_i\,f(z_i)$.
942
943        Similar arrays are available in :code:`gauss20.c` for 20-point
944        quadrature and in :code:`gauss150.c` for 150-point quadrature.
[d0dc9a3]945        The macros :code:`GAUSS_N`, :code:`GAUSS_Z` and :code:`GAUSS_W` are
946        defined so that you can change the order of the integration by
947        selecting an different source without touching the C code.
[990d8df]948
949        :code:`source = ["lib/gauss76.c", ...]`
[870a2f4]950        (`gauss76.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/gauss76.c>`_)
[990d8df]951
952
953
954Problems with C models
955......................
956
957The graphics processor (GPU) in your computer is a specialized computer tuned
958for certain kinds of problems.  This leads to strange restrictions that you
959need to be aware of.  Your code may work fine on some platforms or for some
960models, but then return bad values on other platforms.  Some examples of
961particular problems:
962
963  **(1) Code is too complex, or uses too much memory.** GPU devices only
964  have a limited amount of memory available for each processor. If you run
965  programs which take too much memory, then rather than running multiple
966  values in parallel as it usually does, the GPU may only run a single
967  version of the code at a time, making it slower than running on the CPU.
968  It may fail to run on some platforms, or worse, cause the screen to go
969  blank or the system to reboot.
970
971  **(2) Code takes too long.** Because GPU devices are used for the computer
972  display, the OpenCL drivers are very careful about the amount of time they
973  will allow any code to run. For example, on OS X, the model will stop
974  running after 5 seconds regardless of whether the computation is complete.
975  You may end up with only some of your 2D array defined, with the rest
976  containing random data. Or it may cause the screen to go blank or the
977  system to reboot.
978
979  **(3) Memory is not aligned**. The GPU hardware is specialized to operate
980  on multiple values simultaneously. To keep the GPU simple the values in
981  memory must be aligned with the different GPU compute engines. Not
982  following these rules can lead to unexpected values being loaded into
983  memory, and wrong answers computed. The conclusion from a very long and
984  strange debugging session was that any arrays that you declare in your
985  model should be a multiple of four. For example::
986
987      double Iq(q, p1, p2, ...)
988      {
989          double vector[8];  // Only going to use seven slots, but declare 8
990          ...
991      }
992
993The first step when your model is behaving strangely is to set
994**single=False**. This automatically restricts the model to only run on the
995CPU, or on high-end GPU cards. There can still be problems even on high-end
996cards, so you can force the model off the GPU by setting **opencl=False**.
997This runs the model as a normal C program without any GPU restrictions so
998you know that strange results are probably from your code rather than the
999environment. Once the code is debugged, you can compare your output to the
1000output on the GPU.
1001
1002Although it can be difficult to get your model to work on the GPU, the reward
1003can be a model that runs 1000x faster on a good card.  Even your laptop may
1004show a 50x improvement or more over the equivalent pure python model.
1005
1006
1007.. _Form_Factors:
1008
1009Form Factors
1010............
1011
1012Away from the dilute limit you can estimate scattering including
1013particle-particle interactions using $I(q) = P(q)*S(q)$ where $P(q)$
1014is the form factor and $S(q)$ is the structure factor.  The simplest
1015structure factor is the *hardsphere* interaction, which
1016uses the effective radius of the form factor as an input to the structure
1017factor model.  The effective radius is the average radius of the
1018form averaged over all the polydispersity values.
1019
1020::
1021
1022    def ER(radius, thickness):
1023        """Effective radius of a core-shell sphere."""
1024        return radius + thickness
1025
1026Now consider the *core_shell_sphere*, which has a simple effective radius
1027equal to the radius of the core plus the thickness of the shell, as
1028shown above. Given polydispersity over *(r1, r2, ..., rm)* in radius and
1029*(t1, t2, ..., tn)* in thickness, *ER* is called with a mesh
1030grid covering all possible combinations of radius and thickness.
1031That is, *radius* is *(r1, r2, ..., rm, r1, r2, ..., rm, ...)*
1032and *thickness* is *(t1, t1, ... t1, t2, t2, ..., t2, ...)*.
1033The *ER* function returns one effective radius for each combination.
1034The effective radius calculator weights each of these according to
1035the polydispersity distributions and calls the structure factor
1036with the average *ER*.
1037
1038::
1039
1040    def VR(radius, thickness):
1041        """Sphere and shell volumes for a core-shell sphere."""
1042        whole = 4.0/3.0 * pi * (radius + thickness)**3
1043        core = 4.0/3.0 * pi * radius**3
1044        return whole, whole - core
1045
1046Core-shell type models have an additional volume ratio which scales
1047the structure factor.  The *VR* function returns the volume of
1048the whole sphere and the volume of the shell. Like *ER*, there is
1049one return value for each point in the mesh grid.
1050
1051*NOTE: we may be removing or modifying this feature soon. As of the
1052time of writing, core-shell sphere returns (1., 1.) for VR, giving a volume
1053ratio of 1.0.*
1054
1055Unit Tests
1056..........
1057
1058THESE ARE VERY IMPORTANT. Include at least one test for each model and
1059PLEASE make sure that the answer value is correct (i.e. not a random number).
1060
1061::
1062
1063    tests = [
1064        [{}, 0.2, 0.726362],
1065        [{"scale": 1., "background": 0., "sld": 6., "sld_solvent": 1.,
1066          "radius": 120., "radius_pd": 0.2, "radius_pd_n":45},
1067         0.2, 0.228843],
[304c775]1068        [{"radius": 120., "radius_pd": 0.2, "radius_pd_n":45},
1069         0.1, None, None, 120., None, 1.],  # q, F, F^2, R_eff, V, form:shell
[81751c2]1070        [{"@S": "hardsphere"}, 0.1, None],
[990d8df]1071    ]
1072
1073
[304c775]1074**tests=[[{parameters}, q, Iq], ...]** is a list of lists.
[990d8df]1075Each list is one test and contains, in order:
1076
1077- a dictionary of parameter values. This can be *{}* using the default
1078  parameters, or filled with some parameters that will be different from the
1079  default, such as *{"radius":10.0, "sld":4}*. Unlisted parameters will
1080  be given the default values.
1081- the input $q$ value or tuple of $(q_x, q_y)$ values.
1082- the output $I(q)$ or $I(q_x,q_y)$ expected of the model for the parameters
1083  and input value given.
1084- input and output values can themselves be lists if you have several
1085  $q$ values to test for the same model parameters.
[304c775]1086- for testing effective radius, volume and form:shell volume ratio, use the
1087  extended form of the tests results, with *None, None, R_eff, V, V_r*
1088  instead of *Iq*.  This calls the kernel *Fq* function instead of *Iq*.
1089- for testing F and F^2 (used for beta approximation) do the same as the
1090  effective radius test, but include values for the first two elements,
1091  $<F(q)>$ and $<F^2(q)>$.
[81751c2]1092- for testing interaction between form factor and structure factor, specify
1093  the structure factor name in the parameters as *{"@S": "name", ...}* with
1094  the remaining list of parameters defined by the *P@S* product model.
[990d8df]1095
1096.. _Test_Your_New_Model:
1097
1098Test Your New Model
1099^^^^^^^^^^^^^^^^^^^
1100
1101Minimal Testing
1102...............
1103
1104From SasView either open the Python shell (*Tools* > *Python Shell/Editor*)
1105or the plugin editor (*Fitting* > *Plugin Model Operations* > *Advanced
1106Plugin Editor*), load your model, and then select *Run > Check Model* from
1107the menu bar. An *Info* box will appear with the results of the compilation
1108and a check that the model runs.
1109
1110If you are not using sasmodels from SasView, skip this step.
1111
1112Recommended Testing
1113...................
1114
1115If the model compiles and runs, you can next run the unit tests that
1116you have added using the **test =** values.
1117
1118From SasView, switch to the *Shell* tab and type the following::
1119
1120    from sasmodels.model_test import run_one
1121    run_one("~/.sasview/plugin_models/model.py")
1122
1123This should print::
1124
1125    test_model_python (sasmodels.model_test.ModelTestCase) ... ok
1126
1127To check whether single precision is good enough, type the following::
1128
1129    from sasmodels.compare import main as compare
1130    compare("~/.sasview/plugin_models/model.py")
1131
1132This will pop up a plot showing the difference between single precision
1133and double precision on a range of $q$ values.
1134
1135::
1136
1137  demo = dict(scale=1, background=0,
1138              sld=6, sld_solvent=1,
1139              radius=120,
1140              radius_pd=.2, radius_pd_n=45)
1141
1142**demo={'par': value, ...}** in the model file sets the default values for
1143the comparison. You can include polydispersity parameters such as
1144*radius_pd=0.2, radius_pd_n=45* which would otherwise be zero.
1145
1146These commands can also be run directly in the python interpreter:
1147
1148    $ python -m sasmodels.model_test -v ~/.sasview/plugin_models/model.py
1149    $ python -m sasmodels.compare ~/.sasview/plugin_models/model.py
1150
1151The options to compare are quite extensive; type the following for help::
1152
1153    compare()
1154
1155Options will need to be passed as separate strings.
1156For example to run your model with a random set of parameters::
1157
1158    compare("-random", "-pars", "~/.sasview/plugin_models/model.py")
1159
1160For the random models,
1161
1162- *sld* will be in the range (-0.5,10.5),
1163- angles (*theta, phi, psi*) will be in the range (-180,180),
1164- angular dispersion will be in the range (0,45),
1165- polydispersity will be in the range (0,1)
1166- other values will be in the range (0, 2\ *v*), where *v* is the value
1167  of the parameter in demo.
1168
1169Dispersion parameters *n*\, *sigma* and *type* will be unchanged from
1170demo so that run times are more predictable (polydispersity calculated
1171across multiple parameters can be very slow).
1172
[3048ec6]1173If your model has 2D orientation calculation, then you should also
[990d8df]1174test with::
1175
1176    compare("-2d", "~/.sasview/plugin_models/model.py")
1177
1178Check The Docs
1179^^^^^^^^^^^^^^
1180
1181You can get a rough idea of how the documentation will look using the
1182following::
1183
1184    compare("-help", "~/.sasview/plugin_models/model.py")
1185
1186This does not use the same styling as the rest of the docs, but it will
1187allow you to check that your ReStructuredText and LaTeX formatting.
1188Here are some tools to help with the inevitable syntax errors:
1189
1190- `Sphinx cheat sheet <http://matplotlib.org/sampledoc/cheatsheet.html>`_
1191- `Sphinx Documentation <http://www.sphinx-doc.org/en/stable/>`_
1192- `MathJax <http://www.mathjax.org/>`_
1193- `amsmath <http://www.ams.org/publications/authors/tex/amslatex>`_
1194
1195There is also a neat online WYSIWYG ReStructuredText editor at
1196http://rst.ninjs.org\ .
1197
1198
1199Clean Lint - (Developer Version Only)
1200^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1201
1202**NB: For now we are not providing pylint with the installer version
1203of SasView; so unless you have a SasView build environment available,
1204you can ignore this section!**
1205
1206Run the lint check with::
1207
1208    python -m pylint --rcfile=extra/pylint.rc ~/.sasview/plugin_models/model.py
1209
1210We are not aiming for zero lint just yet, only keeping it to a minimum.
1211For now, don't worry too much about *invalid-name*. If you really want a
1212variable name *Rg* for example because $R_g$ is the right name for the model
1213parameter then ignore the lint errors.  Also, ignore *missing-docstring*
[108e70e]1214for standard model functions *Iq*, *Iqac*, etc.
[990d8df]1215
1216We will have delinting sessions at the SasView Code Camps, where we can
1217decide on standards for model files, parameter names, etc.
1218
1219For now, you can tell pylint to ignore things.  For example, to align your
1220parameters in blocks::
1221
1222    # pylint: disable=bad-whitespace,line-too-long
1223    #   ["name",                  "units", default, [lower, upper], "type", "description"],
1224    parameters = [
1225        ["contrast_factor",       "barns",    10.0,  [-inf, inf], "", "Contrast factor of the polymer"],
1226        ["bjerrum_length",        "Ang",       7.1,  [0, inf],    "", "Bjerrum length"],
1227        ["virial_param",          "1/Ang^2",  12.0,  [-inf, inf], "", "Virial parameter"],
1228        ["monomer_length",        "Ang",      10.0,  [0, inf],    "", "Monomer length"],
1229        ["salt_concentration",    "mol/L",     0.0,  [-inf, inf], "", "Concentration of monovalent salt"],
1230        ["ionization_degree",     "",          0.05, [0, inf],    "", "Degree of ionization"],
1231        ["polymer_concentration", "mol/L",     0.7,  [0, inf],    "", "Polymer molar concentration"],
1232        ]
1233    # pylint: enable=bad-whitespace,line-too-long
1234
1235Don't put in too many pylint statements, though, since they make the code ugly.
1236
1237Share Your Model!
1238^^^^^^^^^^^^^^^^^
1239
1240Once compare and the unit test(s) pass properly and everything is done,
1241consider adding your model to the
1242`Model Marketplace <http://marketplace.sasview.org/>`_ so that others may use it!
1243
1244.. ZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZ
1245
1246*Document History*
1247
1248| 2016-10-25 Steve King
[c654160]1249| 2017-05-07 Paul Kienzle - Moved from sasview to sasmodels docs
Note: See TracBrowser for help on using the repository browser.