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