source: sasmodels/doc/guide/plugin.rst @ 646eeaa

Last change on this file since 646eeaa was aa8c6e0, checked in by Paul Kienzle <pkienzle@…>, 6 years ago

Merge branch 'master' into beta_approx

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