1 | # -*- coding: utf-8 -*- |
---|
2 | #GSASIIobj - data objects for GSAS-II |
---|
3 | ########### SVN repository information ################### |
---|
4 | # $Date: 2013-12-16 16:43:01 +0000 (Mon, 16 Dec 2013) $ |
---|
5 | # $Author: toby $ |
---|
6 | # $Revision: 1168 $ |
---|
7 | # $URL: trunk/GSASIIobj.py $ |
---|
8 | # $Id: GSASIIobj.py 1168 2013-12-16 16:43:01Z toby $ |
---|
9 | ########### SVN repository information ################### |
---|
10 | |
---|
11 | ''' |
---|
12 | *GSASIIobj: Data objects* |
---|
13 | ========================= |
---|
14 | |
---|
15 | This module defines and/or documents the data structures used in GSAS-II, as well |
---|
16 | as provides misc. support routines. |
---|
17 | |
---|
18 | Constraints Tree Item |
---|
19 | ---------------------- |
---|
20 | |
---|
21 | .. _Constraints_table: |
---|
22 | |
---|
23 | .. index:: |
---|
24 | single: Constraints object description |
---|
25 | single: Data object descriptions; Constraints |
---|
26 | |
---|
27 | Constraints are stored in a dict, separated into groups. |
---|
28 | Note that parameter are named in the following pattern, |
---|
29 | p:h:<var>:n, where p is the phase number, h is the histogram number |
---|
30 | <var> is a variable name and n is the parameter number. |
---|
31 | If a parameter does not depend on a histogram or phase or is unnumbered, that |
---|
32 | number is omitted. |
---|
33 | Note that the contents of each dict item is a List where each element in the |
---|
34 | list is a :ref:`constraint definition objects <Constraint_definitions_table>`. |
---|
35 | The constraints in this form are converted in |
---|
36 | :func:`GSASIIstrIO.ProcessConstraints` to the form used in :mod:`GSASIImapvars` |
---|
37 | |
---|
38 | The keys in the Constraints dict are: |
---|
39 | |
---|
40 | .. tabularcolumns:: |l|p{4.5in}| |
---|
41 | |
---|
42 | ========== ==================================================== |
---|
43 | key explanation |
---|
44 | ========== ==================================================== |
---|
45 | Hist This specifies a list of constraints on |
---|
46 | histogram-related parameters, |
---|
47 | which will be of form :h:<var>:n. |
---|
48 | HAP This specifies a list of constraints on parameters |
---|
49 | that are defined for every histogram in each phase |
---|
50 | and are of form p:h:<var>:n. |
---|
51 | Phase This specifies a list of constraints on phase |
---|
52 | parameters, |
---|
53 | which will be of form p::<var>:n. |
---|
54 | Global This specifies a list of constraints on parameters |
---|
55 | that are not tied to a histogram or phase and |
---|
56 | are of form ::<var>:n |
---|
57 | ========== ==================================================== |
---|
58 | |
---|
59 | .. _Constraint_definitions_table: |
---|
60 | |
---|
61 | .. index:: |
---|
62 | single: Constraint definition object description |
---|
63 | single: Data object descriptions; Constraint Definition |
---|
64 | |
---|
65 | Each constraint is defined as an item in a list. Each constraint is of form:: |
---|
66 | |
---|
67 | [[<mult1>, <var1>], [<mult2>, <var2>],..., <fixedval>, <varyflag>, <constype>] |
---|
68 | |
---|
69 | Where the variable pair list item containing two values [<mult>, <var>], where: |
---|
70 | |
---|
71 | * <mult> is a multiplier for the constraint (float) |
---|
72 | * <var> a :class:`G2VarObj` object (previously a str variable name of form |
---|
73 | 'p:h:name[:at]') |
---|
74 | |
---|
75 | Note that the last three items in the list play a special role: |
---|
76 | |
---|
77 | * <fixedval> is the fixed value for a `constant equation` (``constype=c``) |
---|
78 | constraint or is None. For a `New variable` (``constype=f``) constraint, |
---|
79 | a variable name can be specified as a str (used for externally |
---|
80 | generated constraints) |
---|
81 | * <varyflag> is True or False for `New variable` (``constype=f``) constraints |
---|
82 | or is None. This will be implemented in the future to indicate if these variables |
---|
83 | should be refined. |
---|
84 | * <constype> is one of four letters, 'e', 'c', 'h', 'f' that determines the type of constraint: |
---|
85 | |
---|
86 | * 'e' defines a set of equivalent variables. Only the first variable is refined (if the |
---|
87 | appropriate refine flag is set) and and all other equivalent variables in the list |
---|
88 | are generated from that variable, using the appropriate multipliers. |
---|
89 | * 'c' defines a constraint equation of form, |
---|
90 | :math:`m_1 \\times var_1 + m_2 \\times var_2 + ... = c` |
---|
91 | * 'h' defines a variable to hold (not vary). Any variable on this list is not varied, |
---|
92 | even if its refinement flag is set. Only one [mult,var] pair is allowed in a hold |
---|
93 | constraint and the mult value is ignored. |
---|
94 | This is of particular value when needing to hold one or more variables where a |
---|
95 | single flag controls a set of variables such as, coordinates, |
---|
96 | the reciprocal metric tensor or anisotropic displacement parameter. |
---|
97 | * 'f' defines a new variable (function) according to relationship |
---|
98 | :math:`newvar = m_1 \\times var_1 + m_2 \\times var_2 + ...` |
---|
99 | |
---|
100 | Covariance Tree Item |
---|
101 | -------------------- |
---|
102 | |
---|
103 | .. _Covariance_table: |
---|
104 | |
---|
105 | .. index:: |
---|
106 | single: Covariance description |
---|
107 | single: Data object descriptions; Covariance |
---|
108 | |
---|
109 | The Covariance tree item has results from the last least-squares run. They |
---|
110 | are stored in a dict with these keys: |
---|
111 | |
---|
112 | .. tabularcolumns:: |l|l|p{4in}| |
---|
113 | |
---|
114 | ============= =============== ==================================================== |
---|
115 | key sub-key explanation |
---|
116 | ============= =============== ==================================================== |
---|
117 | newCellDict \ dict with lattice parameters computed by |
---|
118 | :func:`GSASIIstrMath.GetNewCellParms` (dict) |
---|
119 | title \ Name of gpx file(?) (str) |
---|
120 | variables \ Values for all N refined variables |
---|
121 | (list of float values, length N, |
---|
122 | ordered to match varyList) |
---|
123 | sig \ Uncertainty values for all N refined variables |
---|
124 | (list of float values, length N, |
---|
125 | ordered to match varyList) |
---|
126 | varyList \ List of directly refined variables |
---|
127 | (list of str values, length N) |
---|
128 | newAtomDict \ dict with atom position values computed in |
---|
129 | :func:`GSASIIstrMath.ApplyXYZshifts` (dict) |
---|
130 | Rvals \ R-factors, GOF, Marquardt value for last |
---|
131 | refinement cycle (dict) |
---|
132 | \ Nobs Number of observed data points (int) |
---|
133 | \ Rwp overall weighted profile R-factor (%, float) |
---|
134 | \ chisq sum[w*(Iobs-Icalc)**2] for all data |
---|
135 | note this is not the reduced chi squared (float) |
---|
136 | \ lamMax Marquardt value applied to Hessian diagonal |
---|
137 | (float) |
---|
138 | \ GOF The goodness-of-fit, aka square root of |
---|
139 | the reduced chi squared. (float) |
---|
140 | covMatrix \ The (NxN) covVariance matrix (np.array) |
---|
141 | ============= =============== ==================================================== |
---|
142 | |
---|
143 | Phase Tree Items |
---|
144 | ---------------- |
---|
145 | |
---|
146 | .. _Phase_table: |
---|
147 | |
---|
148 | .. index:: |
---|
149 | single: Phase object description |
---|
150 | single: Data object descriptions; Phase |
---|
151 | |
---|
152 | Phase information is stored in the GSAS-II data tree as children of the |
---|
153 | Phases item in a dict with keys: |
---|
154 | |
---|
155 | .. tabularcolumns:: |l|l|p{4in}| |
---|
156 | |
---|
157 | ========== =============== ==================================================== |
---|
158 | key sub-key explanation |
---|
159 | ========== =============== ==================================================== |
---|
160 | General \ Overall information for the phase (dict) |
---|
161 | \ AtomPtrs list of four locations to use to pull info |
---|
162 | from the atom records (list) |
---|
163 | \ F000X x-ray F(000) intensity (float) |
---|
164 | \ F000N neutron F(000) intensity (float) |
---|
165 | \ Mydir directory of current .gpx file (str) |
---|
166 | \ MCSA controls Monte Carlo-Simulated Annealing controls (dict) |
---|
167 | \ Cell List with 8 items: cell refinement flag (bool) |
---|
168 | a, b, c, (Angstrom, float) |
---|
169 | alpha, beta & gamma (degrees, float) |
---|
170 | volume (A^3, float) |
---|
171 | \ Type 'nuclear' or 'macromolecular' for now (str) |
---|
172 | \ Map dict of map parameters |
---|
173 | \ SH Texture dict of spherical harmonic preferred orientation |
---|
174 | parameters |
---|
175 | \ Isotope dict of isotopes for each atom type |
---|
176 | \ Isotopes dict of scattering lengths for each isotope |
---|
177 | combination for each element in phase |
---|
178 | \ Name phase name (str) |
---|
179 | \ SGData Space group details as a :ref:`space group (SGData) object <SGData_table>` |
---|
180 | as defined in :func:`GSASIIspc.SpcGroup`. |
---|
181 | \ Pawley neg wt Restraint value for negative Pawley intensities |
---|
182 | (float) |
---|
183 | \ Flip dict of Charge flip controls |
---|
184 | \ Data plot type data plot type ('Mustrain', 'Size' or |
---|
185 | 'Preferred orientation') for powder data (str) |
---|
186 | \ Mass Mass of unit cell contents in g/mol |
---|
187 | \ POhkl March-Dollase preferred orientation direction |
---|
188 | \ Z dict of atomic numbers for each atom type |
---|
189 | \ vdWRadii dict of van der Waals radii for each atom type |
---|
190 | \ Color Colors for atoms (list of (r,b,g) triplets) |
---|
191 | \ AtomTypes List of atom types |
---|
192 | \ AtomMass List of masses for atoms |
---|
193 | \ doPawley Flag for Pawley intensity extraction (bool) |
---|
194 | \ NoAtoms Number of atoms per unit cell of each type (dict) |
---|
195 | \ Pawley dmin maximum Q (as d-space) to use for Pawley |
---|
196 | extraction (float) |
---|
197 | \ BondRadii Default radius for each atom used to compute |
---|
198 | interatomic distances (list of floats) |
---|
199 | \ AngleRadii Default radius for each atom used to compute |
---|
200 | interatomic angles (list of floats) |
---|
201 | \ DisAglCtls Dict with distance/angle search controls, |
---|
202 | which has keys 'Name', 'AtomTypes', |
---|
203 | 'BondRadii', 'AngleRadii' which are as above |
---|
204 | except are possibly edited. Also contains |
---|
205 | 'Factors', which is a 2 element list with |
---|
206 | a multiplier for bond and angle search range |
---|
207 | [typically (0.85,0.85)]. |
---|
208 | ranId \ unique random number Id for phase (int) |
---|
209 | pId \ Phase Id number for current project (int). |
---|
210 | Atoms \ Atoms in phase as a list of lists. The outer list |
---|
211 | is for each atom, the inner list contains varying |
---|
212 | items depending on the type of phase, see |
---|
213 | the :ref:`Atom Records <Atoms_table>` description. |
---|
214 | (list of lists) |
---|
215 | Drawing \ Display parameters (dict) |
---|
216 | \ ballScale Size of spheres in ball-and-stick display (float) |
---|
217 | \ bondList dict with bonds |
---|
218 | \ contourLevel map contour level in e/A^3 (float) |
---|
219 | \ showABC Flag to show view point triplet (bool). True=show. |
---|
220 | \ viewDir cartesian viewing direction (np.array with three |
---|
221 | elements) |
---|
222 | \ Zclip clipping distance in A (float) |
---|
223 | \ backColor background for plot as and R,G,B triplet |
---|
224 | (default = [0, 0, 0], black). |
---|
225 | (list with three atoms) |
---|
226 | \ selectedAtoms List of selected atoms (list of int values) |
---|
227 | \ showRigidBodies Flag to highlight rigid body placement |
---|
228 | \ sizeH Size ratio for H atoms (float) |
---|
229 | \ bondRadius Size of binds in A (float) |
---|
230 | \ atomPtrs positions of x, type, site sym, ADP flag in Draw Atoms (list) |
---|
231 | \ viewPoint list of lists. First item in list is [x,y,z] |
---|
232 | in fractional coordinates for the center of |
---|
233 | the plot. Second item list of previous & current |
---|
234 | atom number viewed (may be [0,0]) |
---|
235 | \ showHydrogen Flag to control plotting of H atoms. |
---|
236 | \ unitCellBox Flag to control display of the unit cell. |
---|
237 | \ ellipseProb Probability limit for display of thermal |
---|
238 | ellipsoids in % (float). |
---|
239 | \ vdwScale Multiplier of van der Waals radius for |
---|
240 | display of vdW spheres. |
---|
241 | \ Atoms A list of lists with an entry for each atom |
---|
242 | that is plotted. |
---|
243 | \ Zstep Step to de/increase Z-clip (float) |
---|
244 | \ Quaternion Viewing quaternion (4 element np.array) |
---|
245 | \ radiusFactor Distance ratio for searching for bonds. ? Bonds |
---|
246 | are located that are within r(Ra+Rb) and (Ra+Rb)/r |
---|
247 | where Ra and Rb are the atomic radii. |
---|
248 | \ oldxy previous view point (list with two floats) |
---|
249 | \ cameraPos Viewing position in A for plot (float) |
---|
250 | \ depthFog True if use depthFog on plot - set currently as False (bool) |
---|
251 | RBModels \ Rigid body assignments (note Rigid body definitions |
---|
252 | are stored in their own main top-level tree entry.) |
---|
253 | Pawley ref \ Pawley reflections |
---|
254 | Histograms \ A dict of dicts. The key for the outer dict is |
---|
255 | the histograms tied to this phase. The inner |
---|
256 | dict contains the combined phase/histogram |
---|
257 | parameters for items such as scale factors, |
---|
258 | size and strain parameters. (dict) |
---|
259 | MCSA \ Monte-Carlo simulated annealing parameters (dict) |
---|
260 | \ |
---|
261 | ========== =============== ==================================================== |
---|
262 | |
---|
263 | Rigid Body Objects |
---|
264 | ------------------ |
---|
265 | |
---|
266 | .. _RBData_table: |
---|
267 | |
---|
268 | .. index:: |
---|
269 | single: Rigid Body Data description |
---|
270 | single: Data object descriptions; Rigid Body Data |
---|
271 | |
---|
272 | Rigid body descriptions are available for two types of rigid bodies: 'Vector' |
---|
273 | and 'Residue'. Vector rigid bodies are developed by a sequence of translations each |
---|
274 | with a refinable magnitude and Residue rigid bodies are described as Cartesian coordinates |
---|
275 | with defined refinable torsion angles. |
---|
276 | |
---|
277 | .. tabularcolumns:: |l|l|p{4in}| |
---|
278 | |
---|
279 | ========== =============== ==================================================== |
---|
280 | key sub-key explanation |
---|
281 | ========== =============== ==================================================== |
---|
282 | Vector RBId vector rigid bodies (dict of dict) |
---|
283 | \ AtInfo Drad, Color: atom drawing radius & color for each atom type (dict) |
---|
284 | \ RBname Name assigned by user to rigid body (str) |
---|
285 | \ VectMag vector magnitudes in A (list) |
---|
286 | \ rbXYZ Cartesian coordinates for Vector rigid body (list of 3 float) |
---|
287 | \ rbRef 3 assigned reference atom nos. in rigid body for origin |
---|
288 | definition, use center of atoms flag (list of 3 int & 1 bool) |
---|
289 | \ VectRef refinement flags for VectMag values (list of bool) |
---|
290 | \ rbTypes Atom types for each atom in rigid body (list of str) |
---|
291 | \ rbVect Cartesian vectors for each translation used to build rigid body (list of lists) |
---|
292 | \ useCount Number of times rigid body is used in any structure (int) |
---|
293 | Residue RBId residue rigid bodies (dict of dict) |
---|
294 | \ AtInfo Drad, Color: atom drawing radius & color for each atom type(dict) |
---|
295 | \ RBname Name assigned by user to rigid body (str) |
---|
296 | \ rbXYZ Cartesian coordinates for Residue rigid body (list of 3 float) |
---|
297 | \ rbTypes Atom types for each atom in rigid body (list of str) |
---|
298 | \ atNames Names of each atom in rigid body (e.g. C1,N2...) (list of str) |
---|
299 | \ rbRef 3 assigned reference atom nos. in rigid body for origin |
---|
300 | definition, use center of atoms flag (list of 3 int & 1 bool) |
---|
301 | \ rbSeq Orig,Piv,angle,Riding (list): definition of internal rigid body |
---|
302 | torsion; origin atom (int), pivot atom (int), torsion angle (float), |
---|
303 | riding atoms (list of int) |
---|
304 | \ SelSeq [int,int] used by SeqSizer to identify objects |
---|
305 | \ useCount Number of times rigid body is used in any structure (int) |
---|
306 | RBIds \ unique Ids generated upon creation of each rigid body (dict) |
---|
307 | \ Vector Ids for each Vector rigid body (list) |
---|
308 | \ Residue Ids for each Residue rigid body (list) |
---|
309 | ========== =============== ==================================================== |
---|
310 | |
---|
311 | Space Group Objects |
---|
312 | ------------------- |
---|
313 | |
---|
314 | .. _SGData_table: |
---|
315 | |
---|
316 | .. index:: |
---|
317 | single: Space Group Data description |
---|
318 | single: Data object descriptions; Space Group Data |
---|
319 | |
---|
320 | Space groups are interpreted by :func:`GSASIIspc.SpcGroup` |
---|
321 | and the information is placed in a SGdata object |
---|
322 | which is a dict with these keys: |
---|
323 | |
---|
324 | .. tabularcolumns:: |l|p{4.5in}| |
---|
325 | |
---|
326 | ========== ==================================================== |
---|
327 | key explanation |
---|
328 | ========== ==================================================== |
---|
329 | SpGrp space group symbol (str) |
---|
330 | Laue one of the following 14 Laue classes: |
---|
331 | -1, 2/m, mmm, 4/m, 4/mmm, 3R, |
---|
332 | 3mR, 3, 3m1, 31m, 6/m, 6/mmm, m3, m3m (str) |
---|
333 | SGInv True if centrosymmetric, False if not (bool) |
---|
334 | SGLatt Lattice centering type. Will be one of |
---|
335 | P, A, B, C, I, F, R (str) |
---|
336 | SGUniq unique axis if monoclinic. Will be |
---|
337 | a, b, or c for monoclinic space groups. |
---|
338 | Will be blank for non-monoclinic. (str) |
---|
339 | SGCen Symmetry cell centering vectors. A (n,3) np.array |
---|
340 | of centers. Will always have at least one row: |
---|
341 | ``np.array([[0, 0, 0]])`` |
---|
342 | SGOps symmetry operations as a list of form |
---|
343 | ``[[M1,T1], [M2,T2],...]`` |
---|
344 | where :math:`M_n` is a 3x3 np.array |
---|
345 | and :math:`T_n` is a length 3 np.array. |
---|
346 | Atom coordinates are transformed where the |
---|
347 | Asymmetric unit coordinates [X is (x,y,z)] |
---|
348 | are transformed using |
---|
349 | :math:`X^\prime = M_n*X+T_n` |
---|
350 | SGSys symmetry unit cell: type one of |
---|
351 | 'triclinic', 'monoclinic', 'orthorhombic', |
---|
352 | 'tetragonal', 'rhombohedral', 'trigonal', |
---|
353 | 'hexagonal', 'cubic' (str) |
---|
354 | SGPolax Axes for space group polarity. Will be one of |
---|
355 | '', 'x', 'y', 'x y', 'z', 'x z', 'y z', |
---|
356 | 'xyz'. In the case where axes are arbitrary |
---|
357 | '111' is used (P 1, and ?). |
---|
358 | ========== ==================================================== |
---|
359 | |
---|
360 | Atom Records |
---|
361 | ------------ |
---|
362 | |
---|
363 | .. _Atoms_table: |
---|
364 | |
---|
365 | .. index:: |
---|
366 | single: Atoms record description |
---|
367 | single: Data object descriptions; Atoms record |
---|
368 | |
---|
369 | |
---|
370 | If ``phasedict`` points to the phase information in the data tree, then |
---|
371 | atoms are contained in a list of atom records (list) in |
---|
372 | ``phasedict['Atoms']``. Also needed to read atom information |
---|
373 | are four pointers, ``cx,ct,cs,cia = phasedict['General']['atomPtrs']``, |
---|
374 | which define locations in the atom record, as shown below. Items shown are |
---|
375 | always present; additional ones for macromolecular phases are marked 'mm' |
---|
376 | |
---|
377 | .. tabularcolumns:: |l|p{4.5in}| |
---|
378 | |
---|
379 | ============== ==================================================== |
---|
380 | location explanation |
---|
381 | ============== ==================================================== |
---|
382 | ct-4 mm - residue number (str) |
---|
383 | ct-3 mm - residue name (e.g. ALA) (str) |
---|
384 | ct-2 mm - chain label (str) |
---|
385 | ct-1 atom label (str) |
---|
386 | ct atom type (str) |
---|
387 | ct+1 refinement flags; combination of 'F', 'X', 'U' (str) |
---|
388 | cx,cx+1,cx+2 the x,y and z coordinates (3 floats) |
---|
389 | cs site symmetry (str) |
---|
390 | cs+1 site multiplicity (int) |
---|
391 | cia ADP flag: Isotropic ('I') or Anisotropic ('A') |
---|
392 | cia+1 Uiso (float) |
---|
393 | cia+2...cia+6 U11, U22, U33, U12, U13, U23 (6 floats) |
---|
394 | atom[-1] unique atom identifier (int) |
---|
395 | ============== ==================================================== |
---|
396 | |
---|
397 | Drawing Atom Records |
---|
398 | -------------------- |
---|
399 | |
---|
400 | .. _Drawing_atoms_table: |
---|
401 | |
---|
402 | .. index:: |
---|
403 | single: Drawing atoms record description |
---|
404 | single: Data object descriptions; Drawing atoms record |
---|
405 | |
---|
406 | |
---|
407 | If ``phasedict`` points to the phase information in the data tree, then |
---|
408 | drawing atoms are contained in a list of drawing atom records (list) in |
---|
409 | ``phasedict['Drawing']['Atoms']``. Also needed to read atom information |
---|
410 | are four pointers, ``cx,ct,cs,ci = phasedict['Drawing']['AtomPtrs']``, |
---|
411 | which define locations in the atom record, as shown below. Items shown are |
---|
412 | always present; additional ones for macromolecular phases are marked 'mm' |
---|
413 | |
---|
414 | .. tabularcolumns:: |l|p{4.5in}| |
---|
415 | |
---|
416 | ============== ==================================================== |
---|
417 | location explanation |
---|
418 | ============== ==================================================== |
---|
419 | ct-4 mm - residue number (str) |
---|
420 | ct-3 mm - residue name (e.g. ALA) (str) |
---|
421 | ct-2 mm - chain label (str) |
---|
422 | ct-1 atom label (str) |
---|
423 | ct atom type (str) |
---|
424 | cx,cx+1,cx+2 the x,y and z coordinates (3 floats) |
---|
425 | cs-1 Sym Op symbol; sym. op number + unit cell id (e.g. '1,0,-1') (str) |
---|
426 | cs atom drawing style; e.g. 'balls & sticks' (str) |
---|
427 | cs+1 atom label style (e.g. 'name') (str) |
---|
428 | cs+2 atom color (RBG triplet) (int) |
---|
429 | cs+3 ADP flag: Isotropic ('I') or Anisotropic ('A') |
---|
430 | cs+4 Uiso (float) |
---|
431 | cs+5...cs+11 U11, U22, U33, U12, U13, U23 (6 floats) |
---|
432 | ci unique atom identifier; matches source atom Id in Atom Records (int) |
---|
433 | ============== ==================================================== |
---|
434 | |
---|
435 | Powder Diffraction Tree Items |
---|
436 | ----------------------------- |
---|
437 | |
---|
438 | .. _Powder_table: |
---|
439 | |
---|
440 | .. index:: |
---|
441 | single: Powder data object description |
---|
442 | single: Data object descriptions; Powder Data |
---|
443 | |
---|
444 | Every powder diffraction histogram is stored in the GSAS-II data tree |
---|
445 | with a top-level entry named beginning with the string "PWDR ". The |
---|
446 | diffraction data for that information are directly associated with |
---|
447 | that tree item and there are a series of children to that item. The |
---|
448 | routines :func:`GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree` |
---|
449 | and :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` will |
---|
450 | load this information into a dictionary where the child tree name is |
---|
451 | used as a key, and the information in the main entry is assigned |
---|
452 | a key of ``Data``, as outlined below. |
---|
453 | |
---|
454 | .. tabularcolumns:: |l|l|p{4in}| |
---|
455 | |
---|
456 | ====================== =============== ==================================================== |
---|
457 | key sub-key explanation |
---|
458 | ====================== =============== ==================================================== |
---|
459 | Limits \ A list of two two element lists, as [[Ld,Hd],[L,H]] |
---|
460 | where L and Ld are the current and default lowest |
---|
461 | two-theta value to be used and |
---|
462 | where H and Hd are the current and default highest |
---|
463 | two-theta value to be used. |
---|
464 | Reflection Lists \ A dict with an entry for each phase in the |
---|
465 | histogram. The contents of each dict item |
---|
466 | is a dict containing reflections, as described in |
---|
467 | the :ref:`Powder Reflections <PowderRefl_table>` |
---|
468 | description. |
---|
469 | Instrument Parameters \ A list containing two dicts where the possible |
---|
470 | keys in each dict are listed below. The value |
---|
471 | for each item is a list containing three values: |
---|
472 | the initial value, the current value and a |
---|
473 | refinement flag which can have a value of |
---|
474 | True, False or 0 where 0 indicates a value that |
---|
475 | cannot be refined. The first and second |
---|
476 | values are floats unless otherwise noted. |
---|
477 | Items in the first dict are noted as [1] |
---|
478 | \ Lam Specifies a wavelength in Angstroms [1] |
---|
479 | \ Lam1 Specifies the primary wavelength in |
---|
480 | Angstrom, when an alpha1, alpha2 |
---|
481 | source is used [1] |
---|
482 | \ Lam2 Specifies the secondary wavelength in |
---|
483 | Angstrom, when an alpha1, alpha2 |
---|
484 | source is used [1] |
---|
485 | I(L2)/I(L1) Ratio of Lam2 to Lam1 [1] |
---|
486 | \ Type Histogram type (str) [1]: |
---|
487 | * 'PXC' for constant wavelength x-ray |
---|
488 | * 'PNC' for constant wavelength neutron |
---|
489 | * 'PNT' for time of flight neutron |
---|
490 | \ Zero Two-theta zero correction in *degrees* [1] |
---|
491 | \ Azimuth Azimuthal setting angle for data recorded |
---|
492 | with differing setting angles [1] |
---|
493 | \ U, V, W Cagliotti profile coefficients |
---|
494 | for Gaussian instrumental broadening, where the |
---|
495 | FWHM goes as |
---|
496 | :math:`U \\tan^2\\theta + V \\tan\\theta + W` [1] |
---|
497 | \ X, Y Cauchy (Lorentzian) instrumental broadening |
---|
498 | coefficients [1] |
---|
499 | \ SH/L Variant of the Finger-Cox-Jephcoat asymmetric |
---|
500 | peak broadening ratio. Note that this is the |
---|
501 | average between S/L and H/L where S is |
---|
502 | sample height, H is the slit height and |
---|
503 | L is the goniometer diameter. [1] |
---|
504 | \ Polariz. Polarization coefficient. [1] |
---|
505 | wtFactor \ A weighting factor to increase or decrease |
---|
506 | the leverage of data in the histogram (float). |
---|
507 | A value of 1.0 weights the data with their |
---|
508 | standard uncertainties and a larger value |
---|
509 | increases the weighting of the data (equivalent |
---|
510 | to decreasing the uncertainties). |
---|
511 | Sample Parameters \ Specifies a dict with parameters that describe how |
---|
512 | the data were collected, as listed |
---|
513 | below. Refinable parameters are a list containing |
---|
514 | a float and a bool, where the second value |
---|
515 | specifies if the value is refined, otherwise |
---|
516 | the value is a float unless otherwise noted. |
---|
517 | \ Scale The histogram scale factor (refinable) |
---|
518 | \ Absorption The sample absorption coefficient as |
---|
519 | :math:`\\mu r` where r is the radius |
---|
520 | (refinable). |
---|
521 | \ DisplaceX, Sample displacement from goniometer center |
---|
522 | DisplaceY where Y is along the beam direction and |
---|
523 | X is perpendicular. Units are :math:`\\mu m` |
---|
524 | (refinable). |
---|
525 | \ Phi, Chi, Goniometer sample setting angles, in degrees. |
---|
526 | Omega |
---|
527 | \ Gonio. radius Radius of the diffractometer in mm |
---|
528 | \ InstrName A name for the instrument, used in preparing |
---|
529 | a CIF (str). |
---|
530 | \ Force, Variables that describe how the measurement |
---|
531 | Temperature, was performed. Not used directly in |
---|
532 | Humidity, any computations. |
---|
533 | Pressure, |
---|
534 | Voltage |
---|
535 | \ ranId The random-number Id for the histogram |
---|
536 | (same value as where top-level key is ranId) |
---|
537 | \ Type Type of diffraction data, may be 'Debye-Scherrer' |
---|
538 | or 'Bragg-Brentano' (str). |
---|
539 | \ Diffuse not in use? |
---|
540 | hId \ The number assigned to the histogram when |
---|
541 | the project is loaded or edited (can change) |
---|
542 | ranId \ A random number id for the histogram |
---|
543 | that does not change |
---|
544 | Background \ The background is stored as a list with where |
---|
545 | the first item in the list is list and the second |
---|
546 | item is a dict. The list contains the background |
---|
547 | function and its coefficients; the dict contains |
---|
548 | Debye diffuse terms and background peaks. |
---|
549 | (TODO: this needs to be expanded.) |
---|
550 | Data \ The data consist of a list of 6 np.arrays |
---|
551 | containing in order: |
---|
552 | |
---|
553 | 1. the x-postions (two-theta in degrees), |
---|
554 | 2. the intensity values (Yobs), |
---|
555 | 3. the weights for each Yobs value |
---|
556 | 4. the computed intensity values (Ycalc) |
---|
557 | 5. the background values |
---|
558 | 6. Yobs-Ycalc |
---|
559 | ====================== =============== ==================================================== |
---|
560 | |
---|
561 | Powder Reflection Data Structure |
---|
562 | -------------------------------- |
---|
563 | |
---|
564 | .. _PowderRefl_table: |
---|
565 | |
---|
566 | .. index:: |
---|
567 | single: Powder reflection object description |
---|
568 | single: Data object descriptions; Powder Reflections |
---|
569 | |
---|
570 | For every phase in a histogram, the ``Reflection Lists`` value is a dict |
---|
571 | one element of which is `'RefList'`, which is a np.array containing |
---|
572 | reflections. The columns in that array are documented below. |
---|
573 | |
---|
574 | ========== ==================================================== |
---|
575 | index explanation |
---|
576 | ========== ==================================================== |
---|
577 | 0,1,2 h,k,l (float) |
---|
578 | 3 multiplicity |
---|
579 | 4 d-space, Angstrom |
---|
580 | 5 pos, two-theta |
---|
581 | 6 sig, Gaussian width |
---|
582 | 7 gam, Lorenzian width |
---|
583 | 8 :math:`F_{obs}^2` |
---|
584 | 9 :math:`F_{calc}^2` |
---|
585 | 10 reflection phase, in degrees |
---|
586 | 11 intensity correction for reflection, this times |
---|
587 | :math:`F_{obs}^2` or :math:`F_{calc}^2` gives Iobs or Icalc |
---|
588 | ========== ==================================================== |
---|
589 | |
---|
590 | Single Crystal Tree Items |
---|
591 | ------------------------- |
---|
592 | |
---|
593 | .. _Xtal_table: |
---|
594 | |
---|
595 | .. index:: |
---|
596 | single: Single Crystal data object description |
---|
597 | single: Data object descriptions; Single crystal data |
---|
598 | |
---|
599 | Every single crystal diffraction histogram is stored in the GSAS-II data tree |
---|
600 | with a top-level entry named beginning with the string "HKLF ". The |
---|
601 | diffraction data for that information are directly associated with |
---|
602 | that tree item and there are a series of children to that item. The |
---|
603 | routines :func:`GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree` |
---|
604 | and :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` will |
---|
605 | load this information into a dictionary where the child tree name is |
---|
606 | used as a key, and the information in the main entry is assigned |
---|
607 | a key of ``Data``, as outlined below. |
---|
608 | |
---|
609 | .. tabularcolumns:: |l|l|p{4in}| |
---|
610 | |
---|
611 | ====================== =============== ==================================================== |
---|
612 | key sub-key explanation |
---|
613 | ====================== =============== ==================================================== |
---|
614 | Data \ A dict that contains the |
---|
615 | reflection table, |
---|
616 | as described in the |
---|
617 | :ref:`Single Crystal Reflections |
---|
618 | <XtalRefl_table>` |
---|
619 | description. |
---|
620 | |
---|
621 | Instrument Parameters \ A list containing two dicts where the possible |
---|
622 | keys in each dict are listed below. The value |
---|
623 | for most items is a list containing two values: |
---|
624 | the initial value, the current value. |
---|
625 | The first and second |
---|
626 | values are floats unless otherwise noted. |
---|
627 | \ Lam Specifies a wavelength in Angstroms (two floats) |
---|
628 | \ Type Histogram type (two str values): |
---|
629 | * 'SXC' for constant wavelength x-ray |
---|
630 | * 'SNC' for constant wavelength neutron |
---|
631 | * 'SNT' for time of flight neutron |
---|
632 | \ InstrName A name for the instrument, used in preparing |
---|
633 | a CIF (str). |
---|
634 | |
---|
635 | wtFactor \ A weighting factor to increase or decrease |
---|
636 | the leverage of data in the histogram (float). |
---|
637 | A value of 1.0 weights the data with their |
---|
638 | standard uncertainties and a larger value |
---|
639 | increases the weighting of the data (equivalent |
---|
640 | to decreasing the uncertainties). |
---|
641 | |
---|
642 | hId \ The number assigned to the histogram when |
---|
643 | the project is loaded or edited (can change) |
---|
644 | ranId \ A random number id for the histogram |
---|
645 | that does not change |
---|
646 | ====================== =============== ==================================================== |
---|
647 | |
---|
648 | Single Crystal Reflection Data Structure |
---|
649 | ---------------------------------------- |
---|
650 | |
---|
651 | .. _XtalRefl_table: |
---|
652 | |
---|
653 | .. index:: |
---|
654 | single: Single Crystal reflection object description |
---|
655 | single: Data object descriptions; Single Crystal Reflections |
---|
656 | |
---|
657 | For every single crystal a histogram, the ``'Data'`` item contains |
---|
658 | the structure factors as an np.array in item `'RefList'`. |
---|
659 | The columns in that array are documented below. |
---|
660 | |
---|
661 | ========== ==================================================== |
---|
662 | index explanation |
---|
663 | ========== ==================================================== |
---|
664 | 0,1,2 h,k,l (float) |
---|
665 | 3 multiplicity |
---|
666 | 4 d-space, Angstrom |
---|
667 | 5 :math:`F_{obs}^2` |
---|
668 | 6 :math:`\sigma(F_{obs}^2)` |
---|
669 | 7 :math:`F_{calc}^2` |
---|
670 | 8 :math:`F_{obs}^2T` |
---|
671 | 9 :math:`F_{calc}^2T` |
---|
672 | 10 reflection phase, in degrees |
---|
673 | 11 intensity correction for reflection, this times |
---|
674 | :math:`F_{obs}^2` or :math:`F_{calc}^2` |
---|
675 | gives Iobs or Icalc |
---|
676 | ========== ==================================================== |
---|
677 | |
---|
678 | |
---|
679 | *Classes and routines* |
---|
680 | ---------------------- |
---|
681 | |
---|
682 | ''' |
---|
683 | import random as ran |
---|
684 | import sys |
---|
685 | import GSASIIpath |
---|
686 | import GSASIImath as G2mth |
---|
687 | |
---|
688 | GSASIIpath.SetVersionNumber("$Revision: 1168 $") |
---|
689 | |
---|
690 | DefaultControls = { |
---|
691 | 'deriv type':'analytic Hessian', #default controls |
---|
692 | 'min dM/M':0.0001,'shift factor':1.,'max cyc':3,'F**2':True, |
---|
693 | 'minF/sig':0, |
---|
694 | 'Author':'no name', |
---|
695 | 'FreeVar1':'Sample humidity (%)', |
---|
696 | 'FreeVar2':'Sample voltage (V)', |
---|
697 | 'FreeVar3':'Applied load (MN)', |
---|
698 | } |
---|
699 | '''Values to be used as defaults for the initial contents of the ``Controls`` |
---|
700 | data tree item. |
---|
701 | ''' |
---|
702 | |
---|
703 | def MakeUniqueLabel(lbl,labellist): |
---|
704 | '''Make sure that every a label is unique against a list by adding |
---|
705 | digits at the end until it is not found in list. |
---|
706 | |
---|
707 | :param str lbl: the input label |
---|
708 | :param list labellist: the labels that have already been encountered |
---|
709 | :returns: lbl if not found in labellist or lbl with ``_1-9` (or |
---|
710 | ``_10-99``, etc.) appended at the end |
---|
711 | ''' |
---|
712 | lbl = lbl.strip() |
---|
713 | if not lbl: # deal with a blank label |
---|
714 | lbl = '_1' |
---|
715 | if lbl not in labellist: |
---|
716 | labellist.append(lbl) |
---|
717 | return lbl |
---|
718 | i = 1 |
---|
719 | prefix = lbl |
---|
720 | if '_' in lbl: |
---|
721 | prefix = lbl[:lbl.rfind('_')] |
---|
722 | suffix = lbl[lbl.rfind('_')+1:] |
---|
723 | try: |
---|
724 | i = int(suffix)+1 |
---|
725 | except: # suffix could not be parsed |
---|
726 | i = 1 |
---|
727 | prefix = lbl |
---|
728 | while prefix+'_'+str(i) in labellist: |
---|
729 | i += 1 |
---|
730 | else: |
---|
731 | lbl = prefix+'_'+str(i) |
---|
732 | labellist.append(lbl) |
---|
733 | return lbl |
---|
734 | |
---|
735 | PhaseIdLookup = {} |
---|
736 | '''dict listing phase name and random Id keyed by sequential phase index as a str; |
---|
737 | best to access this using :func:`LookupPhaseName` |
---|
738 | ''' |
---|
739 | PhaseRanIdLookup = {} |
---|
740 | '''dict listing phase sequential index keyed by phase random Id; |
---|
741 | best to access this using :func:`LookupPhaseId` |
---|
742 | ''' |
---|
743 | HistIdLookup = {} |
---|
744 | '''dict listing histogram name and random Id, keyed by sequential histogram index as a str; |
---|
745 | best to access this using :func:`LookupHistName` |
---|
746 | ''' |
---|
747 | HistRanIdLookup = {} |
---|
748 | '''dict listing histogram sequential index keyed by histogram random Id; |
---|
749 | best to access this using :func:`LookupHistId` |
---|
750 | ''' |
---|
751 | AtomIdLookup = {} |
---|
752 | '''dict listing for each phase index as a str, the atom label and atom random Id, |
---|
753 | keyed by atom sequential index as a str; |
---|
754 | best to access this using :func:`LookupAtomLabel` |
---|
755 | ''' |
---|
756 | AtomRanIdLookup = {} |
---|
757 | '''dict listing for each phase the atom sequential index keyed by atom random Id; |
---|
758 | best to access this using :func:`LookupAtomId` |
---|
759 | ''' |
---|
760 | ShortPhaseNames = {} |
---|
761 | '''a dict containing a possibly shortened and when non-unique numbered |
---|
762 | version of the phase name. Keyed by the phase sequential index. |
---|
763 | ''' |
---|
764 | ShortHistNames = {} |
---|
765 | '''a dict containing a possibly shortened and when non-unique numbered |
---|
766 | version of the histogram name. Keyed by the histogram sequential index. |
---|
767 | ''' |
---|
768 | |
---|
769 | VarDesc = {} |
---|
770 | ''' This dictionary lists descriptions for GSAS-II variables, |
---|
771 | as set in :func:`CompileVarDesc`. See that function for a description |
---|
772 | for how keys and values are written. |
---|
773 | ''' |
---|
774 | |
---|
775 | reVarDesc = {} |
---|
776 | ''' This dictionary lists descriptions for GSAS-II variables with |
---|
777 | the same values as :attr:`VarDesc` except that keys have been compiled as |
---|
778 | regular expressions. Initialized in :func:`CompileVarDesc`. |
---|
779 | ''' |
---|
780 | |
---|
781 | def IndexAllIds(Histograms,Phases): |
---|
782 | '''Scan through the used phases & histograms and create an index |
---|
783 | to the random numbers of phases, histograms and atoms. While doing this, |
---|
784 | confirm that assigned random numbers are unique -- just in case lightning |
---|
785 | strikes twice in the same place. |
---|
786 | |
---|
787 | Note: this code assumes that the atom random Id (ranId) is the last |
---|
788 | element each atom record. |
---|
789 | |
---|
790 | This is called in two places (only) :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` |
---|
791 | (which loads the histograms and phases from a GPX file) and |
---|
792 | :meth:`GSASII.GSASII.GetUsedHistogramsAndPhases` |
---|
793 | (which loads the histograms and phases from the data tree.) |
---|
794 | |
---|
795 | TODO: do we need a lookup for rigid body variables? |
---|
796 | ''' |
---|
797 | # process phases and atoms |
---|
798 | PhaseIdLookup.clear() |
---|
799 | PhaseRanIdLookup.clear() |
---|
800 | AtomIdLookup.clear() |
---|
801 | AtomRanIdLookup.clear() |
---|
802 | ShortPhaseNames.clear() |
---|
803 | for ph in Phases: |
---|
804 | cx,ct,cs,cia = Phases[ph]['General']['AtomPtrs'] |
---|
805 | ranId = Phases[ph]['ranId'] |
---|
806 | while ranId in PhaseRanIdLookup: |
---|
807 | # Found duplicate random Id! note and reassign |
---|
808 | print ("\n\n*** Phase "+str(ph)+" has repeated ranId. Fixing.\n") |
---|
809 | Phases[ph]['ranId'] = ranId = ran.randint(0,sys.maxint) |
---|
810 | pId = str(Phases[ph]['pId']) |
---|
811 | PhaseIdLookup[pId] = (ph,ranId) |
---|
812 | PhaseRanIdLookup[ranId] = pId |
---|
813 | shortname = ph[:10] |
---|
814 | while shortname in ShortPhaseNames.values(): |
---|
815 | shortname = ph[:8] + ' ('+ pId + ')' |
---|
816 | ShortPhaseNames[pId] = shortname |
---|
817 | AtomIdLookup[pId] = {} |
---|
818 | AtomRanIdLookup[pId] = {} |
---|
819 | for iatm,at in enumerate(Phases[ph]['Atoms']): |
---|
820 | ranId = at[-1] |
---|
821 | while ranId in AtomRanIdLookup[pId]: # check for dups |
---|
822 | print ("\n\n*** Phase "+str(ph)+" atom "+str(iatm)+" has repeated ranId. Fixing.\n") |
---|
823 | at[-1] = ranId = ran.randint(0,sys.maxint) |
---|
824 | AtomRanIdLookup[pId][ranId] = str(iatm) |
---|
825 | if Phases[ph]['General']['Type'] == 'macromolecular': |
---|
826 | label = '%s_%s_%s_%s'%(at[ct-1],at[ct-3],at[ct-4],at[ct-2]) |
---|
827 | else: |
---|
828 | label = at[ct-1] |
---|
829 | AtomIdLookup[pId][str(iatm)] = (label,ranId) |
---|
830 | # process histograms |
---|
831 | HistIdLookup.clear() |
---|
832 | HistRanIdLookup.clear() |
---|
833 | ShortHistNames.clear() |
---|
834 | for hist in Histograms: |
---|
835 | ranId = Histograms[hist]['ranId'] |
---|
836 | while ranId in HistRanIdLookup: |
---|
837 | # Found duplicate random Id! note and reassign |
---|
838 | print ("\n\n*** Histogram "+str(hist)+" has repeated ranId. Fixing.\n") |
---|
839 | Histograms[hist]['ranId'] = ranId = ran.randint(0,sys.maxint) |
---|
840 | hId = str(Histograms[hist]['hId']) |
---|
841 | HistIdLookup[hId] = (hist,ranId) |
---|
842 | HistRanIdLookup[ranId] = hId |
---|
843 | shortname = hist[:15] |
---|
844 | while shortname in ShortHistNames.values(): |
---|
845 | shortname = hist[:11] + ' ('+ hId + ')' |
---|
846 | ShortHistNames[hId] = shortname |
---|
847 | |
---|
848 | def LookupAtomId(pId,ranId): |
---|
849 | '''Get the atom number from a phase and atom random Id |
---|
850 | |
---|
851 | :param int/str pId: the sequential number of the phase |
---|
852 | :param int ranId: the random Id assigned to an atom |
---|
853 | |
---|
854 | :returns: the index number of the atom (str) |
---|
855 | ''' |
---|
856 | if not AtomRanIdLookup: |
---|
857 | raise Exception,'Error: LookupAtomId called before IndexAllIds was run' |
---|
858 | if pId is None or pId == '': |
---|
859 | raise KeyError,'Error: phase is invalid (None or blank)' |
---|
860 | pId = str(pId) |
---|
861 | if pId not in AtomRanIdLookup: |
---|
862 | raise KeyError,'Error: LookupAtomId does not have phase '+pId |
---|
863 | if ranId not in AtomRanIdLookup[pId]: |
---|
864 | raise KeyError,'Error: LookupAtomId, ranId '+str(ranId)+' not in AtomRanIdLookup['+pId+']' |
---|
865 | return AtomRanIdLookup[pId][ranId] |
---|
866 | |
---|
867 | def LookupAtomLabel(pId,index): |
---|
868 | '''Get the atom label from a phase and atom index number |
---|
869 | |
---|
870 | :param int/str pId: the sequential number of the phase |
---|
871 | :param int index: the index of the atom in the list of atoms |
---|
872 | |
---|
873 | :returns: the label for the atom (str) and the random Id of the atom (int) |
---|
874 | ''' |
---|
875 | if not AtomIdLookup: |
---|
876 | raise Exception,'Error: LookupAtomLabel called before IndexAllIds was run' |
---|
877 | if pId is None or pId == '': |
---|
878 | raise KeyError,'Error: phase is invalid (None or blank)' |
---|
879 | pId = str(pId) |
---|
880 | if pId not in AtomIdLookup: |
---|
881 | raise KeyError,'Error: LookupAtomLabel does not have phase '+pId |
---|
882 | if index not in AtomIdLookup[pId]: |
---|
883 | raise KeyError,'Error: LookupAtomLabel, ranId '+str(index)+' not in AtomRanIdLookup['+pId+']' |
---|
884 | return AtomIdLookup[pId][index] |
---|
885 | |
---|
886 | def LookupPhaseId(ranId): |
---|
887 | '''Get the phase number and name from a phase random Id |
---|
888 | |
---|
889 | :param int ranId: the random Id assigned to a phase |
---|
890 | :returns: the sequential Id (pId) number for the phase (str) |
---|
891 | ''' |
---|
892 | if not PhaseRanIdLookup: |
---|
893 | raise Exception,'Error: LookupPhaseId called before IndexAllIds was run' |
---|
894 | if ranId not in PhaseRanIdLookup: |
---|
895 | raise KeyError,'Error: LookupPhaseId does not have ranId '+str(ranId) |
---|
896 | return PhaseRanIdLookup[ranId] |
---|
897 | |
---|
898 | def LookupPhaseName(pId): |
---|
899 | '''Get the phase number and name from a phase Id |
---|
900 | |
---|
901 | :param int/str pId: the sequential assigned to a phase |
---|
902 | :returns: (phase,ranId) where phase is the name of the phase (str) |
---|
903 | and ranId is the random # id for the phase (int) |
---|
904 | ''' |
---|
905 | if not PhaseIdLookup: |
---|
906 | raise Exception,'Error: LookupPhaseName called before IndexAllIds was run' |
---|
907 | if pId is None or pId == '': |
---|
908 | raise KeyError,'Error: phase is invalid (None or blank)' |
---|
909 | pId = str(pId) |
---|
910 | if pId not in PhaseIdLookup: |
---|
911 | raise KeyError,'Error: LookupPhaseName does not have index '+pId |
---|
912 | return PhaseIdLookup[pId] |
---|
913 | |
---|
914 | def LookupHistId(ranId): |
---|
915 | '''Get the histogram number and name from a histogram random Id |
---|
916 | |
---|
917 | :param int ranId: the random Id assigned to a histogram |
---|
918 | :returns: the sequential Id (hId) number for the histogram (str) |
---|
919 | ''' |
---|
920 | if not HistRanIdLookup: |
---|
921 | raise Exception,'Error: LookupHistId called before IndexAllIds was run' |
---|
922 | if ranId not in HistRanIdLookup: |
---|
923 | raise KeyError,'Error: LookupHistId does not have ranId '+str(ranId) |
---|
924 | return HistRanIdLookup[ranId] |
---|
925 | |
---|
926 | def LookupHistName(hId): |
---|
927 | '''Get the histogram number and name from a histogram Id |
---|
928 | |
---|
929 | :param int/str hId: the sequential assigned to a histogram |
---|
930 | :returns: (hist,ranId) where hist is the name of the histogram (str) |
---|
931 | and ranId is the random # id for the histogram (int) |
---|
932 | ''' |
---|
933 | if not HistIdLookup: |
---|
934 | raise Exception,'Error: LookupHistName called before IndexAllIds was run' |
---|
935 | if hId is None or hId == '': |
---|
936 | raise KeyError,'Error: histogram is invalid (None or blank)' |
---|
937 | hId = str(hId) |
---|
938 | if hId not in HistIdLookup: |
---|
939 | raise KeyError,'Error: LookupHistName does not have index '+hId |
---|
940 | return HistIdLookup[hId] |
---|
941 | |
---|
942 | def fmtVarDescr(varname): |
---|
943 | '''Return a string with a more complete description for a GSAS-II variable |
---|
944 | |
---|
945 | TODO: This will not handle rigid body parameters yet |
---|
946 | |
---|
947 | :param str name: A full G2 variable name with 2 or 3 |
---|
948 | colons (<p>:<h>:name[:<a>]) |
---|
949 | |
---|
950 | :returns: a string with the description |
---|
951 | ''' |
---|
952 | |
---|
953 | l = getVarDescr(varname) |
---|
954 | if not l: |
---|
955 | return "invalid variable name ("+str(varname)+")!" |
---|
956 | |
---|
957 | if not l[4]: |
---|
958 | l[4] = "(variable needs a definition!)" |
---|
959 | |
---|
960 | s = "" |
---|
961 | if l[0] is not None and l[1] is not None: # HAP: keep short |
---|
962 | lbl = ShortPhaseNames.get(l[0],'? #'+str(l[0])) |
---|
963 | hlbl = ShortHistNames.get(l[1],'? #'+str(l[1])) |
---|
964 | if hlbl[:4] == 'HKLF': |
---|
965 | hlbl = 'Xtl='+hlbl[5:] |
---|
966 | elif hlbl[:4] == 'PWDR': |
---|
967 | hlbl = 'Pwd='+hlbl[5:] |
---|
968 | else: |
---|
969 | hlbl = 'Hist='+hlbl |
---|
970 | s = "Ph="+str(lbl)+" * "+str(hlbl)+": " |
---|
971 | elif l[3] is not None: # atom parameter, |
---|
972 | lbl = ShortPhaseNames.get(l[0],'phase?') |
---|
973 | try: |
---|
974 | albl = LookupAtomLabel(l[0],l[3])[0] |
---|
975 | except KeyError: |
---|
976 | albl = 'Atom?' |
---|
977 | s = "Atom "+str(albl)+" in "+str(lbl)+": " |
---|
978 | elif l[0] is not None: |
---|
979 | lbl = ShortPhaseNames.get(l[0],'phase?') |
---|
980 | s = "Phase "+str(lbl)+": " |
---|
981 | elif l[1] is not None: |
---|
982 | hlbl = ShortHistNames.get(l[1],'? #'+str(l[1])) |
---|
983 | if hlbl[:4] == 'HKLF': |
---|
984 | hlbl = 'Xtl='+hlbl[5:] |
---|
985 | elif hlbl[:4] == 'PWDR': |
---|
986 | hlbl = 'Pwd='+hlbl[5:] |
---|
987 | else: |
---|
988 | hlbl = 'Hist='+hlbl |
---|
989 | s = str(hlbl)+": " |
---|
990 | if not s: |
---|
991 | s = 'Global: ' |
---|
992 | s += l[4] |
---|
993 | return s |
---|
994 | |
---|
995 | def getVarDescr(varname): |
---|
996 | '''Return a short description for a GSAS-II variable |
---|
997 | |
---|
998 | :param str name: A full G2 variable name with 2 or 3 |
---|
999 | colons (<p>:<h>:name[:<a>]) |
---|
1000 | |
---|
1001 | :returns: a five element list as [`p`,`h`,`name`,`a`,`description`], |
---|
1002 | where `p`, `h`, `a` are str values or `None`, for the phase number, |
---|
1003 | the histogram number and the atom number; `name` will always be |
---|
1004 | an str; and `description` is str or `None`. |
---|
1005 | If the variable name is incorrectly formed (for example, wrong |
---|
1006 | number of colons), `None` is returned instead of a list. |
---|
1007 | ''' |
---|
1008 | l = varname.split(':') |
---|
1009 | if len(l) == 3: |
---|
1010 | l += [None] |
---|
1011 | if len(l) != 4: |
---|
1012 | return None |
---|
1013 | for i in (0,1,3): |
---|
1014 | if l[i] == "": |
---|
1015 | l[i] = None |
---|
1016 | l += [getDescr(l[2])] |
---|
1017 | return l |
---|
1018 | |
---|
1019 | def CompileVarDesc(): |
---|
1020 | '''Set the values in the variable description lookup table (:attr:`VarDesc`) |
---|
1021 | into :attr:`reVarDesc`. This is called in :func:`getDescr` so the initialization |
---|
1022 | is always done before use. |
---|
1023 | |
---|
1024 | Note that keys may contain regular expressions, where '[xyz]' |
---|
1025 | matches 'x' 'y' or 'z' (equivalently '[x-z]' describes this as range of values). |
---|
1026 | '.*' matches any string. For example:: |
---|
1027 | |
---|
1028 | 'AUiso':'Atomic isotropic displacement parameter', |
---|
1029 | |
---|
1030 | will match variable ``'p::AUiso:a'``. |
---|
1031 | If parentheses are used in the key, the contents of those parentheses can be |
---|
1032 | used in the value, such as:: |
---|
1033 | |
---|
1034 | 'AU([123][123])':'Atomic anisotropic displacement parameter U\\1', |
---|
1035 | |
---|
1036 | will match ``AU11``, ``AU23``,.. and `U11`, `U23` etc will be displayed |
---|
1037 | in the value when used. |
---|
1038 | |
---|
1039 | ''' |
---|
1040 | import re |
---|
1041 | if reVarDesc: return # already done |
---|
1042 | for key,value in { |
---|
1043 | # Phase vars (p::<var>) |
---|
1044 | 'A([0-5])' : 'Reciprocal metric tensor component \\1', |
---|
1045 | 'Vol' : 'Unit cell volume????', |
---|
1046 | # Atom vars (p::<var>:a) |
---|
1047 | 'dA([xyz])' : 'change to atomic position \\1', |
---|
1048 | 'AUiso':'Atomic isotropic displacement parameter', |
---|
1049 | 'AU([123][123])':'Atomic anisotropic displacement parameter U\\1', |
---|
1050 | 'Afrac': 'Atomic occupancy parameter', |
---|
1051 | # Hist & Phase (HAP) vars (p:h:<var>) |
---|
1052 | 'Bab([AU])': 'Babinet solvent scattering coef. \\1', |
---|
1053 | 'D([123][123])' : 'Anisotropic strain coef. \\1', |
---|
1054 | 'Extinction' : 'Extinction coef.', |
---|
1055 | 'MD' : 'March-Dollase coef.', |
---|
1056 | 'Mustrain;.*' : 'Microstrain coef.', |
---|
1057 | 'Scale' : 'Phase scale factor', |
---|
1058 | 'Size;.*' : 'Crystallite size value', |
---|
1059 | 'eA' : '?', |
---|
1060 | #Histogram vars (:h:<var>) |
---|
1061 | 'Absorption' : 'Absorption coef.', |
---|
1062 | 'Displace([XY])' : 'Debye-Scherrer sample displacement \\1', |
---|
1063 | 'Lam' : 'Wavelength', |
---|
1064 | 'Polariz\.' : 'Polarization correction', |
---|
1065 | 'SH/L' : 'FCJ peak asymmetry correction', |
---|
1066 | 'Scale' : 'Histogram scale factor', |
---|
1067 | '([UVW])' : 'Gaussian instrument broadening \\1', |
---|
1068 | '([XY])' : 'Cauchy instrument broadening \\1', |
---|
1069 | 'Zero' : 'Debye-Scherrer zero correction', |
---|
1070 | 'nDebye' : 'Debye model background corr. terms', |
---|
1071 | 'nPeaks' : 'Fixed peak background corr. terms', |
---|
1072 | # Global vars (::<var>) |
---|
1073 | }.items(): |
---|
1074 | VarDesc[key] = value |
---|
1075 | reVarDesc[re.compile(key)] = value |
---|
1076 | |
---|
1077 | def getDescr(name): |
---|
1078 | '''Return a short description for a GSAS-II variable |
---|
1079 | |
---|
1080 | :param str name: The descriptive part of the variable name without colons (:) |
---|
1081 | |
---|
1082 | :returns: a short description or None if not found |
---|
1083 | ''' |
---|
1084 | |
---|
1085 | CompileVarDesc() # compile the regular expressions, if needed |
---|
1086 | for key in reVarDesc: |
---|
1087 | m = key.match(name) |
---|
1088 | if m: |
---|
1089 | return m.expand(reVarDesc[key]) |
---|
1090 | return None |
---|
1091 | |
---|
1092 | def _lookup(dic,key): |
---|
1093 | '''Lookup a key in a dictionary, where None returns an empty string |
---|
1094 | but an unmatched key returns a question mark. Used in :class:`G2VarObj` |
---|
1095 | ''' |
---|
1096 | if key is None: |
---|
1097 | return "" |
---|
1098 | else: |
---|
1099 | return dic.get(key,'?') |
---|
1100 | |
---|
1101 | class G2VarObj(object): |
---|
1102 | '''Defines a GSAS-II variable either using the phase/atom/histogram |
---|
1103 | unique Id numbers or using a character string that specifies |
---|
1104 | variables by phase/atom/histogram number (which can change). |
---|
1105 | Note that :func:`LoadID` should be used to (re)load the current Ids |
---|
1106 | before creating or later using the G2VarObj object. |
---|
1107 | |
---|
1108 | A :class:`G2VarObj` object can be created with a single parameter: |
---|
1109 | |
---|
1110 | :param str/tuple varname: a single value can be used to create a :class:`G2VarObj` |
---|
1111 | object. If a string, it must be of form "p:h:var" or "p:h:var:a", where |
---|
1112 | |
---|
1113 | * p is the phase number (which may be left blank); |
---|
1114 | * h is the histogram number (which may be left blank); |
---|
1115 | * a is the atom number (which may be left blank in which case the third colon is omitted). |
---|
1116 | |
---|
1117 | Alternately a single tuple of form (Phase,Histogram,VarName,AtomID) can be used, where |
---|
1118 | Phase, Histogram, and AtomID are None or are ranId values and VarName is a string. |
---|
1119 | |
---|
1120 | If four positional arguments are supplied, they are: |
---|
1121 | |
---|
1122 | :param str/int phasenum: The number for the phase |
---|
1123 | :param str/int histnum: The number for the histogram |
---|
1124 | :param str varname: a single value can be used to create a :class:`G2VarObj` |
---|
1125 | :param str/int atomnum: The number for the atom |
---|
1126 | |
---|
1127 | ''' |
---|
1128 | IDdict = {} |
---|
1129 | IDdict['phases'] = {} |
---|
1130 | IDdict['hists'] = {} |
---|
1131 | IDdict['atoms'] = {} |
---|
1132 | def __init__(self,*args): |
---|
1133 | self.phase = None |
---|
1134 | self.histogram = None |
---|
1135 | self.name = '' |
---|
1136 | self.atom = None |
---|
1137 | if len(args) == 1 and (type(args[0]) is list or type(args[0]) is tuple) and len(args[0]) == 4: |
---|
1138 | self.phase,self.histogram,self.name,self.atom = args[0] |
---|
1139 | elif len(args) == 1 and ':' in args[0]: |
---|
1140 | lst = args[0].split(':') |
---|
1141 | self.phase = PhaseIdLookup.get(lst[0],[None,None])[1] |
---|
1142 | self.histogram = HistIdLookup.get(lst[1],[None,None])[1] |
---|
1143 | self.name = lst[2] |
---|
1144 | if len(lst) > 3: |
---|
1145 | self.atom = AtomIdLookup[lst[0]].get(lst[3],[None,None])[1] |
---|
1146 | elif len(args) == 4: |
---|
1147 | self.phase = PhaseIdLookup.get(str(args[0]),[None,None])[1] |
---|
1148 | self.histogram = HistIdLookup.get(str(args[1]),[None,None])[1] |
---|
1149 | self.name = args[2] |
---|
1150 | self.atom = AtomIdLookup[args[0]].get(str(args[3]),[None,None])[1] |
---|
1151 | else: |
---|
1152 | raise Exception,"Incorrectly called GSAS-II parameter name" |
---|
1153 | |
---|
1154 | #print "DEBUG: created ",self.phase,self.histogram,self.name,self.atom |
---|
1155 | |
---|
1156 | def __str__(self): |
---|
1157 | return self.varname() |
---|
1158 | |
---|
1159 | def varname(self): |
---|
1160 | '''Formats the GSAS-II variable name as a "traditional" GSAS-II variable |
---|
1161 | string (p:h:<var>:a) or (p:h:<var>) |
---|
1162 | |
---|
1163 | :returns: the variable name as a str |
---|
1164 | ''' |
---|
1165 | ph = _lookup(PhaseRanIdLookup,self.phase) |
---|
1166 | hist = _lookup(HistRanIdLookup,self.histogram) |
---|
1167 | s = (ph + ":" + hist + ":" + |
---|
1168 | str(self.name)) |
---|
1169 | if self.atom: |
---|
1170 | if ph in AtomRanIdLookup: |
---|
1171 | s += ":" + AtomRanIdLookup[ph].get(self.atom,'?') |
---|
1172 | else: |
---|
1173 | s += ":?" |
---|
1174 | return s |
---|
1175 | |
---|
1176 | def __repr__(self): |
---|
1177 | '''Return the detailed contents of the object |
---|
1178 | ''' |
---|
1179 | s = "<" |
---|
1180 | if self.phase is not None: |
---|
1181 | ph = _lookup(PhaseRanIdLookup,self.phase) |
---|
1182 | s += "Phase: rId=" + str(self.phase) + " (#"+ ph + "); " |
---|
1183 | if self.histogram is not None: |
---|
1184 | hist = _lookup(HistRanIdLookup,self.histogram) |
---|
1185 | s += "Histogram: rId=" + str(self.histogram) + " (#"+ hist + "); " |
---|
1186 | if self.atom is not None: |
---|
1187 | s += "Atom rId=" + str(self.atom) |
---|
1188 | if ph in AtomRanIdLookup: |
---|
1189 | s += " (#" + AtomRanIdLookup[ph].get(self.atom,'?') + "); " |
---|
1190 | else: |
---|
1191 | s += " (#? -- not found!); " |
---|
1192 | s += 'Variable name="' + str(self.name) + '">' |
---|
1193 | return s+"("+self.varname()+")" |
---|
1194 | |
---|
1195 | def __eq__(self, other): |
---|
1196 | if type(other) is type(self): |
---|
1197 | return (self.phase == other.phase and |
---|
1198 | self.histogram == other.histogram and |
---|
1199 | self.name == other.name and |
---|
1200 | self.atom == other.atom) |
---|
1201 | return False |
---|
1202 | |
---|
1203 | def _show(self): |
---|
1204 | 'For testing, shows the current lookup table' |
---|
1205 | print 'phases', self.IDdict['phases'] |
---|
1206 | print 'hists', self.IDdict['hists'] |
---|
1207 | print 'atomDict', self.IDdict['atoms'] |
---|
1208 | |
---|