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