1 | |
---|
2 | |
---|
3 | <!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" |
---|
4 | "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd"> |
---|
5 | |
---|
6 | |
---|
7 | <html xmlns="http://www.w3.org/1999/xhtml"> |
---|
8 | <head> |
---|
9 | <meta http-equiv="Content-Type" content="text/html; charset=utf-8" /> |
---|
10 | |
---|
11 | <title>GSASIIobj: Data objects — GSAS-II 0.2.0 documentation</title> |
---|
12 | |
---|
13 | <link rel="stylesheet" href="_static/default.css" type="text/css" /> |
---|
14 | <link rel="stylesheet" href="_static/pygments.css" type="text/css" /> |
---|
15 | |
---|
16 | <script type="text/javascript"> |
---|
17 | var DOCUMENTATION_OPTIONS = { |
---|
18 | URL_ROOT: '', |
---|
19 | VERSION: '0.2.0', |
---|
20 | COLLAPSE_INDEX: false, |
---|
21 | FILE_SUFFIX: '.html', |
---|
22 | HAS_SOURCE: true |
---|
23 | }; |
---|
24 | </script> |
---|
25 | <script type="text/javascript" src="_static/jquery.js"></script> |
---|
26 | <script type="text/javascript" src="_static/underscore.js"></script> |
---|
27 | <script type="text/javascript" src="_static/doctools.js"></script> |
---|
28 | <script type="text/javascript" src="http://cdn.mathjax.org/mathjax/latest/MathJax.js?config=TeX-AMS-MML_HTMLorMML"></script> |
---|
29 | <link rel="top" title="GSAS-II 0.2.0 documentation" href="index.html" /> |
---|
30 | <link rel="next" title="GSAS-II Utility Modules" href="GSASIIutil.html" /> |
---|
31 | <link rel="prev" title="GSAS-II Main Module" href="GSASII.html" /> |
---|
32 | </head> |
---|
33 | <body> |
---|
34 | <div class="related"> |
---|
35 | <h3>Navigation</h3> |
---|
36 | <ul> |
---|
37 | <li class="right" style="margin-right: 10px"> |
---|
38 | <a href="genindex.html" title="General Index" |
---|
39 | accesskey="I">index</a></li> |
---|
40 | <li class="right" > |
---|
41 | <a href="py-modindex.html" title="Python Module Index" |
---|
42 | >modules</a> |</li> |
---|
43 | <li class="right" > |
---|
44 | <a href="GSASIIutil.html" title="GSAS-II Utility Modules" |
---|
45 | accesskey="N">next</a> |</li> |
---|
46 | <li class="right" > |
---|
47 | <a href="GSASII.html" title="GSAS-II Main Module" |
---|
48 | accesskey="P">previous</a> |</li> |
---|
49 | <li><a href="index.html">GSAS-II 0.2.0 documentation</a> »</li> |
---|
50 | </ul> |
---|
51 | </div> |
---|
52 | |
---|
53 | <div class="document"> |
---|
54 | <div class="documentwrapper"> |
---|
55 | <div class="bodywrapper"> |
---|
56 | <div class="body"> |
---|
57 | |
---|
58 | <span class="target" id="module-GSASIIobj"></span><div class="section" id="gsasiiobj-data-objects"> |
---|
59 | <h1><em>GSASIIobj: Data objects</em><a class="headerlink" href="#gsasiiobj-data-objects" title="Permalink to this headline">¶</a></h1> |
---|
60 | <p>This module defines and/or documents the data structures used in GSAS-II, as well |
---|
61 | as provides misc. support routines.</p> |
---|
62 | <div class="section" id="constraints-tree-item"> |
---|
63 | <h2>Constraints Tree Item<a class="headerlink" href="#constraints-tree-item" title="Permalink to this headline">¶</a></h2> |
---|
64 | <span class="target" id="constraints-table"></span><p id="index-0">Constraints are stored in a dict, separated into groups. |
---|
65 | Note that parameter are named in the following pattern, |
---|
66 | p:h:<var>:n, where p is the phase number, h is the histogram number |
---|
67 | <var> is a variable name and n is the parameter number. |
---|
68 | If a parameter does not depend on a histogram or phase or is unnumbered, that |
---|
69 | number is omitted. |
---|
70 | Note that the contents of each dict item is a List where each element in the |
---|
71 | list is a <a class="reference internal" href="#constraint-definitions-table"><em>constraint definition objects</em></a>. |
---|
72 | The constraints in this form are converted in |
---|
73 | <a class="reference internal" href="GSASIIstruc.html#GSASIIstrIO.ProcessConstraints" title="GSASIIstrIO.ProcessConstraints"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIstrIO.ProcessConstraints()</span></tt></a> to the form used in <a class="reference internal" href="GSASIImapvars.html#module-GSASIImapvars" title="GSASIImapvars"><tt class="xref py py-mod docutils literal"><span class="pre">GSASIImapvars</span></tt></a></p> |
---|
74 | <p>The keys in the Constraints dict are:</p> |
---|
75 | <table border="1" class="docutils"> |
---|
76 | <colgroup> |
---|
77 | <col width="16%" /> |
---|
78 | <col width="84%" /> |
---|
79 | </colgroup> |
---|
80 | <thead valign="bottom"> |
---|
81 | <tr class="row-odd"><th class="head">key</th> |
---|
82 | <th class="head">explanation</th> |
---|
83 | </tr> |
---|
84 | </thead> |
---|
85 | <tbody valign="top"> |
---|
86 | <tr class="row-even"><td>Hist</td> |
---|
87 | <td>This specifies a list of constraints on |
---|
88 | histogram-related parameters, |
---|
89 | which will be of form :h:<var>:n.</td> |
---|
90 | </tr> |
---|
91 | <tr class="row-odd"><td>HAP</td> |
---|
92 | <td>This specifies a list of constraints on parameters |
---|
93 | that are defined for every histogram in each phase |
---|
94 | and are of form p:h:<var>:n.</td> |
---|
95 | </tr> |
---|
96 | <tr class="row-even"><td>Phase</td> |
---|
97 | <td>This specifies a list of constraints on phase |
---|
98 | parameters, |
---|
99 | which will be of form p::<var>:n.</td> |
---|
100 | </tr> |
---|
101 | <tr class="row-odd"><td>Global</td> |
---|
102 | <td>This specifies a list of constraints on parameters |
---|
103 | that are not tied to a histogram or phase and |
---|
104 | are of form ::<var>:n</td> |
---|
105 | </tr> |
---|
106 | </tbody> |
---|
107 | </table> |
---|
108 | <span class="target" id="constraint-definitions-table"></span><p id="index-1">Each constraint is defined as an item in a list. Each constraint is of form:</p> |
---|
109 | <div class="highlight-python"><pre>[[<mult1>, <var1>], [<mult2>, <var2>],..., <fixedval>, <varyflag>, <constype>]</pre> |
---|
110 | </div> |
---|
111 | <p>Where the variable pair list item containing two values [<mult>, <var>], where:</p> |
---|
112 | <blockquote> |
---|
113 | <div><ul> |
---|
114 | <li><p class="first"><mult> is a multiplier for the constraint (float)</p> |
---|
115 | </li> |
---|
116 | <li><dl class="first docutils"> |
---|
117 | <dt><var> a <a class="reference internal" href="#GSASIIobj.G2VarObj" title="GSASIIobj.G2VarObj"><tt class="xref py py-class docutils literal"><span class="pre">G2VarObj</span></tt></a> object (previously a str variable name of form</dt> |
---|
118 | <dd><p class="first last">‘p:h:name[:at]’)</p> |
---|
119 | </dd> |
---|
120 | </dl> |
---|
121 | </li> |
---|
122 | </ul> |
---|
123 | </div></blockquote> |
---|
124 | <p>Note that the last three items in the list play a special role:</p> |
---|
125 | <blockquote> |
---|
126 | <div><ul> |
---|
127 | <li><p class="first"><fixedval> is the fixed value for a <cite>constant equation</cite> (<tt class="docutils literal"><span class="pre">constype=c</span></tt>) |
---|
128 | constraint or is None. For a <cite>New variable</cite> (<tt class="docutils literal"><span class="pre">constype=f</span></tt>) constraint, |
---|
129 | a variable name can be specified as a str (used for externally |
---|
130 | generated constraints)</p> |
---|
131 | </li> |
---|
132 | <li><p class="first"><varyflag> is True or False for <cite>New variable</cite> (<tt class="docutils literal"><span class="pre">constype=f</span></tt>) constraints |
---|
133 | or is None. This will be implemented in the future to indicate if these variables |
---|
134 | should be refined.</p> |
---|
135 | </li> |
---|
136 | <li><p class="first"><constype> is one of four letters, ‘e’, ‘c’, ‘h’, ‘f’ that determines the type of constraint:</p> |
---|
137 | <blockquote> |
---|
138 | <div><ul class="simple"> |
---|
139 | <li>‘e’ defines a set of equivalent variables. Only the first variable is refined (if the |
---|
140 | appropriate refine flag is set) and and all other equivalent variables in the list |
---|
141 | are generated from that variable, using the appropriate multipliers.</li> |
---|
142 | <li>‘c’ defines a constraint equation of form, |
---|
143 | <span class="math">\(m_1 \times var_1 + m_2 \times var_2 + ... = c\)</span></li> |
---|
144 | <li>‘h’ defines a variable to hold (not vary). Any variable on this list is not varied, |
---|
145 | even if its refinement flag is set. Only one [mult,var] pair is allowed in a hold |
---|
146 | constraint and the mult value is ignored. |
---|
147 | This is of particular value when needing to hold one or more variables where a |
---|
148 | single flag controls a set of variables such as, coordinates, |
---|
149 | the reciprocal metric tensor or anisotropic displacement parameter.</li> |
---|
150 | <li>‘f’ defines a new variable (function) according to relationship |
---|
151 | <span class="math">\(newvar = m_1 \times var_1 + m_2 \times var_2 + ...\)</span></li> |
---|
152 | </ul> |
---|
153 | </div></blockquote> |
---|
154 | </li> |
---|
155 | </ul> |
---|
156 | </div></blockquote> |
---|
157 | </div> |
---|
158 | <div class="section" id="covariance-tree-item"> |
---|
159 | <h2>Covariance Tree Item<a class="headerlink" href="#covariance-tree-item" title="Permalink to this headline">¶</a></h2> |
---|
160 | <span class="target" id="covariance-table"></span><p id="index-2">The Covariance tree item has results from the last least-squares run. They |
---|
161 | are stored in a dict with these keys:</p> |
---|
162 | <table border="1" class="docutils"> |
---|
163 | <colgroup> |
---|
164 | <col width="16%" /> |
---|
165 | <col width="19%" /> |
---|
166 | <col width="65%" /> |
---|
167 | </colgroup> |
---|
168 | <thead valign="bottom"> |
---|
169 | <tr class="row-odd"><th class="head">key</th> |
---|
170 | <th class="head">sub-key</th> |
---|
171 | <th class="head">explanation</th> |
---|
172 | </tr> |
---|
173 | </thead> |
---|
174 | <tbody valign="top"> |
---|
175 | <tr class="row-even"><td>newCellDict</td> |
---|
176 | <td></td> |
---|
177 | <td>dict with lattice parameters computed by |
---|
178 | <a class="reference internal" href="GSASIIstruc.html#GSASIIstrMath.GetNewCellParms" title="GSASIIstrMath.GetNewCellParms"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIstrMath.GetNewCellParms()</span></tt></a> (dict)</td> |
---|
179 | </tr> |
---|
180 | <tr class="row-odd"><td>title</td> |
---|
181 | <td></td> |
---|
182 | <td>Name of gpx file(?) (str)</td> |
---|
183 | </tr> |
---|
184 | <tr class="row-even"><td>variables</td> |
---|
185 | <td></td> |
---|
186 | <td>Values for all N refined variables |
---|
187 | (list of float values, length N, |
---|
188 | ordered to match varyList)</td> |
---|
189 | </tr> |
---|
190 | <tr class="row-odd"><td>sig</td> |
---|
191 | <td></td> |
---|
192 | <td>Uncertainty values for all N refined variables |
---|
193 | (list of float values, length N, |
---|
194 | ordered to match varyList)</td> |
---|
195 | </tr> |
---|
196 | <tr class="row-even"><td>varyList</td> |
---|
197 | <td></td> |
---|
198 | <td>List of directly refined variables |
---|
199 | (list of str values, length N)</td> |
---|
200 | </tr> |
---|
201 | <tr class="row-odd"><td>newAtomDict</td> |
---|
202 | <td></td> |
---|
203 | <td>dict with atom position values computed in |
---|
204 | <a class="reference internal" href="GSASIIstruc.html#GSASIIstrMath.ApplyXYZshifts" title="GSASIIstrMath.ApplyXYZshifts"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIstrMath.ApplyXYZshifts()</span></tt></a> (dict)</td> |
---|
205 | </tr> |
---|
206 | <tr class="row-even"><td>Rvals</td> |
---|
207 | <td></td> |
---|
208 | <td>R-factors, GOF, Marquardt value for last |
---|
209 | refinement cycle (dict)</td> |
---|
210 | </tr> |
---|
211 | <tr class="row-odd"><td></td> |
---|
212 | <td>Nobs</td> |
---|
213 | <td>Number of observed data points (int)</td> |
---|
214 | </tr> |
---|
215 | <tr class="row-even"><td></td> |
---|
216 | <td>Rwp</td> |
---|
217 | <td>overall weighted profile R-factor (%, float)</td> |
---|
218 | </tr> |
---|
219 | <tr class="row-odd"><td></td> |
---|
220 | <td>chisq</td> |
---|
221 | <td>sum[w*(Iobs-Icalc)**2] for all data |
---|
222 | note this is not the reduced chi squared (float)</td> |
---|
223 | </tr> |
---|
224 | <tr class="row-even"><td></td> |
---|
225 | <td>lamMax</td> |
---|
226 | <td>Marquardt value applied to Hessian diagonal |
---|
227 | (float)</td> |
---|
228 | </tr> |
---|
229 | <tr class="row-odd"><td></td> |
---|
230 | <td>GOF</td> |
---|
231 | <td>The goodness-of-fit, aka square root of |
---|
232 | the reduced chi squared. (float)</td> |
---|
233 | </tr> |
---|
234 | <tr class="row-even"><td>covMatrix</td> |
---|
235 | <td></td> |
---|
236 | <td>The (NxN) covVariance matrix (np.array)</td> |
---|
237 | </tr> |
---|
238 | </tbody> |
---|
239 | </table> |
---|
240 | </div> |
---|
241 | <div class="section" id="phase-tree-items"> |
---|
242 | <h2>Phase Tree Items<a class="headerlink" href="#phase-tree-items" title="Permalink to this headline">¶</a></h2> |
---|
243 | <span class="target" id="phase-table"></span><p id="index-3">Phase information is stored in the GSAS-II data tree as children of the |
---|
244 | Phases item in a dict with keys:</p> |
---|
245 | <table border="1" class="docutils"> |
---|
246 | <colgroup> |
---|
247 | <col width="10%" /> |
---|
248 | <col width="15%" /> |
---|
249 | <col width="75%" /> |
---|
250 | </colgroup> |
---|
251 | <thead valign="bottom"> |
---|
252 | <tr class="row-odd"><th class="head">key</th> |
---|
253 | <th class="head">sub-key</th> |
---|
254 | <th class="head">explanation</th> |
---|
255 | </tr> |
---|
256 | </thead> |
---|
257 | <tbody valign="top"> |
---|
258 | <tr class="row-even"><td>General</td> |
---|
259 | <td></td> |
---|
260 | <td>Overall information for the phase (dict)</td> |
---|
261 | </tr> |
---|
262 | <tr class="row-odd"><td></td> |
---|
263 | <td>AtomPtrs</td> |
---|
264 | <td>list of four locations to use to pull info |
---|
265 | from the atom records (list)</td> |
---|
266 | </tr> |
---|
267 | <tr class="row-even"><td></td> |
---|
268 | <td>F000X</td> |
---|
269 | <td>x-ray F(000) intensity (float)</td> |
---|
270 | </tr> |
---|
271 | <tr class="row-odd"><td></td> |
---|
272 | <td>F000N</td> |
---|
273 | <td>neutron F(000) intensity (float)</td> |
---|
274 | </tr> |
---|
275 | <tr class="row-even"><td></td> |
---|
276 | <td>Mydir</td> |
---|
277 | <td>directory of current .gpx file (str)</td> |
---|
278 | </tr> |
---|
279 | <tr class="row-odd"><td></td> |
---|
280 | <td>MCSA controls</td> |
---|
281 | <td>Monte Carlo-Simulated Annealing controls (dict)</td> |
---|
282 | </tr> |
---|
283 | <tr class="row-even"><td></td> |
---|
284 | <td>Cell</td> |
---|
285 | <td>List with 8 items: cell refinement flag (bool) |
---|
286 | a, b, c, (Angstrom, float) |
---|
287 | alpha, beta & gamma (degrees, float) |
---|
288 | volume (A^3, float)</td> |
---|
289 | </tr> |
---|
290 | <tr class="row-odd"><td></td> |
---|
291 | <td>Type</td> |
---|
292 | <td>‘nuclear’ or ‘macromolecular’ for now (str)</td> |
---|
293 | </tr> |
---|
294 | <tr class="row-even"><td></td> |
---|
295 | <td>Map</td> |
---|
296 | <td>dict of map parameters</td> |
---|
297 | </tr> |
---|
298 | <tr class="row-odd"><td></td> |
---|
299 | <td>SH Texture</td> |
---|
300 | <td>dict of spherical harmonic preferred orientation |
---|
301 | parameters</td> |
---|
302 | </tr> |
---|
303 | <tr class="row-even"><td></td> |
---|
304 | <td>Isotope</td> |
---|
305 | <td>dict of isotopes for each atom type</td> |
---|
306 | </tr> |
---|
307 | <tr class="row-odd"><td></td> |
---|
308 | <td>Isotopes</td> |
---|
309 | <td>dict of scattering lengths for each isotope |
---|
310 | combination for each element in phase</td> |
---|
311 | </tr> |
---|
312 | <tr class="row-even"><td></td> |
---|
313 | <td>Name</td> |
---|
314 | <td>phase name (str)</td> |
---|
315 | </tr> |
---|
316 | <tr class="row-odd"><td></td> |
---|
317 | <td>SGData</td> |
---|
318 | <td>Space group details as a <a class="reference internal" href="#sgdata-table"><em>space group (SGData) object</em></a> |
---|
319 | as defined in <a class="reference internal" href="GSASIIutil.html#GSASIIspc.SpcGroup" title="GSASIIspc.SpcGroup"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIspc.SpcGroup()</span></tt></a>.</td> |
---|
320 | </tr> |
---|
321 | <tr class="row-even"><td></td> |
---|
322 | <td>Pawley neg wt</td> |
---|
323 | <td>Restraint value for negative Pawley intensities |
---|
324 | (float)</td> |
---|
325 | </tr> |
---|
326 | <tr class="row-odd"><td></td> |
---|
327 | <td>Flip</td> |
---|
328 | <td>dict of Charge flip controls</td> |
---|
329 | </tr> |
---|
330 | <tr class="row-even"><td></td> |
---|
331 | <td>Data plot type</td> |
---|
332 | <td>data plot type (‘Mustrain’, ‘Size’ or |
---|
333 | ‘Preferred orientation’) for powder data (str)</td> |
---|
334 | </tr> |
---|
335 | <tr class="row-odd"><td></td> |
---|
336 | <td>Mass</td> |
---|
337 | <td>Mass of unit cell contents in g/mol</td> |
---|
338 | </tr> |
---|
339 | <tr class="row-even"><td></td> |
---|
340 | <td>POhkl</td> |
---|
341 | <td>March-Dollase preferred orientation direction</td> |
---|
342 | </tr> |
---|
343 | <tr class="row-odd"><td></td> |
---|
344 | <td>Z</td> |
---|
345 | <td>dict of atomic numbers for each atom type</td> |
---|
346 | </tr> |
---|
347 | <tr class="row-even"><td></td> |
---|
348 | <td>vdWRadii</td> |
---|
349 | <td>dict of van der Waals radii for each atom type</td> |
---|
350 | </tr> |
---|
351 | <tr class="row-odd"><td></td> |
---|
352 | <td>Color</td> |
---|
353 | <td>Colors for atoms (list of (r,b,g) triplets)</td> |
---|
354 | </tr> |
---|
355 | <tr class="row-even"><td></td> |
---|
356 | <td>AtomTypes</td> |
---|
357 | <td>List of atom types</td> |
---|
358 | </tr> |
---|
359 | <tr class="row-odd"><td></td> |
---|
360 | <td>AtomMass</td> |
---|
361 | <td>List of masses for atoms</td> |
---|
362 | </tr> |
---|
363 | <tr class="row-even"><td></td> |
---|
364 | <td>doPawley</td> |
---|
365 | <td>Flag for Pawley intensity extraction (bool)</td> |
---|
366 | </tr> |
---|
367 | <tr class="row-odd"><td></td> |
---|
368 | <td>NoAtoms</td> |
---|
369 | <td>Number of atoms per unit cell of each type (dict)</td> |
---|
370 | </tr> |
---|
371 | <tr class="row-even"><td></td> |
---|
372 | <td>Pawley dmin</td> |
---|
373 | <td>maximum Q (as d-space) to use for Pawley |
---|
374 | extraction (float)</td> |
---|
375 | </tr> |
---|
376 | <tr class="row-odd"><td></td> |
---|
377 | <td>BondRadii</td> |
---|
378 | <td>Default radius for each atom used to compute |
---|
379 | interatomic distances (list of floats)</td> |
---|
380 | </tr> |
---|
381 | <tr class="row-even"><td></td> |
---|
382 | <td>AngleRadii</td> |
---|
383 | <td>Default radius for each atom used to compute |
---|
384 | interatomic angles (list of floats)</td> |
---|
385 | </tr> |
---|
386 | <tr class="row-odd"><td></td> |
---|
387 | <td>DisAglCtls</td> |
---|
388 | <td>Dict with distance/angle search controls, |
---|
389 | which has keys ‘Name’, ‘AtomTypes’, |
---|
390 | ‘BondRadii’, ‘AngleRadii’ which are as above |
---|
391 | except are possibly edited. Also contains |
---|
392 | ‘Factors’, which is a 2 element list with |
---|
393 | a multiplier for bond and angle search range |
---|
394 | [typically (0.85,0.85)].</td> |
---|
395 | </tr> |
---|
396 | <tr class="row-even"><td>ranId</td> |
---|
397 | <td></td> |
---|
398 | <td>unique random number Id for phase (int)</td> |
---|
399 | </tr> |
---|
400 | <tr class="row-odd"><td>pId</td> |
---|
401 | <td></td> |
---|
402 | <td>Phase Id number for current project (int).</td> |
---|
403 | </tr> |
---|
404 | <tr class="row-even"><td>Atoms</td> |
---|
405 | <td></td> |
---|
406 | <td>Atoms in phase as a list of lists. The outer list |
---|
407 | is for each atom, the inner list contains varying |
---|
408 | items depending on the type of phase, see |
---|
409 | the <a class="reference internal" href="#atoms-table"><em>Atom Records</em></a> description. |
---|
410 | (list of lists)</td> |
---|
411 | </tr> |
---|
412 | <tr class="row-odd"><td>Drawing</td> |
---|
413 | <td></td> |
---|
414 | <td>Display parameters (dict)</td> |
---|
415 | </tr> |
---|
416 | <tr class="row-even"><td></td> |
---|
417 | <td>ballScale</td> |
---|
418 | <td>Size of spheres in ball-and-stick display (float)</td> |
---|
419 | </tr> |
---|
420 | <tr class="row-odd"><td></td> |
---|
421 | <td>bondList</td> |
---|
422 | <td>dict with bonds</td> |
---|
423 | </tr> |
---|
424 | <tr class="row-even"><td></td> |
---|
425 | <td>contourLevel</td> |
---|
426 | <td>map contour level in e/A^3 (float)</td> |
---|
427 | </tr> |
---|
428 | <tr class="row-odd"><td></td> |
---|
429 | <td>showABC</td> |
---|
430 | <td>Flag to show view point triplet (bool). True=show.</td> |
---|
431 | </tr> |
---|
432 | <tr class="row-even"><td></td> |
---|
433 | <td>viewDir</td> |
---|
434 | <td>cartesian viewing direction (np.array with three |
---|
435 | elements)</td> |
---|
436 | </tr> |
---|
437 | <tr class="row-odd"><td></td> |
---|
438 | <td>Zclip</td> |
---|
439 | <td>clipping distance in A (float)</td> |
---|
440 | </tr> |
---|
441 | <tr class="row-even"><td></td> |
---|
442 | <td>backColor</td> |
---|
443 | <td>background for plot as and R,G,B triplet |
---|
444 | (default = [0, 0, 0], black). |
---|
445 | (list with three atoms)</td> |
---|
446 | </tr> |
---|
447 | <tr class="row-odd"><td></td> |
---|
448 | <td>selectedAtoms</td> |
---|
449 | <td>List of selected atoms (list of int values)</td> |
---|
450 | </tr> |
---|
451 | <tr class="row-even"><td></td> |
---|
452 | <td>showRigidBodies</td> |
---|
453 | <td>Flag to highlight rigid body placement</td> |
---|
454 | </tr> |
---|
455 | <tr class="row-odd"><td></td> |
---|
456 | <td>sizeH</td> |
---|
457 | <td>Size ratio for H atoms (float)</td> |
---|
458 | </tr> |
---|
459 | <tr class="row-even"><td></td> |
---|
460 | <td>bondRadius</td> |
---|
461 | <td>Size of binds in A (float)</td> |
---|
462 | </tr> |
---|
463 | <tr class="row-odd"><td></td> |
---|
464 | <td>atomPtrs</td> |
---|
465 | <td>positions of x, type, site sym, ADP flag in Draw Atoms (list)</td> |
---|
466 | </tr> |
---|
467 | <tr class="row-even"><td></td> |
---|
468 | <td>viewPoint</td> |
---|
469 | <td>list of lists. First item in list is [x,y,z] |
---|
470 | in fractional coordinates for the center of |
---|
471 | the plot. Second item list of previous & current |
---|
472 | atom number viewed (may be [0,0])</td> |
---|
473 | </tr> |
---|
474 | <tr class="row-odd"><td></td> |
---|
475 | <td>showHydrogen</td> |
---|
476 | <td>Flag to control plotting of H atoms.</td> |
---|
477 | </tr> |
---|
478 | <tr class="row-even"><td></td> |
---|
479 | <td>unitCellBox</td> |
---|
480 | <td>Flag to control display of the unit cell.</td> |
---|
481 | </tr> |
---|
482 | <tr class="row-odd"><td></td> |
---|
483 | <td>ellipseProb</td> |
---|
484 | <td>Probability limit for display of thermal |
---|
485 | ellipsoids in % (float).</td> |
---|
486 | </tr> |
---|
487 | <tr class="row-even"><td></td> |
---|
488 | <td>vdwScale</td> |
---|
489 | <td>Multiplier of van der Waals radius for |
---|
490 | display of vdW spheres.</td> |
---|
491 | </tr> |
---|
492 | <tr class="row-odd"><td></td> |
---|
493 | <td>Atoms</td> |
---|
494 | <td>A list of lists with an entry for each atom |
---|
495 | that is plotted.</td> |
---|
496 | </tr> |
---|
497 | <tr class="row-even"><td></td> |
---|
498 | <td>Zstep</td> |
---|
499 | <td>Step to de/increase Z-clip (float)</td> |
---|
500 | </tr> |
---|
501 | <tr class="row-odd"><td></td> |
---|
502 | <td>Quaternion</td> |
---|
503 | <td>Viewing quaternion (4 element np.array)</td> |
---|
504 | </tr> |
---|
505 | <tr class="row-even"><td></td> |
---|
506 | <td>radiusFactor</td> |
---|
507 | <td>Distance ratio for searching for bonds. ? Bonds |
---|
508 | are located that are within r(Ra+Rb) and (Ra+Rb)/r |
---|
509 | where Ra and Rb are the atomic radii.</td> |
---|
510 | </tr> |
---|
511 | <tr class="row-odd"><td></td> |
---|
512 | <td>oldxy</td> |
---|
513 | <td>previous view point (list with two floats)</td> |
---|
514 | </tr> |
---|
515 | <tr class="row-even"><td></td> |
---|
516 | <td>cameraPos</td> |
---|
517 | <td>Viewing position in A for plot (float)</td> |
---|
518 | </tr> |
---|
519 | <tr class="row-odd"><td></td> |
---|
520 | <td>depthFog</td> |
---|
521 | <td>True if use depthFog on plot - set currently as False (bool)</td> |
---|
522 | </tr> |
---|
523 | <tr class="row-even"><td>RBModels</td> |
---|
524 | <td></td> |
---|
525 | <td>Rigid body assignments (note Rigid body definitions |
---|
526 | are stored in their own main top-level tree entry.)</td> |
---|
527 | </tr> |
---|
528 | <tr class="row-odd"><td>Pawley ref</td> |
---|
529 | <td></td> |
---|
530 | <td>Pawley reflections</td> |
---|
531 | </tr> |
---|
532 | <tr class="row-even"><td>Histograms</td> |
---|
533 | <td></td> |
---|
534 | <td>A dict of dicts. The key for the outer dict is |
---|
535 | the histograms tied to this phase. The inner |
---|
536 | dict contains the combined phase/histogram |
---|
537 | parameters for items such as scale factors, |
---|
538 | size and strain parameters. (dict)</td> |
---|
539 | </tr> |
---|
540 | <tr class="row-odd"><td>MCSA</td> |
---|
541 | <td></td> |
---|
542 | <td>Monte-Carlo simulated annealing parameters (dict)</td> |
---|
543 | </tr> |
---|
544 | <tr class="row-even"><td></td> |
---|
545 | <td> </td> |
---|
546 | <td> </td> |
---|
547 | </tr> |
---|
548 | </tbody> |
---|
549 | </table> |
---|
550 | </div> |
---|
551 | <div class="section" id="rigid-body-objects"> |
---|
552 | <h2>Rigid Body Objects<a class="headerlink" href="#rigid-body-objects" title="Permalink to this headline">¶</a></h2> |
---|
553 | <span class="target" id="rbdata-table"></span><p id="index-4">Rigid body descriptions are available for two types of rigid bodies: ‘Vector’ |
---|
554 | and ‘Residue’. Vector rigid bodies are developed by a sequence of translations each |
---|
555 | with a refinable magnitude and Residue rigid bodies are described as Cartesian coordinates |
---|
556 | with defined refinable torsion angles.</p> |
---|
557 | <table border="1" class="docutils"> |
---|
558 | <colgroup> |
---|
559 | <col width="10%" /> |
---|
560 | <col width="14%" /> |
---|
561 | <col width="76%" /> |
---|
562 | </colgroup> |
---|
563 | <thead valign="bottom"> |
---|
564 | <tr class="row-odd"><th class="head">key</th> |
---|
565 | <th class="head">sub-key</th> |
---|
566 | <th class="head">explanation</th> |
---|
567 | </tr> |
---|
568 | </thead> |
---|
569 | <tbody valign="top"> |
---|
570 | <tr class="row-even"><td>Vector</td> |
---|
571 | <td>RBId</td> |
---|
572 | <td>vector rigid bodies (dict of dict)</td> |
---|
573 | </tr> |
---|
574 | <tr class="row-odd"><td></td> |
---|
575 | <td>AtInfo</td> |
---|
576 | <td>Drad, Color: atom drawing radius & color for each atom type (dict)</td> |
---|
577 | </tr> |
---|
578 | <tr class="row-even"><td></td> |
---|
579 | <td>RBname</td> |
---|
580 | <td>Name assigned by user to rigid body (str)</td> |
---|
581 | </tr> |
---|
582 | <tr class="row-odd"><td></td> |
---|
583 | <td>VectMag</td> |
---|
584 | <td>vector magnitudes in A (list)</td> |
---|
585 | </tr> |
---|
586 | <tr class="row-even"><td></td> |
---|
587 | <td>rbXYZ</td> |
---|
588 | <td>Cartesian coordinates for Vector rigid body (list of 3 float)</td> |
---|
589 | </tr> |
---|
590 | <tr class="row-odd"><td></td> |
---|
591 | <td>rbRef</td> |
---|
592 | <td>3 assigned reference atom nos. in rigid body for origin |
---|
593 | definition, use center of atoms flag (list of 3 int & 1 bool)</td> |
---|
594 | </tr> |
---|
595 | <tr class="row-even"><td></td> |
---|
596 | <td>VectRef</td> |
---|
597 | <td>refinement flags for VectMag values (list of bool)</td> |
---|
598 | </tr> |
---|
599 | <tr class="row-odd"><td></td> |
---|
600 | <td>rbTypes</td> |
---|
601 | <td>Atom types for each atom in rigid body (list of str)</td> |
---|
602 | </tr> |
---|
603 | <tr class="row-even"><td></td> |
---|
604 | <td>rbVect</td> |
---|
605 | <td>Cartesian vectors for each translation used to build rigid body (list of lists)</td> |
---|
606 | </tr> |
---|
607 | <tr class="row-odd"><td></td> |
---|
608 | <td>useCount</td> |
---|
609 | <td>Number of times rigid body is used in any structure (int)</td> |
---|
610 | </tr> |
---|
611 | <tr class="row-even"><td>Residue</td> |
---|
612 | <td>RBId</td> |
---|
613 | <td>residue rigid bodies (dict of dict)</td> |
---|
614 | </tr> |
---|
615 | <tr class="row-odd"><td></td> |
---|
616 | <td>AtInfo</td> |
---|
617 | <td>Drad, Color: atom drawing radius & color for each atom type(dict)</td> |
---|
618 | </tr> |
---|
619 | <tr class="row-even"><td></td> |
---|
620 | <td>RBname</td> |
---|
621 | <td>Name assigned by user to rigid body (str)</td> |
---|
622 | </tr> |
---|
623 | <tr class="row-odd"><td></td> |
---|
624 | <td>rbXYZ</td> |
---|
625 | <td>Cartesian coordinates for Residue rigid body (list of 3 float)</td> |
---|
626 | </tr> |
---|
627 | <tr class="row-even"><td></td> |
---|
628 | <td>rbTypes</td> |
---|
629 | <td>Atom types for each atom in rigid body (list of str)</td> |
---|
630 | </tr> |
---|
631 | <tr class="row-odd"><td></td> |
---|
632 | <td>atNames</td> |
---|
633 | <td>Names of each atom in rigid body (e.g. C1,N2...) (list of str)</td> |
---|
634 | </tr> |
---|
635 | <tr class="row-even"><td></td> |
---|
636 | <td>rbRef</td> |
---|
637 | <td>3 assigned reference atom nos. in rigid body for origin |
---|
638 | definition, use center of atoms flag (list of 3 int & 1 bool)</td> |
---|
639 | </tr> |
---|
640 | <tr class="row-odd"><td></td> |
---|
641 | <td>rbSeq</td> |
---|
642 | <td>Orig,Piv,angle,Riding (list): definition of internal rigid body |
---|
643 | torsion; origin atom (int), pivot atom (int), torsion angle (float), |
---|
644 | riding atoms (list of int)</td> |
---|
645 | </tr> |
---|
646 | <tr class="row-even"><td></td> |
---|
647 | <td>SelSeq</td> |
---|
648 | <td>[int,int] used by SeqSizer to identify objects</td> |
---|
649 | </tr> |
---|
650 | <tr class="row-odd"><td></td> |
---|
651 | <td>useCount</td> |
---|
652 | <td>Number of times rigid body is used in any structure (int)</td> |
---|
653 | </tr> |
---|
654 | <tr class="row-even"><td>RBIds</td> |
---|
655 | <td></td> |
---|
656 | <td>unique Ids generated upon creation of each rigid body (dict)</td> |
---|
657 | </tr> |
---|
658 | <tr class="row-odd"><td></td> |
---|
659 | <td>Vector</td> |
---|
660 | <td>Ids for each Vector rigid body (list)</td> |
---|
661 | </tr> |
---|
662 | <tr class="row-even"><td></td> |
---|
663 | <td>Residue</td> |
---|
664 | <td>Ids for each Residue rigid body (list)</td> |
---|
665 | </tr> |
---|
666 | </tbody> |
---|
667 | </table> |
---|
668 | </div> |
---|
669 | <div class="section" id="space-group-objects"> |
---|
670 | <h2>Space Group Objects<a class="headerlink" href="#space-group-objects" title="Permalink to this headline">¶</a></h2> |
---|
671 | <span class="target" id="sgdata-table"></span><p id="index-5">Space groups are interpreted by <a class="reference internal" href="GSASIIutil.html#GSASIIspc.SpcGroup" title="GSASIIspc.SpcGroup"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIspc.SpcGroup()</span></tt></a> |
---|
672 | and the information is placed in a SGdata object |
---|
673 | which is a dict with these keys:</p> |
---|
674 | <table border="1" class="docutils"> |
---|
675 | <colgroup> |
---|
676 | <col width="16%" /> |
---|
677 | <col width="84%" /> |
---|
678 | </colgroup> |
---|
679 | <thead valign="bottom"> |
---|
680 | <tr class="row-odd"><th class="head">key</th> |
---|
681 | <th class="head">explanation</th> |
---|
682 | </tr> |
---|
683 | </thead> |
---|
684 | <tbody valign="top"> |
---|
685 | <tr class="row-even"><td>SpGrp</td> |
---|
686 | <td>space group symbol (str)</td> |
---|
687 | </tr> |
---|
688 | <tr class="row-odd"><td>Laue</td> |
---|
689 | <td>one of the following 14 Laue classes: |
---|
690 | -1, 2/m, mmm, 4/m, 4/mmm, 3R, |
---|
691 | 3mR, 3, 3m1, 31m, 6/m, 6/mmm, m3, m3m (str)</td> |
---|
692 | </tr> |
---|
693 | <tr class="row-even"><td>SGInv</td> |
---|
694 | <td>True if centrosymmetric, False if not (bool)</td> |
---|
695 | </tr> |
---|
696 | <tr class="row-odd"><td>SGLatt</td> |
---|
697 | <td>Lattice centering type. Will be one of |
---|
698 | P, A, B, C, I, F, R (str)</td> |
---|
699 | </tr> |
---|
700 | <tr class="row-even"><td>SGUniq</td> |
---|
701 | <td>unique axis if monoclinic. Will be |
---|
702 | a, b, or c for monoclinic space groups. |
---|
703 | Will be blank for non-monoclinic. (str)</td> |
---|
704 | </tr> |
---|
705 | <tr class="row-odd"><td>SGCen</td> |
---|
706 | <td>Symmetry cell centering vectors. A (n,3) np.array |
---|
707 | of centers. Will always have at least one row: |
---|
708 | <tt class="docutils literal"><span class="pre">np.array([[0,</span> <span class="pre">0,</span> <span class="pre">0]])</span></tt></td> |
---|
709 | </tr> |
---|
710 | <tr class="row-even"><td>SGOps</td> |
---|
711 | <td>symmetry operations as a list of form |
---|
712 | <tt class="docutils literal"><span class="pre">[[M1,T1],</span> <span class="pre">[M2,T2],...]</span></tt> |
---|
713 | where <span class="math">\(M_n\)</span> is a 3x3 np.array |
---|
714 | and <span class="math">\(T_n\)</span> is a length 3 np.array. |
---|
715 | Atom coordinates are transformed where the |
---|
716 | Asymmetric unit coordinates [X is (x,y,z)] |
---|
717 | are transformed using |
---|
718 | <span class="math">\(X^\prime = M_n*X+T_n\)</span></td> |
---|
719 | </tr> |
---|
720 | <tr class="row-odd"><td>SGSys</td> |
---|
721 | <td>symmetry unit cell: type one of |
---|
722 | ‘triclinic’, ‘monoclinic’, ‘orthorhombic’, |
---|
723 | ‘tetragonal’, ‘rhombohedral’, ‘trigonal’, |
---|
724 | ‘hexagonal’, ‘cubic’ (str)</td> |
---|
725 | </tr> |
---|
726 | <tr class="row-even"><td>SGPolax</td> |
---|
727 | <td>Axes for space group polarity. Will be one of |
---|
728 | ‘’, ‘x’, ‘y’, ‘x y’, ‘z’, ‘x z’, ‘y z’, |
---|
729 | ‘xyz’. In the case where axes are arbitrary |
---|
730 | ‘111’ is used (P 1, and ?).</td> |
---|
731 | </tr> |
---|
732 | </tbody> |
---|
733 | </table> |
---|
734 | </div> |
---|
735 | <div class="section" id="atom-records"> |
---|
736 | <h2>Atom Records<a class="headerlink" href="#atom-records" title="Permalink to this headline">¶</a></h2> |
---|
737 | <span class="target" id="atoms-table"></span><p id="index-6">If <tt class="docutils literal"><span class="pre">phasedict</span></tt> points to the phase information in the data tree, then |
---|
738 | atoms are contained in a list of atom records (list) in |
---|
739 | <tt class="docutils literal"><span class="pre">phasedict['Atoms']</span></tt>. Also needed to read atom information |
---|
740 | are four pointers, <tt class="docutils literal"><span class="pre">cx,ct,cs,cia</span> <span class="pre">=</span> <span class="pre">phasedict['General']['atomPtrs']</span></tt>, |
---|
741 | which define locations in the atom record, as shown below. Items shown are |
---|
742 | always present; additional ones for macromolecular phases are marked ‘mm’</p> |
---|
743 | <table border="1" class="docutils"> |
---|
744 | <colgroup> |
---|
745 | <col width="21%" /> |
---|
746 | <col width="79%" /> |
---|
747 | </colgroup> |
---|
748 | <thead valign="bottom"> |
---|
749 | <tr class="row-odd"><th class="head">location</th> |
---|
750 | <th class="head">explanation</th> |
---|
751 | </tr> |
---|
752 | </thead> |
---|
753 | <tbody valign="top"> |
---|
754 | <tr class="row-even"><td>ct-4</td> |
---|
755 | <td>mm - residue number (str)</td> |
---|
756 | </tr> |
---|
757 | <tr class="row-odd"><td>ct-3</td> |
---|
758 | <td>mm - residue name (e.g. ALA) (str)</td> |
---|
759 | </tr> |
---|
760 | <tr class="row-even"><td>ct-2</td> |
---|
761 | <td>mm - chain label (str)</td> |
---|
762 | </tr> |
---|
763 | <tr class="row-odd"><td>ct-1</td> |
---|
764 | <td>atom label (str)</td> |
---|
765 | </tr> |
---|
766 | <tr class="row-even"><td>ct</td> |
---|
767 | <td>atom type (str)</td> |
---|
768 | </tr> |
---|
769 | <tr class="row-odd"><td>ct+1</td> |
---|
770 | <td>refinement flags; combination of ‘F’, ‘X’, ‘U’ (str)</td> |
---|
771 | </tr> |
---|
772 | <tr class="row-even"><td>cx,cx+1,cx+2</td> |
---|
773 | <td>the x,y and z coordinates (3 floats)</td> |
---|
774 | </tr> |
---|
775 | <tr class="row-odd"><td>cs</td> |
---|
776 | <td>site symmetry (str)</td> |
---|
777 | </tr> |
---|
778 | <tr class="row-even"><td>cs+1</td> |
---|
779 | <td>site multiplicity (int)</td> |
---|
780 | </tr> |
---|
781 | <tr class="row-odd"><td>cia</td> |
---|
782 | <td>ADP flag: Isotropic (‘I’) or Anisotropic (‘A’)</td> |
---|
783 | </tr> |
---|
784 | <tr class="row-even"><td>cia+1</td> |
---|
785 | <td>Uiso (float)</td> |
---|
786 | </tr> |
---|
787 | <tr class="row-odd"><td>cia+2...cia+6</td> |
---|
788 | <td>U11, U22, U33, U12, U13, U23 (6 floats)</td> |
---|
789 | </tr> |
---|
790 | <tr class="row-even"><td>atom[-1]</td> |
---|
791 | <td>unique atom identifier (int)</td> |
---|
792 | </tr> |
---|
793 | </tbody> |
---|
794 | </table> |
---|
795 | </div> |
---|
796 | <div class="section" id="drawing-atom-records"> |
---|
797 | <h2>Drawing Atom Records<a class="headerlink" href="#drawing-atom-records" title="Permalink to this headline">¶</a></h2> |
---|
798 | <span class="target" id="drawing-atoms-table"></span><p id="index-7">If <tt class="docutils literal"><span class="pre">phasedict</span></tt> points to the phase information in the data tree, then |
---|
799 | drawing atoms are contained in a list of drawing atom records (list) in |
---|
800 | <tt class="docutils literal"><span class="pre">phasedict['Drawing']['Atoms']</span></tt>. Also needed to read atom information |
---|
801 | are four pointers, <tt class="docutils literal"><span class="pre">cx,ct,cs,ci</span> <span class="pre">=</span> <span class="pre">phasedict['Drawing']['AtomPtrs']</span></tt>, |
---|
802 | which define locations in the atom record, as shown below. Items shown are |
---|
803 | always present; additional ones for macromolecular phases are marked ‘mm’</p> |
---|
804 | <table border="1" class="docutils"> |
---|
805 | <colgroup> |
---|
806 | <col width="17%" /> |
---|
807 | <col width="83%" /> |
---|
808 | </colgroup> |
---|
809 | <thead valign="bottom"> |
---|
810 | <tr class="row-odd"><th class="head">location</th> |
---|
811 | <th class="head">explanation</th> |
---|
812 | </tr> |
---|
813 | </thead> |
---|
814 | <tbody valign="top"> |
---|
815 | <tr class="row-even"><td>ct-4</td> |
---|
816 | <td>mm - residue number (str)</td> |
---|
817 | </tr> |
---|
818 | <tr class="row-odd"><td>ct-3</td> |
---|
819 | <td>mm - residue name (e.g. ALA) (str)</td> |
---|
820 | </tr> |
---|
821 | <tr class="row-even"><td>ct-2</td> |
---|
822 | <td>mm - chain label (str)</td> |
---|
823 | </tr> |
---|
824 | <tr class="row-odd"><td>ct-1</td> |
---|
825 | <td>atom label (str)</td> |
---|
826 | </tr> |
---|
827 | <tr class="row-even"><td>ct</td> |
---|
828 | <td>atom type (str)</td> |
---|
829 | </tr> |
---|
830 | <tr class="row-odd"><td>cx,cx+1,cx+2</td> |
---|
831 | <td>the x,y and z coordinates (3 floats)</td> |
---|
832 | </tr> |
---|
833 | <tr class="row-even"><td>cs-1</td> |
---|
834 | <td>Sym Op symbol; sym. op number + unit cell id (e.g. ‘1,0,-1’) (str)</td> |
---|
835 | </tr> |
---|
836 | <tr class="row-odd"><td>cs</td> |
---|
837 | <td>atom drawing style; e.g. ‘balls & sticks’ (str)</td> |
---|
838 | </tr> |
---|
839 | <tr class="row-even"><td>cs+1</td> |
---|
840 | <td>atom label style (e.g. ‘name’) (str)</td> |
---|
841 | </tr> |
---|
842 | <tr class="row-odd"><td>cs+2</td> |
---|
843 | <td>atom color (RBG triplet) (int)</td> |
---|
844 | </tr> |
---|
845 | <tr class="row-even"><td>cs+3</td> |
---|
846 | <td>ADP flag: Isotropic (‘I’) or Anisotropic (‘A’)</td> |
---|
847 | </tr> |
---|
848 | <tr class="row-odd"><td>cs+4</td> |
---|
849 | <td>Uiso (float)</td> |
---|
850 | </tr> |
---|
851 | <tr class="row-even"><td>cs+5...cs+11</td> |
---|
852 | <td>U11, U22, U33, U12, U13, U23 (6 floats)</td> |
---|
853 | </tr> |
---|
854 | <tr class="row-odd"><td>ci</td> |
---|
855 | <td>unique atom identifier; matches source atom Id in Atom Records (int)</td> |
---|
856 | </tr> |
---|
857 | </tbody> |
---|
858 | </table> |
---|
859 | </div> |
---|
860 | <div class="section" id="powder-diffraction-tree-items"> |
---|
861 | <h2>Powder Diffraction Tree Items<a class="headerlink" href="#powder-diffraction-tree-items" title="Permalink to this headline">¶</a></h2> |
---|
862 | <span class="target" id="powder-table"></span><p id="index-8">Every powder diffraction histogram is stored in the GSAS-II data tree |
---|
863 | with a top-level entry named beginning with the string “PWDR ”. The |
---|
864 | diffraction data for that information are directly associated with |
---|
865 | that tree item and there are a series of children to that item. The |
---|
866 | routines <a class="reference internal" href="GSASII.html#GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree" title="GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree"><tt class="xref py py-func docutils literal"><span class="pre">GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree()</span></tt></a> |
---|
867 | and <a class="reference internal" href="GSASIIstruc.html#GSASIIstrIO.GetUsedHistogramsAndPhases" title="GSASIIstrIO.GetUsedHistogramsAndPhases"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIstrIO.GetUsedHistogramsAndPhases()</span></tt></a> will |
---|
868 | load this information into a dictionary where the child tree name is |
---|
869 | used as a key, and the information in the main entry is assigned |
---|
870 | a key of <tt class="docutils literal"><span class="pre">Data</span></tt>, as outlined below.</p> |
---|
871 | <table border="1" class="docutils"> |
---|
872 | <colgroup> |
---|
873 | <col width="24%" /> |
---|
874 | <col width="17%" /> |
---|
875 | <col width="59%" /> |
---|
876 | </colgroup> |
---|
877 | <thead valign="bottom"> |
---|
878 | <tr class="row-odd"><th class="head">key</th> |
---|
879 | <th class="head">sub-key</th> |
---|
880 | <th class="head">explanation</th> |
---|
881 | </tr> |
---|
882 | </thead> |
---|
883 | <tbody valign="top"> |
---|
884 | <tr class="row-even"><td>Comments</td> |
---|
885 | <td></td> |
---|
886 | <td>Text strings extracted from the original powder |
---|
887 | data header. These cannot be changed by the user; |
---|
888 | it may be empty.</td> |
---|
889 | </tr> |
---|
890 | <tr class="row-odd"><td>Limits</td> |
---|
891 | <td></td> |
---|
892 | <td>A list of two two element lists, as [[Ld,Hd],[L,H]] |
---|
893 | where L and Ld are the current and default lowest |
---|
894 | two-theta value to be used and |
---|
895 | where H and Hd are the current and default highest |
---|
896 | two-theta value to be used.</td> |
---|
897 | </tr> |
---|
898 | <tr class="row-even"><td>Reflection Lists</td> |
---|
899 | <td></td> |
---|
900 | <td>A dict with an entry for each phase in the |
---|
901 | histogram. The contents of each dict item |
---|
902 | is a dict containing reflections, as described in |
---|
903 | the <a class="reference internal" href="#powderrefl-table"><em>Powder Reflections</em></a> |
---|
904 | description.</td> |
---|
905 | </tr> |
---|
906 | <tr class="row-odd"><td>Instrument Parameters</td> |
---|
907 | <td></td> |
---|
908 | <td>A list containing two dicts where the possible |
---|
909 | keys in each dict are listed below. The value |
---|
910 | for each item is a list containing three values: |
---|
911 | the initial value, the current value and a |
---|
912 | refinement flag which can have a value of |
---|
913 | True, False or 0 where 0 indicates a value that |
---|
914 | cannot be refined. The first and second |
---|
915 | values are floats unless otherwise noted. |
---|
916 | Items in the first dict are noted as [1]</td> |
---|
917 | </tr> |
---|
918 | <tr class="row-even"><td></td> |
---|
919 | <td>Lam</td> |
---|
920 | <td>Specifies a wavelength in Angstroms [1]</td> |
---|
921 | </tr> |
---|
922 | <tr class="row-odd"><td></td> |
---|
923 | <td>Lam1</td> |
---|
924 | <td>Specifies the primary wavelength in |
---|
925 | Angstrom, when an alpha1, alpha2 |
---|
926 | source is used [1]</td> |
---|
927 | </tr> |
---|
928 | <tr class="row-even"><td></td> |
---|
929 | <td><p class="first">Lam2</p> |
---|
930 | <p class="last">I(L2)/I(L1)</p> |
---|
931 | </td> |
---|
932 | <td>Specifies the secondary wavelength in |
---|
933 | Angstrom, when an alpha1, alpha2 |
---|
934 | source is used [1] |
---|
935 | Ratio of Lam2 to Lam1 [1]</td> |
---|
936 | </tr> |
---|
937 | <tr class="row-odd"><td></td> |
---|
938 | <td>Type</td> |
---|
939 | <td><dl class="first last docutils"> |
---|
940 | <dt>Histogram type (str) [1]:</dt> |
---|
941 | <dd><ul class="first last simple"> |
---|
942 | <li>‘PXC’ for constant wavelength x-ray</li> |
---|
943 | <li>‘PNC’ for constant wavelength neutron</li> |
---|
944 | <li>‘PNT’ for time of flight neutron</li> |
---|
945 | </ul> |
---|
946 | </dd> |
---|
947 | </dl> |
---|
948 | </td> |
---|
949 | </tr> |
---|
950 | <tr class="row-even"><td></td> |
---|
951 | <td>Zero</td> |
---|
952 | <td>Two-theta zero correction in <em>degrees</em> [1]</td> |
---|
953 | </tr> |
---|
954 | <tr class="row-odd"><td></td> |
---|
955 | <td>Azimuth</td> |
---|
956 | <td>Azimuthal setting angle for data recorded |
---|
957 | with differing setting angles [1]</td> |
---|
958 | </tr> |
---|
959 | <tr class="row-even"><td></td> |
---|
960 | <td>U, V, W</td> |
---|
961 | <td>Cagliotti profile coefficients |
---|
962 | for Gaussian instrumental broadening, where the |
---|
963 | FWHM goes as |
---|
964 | <span class="math">\(U \tan^2\theta + V \tan\theta + W\)</span> [1]</td> |
---|
965 | </tr> |
---|
966 | <tr class="row-odd"><td></td> |
---|
967 | <td>X, Y</td> |
---|
968 | <td>Cauchy (Lorentzian) instrumental broadening |
---|
969 | coefficients [1]</td> |
---|
970 | </tr> |
---|
971 | <tr class="row-even"><td></td> |
---|
972 | <td>SH/L</td> |
---|
973 | <td>Variant of the Finger-Cox-Jephcoat asymmetric |
---|
974 | peak broadening ratio. Note that this is the |
---|
975 | average between S/L and H/L where S is |
---|
976 | sample height, H is the slit height and |
---|
977 | L is the goniometer diameter. [1]</td> |
---|
978 | </tr> |
---|
979 | <tr class="row-odd"><td></td> |
---|
980 | <td>Polariz.</td> |
---|
981 | <td>Polarization coefficient. [1]</td> |
---|
982 | </tr> |
---|
983 | <tr class="row-even"><td>wtFactor</td> |
---|
984 | <td></td> |
---|
985 | <td>A weighting factor to increase or decrease |
---|
986 | the leverage of data in the histogram (float). |
---|
987 | A value of 1.0 weights the data with their |
---|
988 | standard uncertainties and a larger value |
---|
989 | increases the weighting of the data (equivalent |
---|
990 | to decreasing the uncertainties).</td> |
---|
991 | </tr> |
---|
992 | <tr class="row-odd"><td>Sample Parameters</td> |
---|
993 | <td></td> |
---|
994 | <td>Specifies a dict with parameters that describe how |
---|
995 | the data were collected, as listed |
---|
996 | below. Refinable parameters are a list containing |
---|
997 | a float and a bool, where the second value |
---|
998 | specifies if the value is refined, otherwise |
---|
999 | the value is a float unless otherwise noted.</td> |
---|
1000 | </tr> |
---|
1001 | <tr class="row-even"><td></td> |
---|
1002 | <td>Scale</td> |
---|
1003 | <td>The histogram scale factor (refinable)</td> |
---|
1004 | </tr> |
---|
1005 | <tr class="row-odd"><td></td> |
---|
1006 | <td>Absorption</td> |
---|
1007 | <td>The sample absorption coefficient as |
---|
1008 | <span class="math">\(\mu r\)</span> where r is the radius |
---|
1009 | (refinable). Only valid for Debye-Scherrer geometry.</td> |
---|
1010 | </tr> |
---|
1011 | <tr class="row-even"><td></td> |
---|
1012 | <td>SurfaceRoughA</td> |
---|
1013 | <td>Surface roughness parameter A as defined by |
---|
1014 | Surotti,J. Appl. Cryst, 5,325-331, 1972.(refinable - |
---|
1015 | only valid for Bragg-Brentano geometry)</td> |
---|
1016 | </tr> |
---|
1017 | <tr class="row-odd"><td></td> |
---|
1018 | <td>SurfaceRoughB</td> |
---|
1019 | <td>Surface roughness parameter B (refinable - |
---|
1020 | only valid for Bragg-Brentano geometry)</td> |
---|
1021 | </tr> |
---|
1022 | <tr class="row-even"><td></td> |
---|
1023 | <td>DisplaceX, |
---|
1024 | DisplaceY</td> |
---|
1025 | <td>Sample displacement from goniometer center |
---|
1026 | where Y is along the beam direction and |
---|
1027 | X is perpendicular. Units are <span class="math">\(\mu m\)</span> |
---|
1028 | (refinable).</td> |
---|
1029 | </tr> |
---|
1030 | <tr class="row-odd"><td></td> |
---|
1031 | <td>Phi, Chi, |
---|
1032 | Omega</td> |
---|
1033 | <td>Goniometer sample setting angles, in degrees.</td> |
---|
1034 | </tr> |
---|
1035 | <tr class="row-even"><td></td> |
---|
1036 | <td>Gonio. radius</td> |
---|
1037 | <td>Radius of the diffractometer in mm</td> |
---|
1038 | </tr> |
---|
1039 | <tr class="row-odd"><td></td> |
---|
1040 | <td>InstrName</td> |
---|
1041 | <td>A name for the instrument, used in preparing |
---|
1042 | a CIF (str).</td> |
---|
1043 | </tr> |
---|
1044 | <tr class="row-even"><td></td> |
---|
1045 | <td>Force, |
---|
1046 | Temperature, |
---|
1047 | Humidity, |
---|
1048 | Pressure, |
---|
1049 | Voltage</td> |
---|
1050 | <td>Variables that describe how the measurement |
---|
1051 | was performed. Not used directly in |
---|
1052 | any computations.</td> |
---|
1053 | </tr> |
---|
1054 | <tr class="row-odd"><td></td> |
---|
1055 | <td>ranId</td> |
---|
1056 | <td>The random-number Id for the histogram |
---|
1057 | (same value as where top-level key is ranId)</td> |
---|
1058 | </tr> |
---|
1059 | <tr class="row-even"><td></td> |
---|
1060 | <td>Type</td> |
---|
1061 | <td>Type of diffraction data, may be ‘Debye-Scherrer’ |
---|
1062 | or ‘Bragg-Brentano’ (str).</td> |
---|
1063 | </tr> |
---|
1064 | <tr class="row-odd"><td></td> |
---|
1065 | <td>Diffuse</td> |
---|
1066 | <td>not in use?</td> |
---|
1067 | </tr> |
---|
1068 | <tr class="row-even"><td>hId</td> |
---|
1069 | <td></td> |
---|
1070 | <td>The number assigned to the histogram when |
---|
1071 | the project is loaded or edited (can change)</td> |
---|
1072 | </tr> |
---|
1073 | <tr class="row-odd"><td>ranId</td> |
---|
1074 | <td></td> |
---|
1075 | <td>A random number id for the histogram |
---|
1076 | that does not change</td> |
---|
1077 | </tr> |
---|
1078 | <tr class="row-even"><td>Background</td> |
---|
1079 | <td></td> |
---|
1080 | <td>The background is stored as a list with where |
---|
1081 | the first item in the list is list and the second |
---|
1082 | item is a dict. The list contains the background |
---|
1083 | function and its coefficients; the dict contains |
---|
1084 | Debye diffuse terms and background peaks. |
---|
1085 | (TODO: this needs to be expanded.)</td> |
---|
1086 | </tr> |
---|
1087 | <tr class="row-odd"><td>Data</td> |
---|
1088 | <td></td> |
---|
1089 | <td><p class="first">The data consist of a list of 6 np.arrays |
---|
1090 | containing in order:</p> |
---|
1091 | <blockquote class="last"> |
---|
1092 | <div><ol class="arabic simple"> |
---|
1093 | <li>the x-postions (two-theta in degrees),</li> |
---|
1094 | <li>the intensity values (Yobs),</li> |
---|
1095 | <li>the weights for each Yobs value</li> |
---|
1096 | <li>the computed intensity values (Ycalc)</li> |
---|
1097 | <li>the background values</li> |
---|
1098 | <li>Yobs-Ycalc</li> |
---|
1099 | </ol> |
---|
1100 | </div></blockquote> |
---|
1101 | </td> |
---|
1102 | </tr> |
---|
1103 | </tbody> |
---|
1104 | </table> |
---|
1105 | </div> |
---|
1106 | <div class="section" id="powder-reflection-data-structure"> |
---|
1107 | <h2>Powder Reflection Data Structure<a class="headerlink" href="#powder-reflection-data-structure" title="Permalink to this headline">¶</a></h2> |
---|
1108 | <span class="target" id="powderrefl-table"></span><p id="index-9">For every phase in a histogram, the <tt class="docutils literal"><span class="pre">Reflection</span> <span class="pre">Lists</span></tt> value is a dict |
---|
1109 | one element of which is <cite>‘RefList’</cite>, which is a np.array containing |
---|
1110 | reflections. The columns in that array are documented below.</p> |
---|
1111 | <table border="1" class="docutils"> |
---|
1112 | <colgroup> |
---|
1113 | <col width="14%" /> |
---|
1114 | <col width="86%" /> |
---|
1115 | </colgroup> |
---|
1116 | <thead valign="bottom"> |
---|
1117 | <tr class="row-odd"><th class="head">index</th> |
---|
1118 | <th class="head">explanation</th> |
---|
1119 | </tr> |
---|
1120 | </thead> |
---|
1121 | <tbody valign="top"> |
---|
1122 | <tr class="row-even"><td>0,1,2</td> |
---|
1123 | <td>h,k,l (float)</td> |
---|
1124 | </tr> |
---|
1125 | <tr class="row-odd"><td>3</td> |
---|
1126 | <td>multiplicity</td> |
---|
1127 | </tr> |
---|
1128 | <tr class="row-even"><td>4</td> |
---|
1129 | <td>d-space, Angstrom</td> |
---|
1130 | </tr> |
---|
1131 | <tr class="row-odd"><td>5</td> |
---|
1132 | <td>pos, two-theta</td> |
---|
1133 | </tr> |
---|
1134 | <tr class="row-even"><td>6</td> |
---|
1135 | <td>sig, Gaussian width</td> |
---|
1136 | </tr> |
---|
1137 | <tr class="row-odd"><td>7</td> |
---|
1138 | <td>gam, Lorenzian width</td> |
---|
1139 | </tr> |
---|
1140 | <tr class="row-even"><td>8</td> |
---|
1141 | <td><span class="math">\(F_{obs}^2\)</span></td> |
---|
1142 | </tr> |
---|
1143 | <tr class="row-odd"><td>9</td> |
---|
1144 | <td><span class="math">\(F_{calc}^2\)</span></td> |
---|
1145 | </tr> |
---|
1146 | <tr class="row-even"><td>10</td> |
---|
1147 | <td>reflection phase, in degrees</td> |
---|
1148 | </tr> |
---|
1149 | <tr class="row-odd"><td>11</td> |
---|
1150 | <td>intensity correction for reflection, this times |
---|
1151 | <span class="math">\(F_{obs}^2\)</span> or <span class="math">\(F_{calc}^2\)</span> gives Iobs or Icalc</td> |
---|
1152 | </tr> |
---|
1153 | </tbody> |
---|
1154 | </table> |
---|
1155 | </div> |
---|
1156 | <div class="section" id="single-crystal-tree-items"> |
---|
1157 | <h2>Single Crystal Tree Items<a class="headerlink" href="#single-crystal-tree-items" title="Permalink to this headline">¶</a></h2> |
---|
1158 | <span class="target" id="xtal-table"></span><p id="index-10">Every single crystal diffraction histogram is stored in the GSAS-II data tree |
---|
1159 | with a top-level entry named beginning with the string “HKLF ”. The |
---|
1160 | diffraction data for that information are directly associated with |
---|
1161 | that tree item and there are a series of children to that item. The |
---|
1162 | routines <a class="reference internal" href="GSASII.html#GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree" title="GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree"><tt class="xref py py-func docutils literal"><span class="pre">GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree()</span></tt></a> |
---|
1163 | and <a class="reference internal" href="GSASIIstruc.html#GSASIIstrIO.GetUsedHistogramsAndPhases" title="GSASIIstrIO.GetUsedHistogramsAndPhases"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIstrIO.GetUsedHistogramsAndPhases()</span></tt></a> will |
---|
1164 | load this information into a dictionary where the child tree name is |
---|
1165 | used as a key, and the information in the main entry is assigned |
---|
1166 | a key of <tt class="docutils literal"><span class="pre">Data</span></tt>, as outlined below.</p> |
---|
1167 | <table border="1" class="docutils"> |
---|
1168 | <colgroup> |
---|
1169 | <col width="25%" /> |
---|
1170 | <col width="17%" /> |
---|
1171 | <col width="58%" /> |
---|
1172 | </colgroup> |
---|
1173 | <thead valign="bottom"> |
---|
1174 | <tr class="row-odd"><th class="head">key</th> |
---|
1175 | <th class="head">sub-key</th> |
---|
1176 | <th class="head">explanation</th> |
---|
1177 | </tr> |
---|
1178 | </thead> |
---|
1179 | <tbody valign="top"> |
---|
1180 | <tr class="row-even"><td>Data</td> |
---|
1181 | <td></td> |
---|
1182 | <td>A dict that contains the |
---|
1183 | reflection table, |
---|
1184 | as described in the |
---|
1185 | <a class="reference internal" href="#xtalrefl-table"><em>Single Crystal Reflections</em></a> |
---|
1186 | description.</td> |
---|
1187 | </tr> |
---|
1188 | <tr class="row-odd"><td>Instrument Parameters</td> |
---|
1189 | <td></td> |
---|
1190 | <td>A list containing two dicts where the possible |
---|
1191 | keys in each dict are listed below. The value |
---|
1192 | for most items is a list containing two values: |
---|
1193 | the initial value, the current value. |
---|
1194 | The first and second |
---|
1195 | values are floats unless otherwise noted.</td> |
---|
1196 | </tr> |
---|
1197 | <tr class="row-even"><td></td> |
---|
1198 | <td>Lam</td> |
---|
1199 | <td>Specifies a wavelength in Angstroms (two floats)</td> |
---|
1200 | </tr> |
---|
1201 | <tr class="row-odd"><td></td> |
---|
1202 | <td>Type</td> |
---|
1203 | <td><dl class="first last docutils"> |
---|
1204 | <dt>Histogram type (two str values):</dt> |
---|
1205 | <dd><ul class="first last simple"> |
---|
1206 | <li>‘SXC’ for constant wavelength x-ray</li> |
---|
1207 | <li>‘SNC’ for constant wavelength neutron</li> |
---|
1208 | <li>‘SNT’ for time of flight neutron</li> |
---|
1209 | </ul> |
---|
1210 | </dd> |
---|
1211 | </dl> |
---|
1212 | </td> |
---|
1213 | </tr> |
---|
1214 | <tr class="row-even"><td></td> |
---|
1215 | <td>InstrName</td> |
---|
1216 | <td>A name for the instrument, used in preparing |
---|
1217 | a CIF (str).</td> |
---|
1218 | </tr> |
---|
1219 | <tr class="row-odd"><td>wtFactor</td> |
---|
1220 | <td></td> |
---|
1221 | <td>A weighting factor to increase or decrease |
---|
1222 | the leverage of data in the histogram (float). |
---|
1223 | A value of 1.0 weights the data with their |
---|
1224 | standard uncertainties and a larger value |
---|
1225 | increases the weighting of the data (equivalent |
---|
1226 | to decreasing the uncertainties).</td> |
---|
1227 | </tr> |
---|
1228 | <tr class="row-even"><td>hId</td> |
---|
1229 | <td></td> |
---|
1230 | <td>The number assigned to the histogram when |
---|
1231 | the project is loaded or edited (can change)</td> |
---|
1232 | </tr> |
---|
1233 | <tr class="row-odd"><td>ranId</td> |
---|
1234 | <td></td> |
---|
1235 | <td>A random number id for the histogram |
---|
1236 | that does not change</td> |
---|
1237 | </tr> |
---|
1238 | </tbody> |
---|
1239 | </table> |
---|
1240 | </div> |
---|
1241 | <div class="section" id="single-crystal-reflection-data-structure"> |
---|
1242 | <h2>Single Crystal Reflection Data Structure<a class="headerlink" href="#single-crystal-reflection-data-structure" title="Permalink to this headline">¶</a></h2> |
---|
1243 | <span class="target" id="xtalrefl-table"></span><p id="index-11">For every single crystal a histogram, the <tt class="docutils literal"><span class="pre">'Data'</span></tt> item contains |
---|
1244 | the structure factors as an np.array in item <cite>‘RefList’</cite>. |
---|
1245 | The columns in that array are documented below.</p> |
---|
1246 | <table border="1" class="docutils"> |
---|
1247 | <colgroup> |
---|
1248 | <col width="16%" /> |
---|
1249 | <col width="84%" /> |
---|
1250 | </colgroup> |
---|
1251 | <thead valign="bottom"> |
---|
1252 | <tr class="row-odd"><th class="head">index</th> |
---|
1253 | <th class="head">explanation</th> |
---|
1254 | </tr> |
---|
1255 | </thead> |
---|
1256 | <tbody valign="top"> |
---|
1257 | <tr class="row-even"><td>0,1,2</td> |
---|
1258 | <td>h,k,l (float)</td> |
---|
1259 | </tr> |
---|
1260 | <tr class="row-odd"><td>3</td> |
---|
1261 | <td>multiplicity</td> |
---|
1262 | </tr> |
---|
1263 | <tr class="row-even"><td>4</td> |
---|
1264 | <td>d-space, Angstrom</td> |
---|
1265 | </tr> |
---|
1266 | <tr class="row-odd"><td>5</td> |
---|
1267 | <td><span class="math">\(F_{obs}^2\)</span></td> |
---|
1268 | </tr> |
---|
1269 | <tr class="row-even"><td>6</td> |
---|
1270 | <td><span class="math">\(\sigma(F_{obs}^2)\)</span></td> |
---|
1271 | </tr> |
---|
1272 | <tr class="row-odd"><td>7</td> |
---|
1273 | <td><span class="math">\(F_{calc}^2\)</span></td> |
---|
1274 | </tr> |
---|
1275 | <tr class="row-even"><td>8</td> |
---|
1276 | <td><span class="math">\(F_{obs}^2T\)</span></td> |
---|
1277 | </tr> |
---|
1278 | <tr class="row-odd"><td>9</td> |
---|
1279 | <td><span class="math">\(F_{calc}^2T\)</span></td> |
---|
1280 | </tr> |
---|
1281 | <tr class="row-even"><td>10</td> |
---|
1282 | <td>reflection phase, in degrees</td> |
---|
1283 | </tr> |
---|
1284 | <tr class="row-odd"><td>11</td> |
---|
1285 | <td>intensity correction for reflection, this times |
---|
1286 | <span class="math">\(F_{obs}^2\)</span> or <span class="math">\(F_{calc}^2\)</span> |
---|
1287 | gives Iobs or Icalc</td> |
---|
1288 | </tr> |
---|
1289 | </tbody> |
---|
1290 | </table> |
---|
1291 | </div> |
---|
1292 | <div class="section" id="image-data-structure"> |
---|
1293 | <h2>Image Data Structure<a class="headerlink" href="#image-data-structure" title="Permalink to this headline">¶</a></h2> |
---|
1294 | <span class="target" id="image-table"></span><p id="index-12">Every 2-dimensional image is stored in the GSAS-II data tree |
---|
1295 | with a top-level entry named beginning with the string “IMG ”. The |
---|
1296 | image data are directly associated with that tree item and there |
---|
1297 | are a series of children to that item. The routines <a class="reference internal" href="GSASII.html#GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree" title="GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree"><tt class="xref py py-func docutils literal"><span class="pre">GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree()</span></tt></a> |
---|
1298 | and <a class="reference internal" href="GSASIIstruc.html#GSASIIstrIO.GetUsedHistogramsAndPhases" title="GSASIIstrIO.GetUsedHistogramsAndPhases"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIstrIO.GetUsedHistogramsAndPhases()</span></tt></a> will |
---|
1299 | load this information into a dictionary where the child tree name is |
---|
1300 | used as a key, and the information in the main entry is assigned |
---|
1301 | a key of <tt class="docutils literal"><span class="pre">Data</span></tt>, as outlined below.</p> |
---|
1302 | <table border="1" class="docutils"> |
---|
1303 | <colgroup> |
---|
1304 | <col width="12%" /> |
---|
1305 | <col width="12%" /> |
---|
1306 | <col width="76%" /> |
---|
1307 | </colgroup> |
---|
1308 | <thead valign="bottom"> |
---|
1309 | <tr class="row-odd"><th class="head">key</th> |
---|
1310 | <th class="head">sub-key</th> |
---|
1311 | <th class="head">explanation</th> |
---|
1312 | </tr> |
---|
1313 | </thead> |
---|
1314 | <tbody valign="top"> |
---|
1315 | <tr class="row-even"><td>Comments</td> |
---|
1316 | <td></td> |
---|
1317 | <td>Text strings extracted from the original image data |
---|
1318 | header or a metafile. These cannot be changed by |
---|
1319 | the user; it may be empty.</td> |
---|
1320 | </tr> |
---|
1321 | <tr class="row-odd"><td>Image Controls</td> |
---|
1322 | <td>azmthOff</td> |
---|
1323 | <td>(float) The offset to be applied to an azimuthal |
---|
1324 | value. Accomodates |
---|
1325 | detector orientations other than with the detector |
---|
1326 | X-axis |
---|
1327 | horizontal.</td> |
---|
1328 | </tr> |
---|
1329 | <tr class="row-even"><td></td> |
---|
1330 | <td>background image</td> |
---|
1331 | <td>(list:str,float) The name of a tree item (“IMG ...”) that is to be subtracted |
---|
1332 | during image integration multiplied by value. It must have the same size/shape as |
---|
1333 | the integrated image. NB: value < 0 for subtraction.</td> |
---|
1334 | </tr> |
---|
1335 | <tr class="row-odd"><td></td> |
---|
1336 | <td>calibrant</td> |
---|
1337 | <td>(str) The material used for determining the position/orientation |
---|
1338 | of the image. The data is obtained from <a class="reference internal" href="GSASIIutil.html#module-ImageCalibrants" title="ImageCalibrants"><tt class="xref py py-func docutils literal"><span class="pre">ImageCalibrants()</span></tt></a> |
---|
1339 | and UserCalibrants.py (supplied by user).</td> |
---|
1340 | </tr> |
---|
1341 | <tr class="row-even"><td></td> |
---|
1342 | <td>calibdmin</td> |
---|
1343 | <td>(float) The minimum d-spacing used during the last calibration run.</td> |
---|
1344 | </tr> |
---|
1345 | <tr class="row-odd"><td></td> |
---|
1346 | <td>calibskip</td> |
---|
1347 | <td>(int) The number of expected diffraction lines skipped during the last |
---|
1348 | calibration run.</td> |
---|
1349 | </tr> |
---|
1350 | <tr class="row-even"><td></td> |
---|
1351 | <td>center</td> |
---|
1352 | <td>(list:floats) The [X,Y] point in detector coordinates (mm) where the direct beam |
---|
1353 | strikes the detector plane as determined by calibration. This point |
---|
1354 | does not have to be within the limits of the detector boundaries.</td> |
---|
1355 | </tr> |
---|
1356 | <tr class="row-odd"><td></td> |
---|
1357 | <td>centerAzm</td> |
---|
1358 | <td>(bool) If True then the azimuth reported for the integrated slice |
---|
1359 | of the image is at the center line otherwise it is at the leading edge.</td> |
---|
1360 | </tr> |
---|
1361 | <tr class="row-even"><td></td> |
---|
1362 | <td>color</td> |
---|
1363 | <td>(str) The name of the colormap used to display the image. Default = ‘Paired’.</td> |
---|
1364 | </tr> |
---|
1365 | <tr class="row-odd"><td></td> |
---|
1366 | <td>cutoff</td> |
---|
1367 | <td>(float) The minimum value of I/Ib for a point selected in a diffraction ring for |
---|
1368 | calibration calculations. See pixLimit for details as how point is found.</td> |
---|
1369 | </tr> |
---|
1370 | <tr class="row-even"><td></td> |
---|
1371 | <td>DetDepth</td> |
---|
1372 | <td>(float) Coefficient for penetration correction to distance; accounts for diffraction |
---|
1373 | ring offset at higher angles. Optionally determined by calibration.</td> |
---|
1374 | </tr> |
---|
1375 | <tr class="row-odd"><td></td> |
---|
1376 | <td>DetDepthRef</td> |
---|
1377 | <td>(bool) If True then refine DetDepth during calibration/recalibration calculation.</td> |
---|
1378 | </tr> |
---|
1379 | <tr class="row-even"><td></td> |
---|
1380 | <td>distance</td> |
---|
1381 | <td>(float) The distance (mm) from sample to detector plane.</td> |
---|
1382 | </tr> |
---|
1383 | <tr class="row-odd"><td></td> |
---|
1384 | <td>ellipses</td> |
---|
1385 | <td>(list:lists) Each object in ellipses is a list [center,phi,radii,color] where |
---|
1386 | center (list) is location (mm) of the ellipse center on the detector plane, phi is the |
---|
1387 | rotation of the ellipse minor axis from the x-axis, and radii are the minor & major |
---|
1388 | radii of the ellipse. If radii[0] is negative then parameters describe a hyperbola. Color |
---|
1389 | is the selected drawing color (one of ‘b’, ‘g’ ,’r’) for the ellipse/hyperbola.</td> |
---|
1390 | </tr> |
---|
1391 | <tr class="row-even"><td></td> |
---|
1392 | <td>edgemin</td> |
---|
1393 | <td>(float) Not used; parameter in EdgeFinder code.</td> |
---|
1394 | </tr> |
---|
1395 | <tr class="row-odd"><td></td> |
---|
1396 | <td>fullIntegrate</td> |
---|
1397 | <td>(bool) If True then integrate over full 360 deg azimuthal range.</td> |
---|
1398 | </tr> |
---|
1399 | <tr class="row-even"><td></td> |
---|
1400 | <td>GonioAngles</td> |
---|
1401 | <td>(list:floats) The ‘Omega’,’Chi’,’Phi’ goniometer angles used for this image. |
---|
1402 | Required for texture calculations.</td> |
---|
1403 | </tr> |
---|
1404 | <tr class="row-odd"><td></td> |
---|
1405 | <td>invert_x</td> |
---|
1406 | <td>(bool) If True display the image with the x-axis inverted.</td> |
---|
1407 | </tr> |
---|
1408 | <tr class="row-even"><td></td> |
---|
1409 | <td>invert_y</td> |
---|
1410 | <td>(bool) If True display the image with the y-axis inverted.</td> |
---|
1411 | </tr> |
---|
1412 | <tr class="row-odd"><td></td> |
---|
1413 | <td>IOtth</td> |
---|
1414 | <td>(list:floats) The minimum and maximum 2-theta values to be used for integration.</td> |
---|
1415 | </tr> |
---|
1416 | <tr class="row-even"><td></td> |
---|
1417 | <td>LRazimuth</td> |
---|
1418 | <td>(list:floats) The minimum and maximum azimuth values to be used for integration.</td> |
---|
1419 | </tr> |
---|
1420 | <tr class="row-odd"><td></td> |
---|
1421 | <td>Oblique</td> |
---|
1422 | <td>(list:float,bool) If True apply a detector absorption correction using the value to the |
---|
1423 | intensities obtained during integration.</td> |
---|
1424 | </tr> |
---|
1425 | <tr class="row-even"><td></td> |
---|
1426 | <td>outAzimuths</td> |
---|
1427 | <td>(int) The number of azimuth pie slices.</td> |
---|
1428 | </tr> |
---|
1429 | <tr class="row-odd"><td></td> |
---|
1430 | <td>outChannels</td> |
---|
1431 | <td>(int) The number of 2-theta steps.</td> |
---|
1432 | </tr> |
---|
1433 | <tr class="row-even"><td></td> |
---|
1434 | <td>pixelSize</td> |
---|
1435 | <td>(list:ints) The X,Y dimensions (microns) of each pixel.</td> |
---|
1436 | </tr> |
---|
1437 | <tr class="row-odd"><td></td> |
---|
1438 | <td>pixLimit</td> |
---|
1439 | <td>(int) A box in the image with 2*pixLimit+1 edges is searched to find the maximum. |
---|
1440 | This value (I) along with the minimum (Ib) in the box is reported by <a class="reference internal" href="GSASIIimage.html#GSASIIimage.ImageLocalMax" title="GSASIIimage.ImageLocalMax"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIimage.ImageLocalMax()</span></tt></a> |
---|
1441 | and subject to cutoff in <a class="reference internal" href="GSASIIimage.html#GSASIIimage.makeRing" title="GSASIIimage.makeRing"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIimage.makeRing()</span></tt></a>. |
---|
1442 | Locations are used to construct rings of points for calibration calcualtions.</td> |
---|
1443 | </tr> |
---|
1444 | <tr class="row-even"><td></td> |
---|
1445 | <td>PolaVal</td> |
---|
1446 | <td>(list:float,bool) If type=’SASD’ and if True, apply polarization correction to intensities from |
---|
1447 | integration using value.</td> |
---|
1448 | </tr> |
---|
1449 | <tr class="row-odd"><td></td> |
---|
1450 | <td>rings</td> |
---|
1451 | <td>(list:lists) Each entry is [X,Y,dsp] where X & Y are lists of x,y coordinates around a |
---|
1452 | diffraction ring with the same d-spacing (dsp)</td> |
---|
1453 | </tr> |
---|
1454 | <tr class="row-even"><td></td> |
---|
1455 | <td>ring</td> |
---|
1456 | <td>(list) The x,y coordinates of the >5 points on an inner ring |
---|
1457 | selected by the user,</td> |
---|
1458 | </tr> |
---|
1459 | <tr class="row-odd"><td></td> |
---|
1460 | <td>Range</td> |
---|
1461 | <td>(list) The minimum & maximum values of the image</td> |
---|
1462 | </tr> |
---|
1463 | <tr class="row-even"><td></td> |
---|
1464 | <td>rotation</td> |
---|
1465 | <td>(float) The angle between the x-axis and the vector about which the |
---|
1466 | detector is tilted. Constrained to -180 to 180 deg.</td> |
---|
1467 | </tr> |
---|
1468 | <tr class="row-odd"><td></td> |
---|
1469 | <td>SampleShape</td> |
---|
1470 | <td>(str) Currently only ‘Cylinder’. Sample shape for Debye-Scherrer experiments; used for absorption |
---|
1471 | calculations.</td> |
---|
1472 | </tr> |
---|
1473 | <tr class="row-even"><td></td> |
---|
1474 | <td>SampleAbs</td> |
---|
1475 | <td>(list: float,bool) Value of absorption coefficient for Debye-Scherrer experimnents, flag if True |
---|
1476 | to cause correction to be applied.</td> |
---|
1477 | </tr> |
---|
1478 | <tr class="row-odd"><td></td> |
---|
1479 | <td>setDefault</td> |
---|
1480 | <td>(bool) If True the use the image controls values for all new images to be read. (might be removed)</td> |
---|
1481 | </tr> |
---|
1482 | <tr class="row-even"><td></td> |
---|
1483 | <td>setRings</td> |
---|
1484 | <td>(bool) If True then display all the selected x,y ring positions (vida supra rings) used in the calibration.</td> |
---|
1485 | </tr> |
---|
1486 | <tr class="row-odd"><td></td> |
---|
1487 | <td>showLines</td> |
---|
1488 | <td>(bool) If True then isplay the integration limits to be used.</td> |
---|
1489 | </tr> |
---|
1490 | <tr class="row-even"><td></td> |
---|
1491 | <td>size</td> |
---|
1492 | <td>(list:int) The number of pixels on the image x & y axes</td> |
---|
1493 | </tr> |
---|
1494 | <tr class="row-odd"><td></td> |
---|
1495 | <td>type</td> |
---|
1496 | <td>(str) One of ‘PWDR’, ‘SASD’ or ‘REFL’ for powder, small angle or reflectometry data, respectively.</td> |
---|
1497 | </tr> |
---|
1498 | <tr class="row-even"><td></td> |
---|
1499 | <td>tilt</td> |
---|
1500 | <td>(float) The angle the detector normal makes with the incident beam; range -90 to 90.</td> |
---|
1501 | </tr> |
---|
1502 | <tr class="row-odd"><td></td> |
---|
1503 | <td>wavelength</td> |
---|
1504 | <td>(float) Tha radiation wavelength (Angstroms) as entered by the user (or someday obtained from the image header).</td> |
---|
1505 | </tr> |
---|
1506 | <tr class="row-even"><td>Masks</td> |
---|
1507 | <td>Arcs</td> |
---|
1508 | <td>(list: lists) Each entry [2-theta,[azimuth[0],azimuth[1]],thickness] describes an arc mask |
---|
1509 | to be excluded from integration</td> |
---|
1510 | </tr> |
---|
1511 | <tr class="row-odd"><td></td> |
---|
1512 | <td>Frames</td> |
---|
1513 | <td>(list:lists) Each entry describes the x,y points (3 or more - mm) that describe a frame outside |
---|
1514 | of which is excluded from recalibration and integration. Only one frame is allowed.</td> |
---|
1515 | </tr> |
---|
1516 | <tr class="row-even"><td></td> |
---|
1517 | <td>Points</td> |
---|
1518 | <td>(list:lists) Each entry [x,y,radius] (mm) describes an excluded spot on the image to be excluded |
---|
1519 | from integration.</td> |
---|
1520 | </tr> |
---|
1521 | <tr class="row-odd"><td></td> |
---|
1522 | <td>Polygons</td> |
---|
1523 | <td>(list:lists) Each entry is a list of 3+ [x,y] points (mm) that describe a polygon on the image |
---|
1524 | to be excluded from integration.</td> |
---|
1525 | </tr> |
---|
1526 | <tr class="row-even"><td></td> |
---|
1527 | <td>Rings</td> |
---|
1528 | <td>(list: lists) Each entry [2-theta,thickness] describes a ring mask |
---|
1529 | to be excluded from integration.</td> |
---|
1530 | </tr> |
---|
1531 | <tr class="row-odd"><td></td> |
---|
1532 | <td>Thresholds</td> |
---|
1533 | <td>(list:[tuple,list]) [(Imin,Imax),[Imin,Imax]] This gives lower and upper limits for points on the image to be included |
---|
1534 | in integrsation. The tuple is the image intensity limits and the list are those set by the user.</td> |
---|
1535 | </tr> |
---|
1536 | <tr class="row-even"><td>Stress/Strain</td> |
---|
1537 | <td>Sample phi</td> |
---|
1538 | <td>(float) Sample rotation about vertical axis.</td> |
---|
1539 | </tr> |
---|
1540 | <tr class="row-odd"><td></td> |
---|
1541 | <td>Sample z</td> |
---|
1542 | <td>(float) Sample translation from the calibration sample position (for Sample phi = 0)</td> |
---|
1543 | </tr> |
---|
1544 | <tr class="row-even"><td></td> |
---|
1545 | <td>strain</td> |
---|
1546 | <td>(list: 3x3 array of float) The strain tensor coefficients [[‘ e11’,’e12’,’e13’],[‘ e21’,’e22’,’e23’],[‘ e31’,’e32’,’e33’]]. |
---|
1547 | These will be restricted by space group symmetry; result of strain fit refinement.</td> |
---|
1548 | </tr> |
---|
1549 | <tr class="row-odd"><td></td> |
---|
1550 | <td>Type</td> |
---|
1551 | <td>(str) ‘True’ or ‘Conventional’: The strain model used for the calculation.</td> |
---|
1552 | </tr> |
---|
1553 | <tr class="row-even"><td></td> |
---|
1554 | <td>d-zero</td> |
---|
1555 | <td>(list:dict) Each item is for a diffraction ring on the image; all items are from the same phase and are used to determine the strain tensor. |
---|
1556 | The dictionary items are: |
---|
1557 | ‘Dset’: (float) True d-spacing for the diffraction ring; entered by the user. |
---|
1558 | ‘Dcalc’: (float) d-spacing... |
---|
1559 | ‘pixLimit’: (int) Search range to find highest point on ring for each data point |
---|
1560 | ‘cutoff’: (float) I/Ib cutoff for searching. |
---|
1561 | ‘ImxyObs’: (list:lists) [[X],[Y]] observed points to be used for strain calculations. |
---|
1562 | ‘ImxyCalc’:(list:lists) [[X],[Y]] calculated points based on refined strain.</td> |
---|
1563 | </tr> |
---|
1564 | </tbody> |
---|
1565 | </table> |
---|
1566 | </div> |
---|
1567 | <div class="section" id="parameter-dictionary"> |
---|
1568 | <h2>Parameter Dictionary<a class="headerlink" href="#parameter-dictionary" title="Permalink to this headline">¶</a></h2> |
---|
1569 | <span class="target" id="parmdict-table"></span><p id="index-13">The parameter dictionary contains all of the variable parameters for the refinement. |
---|
1570 | The dictionary keys are the name of the parameter (<phase>:<hist>:<name>:<atom>). |
---|
1571 | It is prepared in two ways. When loaded from the tree |
---|
1572 | (in <a class="reference internal" href="GSASII.html#GSASII.GSASII.MakeLSParmDict" title="GSASII.GSASII.MakeLSParmDict"><tt class="xref py py-meth docutils literal"><span class="pre">GSASII.GSASII.MakeLSParmDict()</span></tt></a> and |
---|
1573 | <a class="reference internal" href="GSASIIGUIr.html#GSASIIIO.ExportBaseclass.loadParmDict" title="GSASIIIO.ExportBaseclass.loadParmDict"><tt class="xref py py-meth docutils literal"><span class="pre">GSASIIIO.ExportBaseclass.loadParmDict()</span></tt></a>), |
---|
1574 | the values are lists with two elements: <tt class="docutils literal"><span class="pre">[value,</span> <span class="pre">refine</span> <span class="pre">flag]</span></tt></p> |
---|
1575 | <p>When loaded from the GPX file (in |
---|
1576 | <a class="reference internal" href="GSASIIstruc.html#GSASIIstrMain.Refine" title="GSASIIstrMain.Refine"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIstrMain.Refine()</span></tt></a> and <a class="reference internal" href="GSASIIstruc.html#GSASIIstrMain.SeqRefine" title="GSASIIstrMain.SeqRefine"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIstrMain.SeqRefine()</span></tt></a>), the value in the |
---|
1577 | dict is the actual parameter value (usually a float, but sometimes a |
---|
1578 | letter or string flag value (such as I or A for iso/anisotropic).</p> |
---|
1579 | </div> |
---|
1580 | <div class="section" id="classes-and-routines"> |
---|
1581 | <h2><em>Classes and routines</em><a class="headerlink" href="#classes-and-routines" title="Permalink to this headline">¶</a></h2> |
---|
1582 | <dl class="data"> |
---|
1583 | <dt id="GSASIIobj.AtomIdLookup"> |
---|
1584 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">AtomIdLookup</tt><em class="property"> = {}</em><a class="headerlink" href="#GSASIIobj.AtomIdLookup" title="Permalink to this definition">¶</a></dt> |
---|
1585 | <dd><p>dict listing for each phase index as a str, the atom label and atom random Id, |
---|
1586 | keyed by atom sequential index as a str; |
---|
1587 | best to access this using <a class="reference internal" href="#GSASIIobj.LookupAtomLabel" title="GSASIIobj.LookupAtomLabel"><tt class="xref py py-func docutils literal"><span class="pre">LookupAtomLabel()</span></tt></a></p> |
---|
1588 | </dd></dl> |
---|
1589 | |
---|
1590 | <dl class="data"> |
---|
1591 | <dt id="GSASIIobj.AtomRanIdLookup"> |
---|
1592 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">AtomRanIdLookup</tt><em class="property"> = {}</em><a class="headerlink" href="#GSASIIobj.AtomRanIdLookup" title="Permalink to this definition">¶</a></dt> |
---|
1593 | <dd><p>dict listing for each phase the atom sequential index keyed by atom random Id; |
---|
1594 | best to access this using <a class="reference internal" href="#GSASIIobj.LookupAtomId" title="GSASIIobj.LookupAtomId"><tt class="xref py py-func docutils literal"><span class="pre">LookupAtomId()</span></tt></a></p> |
---|
1595 | </dd></dl> |
---|
1596 | |
---|
1597 | <dl class="function"> |
---|
1598 | <dt id="GSASIIobj.CompileVarDesc"> |
---|
1599 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">CompileVarDesc</tt><big>(</big><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#CompileVarDesc"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.CompileVarDesc" title="Permalink to this definition">¶</a></dt> |
---|
1600 | <dd><p>Set the values in the variable description lookup table (<a class="reference internal" href="#GSASIIobj.VarDesc" title="GSASIIobj.VarDesc"><tt class="xref py py-attr docutils literal"><span class="pre">VarDesc</span></tt></a>) |
---|
1601 | into <a class="reference internal" href="#GSASIIobj.reVarDesc" title="GSASIIobj.reVarDesc"><tt class="xref py py-attr docutils literal"><span class="pre">reVarDesc</span></tt></a>. This is called in <a class="reference internal" href="#GSASIIobj.getDescr" title="GSASIIobj.getDescr"><tt class="xref py py-func docutils literal"><span class="pre">getDescr()</span></tt></a> so the initialization |
---|
1602 | is always done before use.</p> |
---|
1603 | <p>Note that keys may contain regular expressions, where ‘[xyz]’ |
---|
1604 | matches ‘x’ ‘y’ or ‘z’ (equivalently ‘[x-z]’ describes this as range of values). |
---|
1605 | ‘.*’ matches any string. For example:</p> |
---|
1606 | <div class="highlight-python"><pre>'AUiso':'Atomic isotropic displacement parameter',</pre> |
---|
1607 | </div> |
---|
1608 | <p>will match variable <tt class="docutils literal"><span class="pre">'p::AUiso:a'</span></tt>. |
---|
1609 | If parentheses are used in the key, the contents of those parentheses can be |
---|
1610 | used in the value, such as:</p> |
---|
1611 | <div class="highlight-python"><pre>'AU([123][123])':'Atomic anisotropic displacement parameter U\1',</pre> |
---|
1612 | </div> |
---|
1613 | <p>will match <tt class="docutils literal"><span class="pre">AU11</span></tt>, <tt class="docutils literal"><span class="pre">AU23</span></tt>,.. and <cite>U11</cite>, <cite>U23</cite> etc will be displayed |
---|
1614 | in the value when used.</p> |
---|
1615 | </dd></dl> |
---|
1616 | |
---|
1617 | <dl class="data"> |
---|
1618 | <dt id="GSASIIobj.DefaultControls"> |
---|
1619 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">DefaultControls</tt><em class="property"> = {'F**2': True, 'shift factor': 1.0, 'deriv type': 'analytic Hessian', 'max cyc': 3, 'minF/sig': 0, 'FreeVar1': 'Sample humidity (%)', 'Author': 'no name', 'FreeVar2': 'Sample voltage (V)', 'FreeVar3': 'Applied load (MN)', 'min dM/M': 0.0001}</em><a class="headerlink" href="#GSASIIobj.DefaultControls" title="Permalink to this definition">¶</a></dt> |
---|
1620 | <dd><p>Values to be used as defaults for the initial contents of the <tt class="docutils literal"><span class="pre">Controls</span></tt> |
---|
1621 | data tree item.</p> |
---|
1622 | </dd></dl> |
---|
1623 | |
---|
1624 | <dl class="class"> |
---|
1625 | <dt id="GSASIIobj.G2VarObj"> |
---|
1626 | <em class="property">class </em><tt class="descclassname">GSASIIobj.</tt><tt class="descname">G2VarObj</tt><big>(</big><em>*args</em><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#G2VarObj"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.G2VarObj" title="Permalink to this definition">¶</a></dt> |
---|
1627 | <dd><p>Defines a GSAS-II variable either using the phase/atom/histogram |
---|
1628 | unique Id numbers or using a character string that specifies |
---|
1629 | variables by phase/atom/histogram number (which can change). |
---|
1630 | Note that <tt class="xref py py-func docutils literal"><span class="pre">LoadID()</span></tt> should be used to (re)load the current Ids |
---|
1631 | before creating or later using the G2VarObj object.</p> |
---|
1632 | <p>A <a class="reference internal" href="#GSASIIobj.G2VarObj" title="GSASIIobj.G2VarObj"><tt class="xref py py-class docutils literal"><span class="pre">G2VarObj</span></tt></a> object can be created with a single parameter:</p> |
---|
1633 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1634 | <col class="field-name" /> |
---|
1635 | <col class="field-body" /> |
---|
1636 | <tbody valign="top"> |
---|
1637 | <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>varname</strong> (<em>str/tuple</em>) – <dl class="docutils"> |
---|
1638 | <dt>a single value can be used to create a <a class="reference internal" href="#GSASIIobj.G2VarObj" title="GSASIIobj.G2VarObj"><tt class="xref py py-class docutils literal"><span class="pre">G2VarObj</span></tt></a></dt> |
---|
1639 | <dd>object. If a string, it must be of form “p:h:var” or “p:h:var:a”, where</dd> |
---|
1640 | </dl> |
---|
1641 | <ul class="simple"> |
---|
1642 | <li>p is the phase number (which may be left blank);</li> |
---|
1643 | <li>h is the histogram number (which may be left blank);</li> |
---|
1644 | <li>a is the atom number (which may be left blank in which case the third colon is omitted).</li> |
---|
1645 | </ul> |
---|
1646 | <blockquote> |
---|
1647 | <div>Alternately a single tuple of form (Phase,Histogram,VarName,AtomID) can be used, where |
---|
1648 | Phase, Histogram, and AtomID are None or are ranId values and VarName is a string.</div></blockquote> |
---|
1649 | </td> |
---|
1650 | </tr> |
---|
1651 | </tbody> |
---|
1652 | </table> |
---|
1653 | <p>If four positional arguments are supplied, they are:</p> |
---|
1654 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1655 | <col class="field-name" /> |
---|
1656 | <col class="field-body" /> |
---|
1657 | <tbody valign="top"> |
---|
1658 | <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><ul class="first last simple"> |
---|
1659 | <li><strong>phasenum</strong> (<em>str/int</em>) – The number for the phase</li> |
---|
1660 | <li><strong>histnum</strong> (<em>str/int</em>) – The number for the histogram</li> |
---|
1661 | <li><strong>varname</strong> (<em>str</em>) – a single value can be used to create a <a class="reference internal" href="#GSASIIobj.G2VarObj" title="GSASIIobj.G2VarObj"><tt class="xref py py-class docutils literal"><span class="pre">G2VarObj</span></tt></a></li> |
---|
1662 | <li><strong>atomnum</strong> (<em>str/int</em>) – The number for the atom</li> |
---|
1663 | </ul> |
---|
1664 | </td> |
---|
1665 | </tr> |
---|
1666 | </tbody> |
---|
1667 | </table> |
---|
1668 | <dl class="method"> |
---|
1669 | <dt id="GSASIIobj.G2VarObj.varname"> |
---|
1670 | <tt class="descname">varname</tt><big>(</big><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#G2VarObj.varname"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.G2VarObj.varname" title="Permalink to this definition">¶</a></dt> |
---|
1671 | <dd><p>Formats the GSAS-II variable name as a “traditional” GSAS-II variable |
---|
1672 | string (p:h:<var>:a) or (p:h:<var>)</p> |
---|
1673 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1674 | <col class="field-name" /> |
---|
1675 | <col class="field-body" /> |
---|
1676 | <tbody valign="top"> |
---|
1677 | <tr class="field-odd field"><th class="field-name">Returns:</th><td class="field-body">the variable name as a str</td> |
---|
1678 | </tr> |
---|
1679 | </tbody> |
---|
1680 | </table> |
---|
1681 | </dd></dl> |
---|
1682 | |
---|
1683 | </dd></dl> |
---|
1684 | |
---|
1685 | <dl class="data"> |
---|
1686 | <dt id="GSASIIobj.HistIdLookup"> |
---|
1687 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">HistIdLookup</tt><em class="property"> = {}</em><a class="headerlink" href="#GSASIIobj.HistIdLookup" title="Permalink to this definition">¶</a></dt> |
---|
1688 | <dd><p>dict listing histogram name and random Id, keyed by sequential histogram index as a str; |
---|
1689 | best to access this using <a class="reference internal" href="#GSASIIobj.LookupHistName" title="GSASIIobj.LookupHistName"><tt class="xref py py-func docutils literal"><span class="pre">LookupHistName()</span></tt></a></p> |
---|
1690 | </dd></dl> |
---|
1691 | |
---|
1692 | <dl class="data"> |
---|
1693 | <dt id="GSASIIobj.HistRanIdLookup"> |
---|
1694 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">HistRanIdLookup</tt><em class="property"> = {}</em><a class="headerlink" href="#GSASIIobj.HistRanIdLookup" title="Permalink to this definition">¶</a></dt> |
---|
1695 | <dd><p>dict listing histogram sequential index keyed by histogram random Id; |
---|
1696 | best to access this using <a class="reference internal" href="#GSASIIobj.LookupHistId" title="GSASIIobj.LookupHistId"><tt class="xref py py-func docutils literal"><span class="pre">LookupHistId()</span></tt></a></p> |
---|
1697 | </dd></dl> |
---|
1698 | |
---|
1699 | <dl class="function"> |
---|
1700 | <dt id="GSASIIobj.IndexAllIds"> |
---|
1701 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">IndexAllIds</tt><big>(</big><em>Histograms</em>, <em>Phases</em><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#IndexAllIds"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.IndexAllIds" title="Permalink to this definition">¶</a></dt> |
---|
1702 | <dd><p>Scan through the used phases & histograms and create an index |
---|
1703 | to the random numbers of phases, histograms and atoms. While doing this, |
---|
1704 | confirm that assigned random numbers are unique – just in case lightning |
---|
1705 | strikes twice in the same place.</p> |
---|
1706 | <p>Note: this code assumes that the atom random Id (ranId) is the last |
---|
1707 | element each atom record.</p> |
---|
1708 | <p>This is called in two places (only) <a class="reference internal" href="GSASIIstruc.html#GSASIIstrIO.GetUsedHistogramsAndPhases" title="GSASIIstrIO.GetUsedHistogramsAndPhases"><tt class="xref py py-func docutils literal"><span class="pre">GSASIIstrIO.GetUsedHistogramsAndPhases()</span></tt></a> |
---|
1709 | (which loads the histograms and phases from a GPX file) and |
---|
1710 | <a class="reference internal" href="GSASII.html#GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree" title="GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree"><tt class="xref py py-meth docutils literal"><span class="pre">GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree()</span></tt></a> |
---|
1711 | (which loads the histograms and phases from the data tree.)</p> |
---|
1712 | <p>TODO: do we need a lookup for rigid body variables?</p> |
---|
1713 | </dd></dl> |
---|
1714 | |
---|
1715 | <dl class="function"> |
---|
1716 | <dt id="GSASIIobj.LookupAtomId"> |
---|
1717 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">LookupAtomId</tt><big>(</big><em>pId</em>, <em>ranId</em><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#LookupAtomId"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.LookupAtomId" title="Permalink to this definition">¶</a></dt> |
---|
1718 | <dd><p>Get the atom number from a phase and atom random Id</p> |
---|
1719 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1720 | <col class="field-name" /> |
---|
1721 | <col class="field-body" /> |
---|
1722 | <tbody valign="top"> |
---|
1723 | <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><ul class="first simple"> |
---|
1724 | <li><strong>pId</strong> (<em>int/str</em>) – the sequential number of the phase</li> |
---|
1725 | <li><strong>ranId</strong> (<em>int</em>) – the random Id assigned to an atom</li> |
---|
1726 | </ul> |
---|
1727 | </td> |
---|
1728 | </tr> |
---|
1729 | <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body"><p class="first last">the index number of the atom (str)</p> |
---|
1730 | </td> |
---|
1731 | </tr> |
---|
1732 | </tbody> |
---|
1733 | </table> |
---|
1734 | </dd></dl> |
---|
1735 | |
---|
1736 | <dl class="function"> |
---|
1737 | <dt id="GSASIIobj.LookupAtomLabel"> |
---|
1738 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">LookupAtomLabel</tt><big>(</big><em>pId</em>, <em>index</em><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#LookupAtomLabel"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.LookupAtomLabel" title="Permalink to this definition">¶</a></dt> |
---|
1739 | <dd><p>Get the atom label from a phase and atom index number</p> |
---|
1740 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1741 | <col class="field-name" /> |
---|
1742 | <col class="field-body" /> |
---|
1743 | <tbody valign="top"> |
---|
1744 | <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><ul class="first simple"> |
---|
1745 | <li><strong>pId</strong> (<em>int/str</em>) – the sequential number of the phase</li> |
---|
1746 | <li><strong>index</strong> (<em>int</em>) – the index of the atom in the list of atoms</li> |
---|
1747 | </ul> |
---|
1748 | </td> |
---|
1749 | </tr> |
---|
1750 | <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body"><p class="first last">the label for the atom (str) and the random Id of the atom (int)</p> |
---|
1751 | </td> |
---|
1752 | </tr> |
---|
1753 | </tbody> |
---|
1754 | </table> |
---|
1755 | </dd></dl> |
---|
1756 | |
---|
1757 | <dl class="function"> |
---|
1758 | <dt id="GSASIIobj.LookupHistId"> |
---|
1759 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">LookupHistId</tt><big>(</big><em>ranId</em><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#LookupHistId"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.LookupHistId" title="Permalink to this definition">¶</a></dt> |
---|
1760 | <dd><p>Get the histogram number and name from a histogram random Id</p> |
---|
1761 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1762 | <col class="field-name" /> |
---|
1763 | <col class="field-body" /> |
---|
1764 | <tbody valign="top"> |
---|
1765 | <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>ranId</strong> (<em>int</em>) – the random Id assigned to a histogram</td> |
---|
1766 | </tr> |
---|
1767 | <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body">the sequential Id (hId) number for the histogram (str)</td> |
---|
1768 | </tr> |
---|
1769 | </tbody> |
---|
1770 | </table> |
---|
1771 | </dd></dl> |
---|
1772 | |
---|
1773 | <dl class="function"> |
---|
1774 | <dt id="GSASIIobj.LookupHistName"> |
---|
1775 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">LookupHistName</tt><big>(</big><em>hId</em><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#LookupHistName"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.LookupHistName" title="Permalink to this definition">¶</a></dt> |
---|
1776 | <dd><p>Get the histogram number and name from a histogram Id</p> |
---|
1777 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1778 | <col class="field-name" /> |
---|
1779 | <col class="field-body" /> |
---|
1780 | <tbody valign="top"> |
---|
1781 | <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>hId</strong> (<em>int/str</em>) – the sequential assigned to a histogram</td> |
---|
1782 | </tr> |
---|
1783 | <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body">(hist,ranId) where hist is the name of the histogram (str) |
---|
1784 | and ranId is the random # id for the histogram (int)</td> |
---|
1785 | </tr> |
---|
1786 | </tbody> |
---|
1787 | </table> |
---|
1788 | </dd></dl> |
---|
1789 | |
---|
1790 | <dl class="function"> |
---|
1791 | <dt id="GSASIIobj.LookupPhaseId"> |
---|
1792 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">LookupPhaseId</tt><big>(</big><em>ranId</em><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#LookupPhaseId"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.LookupPhaseId" title="Permalink to this definition">¶</a></dt> |
---|
1793 | <dd><p>Get the phase number and name from a phase random Id</p> |
---|
1794 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1795 | <col class="field-name" /> |
---|
1796 | <col class="field-body" /> |
---|
1797 | <tbody valign="top"> |
---|
1798 | <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>ranId</strong> (<em>int</em>) – the random Id assigned to a phase</td> |
---|
1799 | </tr> |
---|
1800 | <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body">the sequential Id (pId) number for the phase (str)</td> |
---|
1801 | </tr> |
---|
1802 | </tbody> |
---|
1803 | </table> |
---|
1804 | </dd></dl> |
---|
1805 | |
---|
1806 | <dl class="function"> |
---|
1807 | <dt id="GSASIIobj.LookupPhaseName"> |
---|
1808 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">LookupPhaseName</tt><big>(</big><em>pId</em><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#LookupPhaseName"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.LookupPhaseName" title="Permalink to this definition">¶</a></dt> |
---|
1809 | <dd><p>Get the phase number and name from a phase Id</p> |
---|
1810 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1811 | <col class="field-name" /> |
---|
1812 | <col class="field-body" /> |
---|
1813 | <tbody valign="top"> |
---|
1814 | <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>pId</strong> (<em>int/str</em>) – the sequential assigned to a phase</td> |
---|
1815 | </tr> |
---|
1816 | <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body">(phase,ranId) where phase is the name of the phase (str) |
---|
1817 | and ranId is the random # id for the phase (int)</td> |
---|
1818 | </tr> |
---|
1819 | </tbody> |
---|
1820 | </table> |
---|
1821 | </dd></dl> |
---|
1822 | |
---|
1823 | <dl class="function"> |
---|
1824 | <dt id="GSASIIobj.MakeUniqueLabel"> |
---|
1825 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">MakeUniqueLabel</tt><big>(</big><em>lbl</em>, <em>labellist</em><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#MakeUniqueLabel"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.MakeUniqueLabel" title="Permalink to this definition">¶</a></dt> |
---|
1826 | <dd><p>Make sure that every a label is unique against a list by adding |
---|
1827 | digits at the end until it is not found in list.</p> |
---|
1828 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1829 | <col class="field-name" /> |
---|
1830 | <col class="field-body" /> |
---|
1831 | <tbody valign="top"> |
---|
1832 | <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><ul class="first simple"> |
---|
1833 | <li><strong>lbl</strong> (<em>str</em>) – the input label</li> |
---|
1834 | <li><strong>labellist</strong> (<em>list</em>) – the labels that have already been encountered</li> |
---|
1835 | </ul> |
---|
1836 | </td> |
---|
1837 | </tr> |
---|
1838 | <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body"><p class="first last">lbl if not found in labellist or lbl with <tt class="docutils literal"><span class="pre">_1-9</span></tt> (or |
---|
1839 | <tt class="docutils literal"><span class="pre">_10-99</span></tt>, etc.) appended at the end</p> |
---|
1840 | </td> |
---|
1841 | </tr> |
---|
1842 | </tbody> |
---|
1843 | </table> |
---|
1844 | </dd></dl> |
---|
1845 | |
---|
1846 | <dl class="data"> |
---|
1847 | <dt id="GSASIIobj.PhaseIdLookup"> |
---|
1848 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">PhaseIdLookup</tt><em class="property"> = {}</em><a class="headerlink" href="#GSASIIobj.PhaseIdLookup" title="Permalink to this definition">¶</a></dt> |
---|
1849 | <dd><p>dict listing phase name and random Id keyed by sequential phase index as a str; |
---|
1850 | best to access this using <a class="reference internal" href="#GSASIIobj.LookupPhaseName" title="GSASIIobj.LookupPhaseName"><tt class="xref py py-func docutils literal"><span class="pre">LookupPhaseName()</span></tt></a></p> |
---|
1851 | </dd></dl> |
---|
1852 | |
---|
1853 | <dl class="data"> |
---|
1854 | <dt id="GSASIIobj.PhaseRanIdLookup"> |
---|
1855 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">PhaseRanIdLookup</tt><em class="property"> = {}</em><a class="headerlink" href="#GSASIIobj.PhaseRanIdLookup" title="Permalink to this definition">¶</a></dt> |
---|
1856 | <dd><p>dict listing phase sequential index keyed by phase random Id; |
---|
1857 | best to access this using <a class="reference internal" href="#GSASIIobj.LookupPhaseId" title="GSASIIobj.LookupPhaseId"><tt class="xref py py-func docutils literal"><span class="pre">LookupPhaseId()</span></tt></a></p> |
---|
1858 | </dd></dl> |
---|
1859 | |
---|
1860 | <dl class="data"> |
---|
1861 | <dt id="GSASIIobj.ShortHistNames"> |
---|
1862 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">ShortHistNames</tt><em class="property"> = {}</em><a class="headerlink" href="#GSASIIobj.ShortHistNames" title="Permalink to this definition">¶</a></dt> |
---|
1863 | <dd><p>a dict containing a possibly shortened and when non-unique numbered |
---|
1864 | version of the histogram name. Keyed by the histogram sequential index.</p> |
---|
1865 | </dd></dl> |
---|
1866 | |
---|
1867 | <dl class="data"> |
---|
1868 | <dt id="GSASIIobj.ShortPhaseNames"> |
---|
1869 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">ShortPhaseNames</tt><em class="property"> = {}</em><a class="headerlink" href="#GSASIIobj.ShortPhaseNames" title="Permalink to this definition">¶</a></dt> |
---|
1870 | <dd><p>a dict containing a possibly shortened and when non-unique numbered |
---|
1871 | version of the phase name. Keyed by the phase sequential index.</p> |
---|
1872 | </dd></dl> |
---|
1873 | |
---|
1874 | <dl class="data"> |
---|
1875 | <dt id="GSASIIobj.VarDesc"> |
---|
1876 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">VarDesc</tt><em class="property"> = {}</em><a class="headerlink" href="#GSASIIobj.VarDesc" title="Permalink to this definition">¶</a></dt> |
---|
1877 | <dd><p>This dictionary lists descriptions for GSAS-II variables, |
---|
1878 | as set in <a class="reference internal" href="#GSASIIobj.CompileVarDesc" title="GSASIIobj.CompileVarDesc"><tt class="xref py py-func docutils literal"><span class="pre">CompileVarDesc()</span></tt></a>. See that function for a description |
---|
1879 | for how keys and values are written.</p> |
---|
1880 | </dd></dl> |
---|
1881 | |
---|
1882 | <dl class="function"> |
---|
1883 | <dt id="GSASIIobj.fmtVarDescr"> |
---|
1884 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">fmtVarDescr</tt><big>(</big><em>varname</em><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#fmtVarDescr"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.fmtVarDescr" title="Permalink to this definition">¶</a></dt> |
---|
1885 | <dd><p>Return a string with a more complete description for a GSAS-II variable</p> |
---|
1886 | <p>TODO: This will not handle rigid body parameters yet</p> |
---|
1887 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1888 | <col class="field-name" /> |
---|
1889 | <col class="field-body" /> |
---|
1890 | <tbody valign="top"> |
---|
1891 | <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>name</strong> (<em>str</em>) – A full G2 variable name with 2 or 3 |
---|
1892 | colons (<p>:<h>:name[:<a>])</td> |
---|
1893 | </tr> |
---|
1894 | <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body">a string with the description</td> |
---|
1895 | </tr> |
---|
1896 | </tbody> |
---|
1897 | </table> |
---|
1898 | </dd></dl> |
---|
1899 | |
---|
1900 | <dl class="function"> |
---|
1901 | <dt id="GSASIIobj.getDescr"> |
---|
1902 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">getDescr</tt><big>(</big><em>name</em><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#getDescr"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.getDescr" title="Permalink to this definition">¶</a></dt> |
---|
1903 | <dd><p>Return a short description for a GSAS-II variable</p> |
---|
1904 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1905 | <col class="field-name" /> |
---|
1906 | <col class="field-body" /> |
---|
1907 | <tbody valign="top"> |
---|
1908 | <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>name</strong> (<em>str</em>) – The descriptive part of the variable name without colons (:)</td> |
---|
1909 | </tr> |
---|
1910 | <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body">a short description or None if not found</td> |
---|
1911 | </tr> |
---|
1912 | </tbody> |
---|
1913 | </table> |
---|
1914 | </dd></dl> |
---|
1915 | |
---|
1916 | <dl class="function"> |
---|
1917 | <dt id="GSASIIobj.getVarDescr"> |
---|
1918 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">getVarDescr</tt><big>(</big><em>varname</em><big>)</big><a class="reference internal" href="_modules/GSASIIobj.html#getVarDescr"><span class="viewcode-link">[source]</span></a><a class="headerlink" href="#GSASIIobj.getVarDescr" title="Permalink to this definition">¶</a></dt> |
---|
1919 | <dd><p>Return a short description for a GSAS-II variable</p> |
---|
1920 | <table class="docutils field-list" frame="void" rules="none"> |
---|
1921 | <col class="field-name" /> |
---|
1922 | <col class="field-body" /> |
---|
1923 | <tbody valign="top"> |
---|
1924 | <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>name</strong> (<em>str</em>) – A full G2 variable name with 2 or 3 |
---|
1925 | colons (<p>:<h>:name[:<a>])</td> |
---|
1926 | </tr> |
---|
1927 | <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body">a five element list as [<cite>p</cite>,`h`,`name`,`a`,`description`], |
---|
1928 | where <cite>p</cite>, <cite>h</cite>, <cite>a</cite> are str values or <cite>None</cite>, for the phase number, |
---|
1929 | the histogram number and the atom number; <cite>name</cite> will always be |
---|
1930 | an str; and <cite>description</cite> is str or <cite>None</cite>. |
---|
1931 | If the variable name is incorrectly formed (for example, wrong |
---|
1932 | number of colons), <cite>None</cite> is returned instead of a list.</td> |
---|
1933 | </tr> |
---|
1934 | </tbody> |
---|
1935 | </table> |
---|
1936 | </dd></dl> |
---|
1937 | |
---|
1938 | <dl class="data"> |
---|
1939 | <dt id="GSASIIobj.reVarDesc"> |
---|
1940 | <tt class="descclassname">GSASIIobj.</tt><tt class="descname">reVarDesc</tt><em class="property"> = {}</em><a class="headerlink" href="#GSASIIobj.reVarDesc" title="Permalink to this definition">¶</a></dt> |
---|
1941 | <dd><p>This dictionary lists descriptions for GSAS-II variables with |
---|
1942 | the same values as <a class="reference internal" href="#GSASIIobj.VarDesc" title="GSASIIobj.VarDesc"><tt class="xref py py-attr docutils literal"><span class="pre">VarDesc</span></tt></a> except that keys have been compiled as |
---|
1943 | regular expressions. Initialized in <a class="reference internal" href="#GSASIIobj.CompileVarDesc" title="GSASIIobj.CompileVarDesc"><tt class="xref py py-func docutils literal"><span class="pre">CompileVarDesc()</span></tt></a>.</p> |
---|
1944 | </dd></dl> |
---|
1945 | |
---|
1946 | </div> |
---|
1947 | </div> |
---|
1948 | |
---|
1949 | |
---|
1950 | </div> |
---|
1951 | </div> |
---|
1952 | </div> |
---|
1953 | <div class="sphinxsidebar"> |
---|
1954 | <div class="sphinxsidebarwrapper"> |
---|
1955 | <p class="logo"><a href="index.html"> |
---|
1956 | <img class="logo" src="_static/G2_html_logo.png" alt="Logo"/> |
---|
1957 | </a></p> |
---|
1958 | <h3><a href="index.html">Table Of Contents</a></h3> |
---|
1959 | <ul> |
---|
1960 | <li><a class="reference internal" href="#"><em>GSASIIobj: Data objects</em></a><ul> |
---|
1961 | <li><a class="reference internal" href="#constraints-tree-item">Constraints Tree Item</a></li> |
---|
1962 | <li><a class="reference internal" href="#covariance-tree-item">Covariance Tree Item</a></li> |
---|
1963 | <li><a class="reference internal" href="#phase-tree-items">Phase Tree Items</a></li> |
---|
1964 | <li><a class="reference internal" href="#rigid-body-objects">Rigid Body Objects</a></li> |
---|
1965 | <li><a class="reference internal" href="#space-group-objects">Space Group Objects</a></li> |
---|
1966 | <li><a class="reference internal" href="#atom-records">Atom Records</a></li> |
---|
1967 | <li><a class="reference internal" href="#drawing-atom-records">Drawing Atom Records</a></li> |
---|
1968 | <li><a class="reference internal" href="#powder-diffraction-tree-items">Powder Diffraction Tree Items</a></li> |
---|
1969 | <li><a class="reference internal" href="#powder-reflection-data-structure">Powder Reflection Data Structure</a></li> |
---|
1970 | <li><a class="reference internal" href="#single-crystal-tree-items">Single Crystal Tree Items</a></li> |
---|
1971 | <li><a class="reference internal" href="#single-crystal-reflection-data-structure">Single Crystal Reflection Data Structure</a></li> |
---|
1972 | <li><a class="reference internal" href="#image-data-structure">Image Data Structure</a></li> |
---|
1973 | <li><a class="reference internal" href="#parameter-dictionary">Parameter Dictionary</a></li> |
---|
1974 | <li><a class="reference internal" href="#classes-and-routines"><em>Classes and routines</em></a></li> |
---|
1975 | </ul> |
---|
1976 | </li> |
---|
1977 | </ul> |
---|
1978 | |
---|
1979 | <h4>Previous topic</h4> |
---|
1980 | <p class="topless"><a href="GSASII.html" |
---|
1981 | title="previous chapter"><em>GSAS-II Main Module</em></a></p> |
---|
1982 | <h4>Next topic</h4> |
---|
1983 | <p class="topless"><a href="GSASIIutil.html" |
---|
1984 | title="next chapter"><em>GSAS-II Utility Modules</em></a></p> |
---|
1985 | <h3>This Page</h3> |
---|
1986 | <ul class="this-page-menu"> |
---|
1987 | <li><a href="_sources/GSASIIobj.txt" |
---|
1988 | rel="nofollow">Show Source</a></li> |
---|
1989 | </ul> |
---|
1990 | <div id="searchbox" style="display: none"> |
---|
1991 | <h3>Quick search</h3> |
---|
1992 | <form class="search" action="search.html" method="get"> |
---|
1993 | <input type="text" name="q" /> |
---|
1994 | <input type="submit" value="Go" /> |
---|
1995 | <input type="hidden" name="check_keywords" value="yes" /> |
---|
1996 | <input type="hidden" name="area" value="default" /> |
---|
1997 | </form> |
---|
1998 | <p class="searchtip" style="font-size: 90%"> |
---|
1999 | Enter search terms or a module, class or function name. |
---|
2000 | </p> |
---|
2001 | </div> |
---|
2002 | <script type="text/javascript">$('#searchbox').show(0);</script> |
---|
2003 | </div> |
---|
2004 | </div> |
---|
2005 | <div class="clearer"></div> |
---|
2006 | </div> |
---|
2007 | <div class="related"> |
---|
2008 | <h3>Navigation</h3> |
---|
2009 | <ul> |
---|
2010 | <li class="right" style="margin-right: 10px"> |
---|
2011 | <a href="genindex.html" title="General Index" |
---|
2012 | >index</a></li> |
---|
2013 | <li class="right" > |
---|
2014 | <a href="py-modindex.html" title="Python Module Index" |
---|
2015 | >modules</a> |</li> |
---|
2016 | <li class="right" > |
---|
2017 | <a href="GSASIIutil.html" title="GSAS-II Utility Modules" |
---|
2018 | >next</a> |</li> |
---|
2019 | <li class="right" > |
---|
2020 | <a href="GSASII.html" title="GSAS-II Main Module" |
---|
2021 | >previous</a> |</li> |
---|
2022 | <li><a href="index.html">GSAS-II 0.2.0 documentation</a> »</li> |
---|
2023 | </ul> |
---|
2024 | </div> |
---|
2025 | <div class="footer"> |
---|
2026 | © Copyright 2013, Von Dreele and Toby for Argonne National Laboratory. |
---|
2027 | Created using <a href="http://sphinx.pocoo.org/">Sphinx</a> 1.1.2. |
---|
2028 | </div> |
---|
2029 | </body> |
---|
2030 | </html> |
---|