source: trunk/GSASIIobj.py

Last change on this file was 5025, checked in by vondreele, 7 days ago

Add shell plots to 3D reflection plots.
Fix escape character issues in CopyRietveld2Origin
Revise LeBail? operation - now works in sequential refinements. Remove "Fit LeBail?" from Calculate menu.
OnLeBail? does 10 cycles instead of 3 in the "zero cycle" step.
Initial FOSQ set to Scale*Phase fr. instead of 1.0 - improves initial LeBaill? fit.
Put newLeBail flag in Controls - as more global.
Allow LeBail? refinements in sequential refinements.
Remove a wx.Yield() - obsolete
Disable deletion of selected parameters in a LeBail? refinement & let the LS handle any singularities that may ensue
Derivatives no longer controlled by LeBail? flag.

  • Property svn:eol-style set to native
  • Property svn:keywords set to Date Author Revision URL Id
File size: 164.4 KB
Line 
1# -*- coding: utf-8 -*-
2#GSASIIobj - data objects for GSAS-II
3########### SVN repository information ###################
4# $Date: 2021-09-14 16:55:59 +0000 (Tue, 14 Sep 2021) $
5# $Author: vondreele $
6# $Revision: 5025 $
7# $URL: trunk/GSASIIobj.py $
8# $Id: GSASIIobj.py 5025 2021-09-14 16:55:59Z vondreele $
9########### SVN repository information ###################
10
11'''
12*GSASIIobj: Data objects*
13=========================
14
15This module defines and/or documents the data structures used in GSAS-II, as well
16as provides misc. support routines.
17
18.. Next command allows \\AA to be used in HTML
19
20.. only:: html
21
22   :math:`\\require{mediawiki-texvc}`
23
24.. _Constraints_table:
25
26.. index::
27   single: Constraints object description
28   single: Data object descriptions; Constraints
29
30Constraints Tree Item
31----------------------
32
33Constraints are stored in a dict, separated into groups.
34Note that parameter are named in the following pattern,
35p:h:<var>:n, where p is the phase number, h is the histogram number
36<var> is a variable name and n is the parameter number.
37If a parameter does not depend on a histogram or phase or is unnumbered, that
38number is omitted.
39Note that the contents of each dict item is a List where each element in the
40list is a :ref:`constraint definition objects <Constraint_definitions_table>`.
41The constraints in this form are converted in
42:func:`GSASIIstrIO.ProcessConstraints` to the form used in :mod:`GSASIImapvars`
43
44The keys in the Constraints dict are:
45
46.. tabularcolumns:: |l|p{4.5in}|
47
48==========  ====================================================
49  key         explanation
50==========  ====================================================
51Hist        This specifies a list of constraints on
52            histogram-related parameters,
53            which will be of form :h:<var>:n.
54HAP         This specifies a list of constraints on parameters
55            that are defined for every histogram in each phase
56            and are of form p:h:<var>:n.
57Phase       This specifies a list of constraints on phase
58            parameters,
59            which will be of form p::<var>:n.
60Global      This specifies a list of constraints on parameters
61            that are not tied to a histogram or phase and
62            are of form ::<var>:n
63==========  ====================================================
64
65.. _Constraint_definitions_table:
66
67.. index::
68   single: Constraint definition object description
69   single: Data object descriptions; Constraint Definition
70
71Each constraint is defined as an item in a list. Each constraint is of form::
72
73[[<mult1>, <var1>], [<mult2>, <var2>],..., <fixedval>, <varyflag>, <constype>]
74
75Where the variable pair list item containing two values [<mult>, <var>], where:
76
77  * <mult> is a multiplier for the constraint (float)
78  * <var> a :class:`G2VarObj` object (previously a str variable name of form
79      'p:h:name[:at]')
80
81Note that the last three items in the list play a special role:
82
83 * <fixedval> is the fixed value for a `constant equation` (``constype=c``)
84   constraint or is None. For a `New variable` (``constype=f``) constraint,
85   a variable name can be specified as a str (used for externally
86   generated constraints)
87 * <varyflag> is True or False for `New variable` (``constype=f``) constraints
88   or is None. This indicates if this variable should be refined.
89 * <constype> is one of four letters, 'e', 'c', 'h', 'f' that determines the type of constraint:
90
91    * 'e' defines a set of equivalent variables. Only the first variable is refined (if the
92      appropriate refine flag is set) and and all other equivalent variables in the list
93      are generated from that variable, using the appropriate multipliers.
94    * 'c' defines a constraint equation of form,
95      :math:`m_1 \\times var_1 + m_2 \\times var_2 + ... = c`
96    * 'h' defines a variable to hold (not vary). Any variable on this list is not varied,
97      even if its refinement flag is set. Only one [mult,var] pair is allowed in a hold
98      constraint and the mult value is ignored.
99      This is of particular value when needing to hold one or more variables where a
100      single flag controls a set of variables such as, coordinates,
101      the reciprocal metric tensor or anisotropic displacement parameter.
102    * 'f' defines a new variable (function) according to relationship
103      :math:`newvar = m_1 \\times var_1 + m_2 \\times var_2 + ...`
104
105.. _Covariance_table:
106
107.. index::
108   single: Covariance description
109   single: Data object descriptions; Covariance
110
111Covariance Tree Item
112--------------------
113
114The Covariance tree item has results from the last least-squares run. They
115are stored in a dict with these keys:
116
117.. tabularcolumns:: |l|l|p{4in}|
118
119=============  ===============  ====================================================
120  key            sub-key        explanation
121=============  ===============  ====================================================
122newCellDict    \\                (dict) ith lattice parameters computed by
123                                :func:`GSASIIstrMath.GetNewCellParms`
124title          \\                (str) Name of gpx file(?)
125variables      \\                (list) Values for all N refined variables
126                                (list of float values, length N,
127                                ordered to match varyList)
128sig            \\                (list) Uncertainty values for all N refined variables
129                                (list of float values, length N,
130                                ordered to match varyList)
131varyList       \\                (list of str values, length N) List of directly refined variables
132                               
133newAtomDict    \\                (dict) atom position values computed in
134                                :func:`GSASIIstrMath.ApplyXYZshifts`
135Rvals          \\                (dict) R-factors, GOF, Marquardt value for last
136                                refinement cycle
137\\              Nobs             (int) Number of observed data points
138\\              Rwp              (float) overall weighted profile R-factor (%)
139\\              chisq            (float) :math:`\\sum w*(I_{obs}-I_{calc})^2`                               
140                                for all data.
141                                Note: this is not the reduced :math:`\\chi^2`.
142\\              lamMax           (float) Marquardt value applied to Hessian diagonal
143\\              GOF              (float) The goodness-of-fit, aka square root of
144                                the reduced chi squared.
145covMatrix      \\                (np.array) The (NxN) covVariance matrix
146=============  ===============  ====================================================
147
148.. _Phase_table:
149
150.. index::
151   single: Phase object description
152   single: Data object descriptions; Phase
153
154Phase Tree Items
155----------------
156
157Phase information is stored in the GSAS-II data tree as children of the
158Phases item in a dict with keys:
159
160.. tabularcolumns:: |l|l|p{4in}|
161
162==========  ===============     =====================================================================================================
163  key         sub-key           explanation
164==========  ===============     =====================================================================================================
165General         \\               (dict) Overall information for the phase
166  \\         3Dproj              (list of str) projections for 3D pole distribution plots
167  \\         AngleRadii          (list of floats) Default radius for each atom used to compute
168                                interatomic angles
169  \\         AtomMass            (list of floats) Masses for atoms
170  \\         AtomPtrs            (list of int) four locations (cx,ct,cs & cu) to use to pull info
171                                from the atom records
172  \\         AtomTypes           (llist of str) Atom types
173  \\         BondRadii           (list of floats) Default radius for each atom used to compute
174                                interatomic distances
175  \\         Cell                Unit cell parameters & ref. flag
176                                (list with 8 items. All but first item are float.)
177
178                                 0: cell refinement flag (True/False),
179                                 1-3: a, b, c, (:math:`\\AA`)
180                                 4-6: alpha, beta & gamma, (degrees)
181                                 7: volume (:math:`\\AA^3`)
182  \\         Color               (list of (r,b,g) triplets) Colors for atoms
183  \\         Compare             (dict) Polygon comparison parameters
184  \\         Data plot type      (str) data plot type ('Mustrain', 'Size' or
185                                'Preferred orientation') for powder data
186  \\         DisAglCtls          (dDict) with distance/angle search controls,
187                                which has keys 'Name', 'AtomTypes',
188                                'BondRadii', 'AngleRadii' which are as above
189                                except are possibly edited. Also contains
190                                'Factors', which is a 2 element list with
191                                a multiplier for bond and angle search range
192                                [typically (0.85,0.85)].
193  \\         F000X               (float) x-ray F(000) intensity
194  \\         F000N               (float) neutron F(000) intensity
195  \\         Flip                (dict) Charge flip controls
196  \\         HydIds              (dict) geometrically generated hydrogen atoms
197  \\         Isotope             (dict) Isotopes for each atom type
198  \\         Isotopes            (dict) Scattering lengths for each isotope
199                                combination for each element in phase
200  \\         MCSA controls       (dict) Monte Carlo-Simulated Annealing controls
201  \\         Map                 (dict) Map parameters
202  \\         Mass                (float) Mass of unit cell contents in g/mol
203  \\         Modulated           (bool) True if phase modulated
204  \\         Mydir               (str) Directory of current .gpx file
205  \\         Name                (str) Phase name
206  \\         NoAtoms             (dict) Number of atoms per unit cell of each type
207  \\         POhkl               (list) March-Dollase preferred orientation direction
208  \\         Pawley dmin         (float) maximum Q (as d-space) to use for Pawley extraction
209  \\         Pawley dmax         (float) minimum Q (as d-space) to use for Pawley extraction
210  \\         Pawley neg wt       (float) Restraint value for negative Pawley intensities
211  \\         SGData              (object) Space group details as a
212                                :ref:`space group (SGData) <SGData_table>`
213                                object, as defined in :func:`GSASIIspc.SpcGroup`.
214  \\         SH Texture          (dict) Spherical harmonic preferred orientation parameters
215  \\         Super               (int) dimension of super group (0,1 only)
216  \\         Type                (str) phase type (e.g. 'nuclear')
217  \\         Z                   (dict) Atomic numbers for each atom type
218  \\         doDysnomia          (bool) flag for max ent map modification via Dysnomia
219  \\         doPawley            (bool) Flag for Pawley intensity extraction
220  \\         vdWRadii            (dict) Van der Waals radii for each atom type
221ranId           \\               (int) unique random number Id for phase
222pId             \\               (int) Phase Id number for current project.
223Atoms           \\               (list of lists) Atoms in phase as a list of lists. The outer list
224                                is for each atom, the inner list contains varying
225                                items depending on the type of phase, see
226                                the :ref:`Atom Records <Atoms_table>` description.
227Drawing         \\               (dict) Display parameters
228\\           Atoms               (list of lists) with an entry for each atom that is drawn
229\\           Plane               (list) Controls for contour density plane display
230\\           Quaternion          (4 element np.array) Viewing quaternion
231\\           Zclip               (float) clipping distance in :math:`\\AA`
232\\           Zstep               (float) Step to de/increase Z-clip
233\\           atomPtrs            (list) positions of x, type, site sym, ADP flag in Draw Atoms
234\\           backColor           (list) background for plot as and R,G,B triplet
235                                (default = [0, 0, 0], black).
236\\           ballScale           (float) Radius of spheres in ball-and-stick display
237\\           bondList            (dict) Bonds
238\\           bondRadius          (float) Radius of binds in :math:`\\AA`
239\\           cameraPos           (float) Viewing position in :math:`\\AA` for plot
240\\           contourLevel        (float) map contour level in :math:`e/\\AA^3`
241\\           contourMax          (float) map contour maximum
242\\           depthFog            (bool) True if use depthFog on plot - set currently as False
243\\           ellipseProb         (float) Probability limit for display of thermal
244                                ellipsoids in % .
245\\           magMult             (float) multiplier for magnetic moment arrows
246\\           mapSize             (float) x & y dimensions of contourmap (fixed internally)
247\\           modelView           (4,4 array) from openGL drawing transofmation matrix
248\\           oldxy               (list with two floats) previous view point
249\\           radiusFactor        (float) Distance ratio for searching for bonds. Bonds
250                                are located that are within r(Ra+Rb) and (Ra+Rb)/r
251                                where Ra and Rb are the atomic radii.
252\\           selectedAtoms       (list of int values) List of selected atoms
253\\           showABC             (bool) Flag to show view point triplet. True=show.
254\\           showHydrogen        (bool) Flag to control plotting of H atoms.
255\\           showRigidBodies     (bool) Flag to highlight rigid body placement
256\\           showSlice           (bool) flag to show contour map
257\\           sizeH               (float) Size ratio for H atoms
258\\           unitCellBox         (bool) Flag to control display of the unit cell.
259\\           vdwScale            (float) Multiplier of van der Waals radius for display of vdW spheres.
260\\           viewDir             (np.array with three floats) cartesian viewing direction
261\\           viewPoint           (list of lists) First item in list is [x,y,z]
262                                in fractional coordinates for the center of
263                                the plot. Second item list of previous & current
264                                atom number viewed (may be [0,0])
265RBModels        \\               Rigid body assignments (note Rigid body definitions
266                                are stored in their own main top-level tree entry.)
267RMC             \\               (dict) RMCProfile & rmcfull controls
268Pawley ref      \\               (list) Pawley reflections
269Histograms      \\               (dict of dicts) The key for the outer dict is
270                                the histograms tied to this phase. The inner
271                                dict contains the combined phase/histogram
272                                parameters for items such as scale factors,
273                                size and strain parameters. The following are the
274                                keys to the inner dict. (dict)
275\\           Babinet             (dict) For protein crystallography. Dictionary with two
276                                entries, 'BabA', 'BabU'
277\\           Extinction          (list of float, bool) Extinction parameter
278\\           Flack               (list of [float, bool]) Flack parameter & refine flag
279\\           HStrain             (list of two lists) Hydrostatic strain. The first is
280                                a list of the HStrain parameters (1, 2, 3, 4, or 6
281                                depending on unit cell), the second is a list of boolean
282                                refinement parameters (same length)
283\\           Histogram           (str) The name of the associated histogram
284\\           Layer Disp          (list of [float, bool]) Layer displacement in beam direction & refine flag
285\\           LeBail              (bool) Flag for LeBail extraction
286\\           Mustrain            (list) Microstrain parameters, in order:
287   
288                                0. Type, one of  u'isotropic', u'uniaxial', u'generalized'
289                                1. Isotropic/uniaxial parameters - list of 3 floats
290                                2. Refinement flags - list of 3 bools
291                                3. Microstrain axis - list of 3 ints, [h, k, l]
292                                4. Generalized mustrain parameters - list of 2-6 floats, depending on space group
293                                5. Generalized refinement flags - list of bools, corresponding to the parameters of (4)
294\\           Pref.Ori.           (list) Preferred Orientation. List of eight parameters.
295                                Items marked SH are only used for Spherical Harmonics.
296                               
297                                0. (str) Type, 'MD' for March-Dollase or 'SH' for Spherical Harmonics
298                                1. (float) Value
299                                2. (bool) Refinement flag
300                                3. (list) Preferred direction, list of ints, [h, k, l]
301                                4. (int) SH - number of terms
302                                5. (dict) SH - 
303                                6. (list) SH
304                                7. (float) SH
305\\           Scale               (list of [float, bool]) Phase fraction & refine flag
306\\           Size                List of crystallite size parameters, in order:
307
308                                0. (str) Type, one of  u'isotropic', u'uniaxial', u'ellipsoidal'
309                                1. (list) Isotropic/uniaxial parameters - list of 3 floats
310                                2. (list) Refinement flags - list of 3 bools
311                                3. (list) Size axis - list of 3 ints, [h, k, l]
312                                4. (list) Ellipsoidal size parameters - list of 6 floats
313                                5. (list) Ellipsoidal refinement flags - list of bools, corresponding to the parameters of (4)
314\\           Use                 (bool) True if this histogram is to be used in refinement
315MCSA            \\               (dict) Monte-Carlo simulated annealing parameters
316==========  ===============     =====================================================================================================
317
318.. _RBData_table:
319
320.. index::
321   single: Rigid Body Data description
322   single: Data object descriptions; Rigid Body Data
323
324Rigid Body Objects
325------------------
326
327Rigid body descriptions are available for two types of rigid bodies: 'Vector'
328and 'Residue'. Vector rigid bodies are developed by a sequence of translations each
329with a refinable magnitude and Residue rigid bodies are described as Cartesian coordinates
330with defined refinable torsion angles.
331
332.. tabularcolumns:: |l|l|p{4in}|
333
334==========  ===============     ====================================================
335  key         sub-key           explanation
336==========  ===============     ====================================================
337Vector      RBId                (dict of dict) vector rigid bodies
338\\           AtInfo              (dict) Drad, Color: atom drawing radius & color for each atom type
339\\           RBname              (str) Name assigned by user to rigid body
340\\           VectMag             (list) vector magnitudes in :math:`\\AA`
341\\           rbXYZ               (list of 3 float Cartesian coordinates for Vector rigid body )
342\\           rbRef               (list of 3 int & 1 bool) 3 assigned reference atom nos. in rigid body for origin
343                                definition, use center of atoms flag
344\\           VectRef             (list of bool refinement flags for VectMag values )
345\\           rbTypes             (list of str) Atom types for each atom in rigid body
346\\           rbVect              (list of lists) Cartesian vectors for each translation used to build rigid body
347\\           useCount            (int) Number of times rigid body is used in any structure
348Residue     RBId                (dict of dict) residue rigid bodies
349\\           AtInfo              (dict) Drad, Color: atom drawing radius & color for each atom type
350\\           RBname              (str) Name assigned by user to rigid body
351\\           rbXYZ               (list of 3 float) Cartesian coordinates for Residue rigid body
352\\           rbTypes             (list of str) Atom types for each atom in rigid body
353\\           atNames             (list of str) Names of each atom in rigid body (e.g. C1,N2...)
354\\           rbRef               (list of 3 int & 1 bool) 3 assigned reference atom nos. in rigid body for origin
355                                definition, use center of atoms flag
356\\           rbSeq               (list) Orig,Piv,angle,Riding : definition of internal rigid body
357                                torsion; origin atom (int), pivot atom (int), torsion angle (float),
358                                riding atoms (list of int)
359\\           SelSeq              (int,int) used by SeqSizer to identify objects
360\\           useCount            (int)Number of times rigid body is used in any structure
361RBIds           \\               (dict) unique Ids generated upon creation of each rigid body
362\\           Vector              (list) Ids for each Vector rigid body
363\\           Residue             (list) Ids for each Residue rigid body
364==========  ===============     ====================================================
365
366.. _SGData_table:
367
368.. index::
369   single: Space Group Data description
370   single: Data object descriptions; Space Group Data
371
372Space Group Objects
373-------------------
374
375Space groups are interpreted by :func:`GSASIIspc.SpcGroup`
376and the information is placed in a SGdata object
377which is a dict with these keys. Magnetic ones are marked "mag"
378
379.. tabularcolumns:: |l|p{4.5in}|
380
381==========  ========================================================================================
382  key         explanation
383==========  ========================================================================================
384BNSlattsym  mag - (str) BNS magnetic space group symbol and centering vector
385GenFlg      mag - (list) symmetry generators indices
386GenSym      mag - (list) names for each generator
387MagMom      mag - (list) "time reversals" for each magnetic operator
388MagPtGp     mag - (str) Magnetic point group symbol
389MagSpGrp    mag - (str) Magnetic space group symbol
390OprNames    mag - (list) names for each space group operation
391SGCen       (np.array) Symmetry cell centering vectors. A (n,3) np.array
392            of centers. Will always have at least one row: ``np.array([[0, 0, 0]])``
393SGFixed     (bool) Only True if phase mported from a magnetic cif file
394            then the space group can not be changed by the user because
395            operator set from cif may be nonstandard
396SGGen       (list) generators
397SGGray      (bool) True if space group is a gray group (incommensurate magnetic structures)
398SGInv       (bool) True if centrosymmetric, False if not
399SGLatt      (str)Lattice centering type. Will be one of
400            P, A, B, C, I, F, R
401SGLaue      (str) one of the following 14 Laue classes:
402            -1, 2/m, mmm, 4/m, 4/mmm, 3R,
403            3mR, 3, 3m1, 31m, 6/m, 6/mmm, m3, m3m
404SGOps       (list) symmetry operations as a list of form
405            ``[[M1,T1], [M2,T2],...]``
406            where :math:`M_n` is a 3x3 np.array
407            and :math:`T_n` is a length 3 np.array.
408            Atom coordinates are transformed where the
409            Asymmetric unit coordinates [X is (x,y,z)]
410            are transformed using
411            :math:`X^\\prime = M_n*X+T_n`
412SGPolax     (str) Axes for space group polarity. Will be one of
413            '', 'x', 'y', 'x y', 'z', 'x z', 'y z',
414            'xyz'. In the case where axes are arbitrary
415            '111' is used (P 1, and ?).
416SGPtGrp     (str) Point group of the space group
417SGUniq      unique axis if monoclinic. Will be
418            a, b, or c for monoclinic space groups.
419            Will be blank for non-monoclinic.
420SGSpin      mag - (list) of spin flip operatiors (+1 or -1) for the space group operations
421SGSys       (str) symmetry unit cell: type one of
422            'triclinic', 'monoclinic', 'orthorhombic',
423            'tetragonal', 'rhombohedral', 'trigonal',
424            'hexagonal', 'cubic'
425SSGK1       (list) Superspace multipliers
426SpGrp       (str) space group symbol
427SpnFlp      mag - (list) Magnetic spin flips for every magnetic space group operator
428==========  ========================================================================================
429
430.. _SSGData_table:
431
432.. index::
433   single: Superspace Group Data description
434   single: Data object descriptions; Superspace Group Data
435
436Superspace groups [3+1] are interpreted by :func:`GSASIIspc.SSpcGroup`
437and the information is placed in a SSGdata object
438which is a dict with these keys:
439
440.. tabularcolumns:: |l|p{4.5in}|
441
442==========  ====================================================
443  key         explanation
444==========  ====================================================
445SSGCen      (list) 4D cell centering vectors [0,0,0,0] at least
446SSGK1       (list) Superspace multipliers
447SSGOps      (list) 4D symmetry operations as [M,T] so that M*x+T = x'
448SSpGrp      (str) superspace group symbol extension to space group
449            symbol, accidental spaces removed
450modQ        (list) modulation/propagation vector
451modSymb     (list of str) Modulation symbols
452==========  ====================================================
453
454
455Phase Information
456--------------------
457
458.. index::
459   single: Phase information record description
460
461Phase information is placed in one of the following keys:
462
463.. tabularcolumns:: |l|p{4.5in}|
464
465==========  ==============================================================
466  key         explanation
467==========  ==============================================================
468General       Overall information about a phase
469Histograms    Information about each histogram linked to the
470              current phase as well as parameters that
471              are defined for each histogram and phase
472              (such as sample peak widths and preferred
473              orientation parameters.
474Atoms         Contains a list of atoms, as described in the
475              :ref:`Atom Records <Atoms_table>` description.
476Drawing       Parameters that determine how the phase is
477              displayed, including a list of atoms to be
478              included, as described in the
479              :ref:`Drawing Atom Records <Drawing_atoms_table>`
480              description
481MCSA          Monte-Carlo simulated annealing parameters
482pId           The index of each phase in the project, numbered
483              starting at 0
484ranId         An int value with a unique value for each phase
485RBModels      A list of dicts with parameters for each
486              rigid body inserted into the current phase,
487              as defined in the
488              :ref:`Rigid Body Insertions <Rigid_Body_Insertions>`.
489              Note that the rigid bodies are defined as
490              :ref:`Rigid Body Objects <RBData_table>`
491RMC           PDF modeling parameters
492Pawley ref    Pawley refinement parameters
493
494==========  ==============================================================
495
496.. _Atoms_table:
497
498.. index::
499   single: Atoms record description
500   single: Data object descriptions; Atoms record
501
502--------------------
503Atom Records
504--------------------
505
506If ``phasedict`` points to the phase information in the data tree, then
507atoms are contained in a list of atom records (list) in
508``phasedict['Atoms']``. Also needed to read atom information
509are four pointers, ``cx,ct,cs,cia = phasedict['General']['AtomPtrs']``,
510which define locations in the atom record, as shown below. Items shown are
511always present; additional ones for macromolecular phases are marked 'mm',
512and those for magnetic structures are marked 'mg'
513
514.. tabularcolumns:: |l|p{4.5in}|
515
516==============      ====================================================
517location            explanation
518==============      ====================================================
519ct-4                mm - (str) residue number
520ct-3                mm - (str) residue name (e.g. ALA)
521ct-2                mm - (str) chain label
522ct-1                (str) atom label
523ct                  (str) atom type
524ct+1                (str) refinement flags; combination of 'F', 'X', 'U', 'M'
525cx,cx+1,cx+2        (3 floats) the x,y and z coordinates
526cx+3                (float) site occupancy
527cx+4,cx+5,cx+6      mg - (list) atom magnetic moment along a,b,c in Bohr magnetons
528cs                  (str) site symmetry
529cs+1                (int) site multiplicity
530cia                 (str) ADP flag: Isotropic ('I') or Anisotropic ('A')
531cia+1               (float) Uiso
532cia+2...cia+7       (6 floats) U11, U22, U33, U12, U13, U23
533atom[cia+8]         (int) unique atom identifier
534
535==============      ====================================================
536
537.. _Drawing_atoms_table:
538
539.. index::
540   single: Drawing atoms record description
541   single: Data object descriptions; Drawing atoms record
542
543----------------------------
544Drawing Atom Records
545----------------------------
546
547If ``phasedict`` points to the phase information in the data tree, then
548drawing atoms are contained in a list of drawing atom records (list) in
549``phasedict['Drawing']['Atoms']``. Also needed to read atom information
550are four pointers, ``cx,ct,cs,ci = phasedict['Drawing']['AtomPtrs']``,
551which define locations in the atom record, as shown below. Items shown are
552always present; additional ones for macromolecular phases are marked 'mm',
553and those for magnetic structures are marked 'mg'
554
555.. tabularcolumns:: |l|p{4.5in}|
556
557==============   ===================================================================================
558location            explanation
559==============   ===================================================================================
560ct-4                mm - (str) residue number
561ct-3                mm - (str) residue name (e.g. ALA)
562ct-2                mm - (str) chain label
563ct-1                (str) atom label
564ct                  (str) atom type
565cx,cx+1,cx+2        (3 floats) the x,y and z coordinates
566cx+3,cx+4,cx+5      mg - (3 floats) atom magnetic moment along a,b,c in Bohr magnetons
567cs-1                (str) Sym Op symbol; sym. op number + unit cell id (e.g. '1,0,-1')
568cs                  (str) atom drawing style; e.g. 'balls & sticks'
569cs+1                (str) atom label style (e.g. 'name')
570cs+2                (int) atom color (RBG triplet)
571cs+3                (str) ADP flag: Isotropic ('I') or Anisotropic ('A')
572cs+4                (float) Uiso
573cs+5...cs+11        (6 floats) U11, U22, U33, U12, U13, U23
574ci                  (int) unique atom identifier; matches source atom Id in Atom Records
575==============   ===================================================================================
576
577.. _Rigid_Body_Insertions:
578
579----------------------------
580Rigid Body Insertions
581----------------------------
582
583If ``phasedict`` points to the phase information in the data tree, then
584rigid body information is contained in list(s) in
585``phasedict['RBModels']['Residue']`` and/or ``phasedict['RBModels']['Vector']``
586for each rigid body inserted into the current phase.
587
588.. tabularcolumns:: |l|p{4.5in}|
589
590==============   ===================================================================================
591key              explanation
592==============   ===================================================================================
593fixOrig           Should the origin be fixed (when editing, not the refinement flag)
594Ids               Ids for assignment of atoms in the rigid body
595numChain          Chain number for macromolecular fits
596Orient            Orientation of the RB as a quaternion and a refinement flag (' ', 'A' or 'AV')
597OrientVec         Orientation of the RB expressed as a vector and azimuthal rotation angle
598Orig              Origin of the RB in fractional coordinates and refinement flag (bool)
599RBId              References the unique ID of a rigid body in the
600                  :ref:`Rigid Body Objects <RBData_table>`
601RBname            The name for the rigid body (str)
602AtomFrac          The atom fractions for the rigid body
603ThermalMotion     The thermal motion description for the rigid body, which includes a choice for
604                  the model and can include TLS parameters or an overall Uiso value.
605Torsions          Defines the torsion angle and refinement flag for each torsion defined in
606                  the :ref:`Rigid Body Object <RBData_table>`
607==============   ===================================================================================
608
609.. _Powder_table:
610
611.. index::
612   single: Powder data object description
613   single: Data object descriptions; Powder Data
614
615Powder Diffraction Tree Items
616-----------------------------
617
618Every powder diffraction histogram is stored in the GSAS-II data tree
619with a top-level entry named beginning with the string "PWDR ". The
620diffraction data for that information are directly associated with
621that tree item and there are a series of children to that item. The
622routines :func:`GSASIIdataGUI.GSASII.GetUsedHistogramsAndPhasesfromTree`
623and :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` will
624load this information into a dictionary where the child tree name is
625used as a key, and the information in the main entry is assigned
626a key of ``Data``, as outlined below.
627
628.. tabularcolumns:: |p{1in}|p{1in}|p{4in}|
629
630======================     ===============  ===========================================================
631  key                       sub-key          explanation
632======================     ===============  ===========================================================
633Comments                    \\               (list of str) Text strings extracted from the original powder
634                                            data header. These cannot be changed by the user;
635                                            it may be empty.
636Limits                      \\               (list) two two element lists, as [[Ld,Hd],[L,H]]
637                                            where L and Ld are the current and default lowest
638                                            two-theta value to be used and
639                                            where H and Hd are the current and default highest
640                                            two-theta value to be used.
641Reflection Lists            \\               (dict of dicts) with an entry for each phase in the
642                                            histogram. The contents of each dict item
643                                            is a dict containing reflections, as described in
644                                            the :ref:`Powder Reflections <PowderRefl_table>`
645                                            description.
646Instrument Parameters       \\               (dict) The instrument parameters uses different dicts
647                                            for the constant wavelength (CW) and time-of-flight (TOF)
648                                            cases. See below for the descriptions of each.
649wtFactor                    \\               (float) A weighting factor to increase or decrease
650                                            the leverage of data in the histogram .
651                                            A value of 1.0 weights the data with their
652                                            standard uncertainties and a larger value
653                                            increases the weighting of the data (equivalent
654                                            to decreasing the uncertainties).
655Sample Parameters           \\               (dict) Parameters that describe how
656                                            the data were collected, as listed
657                                            below. Refinable parameters are a list containing
658                                            a float and a bool, where the second value
659                                            specifies if the value is refined, otherwise
660                                            the value is a float unless otherwise noted.
661\\                           Scale           The histogram scale factor (refinable)
662\\                           Absorption      The sample absorption coefficient as
663                                            :math:`\\mu r` where r is the radius
664                                            (refinable). Only valid for Debye-Scherrer geometry.
665\\                           SurfaceRoughA   Surface roughness parameter A as defined by
666                                            Surotti, *J. Appl. Cryst*, **5**, 325-331, 1972.
667                                            (refinable - only valid for Bragg-Brentano geometry)
668\\                           SurfaceRoughB   Surface roughness parameter B (refinable -
669                                            only valid for Bragg-Brentano geometry)
670\\                           DisplaceX,      Sample displacement from goniometer center
671                            DisplaceY       where Y is along the beam direction and
672                                            X is perpendicular. Units are :math:`\\mu m`
673                                            (refinable).
674\\                           Phi, Chi,       Goniometer sample setting angles, in degrees.
675                            Omega
676\\                           Gonio. radius   Radius of the diffractometer in mm
677\\                           InstrName       (str) A name for the instrument, used in preparing
678                                            a CIF .
679\\                           Force,          Variables that describe how the measurement
680                            Temperature,    was performed. Not used directly in
681                            Humidity,       any computations.
682                            Pressure,
683                            Voltage
684\\                           ranId           (int) The random-number Id for the histogram
685                                            (same value as where top-level key is ranId)
686\\                           Type            (str) Type of diffraction data, may be 'Debye-Scherrer'
687                                            or 'Bragg-Brentano' .
688hId                         \\               (int) The number assigned to the histogram when
689                                            the project is loaded or edited (can change)
690ranId                       \\               (int) A random number id for the histogram
691                                            that does not change
692Background                  \\               (list) The background is stored as a list with where
693                                            the first item in the list is list and the second
694                                            item is a dict. The list contains the background
695                                            function and its coefficients; the dict contains
696                                            Debye diffuse terms and background peaks.
697                                            (TODO: this needs to be expanded.)
698Data                        \\               (list) The data consist of a list of 6 np.arrays
699                                            containing in order:
700
701                                            0. the x-postions (two-theta in degrees),
702                                            1. the intensity values (Yobs),
703                                            2. the weights for each Yobs value
704                                            3. the computed intensity values (Ycalc)
705                                            4. the background values
706                                            5. Yobs-Ycalc
707======================     ===============  ===========================================================
708
709.. _CWPowder_table:
710
711.. index::
712   single: Powder data CW Instrument Parameters
713
714-----------------------------
715CW Instrument Parameters
716-----------------------------
717
718Instrument Parameters are placed in a list of two dicts,
719where the keys in the first dict are listed below. Note that the dict contents are different for
720constant wavelength (CW) vs. time-of-flight (TOF) histograms.
721The value for each item is a list containing three values: the initial value, the current value
722and a refinement flag which can have a value of True, False or 0 where 0 indicates a value that
723cannot be refined. The first and second values are floats unless otherwise noted.
724Items not refined are noted as [*]
725
726.. tabularcolumns:: |l|p{1in}|p{4in}|
727
728========================    ===============  ===========================================================
729  key                       sub-key           explanation
730========================    ===============  ===========================================================
731Instrument Parameters[0]    Type [*]            (str) Histogram type:
732                                                * 'PXC' for constant wavelength x-ray
733                                                * 'PNC' for constant wavelength neutron
734\\                           Bank [*]            (int) Data set number in a multidata file (usually 1)
735\\                           Lam                 (float) Specifies a wavelength in :math:`\\AA`
736\\                           Lam1 [*]            (float) Specifies the primary wavelength in
737                                                :math:`\\AA`, used in place of Lam
738                                                when an :math:`\\alpha_1, \\alpha_2`
739                                                source is used.
740\\                           Lam2 [*]            (float) Specifies the secondary wavelength in
741                                                :math:`\\AA`, used with Lam1
742\\                           I(L2)/I(L1)         (float) Ratio of Lam2 to Lam1, used with Lam1
743\\                           Zero                (float) Two-theta zero correction in *degrees*
744\\                           Azimuth [*]         (float) Azimuthal setting angle for data recorded with differing setting angles
745\\                           U, V, W             (float) Cagliotti profile coefficients
746                                                for Gaussian instrumental broadening, where the
747                                                FWHM goes as
748                                                :math:`U \\tan^2\\theta + V \\tan\\theta + W`
749\\                           X, Y, Z             (float) Cauchy (Lorentzian) instrumental broadening coefficients
750\\                           SH/L                (float) Variant of the Finger-Cox-Jephcoat asymmetric
751                                                peak broadening ratio. Note that this is the
752                                                sum of S/L and H/L where S is
753                                                sample height, H is the slit height and
754                                                L is the goniometer diameter.
755\\                           Polariz.            (float) Polarization coefficient.
756Instrument Parameters[1]                        (empty dict)
757========================    ===============  ===========================================================
758
759.. _TOFPowder_table:
760
761.. index::
762   single: Powder data TOF Instrument Parameters
763
764-----------------------------
765TOF Instrument Parameters
766-----------------------------
767
768Instrument Parameters are also placed in a list of two dicts,
769where the keys in each dict listed below, but here for
770time-of-flight (TOF) histograms.
771The value for each item is a list containing three values: the initial value, the current value
772and a refinement flag which can have a value of True, False or 0 where 0 indicates a value that
773cannot be refined. The first and second values are floats unless otherwise noted.
774Items not refined are noted as [*]
775
776.. tabularcolumns:: |l|p{1.5in}|p{4in}|
777
778========================    ===============  ===========================================================
779  key                        sub-key          explanation
780========================    ===============  ===========================================================
781Instrument Parameters[0]    Type [*]            (str) Histogram type:
782                                                * 'PNT' for time of flight neutron
783\\                           Bank                (int) Data set number in a multidata file
784\\                           2-theta [*]         (float) Nominal scattering angle for the detector
785\\                           fltPath [*]         (float) Total flight path source-sample-detector
786\\                           Azimuth [*]         (float) Azimuth angle for detector right hand rotation
787                                                from horizontal away from source
788\\                           difC,difA,          (float) Diffractometer constants for conversion of d-spacing to TOF
789                            difB                in microseconds
790\\                           Zero                (float) Zero point offset (microseconds)
791\\                           alpha               (float) Exponential rise profile coefficients
792\\                           beta-0              (float) Exponential decay profile coefficients
793                            beta-1
794                            beta-q
795\\                           sig-0               (float) Gaussian profile coefficients
796                            sig-1
797                            sig-2
798                            sig-q   
799\\                           X,Y,Z               (float) Lorentzian profile coefficients
800Instrument Parameters[1]    Pdabc               (list of 4 float lists) Originally created for use in gsas as optional tables
801                                                of d, alp, bet, d-true; for a reflection alpha & beta are obtained via interpolation
802                                                from the d-spacing and these tables. The d-true column is apparently unused.
803========================    ===============  ===========================================================
804
805
806.. _PowderRefl_table:
807
808.. index::
809   single: Powder reflection object description
810   single: Data object descriptions; Powder Reflections
811
812Powder Reflection Data Structure
813--------------------------------
814
815For every phase in a histogram, the ``Reflection Lists`` value is a dict
816one element of which is `'RefList'`, which is a np.array containing
817reflections. The columns in that array are documented below.
818
819==========  ====================================================
820  index         explanation
821==========  ====================================================
822 0,1,2          h,k,l (float)
823 3              (int) multiplicity
824 4              (float) d-space, :math:`\\AA`
825 5              (float) pos, two-theta
826 6              (float) sig, Gaussian width
827 7              (float) gam, Lorenzian width
828 8              (float) :math:`F_{obs}^2`
829 9              (float) :math:`F_{calc}^2`
830 10             (float) reflection phase, in degrees
831 11             (float) intensity correction for reflection, this times
832                :math:`F_{obs}^2` or :math:`F_{calc}^2` gives Iobs or Icalc
833 12             (float) Preferred orientation correction
834 13             (float) Transmission (absorption correction)
835 14             (float) Extinction correction
836==========  ====================================================
837
838.. _Xtal_table:
839
840.. index::
841   single: Single Crystal data object description
842   single: Data object descriptions; Single crystal data
843
844Single Crystal Tree Items
845-------------------------
846
847Every single crystal diffraction histogram is stored in the GSAS-II data tree
848with a top-level entry named beginning with the string "HKLF ". The
849diffraction data for that information are directly associated with
850that tree item and there are a series of children to that item. The
851routines :func:`GSASIIdataGUI.GSASII.GetUsedHistogramsAndPhasesfromTree`
852and :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` will
853load this information into a dictionary where the child tree name is
854used as a key, and the information in the main entry is assigned
855a key of ``Data``, as outlined below.
856
857.. tabularcolumns:: |l|l|p{4in}|
858
859======================  ===============     ====================================================
860  key                      sub-key          explanation
861======================  ===============     ====================================================
862Data                        \\               (dict) that contains the
863                                            reflection table,
864                                            as described in the
865                                            :ref:`Single Crystal Reflections
866                                            <XtalRefl_table>`
867                                            description.
868
869Instrument Parameters       \\               (list) containing two dicts where the possible
870                                            keys in each dict are listed below. The value
871                                            for most items is a list containing two values:
872                                            the initial value, the current value.
873                                            The first and second
874                                            values are floats unless otherwise noted.
875\\                           Lam             (two floats) Specifies a wavelength in :math:`\\AA`
876\\                           Type            (two str values) Histogram type :
877                                            * 'SXC' for constant wavelength x-ray
878                                            * 'SNC' for constant wavelength neutron
879                                            * 'SNT' for time of flight neutron
880\\                           InstrName       (str) A name for the instrument, used in preparing a CIF
881wtFactor                    \\               (float) A weighting factor to increase or decrease
882                                            the leverage of data in the histogram.
883                                            A value of 1.0 weights the data with their
884                                            standard uncertainties and a larger value
885                                            increases the weighting of the data (equivalent
886                                            to decreasing the uncertainties).
887
888hId                         \\               (int) The number assigned to the histogram when
889                                            the project is loaded or edited (can change)
890ranId                       \\               (int) A random number id for the histogram
891                                            that does not change
892======================  ===============     ====================================================
893
894.. _XtalRefl_table:
895
896.. index::
897   single: Single Crystal reflection object description
898   single: Data object descriptions; Single Crystal Reflections
899
900Single Crystal Reflection Data Structure
901----------------------------------------
902
903For every single crystal a histogram, the ``'Data'`` item contains
904the structure factors as an np.array in item `'RefList'`.
905The columns in that array are documented below.
906
907.. tabularcolumns:: |l|p{4in}|
908
909==========  ====================================================
910  index         explanation
911==========  ====================================================
912 0,1,2          (float) h,k,l
913 3              (int) multiplicity
914 4              (float) d-space, :math:`\\AA`
915 5              (float) :math:`F_{obs}^2`
916 6              (float) :math:`\\sigma(F_{obs}^2)`
917 7              (float) :math:`F_{calc}^2`
918 8              (float) :math:`F_{obs}^2T`
919 9              (float) :math:`F_{calc}^2T`
920 10             (float) reflection phase, in degrees
921 11             (float) intensity correction for reflection, this times
922                :math:`F_{obs}^2` or :math:`F_{calc}^2`
923                gives Iobs or Icalc
924==========  ====================================================
925
926.. _Image_table:
927
928.. index::
929   image: Image data object description
930   image: Image object descriptions
931
932Image Data Structure
933--------------------
934
935Every 2-dimensional image is stored in the GSAS-II data tree
936with a top-level entry named beginning with the string "IMG ". The
937image data are directly associated with that tree item and there
938are a series of children to that item. The routines :func:`GSASIIdataGUI.GSASII.GetUsedHistogramsAndPhasesfromTree`
939and :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` will
940load this information into a dictionary where the child tree name is
941used as a key, and the information in the main entry is assigned
942a key of ``Data``, as outlined below.
943
944.. tabularcolumns:: |l|l|p{4in}|
945
946======================  ======================  ====================================================
947  key                      sub-key              explanation
948======================  ======================  ====================================================
949Comments                    \\                   (list of str) Text strings extracted from the original image data
950                                                header or a metafile. These cannot be changed by
951                                                the user; it may be empty.
952Image Controls              azmthOff            (float) The offset to be applied to an azimuthal
953                                                value. Accomodates
954                                                detector orientations other than with the detector
955                                                X-axis
956                                                horizontal.
957\\                           background image    (list:str,float) The name of a tree item ("IMG ...") that is to be subtracted
958                                                during image integration multiplied by value. It must have the same size/shape as
959                                                the integrated image. NB: value < 0 for subtraction.
960\\                           calibrant           (str) The material used for determining the position/orientation
961                                                of the image. The data is obtained from :func:`ImageCalibrants`
962                                                and UserCalibrants.py (supplied by user).
963\\                           calibdmin           (float) The minimum d-spacing used during the last calibration run.
964\\                           calibskip           (int) The number of expected diffraction lines skipped during the last
965                                                calibration run.
966\\                           center              (list:floats) The [X,Y] point in detector coordinates (mm) where the direct beam
967                                                strikes the detector plane as determined by calibration. This point
968                                                does not have to be within the limits of the detector boundaries.
969\\                           centerAzm           (bool) If True then the azimuth reported for the integrated slice
970                                                of the image is at the center line otherwise it is at the leading edge.
971\\                           color               (str) The name of the colormap used to display the image. Default = 'Paired'.
972\\                           cutoff              (float) The minimum value of I/Ib for a point selected in a diffraction ring for
973                                                calibration calculations. See pixLimit for details as how point is found.
974\\                           DetDepth            (float) Coefficient for penetration correction to distance; accounts for diffraction
975                                                ring offset at higher angles. Optionally determined by calibration.
976\\                           DetDepthRef         (bool) If True then refine DetDepth during calibration/recalibration calculation.
977\\                           distance            (float) The distance (mm) from sample to detector plane.
978\\                           ellipses            (list:lists) Each object in ellipses is a list [center,phi,radii,color] where
979                                                center (list) is location (mm) of the ellipse center on the detector plane, phi is the
980                                                rotation of the ellipse minor axis from the x-axis, and radii are the minor & major
981                                                radii of the ellipse. If radii[0] is negative then parameters describe a hyperbola. Color
982                                                is the selected drawing color (one of 'b', 'g' ,'r') for the ellipse/hyperbola.
983\\                           edgemin             (float) Not used;  parameter in EdgeFinder code.
984\\                           fullIntegrate       (bool) If True then integrate over full 360 deg azimuthal range.
985\\                           GonioAngles         (list:floats) The 'Omega','Chi','Phi' goniometer angles used for this image.
986                                                Required for texture calculations.
987\\                           invert_x            (bool) If True display the image with the x-axis inverted.
988\\                           invert_y            (bool) If True display the image with the y-axis inverted.
989\\                           IOtth               (list:floats) The minimum and maximum 2-theta values to be used for integration.
990\\                           LRazimuth           (list:floats) The minimum and maximum azimuth values to be used for integration.
991\\                           Oblique             (list:float,bool) If True apply a detector absorption correction using the value to the
992                                                intensities obtained during integration.
993\\                           outAzimuths         (int) The number of azimuth pie slices.
994\\                           outChannels         (int) The number of 2-theta steps.
995\\                           pixelSize           (list:ints) The X,Y dimensions (microns) of each pixel.
996\\                           pixLimit            (int) A box in the image with 2*pixLimit+1 edges is searched to find the maximum.
997                                                This value (I) along with the minimum (Ib) in the box is reported by :func:`GSASIIimage.ImageLocalMax`
998                                                and subject to cutoff in :func:`GSASIIimage.makeRing`.
999                                                Locations are used to construct rings of points for calibration calcualtions.
1000\\                           PolaVal             (list:float,bool) If type='SASD' and if True, apply polarization correction to intensities from
1001                                                integration using value.
1002\\                           rings               (list:lists) Each entry is [X,Y,dsp] where X & Y are lists of x,y coordinates around a
1003                                                diffraction ring with the same d-spacing (dsp)
1004\\                           ring                (list) The x,y coordinates of the >5 points on an inner ring
1005                                                selected by the user,
1006\\                           Range               (list) The minimum & maximum values of the image
1007\\                           rotation            (float) The angle between the x-axis and the vector about which the
1008                                                detector is tilted. Constrained to -180 to 180 deg.
1009\\                           SampleShape         (str) Currently only 'Cylinder'. Sample shape for Debye-Scherrer experiments; used for absorption
1010                                                calculations.
1011\\                           SampleAbs           (list: float,bool) Value of absorption coefficient for Debye-Scherrer experimnents, flag if True
1012                                                to cause correction to be applied.
1013\\                           setDefault          (bool) If True the use the image controls values for all new images to be read. (might be removed)
1014\\                           setRings            (bool) If True then display all the selected x,y ring positions (vida supra rings) used in the calibration.
1015\\                           showLines           (bool) If True then isplay the integration limits to be used.
1016\\                           size                (list:int) The number of pixels on the image x & y axes
1017\\                           type                (str) One of 'PWDR', 'SASD' or 'REFL' for powder, small angle or reflectometry data, respectively.
1018\\                           tilt                (float) The angle the detector normal makes with the incident beam; range -90 to 90.
1019\\                           wavelength          (float) The radiation wavelength (:math:`\\AA`) as entered by the user
1020                                                (or someday obtained from the image header).
1021Masks                       Arcs                (list: lists) Each entry [2-theta,[azimuth[0],azimuth[1]],thickness] describes an arc mask
1022                                                to be excluded from integration
1023\\                           Frames              (list:lists) Each entry describes the x,y points (3 or more - mm) that describe a frame outside
1024                                                of which is excluded from recalibration and integration. Only one frame is allowed.
1025\\                           Points              (list:lists) Each entry [x,y,radius] (mm) describes an excluded spot on the image to be excluded
1026                                                from integration.
1027\\                           Polygons            (list:lists) Each entry is a list of 3+ [x,y] points (mm) that describe a polygon on the image
1028                                                to be excluded from integration.
1029\\                           Rings               (list: lists) Each entry [2-theta,thickness] describes a ring mask
1030                                                to be excluded from integration.
1031\\                           Thresholds          (list:[tuple,list]) [(Imin,Imax),[Imin,Imax]] This gives lower and upper limits for points on the image to be included
1032                                                in integrsation. The tuple is the image intensity limits and the list are those set by the user.
1033\\                           SpotMask            (dict: int & array)
1034                                                'esdMul'(int) number of standard deviations above mean ring intensity to mask
1035                                                'spotMask' (bool array) the spot mask for every pixel in image         
1036
1037Stress/Strain               Sample phi          (float) Sample rotation about vertical axis.
1038\\                           Sample z            (float) Sample translation from the calibration sample position (for Sample phi = 0)
1039                                                These will be restricted by space group symmetry; result of strain fit refinement.
1040\\                           Type                (str) 'True' or 'Conventional': The strain model used for the calculation.
1041\\                           d-zero              (list:dict) Each item is for a diffraction ring on the image; all items are from the same phase
1042                                                and are used to determine the strain tensor.
1043                                                The dictionary items are:
1044                                                'Dset': (float) True d-spacing for the diffraction ring; entered by the user.
1045                                                'Dcalc': (float) Average calculated d-spacing determined from strain coeff.
1046                                                'Emat': (list: float) The strain tensor elements e11, e12 & e22 (e21=e12, rest are 0)
1047                                                'Esig': (list: float) Esds for Emat from fitting.
1048                                                'pixLimit': (int) Search range to find highest point on ring for each data point
1049                                                'cutoff': (float) I/Ib cutoff for searching.
1050                                                'ImxyObs': (list: lists) [[X],[Y]] observed points to be used for strain calculations.
1051                                                'ImtaObs': (list: lists) [[d],[azm]] transformed via detector calibration from ImxyObs.
1052                                                'ImtaCalc': (list: lists [[d],[azm]] calculated d-spacing & azimuth from fit.
1053
1054======================  ======================  ====================================================
1055
1056.. _parmDict_table:
1057
1058.. index::
1059   single: Parameter dictionary
1060
1061Parameter Dictionary
1062-------------------------
1063
1064The parameter dictionary contains all of the variable parameters for the refinement.
1065The dictionary keys are the name of the parameter (<phase>:<hist>:<name>:<atom>).
1066It is prepared in two ways. When loaded from the tree
1067(in :meth:`GSASIIdataGUI.GSASII.MakeLSParmDict` and
1068:meth:`GSASIIIO.ExportBaseclass.loadParmDict`),
1069the values are lists with two elements: ``[value, refine flag]``
1070
1071When loaded from the GPX file (in
1072:func:`GSASIIstrMain.Refine` and :func:`GSASIIstrMain.SeqRefine`), the value in the
1073dict is the actual parameter value (usually a float, but sometimes a
1074letter or string flag value (such as I or A for iso/anisotropic).
1075
1076Texture implementation
1077------------------------------
1078
1079There are two different places where texture can be treated in GSAS-II.
1080One is for mitigating the effects of texture in a structural refinement.
1081The other is for texture characterization.
1082
1083For reducing the effect of texture in a structural refinement
1084there are entries labeled preferred orientation in each phase's
1085data tab. Two different approaches can be used for this, the March-Dollase
1086model and spherical harmonics.
1087
1088For the March-Dollase model, one axis in reciprocal space is designated as
1089unique (defaulting to the 001 axis) and reflections are corrected
1090according to the angle they make with this axis depending on
1091the March-Dollase ratio. (If unity, no correction is made).
1092The ratio can be greater than one or less than one depending on if
1093crystallites oriented along the designated axis are
1094overrepresented or underrepresented. For most crystal systems there is an
1095obvious choice for the direction of the unique axis and then only a single
1096term needs to be refined. If the number is close to 1, then the correction
1097is not needed.
1098
1099The second method for reducing the effect of texture in a structural
1100refinement is to create a crystallite orientation probability surface as an
1101expansion in terms spherical harmonic functions. Only functions consistent with
1102cylindrical diffraction suymmetry and having texture symmetry
1103consistent with the Laue class of phase are used and are allowed,
1104so the higher the symmetry the fewer terms that are available for a given spherical harmonics order.
1105To use this correction, select the lowest order that provides
1106refinable terms and perform a refinement. If the texture index remains close to
1107one, then the correction is not needed. If a significant improvement is
1108noted in the profile Rwp, one may wish to see if a higher order expansion
1109gives an even larger improvement.
1110
1111To characterize texture in a material, generally one needs data collected with the
1112sample at multiple orientations or, for TOF, with detectors at multiple
1113locations around the sample. In this case the detector orientation is given in
1114each histogram's Sample Parameters and the sample's orientation is described
1115with the Euler angles specifed on the phase's Texture tab, which is also
1116where the texture type (cylindrical, rolling,...) and the spherical
1117harmonic order is selected. This should not be used with a single dataset and
1118should not be used if the preferred orientations corrections are used.
1119
1120The coordinate system used for texture characterization is defined where
1121the sample coordinates (Psi, gamma) are defined with an instrument coordinate
1122system (I, J, K) such that K is normal to the diffraction plane and J is coincident with the
1123direction of the incident radiation beam toward the source. We further define
1124a standard set of right-handed goniometer eulerian angles (Omega, Chi, Phi) so that Omega and Phi are
1125rotations about K and Chi is a rotation about J when Omega = 0. Finally, as the sample
1126may be mounted so that the sample coordinate system (Is, Js, Ks) does not coincide with
1127the instrument coordinate system (I, J, K), we define three eulerian sample rotation angles
1128(Omega-s, Chi-s, Phi-s) that describe the rotation from (Is, Js, Ks) to (I, J, K). The sample rotation
1129angles are defined so that with the goniometer angles at zero Omega-s and Phi-s are rotations
1130about K and Chi-s is a rotation about J.
1131
1132Three typical examples:
1133    1) Bragg-Brentano laboratory diffractometer: Chi=0
1134    2) Debye-Scherrer counter detector; sample capillary axis perpendicular to diffraction plane: Chi=90
1135    3) Debye-Scherrer 2D area detector positioned directly behind sample; sample capillary axis horizontal; Chi=0
1136            NB: The area detector azimuthal angle = 0 in horizontal plane to right as viewed from x-ray source & 90 at vertical "up" direction
1137           
1138ISODISTORT implementation
1139------------------------------
1140
1141CIFs prepared with the ISODISTORT web site
1142https://stokes.byu.edu/iso/isodistort_version5.6.1/isodistort.php
1143[B. J. Campbell, H. T. Stokes, D. E. Tanner, and D. M. Hatch, "ISODISPLACE: An Internet Tool for Exploring Structural Distortions." J. Appl. Cryst. 39, 607-614 (2006).]
1144can be read into GSAS-II using import CIF. This will cause constraints to be established for structural distortion modes read from the CIF.
1145At present, of the five types of modes  only displacive(``_iso_displacivemode``...) and occupancy (``_iso_occupancymode``...) are processed. Not yet processed: ``_iso_magneticmode``..., ``_iso_rotationalmode``... & ``_iso_strainmode``...
1146
1147The CIF importer :mod:`G2phase_CIF` implements class :class:`G2phase_CIF.CIFPhaseReader` which offers two methods associated with ISODISTORT (ID) input. Method :meth:`G2phase_CIF.CIFPhaseReader.ISODISTORT_test`
1148checks to see if a CIF block contains the loops with ``_iso_displacivemode_label`` or  ``_iso_occupancymode_label`` items.
1149If so, method :meth:`G2phase_CIF.CIFPhaseReader.ISODISTORT_proc` is called to read and interpret them.
1150The results are placed into the reader object's ``.Phase`` class variable as a dict item with key ``'ISODISTORT'``.
1151
1152Note that each mode ID has a long label with a name such as  Pm-3m[1/2,1/2,1/2]R5+(a,a,0)[La:b:dsp]T1u(a). Function :func:`G2phase_CIF.ISODISTORT_shortLbl` is used to create a short name for this,
1153such as R5_T1u(a) which is made unique by addition of _n if the short name is duplicated.
1154As each mode is processed, a constraint corresponding to that mode is
1155created and is added to list in the reader object's ``.Constraints`` class variable. Items placed into that list can either be a list, which corresponds to a function (new var) type :ref:`constraint definition <Constraints_table>` entry, or an item can be a dict, which provides help information for each constraint.
1156
1157------------------------------
1158Displacive modes
1159------------------------------
1160
1161The coordinate variables, as named by ISODISTORT, are placed in
1162``.Phase['ISODISTORT']['IsoVarList']`` and the corresponding :class:`GSASIIobj.G2VarObj` objects for each are placed in ``.Phase['ISODISTORT']['G2VarList']``. The mode variables,
1163as named by ISODISTORT, are placed in ``.Phase['ISODISTORT']['IsoModeList']`` and the corresponding :class:`GSASIIobj.G2VarObj` objects for each are placed in ``.Phase['ISODISTORT']['G2ModeList']``.
1164[Use ``str(G2VarObj)`` to get the variable name from the G2VarObj object, but note that the phase number, *n*, for the prefix "*n*::" cannot be determined as the phase number is not yet assigned.]
1165
1166Displacive modes are a bit complex in that they relate to delta displacements, relative to an offset value for each coordinate, and because the modes are normalized, while GSAS-II also uses displacements, but these are added to the coordinates after each refinement cycle and the delta values are set to zero. ISODISTORT uses fixed offsets (subtracted from the actual position to obtain the delta values) that are taken from ``_iso_coordinate_formula`` and these are placed in ``.Phase['ISODISTORT']['ParentStructure]`` (keyed by atom label). The normalization factors (which the delta values are divided by) are taken from ``_iso_displacivemodenorm_value`` and are placed in ``.Phase['ISODISTORT']['NormList']`` in the same order as as ``...['IsoModeList']`` and
1167``...['G2ModeList']``.
1168
1169The CIF contains a sparse matrix, from the ``loop_`` containing
1170``_iso_displacivemodematrix_value`` which provides the equations for determining the mode values from the coordinates, that matrix is placed in ``.Phase['ISODISTORT']['Mode2VarMatrix']``. The matrix is inverted to produce
1171``.Phase['ISODISTORT']['Var2ModeMatrix']``, which determines how to compute the
1172mode values from the delta coordinate values. These values are used for the
1173in :func:`GSASIIconstrGUI.ShowIsoDistortCalc` which shows coordinate and mode
1174values, the latter with s.u. values.
1175
1176
1177------------------------------
1178Occupancy modes
1179------------------------------
1180
1181
1182The delta occupancy variables, as named by ISODISTORT, are placed in
1183``.Phase['ISODISTORT']['OccVarList']`` and the corresponding :class:`GSASIIobj.G2VarObj` objects for each are placed in ``.Phase['ISODISTORT']['G2OccVarList']``. The mode variables, as named by ISODISTORT, are placed in ``.Phase['ISODISTORT']['OccModeList']`` and the corresponding :class:`GSASIIobj.G2VarObj` objects for each are placed in ``.Phase['ISODISTORT']['G2OccModeList']``.
1184
1185Occupancy modes, like Displacive modes, are also refined as delta values.  However, GSAS-II directly refines the fractional occupancies.
1186Offset values for each atom, are taken from ``_iso_occupancy_formula`` and are placed in
1187``.Phase['ISODISTORT']['ParentOcc]``.
1188(Offset values are subtracted from the actual position to obtain the delta values.)
1189Modes are normalized (where the mode values are divided by the normalization factor)
1190are taken from ``_iso_occupancymodenorm_value`` and are placed in ``.Phase['ISODISTORT']['OccNormList']`` in the same order as as ``...['OccModeList']`` and
1191``...['G2OccModeList']``.
1192
1193The CIF contains a sparse matrix, from the ``loop_`` containing
1194``_iso_occupancymodematrix_value``, which provides the equations for determining the mode values from the coordinates. That matrix is placed in ``.Phase['ISODISTORT']['Occ2VarMatrix']``. The matrix is inverted to produce
1195``.Phase['ISODISTORT']['Var2OccMatrix']``, which determines how to compute the
1196mode values from the delta coordinate values.
1197
1198
1199------------------------------
1200Mode Computations
1201------------------------------
1202
1203Constraints are processed after the CIF has been read in
1204:meth:`GSASIIdataGUI.GSASII.OnImportPhase` or 
1205:meth:`GSASIIscriptable.G2Project.add_phase` by moving them from the reader
1206object's ``.Constraints`` class variable to
1207the Constraints tree entry's ['Phase'] list (for list items defining constraints) or
1208the Constraints tree entry's ['_Explain'] dict (for dict items defining constraint help information)
1209
1210The information in ``.Phase['ISODISTORT']`` is used
1211in :func:`GSASIIconstrGUI.ShowIsoDistortCalc` which shows coordinate and mode
1212values, the latter with s.u. values. This can be called from the Constraints and Phase/Atoms tree items.
1213
1214Before each refinement, constraints are processed as
1215:ref:`described elsewhere <Constraints_processing>`. After a refinement
1216is complete, :func:`GSASIImapvars.PrintIndependentVars` shows the
1217shifts and s.u.'s on the refined modes, using GSAS-II values, but
1218:func:`GSASIIstrIO.PrintISOmodes` prints the ISODISTORT modes as computed
1219in the web site.
1220
1221
1222.. _ParameterLimits:
1223
1224.. index::
1225   single: Parameter limits
1226
1227Parameter Limits
1228------------------------------
1229
1230One of the most often requested "enhancements" for GSAS-II would be the inclusion
1231of constraints to force parameters such as occupancies or Uiso values to stay within
1232expected ranges. While it is possible for users to supply their own restraints that would
1233perform this by supplying an appropriate expression with the "General" restraints, the
1234GSAS-II authors do not feel that use of restraints or constraints are a good solution for
1235this common problem where parameters refine to non-physical values. This is because when
1236this occurs, most likely one of the following cases is occurring:
1237
1238#. there is a significant problem
1239   with the model, for example for an x-ray fit if an O atom is placed where a S is actually
1240   present, the Uiso will refine artificially small or the occupancy much larger than unity
1241   to try to compensate for the missing electrons; or
1242 
1243#. the data are simply insensitive
1244   to the parameter or combination of parameters, for example unless very high-Q data
1245   are included, the effects of a occupancy and Uiso value can have compensating effects,
1246   so an assumption must be made; likewise, with neutron data natural-abundance V atoms
1247   are nearly invisible due to weak coherent scattering. No parameters can be fit for a
1248   V atom with neutrons.
1249
1250#. the parameter is non-physical (such as a negative Uiso value) but within
1251   two sigma (sigma = standard uncertainty, aka e.s.d.) of a reasonable value,
1252   in which case the
1253   value is not problematic as it is experimentally indistinguishable from an
1254   expected value. 
1255
1256#. there is a systematic problem with the data (experimental error)
1257
1258In all these cases, this situation needs to be reviewed by a crystallographer to decide
1259how to best determine a structural model for these data. An implementation with a constraint
1260or restraint is likely to simply hide the problem from the user, making it more probable
1261that a poor model choice is obtained.
1262
1263What GSAS-II does implement is to allow users to specify ranges for parameters
1264that works by disabling
1265refinement of parameters that refine beyond either a lower limit or an upper limit, where
1266either or both may be optionally specified. Parameters limits are specified in the Controls
1267tree entry in dicts named as ``Controls['parmMaxDict']`` and ``Controls['parmMinDict']``, where
1268the keys are :class:`G2VarObj` objects corresponding to standard GSAS-II variable
1269(see :func:`getVarDescr` and :func:`CompileVarDesc`) names, where a
1270wildcard ('*') may optionally be used for histogram number or atom number
1271(phase number is intentionally not  allowed as a wildcard as it makes little sense
1272to group the same parameter together different phases). Note
1273that :func:`prmLookup` is used to see if a name matches a wildcard. The upper or lower limit
1274is placed into these dicts as a float value. These values can be edited using the window
1275created by the Calculate/"View LS parms" menu command or in scripting with the
1276:meth:`GSASIIscriptable.G2Project.set_Controls` function.
1277In the GUI, a checkbox labeled "match all histograms/atoms" is used to insert a wildcard
1278into the appropriate part of the variable name.
1279
1280When a refinement is conducted, routine :func:`GSASIIstrMain.dropOOBvars` is used to
1281find parameters that have refined to values outside their limits. If this occurs, the parameter
1282is set to the limiting value and the variable name is added to a list of frozen variables
1283(as a :class:`G2VarObj` objects) kept in a list in the
1284``Controls['parmFrozen']`` dict. In a sequential refinement, this is kept separate for
1285each histogram as a list in
1286``Controls['parmFrozen'][histogram]`` (where the key is the histogram name) or as a list in
1287``Controls['parmFrozen']['FrozenList']`` for a non-sequential fit.
1288This allows different variables
1289to be frozen in each section of a sequential fit.
1290Frozen parameters are not included in refinements through removal from the
1291list of parameters to be refined (``varyList``) in :func:`GSASIIstrMain.Refine` or
1292:func:`GSASIIstrMain.SeqRefine`.
1293The data window for the Controls tree item shows the number of Frozen variables and
1294the individual variables can be viewed with the Calculate/"View LS parms" menu window or
1295obtained with :meth:`GSASIIscriptable.G2Project.get_Frozen`.
1296Once a variable is frozen, it will not be refined in any
1297future refinements unless the the variable is removed (manually) from the list. This can also
1298be done with the Calculate/"View LS parms" menu window or
1299:meth:`GSASIIscriptable.G2Project.set_Frozen`.
1300
1301
1302.. seealso::
1303  :class:`G2VarObj`
1304  :func:`getVarDescr`
1305  :func:`CompileVarDesc`
1306  :func:`prmLookup`
1307  :class:`GSASIIctrlGUI.ShowLSParms`
1308  :class:`GSASIIctrlGUI.VirtualVarBox`
1309  :func:`GSASIIstrIO.SetUsedHistogramsAndPhases`
1310  :func:`GSASIIstrIO.SaveUpdatedHistogramsAndPhases`
1311  :func:`GSASIIstrIO.SetSeqResult`
1312  :func:`GSASIIstrMain.dropOOBvars`
1313  :meth:`GSASIIscriptable.G2Project.set_Controls`
1314  :meth:`GSASIIscriptable.G2Project.get_Frozen`
1315  :meth:`GSASIIscriptable.G2Project.set_Frozen`
1316
1317*Classes and routines*
1318----------------------
1319
1320'''
1321from __future__ import division, print_function
1322import platform
1323import re
1324import random as ran
1325import sys
1326import os.path as ospath
1327if '2' in platform.python_version_tuple()[0]:
1328    import cPickle
1329else:
1330    import pickle as cPickle
1331import GSASIIpath
1332import GSASIImath as G2mth
1333import GSASIIspc as G2spc
1334import numpy as np
1335
1336GSASIIpath.SetVersionNumber("$Revision: 5025 $")
1337
1338DefaultControls = {
1339    'deriv type':'analytic Hessian',
1340    'min dM/M':0.001,'shift factor':1.,'max cyc':3,'F**2':False,'SVDtol':1.e-6,
1341    'UsrReject':{'minF/sig':0.,'MinExt':0.01,'MaxDF/F':100.,'MaxD':500.,'MinD':0.05},
1342    'Copy2Next':False,'Reverse Seq':False,'HatomFix':False,
1343    'Author':'no name','newLeBail':False,
1344    'FreePrm1':'Sample humidity (%)',
1345    'FreePrm2':'Sample voltage (V)',
1346    'FreePrm3':'Applied load (MN)',
1347    'ShowCell':False,
1348    }
1349'''Values to be used as defaults for the initial contents of the ``Controls``
1350data tree item.
1351'''
1352def StripUnicode(string,subs='.'):
1353    '''Strip non-ASCII characters from strings
1354
1355    :param str string: string to strip Unicode characters from
1356    :param str subs: character(s) to place into string in place of each
1357      Unicode character. Defaults to '.'
1358
1359    :returns: a new string with only ASCII characters
1360    '''
1361    s = ''
1362    for c in string:
1363        if ord(c) < 128:
1364            s += c
1365        else:
1366            s += subs
1367    return s
1368#    return s.encode('ascii','replace')
1369
1370def MakeUniqueLabel(lbl,labellist):
1371    '''Make sure that every a label is unique against a list by adding
1372    digits at the end until it is not found in list.
1373
1374    :param str lbl: the input label
1375    :param list labellist: the labels that have already been encountered
1376    :returns: lbl if not found in labellist or lbl with ``_1-9`` (or
1377      ``_10-99``, etc.) appended at the end
1378    '''
1379    lbl = StripUnicode(lbl.strip(),'_')
1380    if not lbl: # deal with a blank label
1381        lbl = '_1'
1382    if lbl not in labellist:
1383        labellist.append(lbl)
1384        return lbl
1385    i = 1
1386    prefix = lbl
1387    if '_' in lbl:
1388        prefix = lbl[:lbl.rfind('_')]
1389        suffix = lbl[lbl.rfind('_')+1:]
1390        try:
1391            i = int(suffix)+1
1392        except: # suffix could not be parsed
1393            i = 1
1394            prefix = lbl
1395    while prefix+'_'+str(i) in labellist:
1396        i += 1
1397    else:
1398        lbl = prefix+'_'+str(i)
1399        labellist.append(lbl)
1400    return lbl
1401
1402PhaseIdLookup = {}
1403'''dict listing phase name and random Id keyed by sequential phase index as a str;
1404best to access this using :func:`LookupPhaseName`
1405'''
1406PhaseRanIdLookup = {}
1407'''dict listing phase sequential index keyed by phase random Id;
1408best to access this using :func:`LookupPhaseId`
1409'''
1410HistIdLookup = {}
1411'''dict listing histogram name and random Id, keyed by sequential histogram index as a str;
1412best to access this using :func:`LookupHistName`
1413'''
1414HistRanIdLookup = {}
1415'''dict listing histogram sequential index keyed by histogram random Id;
1416best to access this using :func:`LookupHistId`
1417'''
1418AtomIdLookup = {}
1419'''dict listing for each phase index as a str, the atom label and atom random Id,
1420keyed by atom sequential index as a str;
1421best to access this using :func:`LookupAtomLabel`
1422'''
1423AtomRanIdLookup = {}
1424'''dict listing for each phase the atom sequential index keyed by atom random Id;
1425best to access this using :func:`LookupAtomId`
1426'''
1427ShortPhaseNames = {}
1428'''a dict containing a possibly shortened and when non-unique numbered
1429version of the phase name. Keyed by the phase sequential index.
1430'''
1431ShortHistNames = {}
1432'''a dict containing a possibly shortened and when non-unique numbered
1433version of the histogram name. Keyed by the histogram sequential index.
1434'''
1435
1436#VarDesc = {}  # removed 1/30/19 BHT as no longer needed (I think)
1437#''' This dictionary lists descriptions for GSAS-II variables,
1438#as set in :func:`CompileVarDesc`. See that function for a description
1439#for how keys and values are written.
1440#'''
1441
1442reVarDesc = {}
1443''' This dictionary lists descriptions for GSAS-II variables where
1444keys are compiled regular expressions that will match the name portion
1445of a parameter name. Initialized in :func:`CompileVarDesc`.
1446'''
1447
1448reVarStep = {}
1449''' This dictionary lists the preferred step size for numerical
1450derivative computation w/r to a GSAS-II variable. Keys are compiled
1451regular expressions and values are the step size for that parameter.
1452Initialized in :func:`CompileVarDesc`.
1453'''
1454# create a default space group object for P1; N.B. fails when building documentation
1455try:
1456    P1SGData = G2spc.SpcGroup('P 1')[1] # data structure for default space group
1457except:
1458    pass
1459
1460def GetPhaseNames(fl):
1461    ''' Returns a list of phase names found under 'Phases' in GSASII gpx file
1462    NB: there is another one of these in GSASIIstrIO.py that uses the gpx filename
1463
1464    :param file fl: opened .gpx file
1465    :return: list of phase names
1466    '''
1467    PhaseNames = []
1468    while True:
1469        try:
1470            data = cPickle.load(fl)
1471        except EOFError:
1472            break
1473        datum = data[0]
1474        if 'Phases' == datum[0]:
1475            for datus in data[1:]:
1476                PhaseNames.append(datus[0])
1477    fl.seek(0)          #reposition file
1478    return PhaseNames
1479
1480def SetNewPhase(Name='New Phase',SGData=None,cell=None,Super=None):
1481    '''Create a new phase dict with default values for various parameters
1482
1483    :param str Name: Name for new Phase
1484
1485    :param dict SGData: space group data from :func:`GSASIIspc:SpcGroup`;
1486      defaults to data for P 1
1487
1488    :param list cell: unit cell parameter list; defaults to
1489      [1.0,1.0,1.0,90.,90,90.,1.]
1490
1491    '''
1492    if SGData is None: SGData = P1SGData
1493    if cell is None: cell=[1.0,1.0,1.0,90.,90.,90.,1.]
1494    phaseData = {
1495        'ranId':ran.randint(0,sys.maxsize),
1496        'General':{
1497            'Name':Name,
1498            'Type':'nuclear',
1499            'Modulated':False,
1500            'AtomPtrs':[3,1,7,9],
1501            'SGData':SGData,
1502            'Cell':[False,]+cell,
1503            'Pawley dmin':1.0,
1504            'Data plot type':'None',
1505            'SH Texture':{
1506                'Order':0,
1507                'Model':'cylindrical',
1508                'Sample omega':[False,0.0],
1509                'Sample chi':[False,0.0],
1510                'Sample phi':[False,0.0],
1511                'SH Coeff':[False,{}],
1512                'SHShow':False,
1513                'PFhkl':[0,0,1],
1514                'PFxyz':[0,0,1],
1515                'PlotType':'Pole figure',
1516                'Penalty':[['',],0.1,False,1.0]}},
1517        'Atoms':[],
1518        'Drawing':{},
1519        'Histograms':{},
1520        'Pawley ref':[],
1521        'RBModels':{},
1522        }
1523    if Super and Super.get('Use',False):
1524        phaseData['General'].update({'Modulated':True,'Super':True,'SuperSg':Super['ssSymb']})
1525        phaseData['General']['SSGData'] = G2spc.SSpcGroup(SGData,Super['ssSymb'])[1]
1526        phaseData['General']['SuperVec'] = [Super['ModVec'],False,Super['maxH']]
1527
1528    return phaseData
1529
1530def ReadCIF(URLorFile):
1531    '''Open a CIF, which may be specified as a file name or as a URL using PyCifRW
1532    (from James Hester).
1533    The open routine gets confused with DOS names that begin with a letter and colon
1534    "C:\\dir\" so this routine will try to open the passed name as a file and if that
1535    fails, try it as a URL
1536
1537    :param str URLorFile: string containing a URL or a file name. Code will try first
1538      to open it as a file and then as a URL.
1539
1540    :returns: a PyCifRW CIF object.
1541    '''
1542    import CifFile as cif # PyCifRW from James Hester
1543
1544    # alternate approach:
1545    #import urllib
1546    #ciffile = 'file:'+urllib.pathname2url(filename)
1547
1548    try:
1549        fp = open(URLorFile,'r')
1550        cf = cif.ReadCif(fp)
1551        fp.close()
1552        return cf
1553    except IOError:
1554        return cif.ReadCif(URLorFile)
1555
1556def TestIndexAll():
1557    '''Test if :func:`IndexAllIds` has been called to index all phases and
1558    histograms (this is needed before :func:`G2VarObj` can be used.
1559
1560    :returns: Returns True if indexing is needed.
1561    '''
1562    if PhaseIdLookup or AtomIdLookup or HistIdLookup:
1563        return False
1564    return True
1565       
1566def IndexAllIds(Histograms,Phases):
1567    '''Scan through the used phases & histograms and create an index
1568    to the random numbers of phases, histograms and atoms. While doing this,
1569    confirm that assigned random numbers are unique -- just in case lightning
1570    strikes twice in the same place.
1571
1572    Note: this code assumes that the atom random Id (ranId) is the last
1573    element each atom record.
1574
1575    This is called in three places (only): :func:`GSASIIstrIO.GetUsedHistogramsAndPhases`
1576    (which loads the histograms and phases from a GPX file),
1577    :meth:`~GSASIIdataGUI.GSASII.GetUsedHistogramsAndPhasesfromTree`
1578    (which loads the histograms and phases from the data tree.) and
1579    :meth:`GSASIIconstrGUI.UpdateConstraints`
1580    (which displays & edits the constraints in a GUI)
1581
1582    TODO: do we need a lookup for rigid body variables?
1583    '''
1584    # process phases and atoms
1585    PhaseIdLookup.clear()
1586    PhaseRanIdLookup.clear()
1587    AtomIdLookup.clear()
1588    AtomRanIdLookup.clear()
1589    ShortPhaseNames.clear()
1590    for ph in Phases:
1591        cx,ct,cs,cia = Phases[ph]['General']['AtomPtrs']
1592        ranId = Phases[ph]['ranId']
1593        while ranId in PhaseRanIdLookup:
1594            # Found duplicate random Id! note and reassign
1595            print ("\n\n*** Phase "+str(ph)+" has repeated ranId. Fixing.\n")
1596            Phases[ph]['ranId'] = ranId = ran.randint(0,sys.maxsize)
1597        pId = str(Phases[ph]['pId'])
1598        PhaseIdLookup[pId] = (ph,ranId)
1599        PhaseRanIdLookup[ranId] = pId
1600        shortname = ph  #[:10]
1601        while shortname in ShortPhaseNames.values():
1602            shortname = ph[:8] + ' ('+ pId + ')'
1603        ShortPhaseNames[pId] = shortname
1604        AtomIdLookup[pId] = {}
1605        AtomRanIdLookup[pId] = {}
1606        for iatm,at in enumerate(Phases[ph]['Atoms']):
1607            ranId = at[cia+8]
1608            while ranId in AtomRanIdLookup[pId]: # check for dups
1609                print ("\n\n*** Phase "+str(ph)+" atom "+str(iatm)+" has repeated ranId. Fixing.\n")
1610                at[cia+8] = ranId = ran.randint(0,sys.maxsize)
1611            AtomRanIdLookup[pId][ranId] = str(iatm)
1612            if Phases[ph]['General']['Type'] == 'macromolecular':
1613                label = '%s_%s_%s_%s'%(at[ct-1],at[ct-3],at[ct-4],at[ct-2])
1614            else:
1615                label = at[ct-1]
1616            AtomIdLookup[pId][str(iatm)] = (label,ranId)
1617    # process histograms
1618    HistIdLookup.clear()
1619    HistRanIdLookup.clear()
1620    ShortHistNames.clear()
1621    for hist in Histograms:
1622        ranId = Histograms[hist]['ranId']
1623        while ranId in HistRanIdLookup:
1624            # Found duplicate random Id! note and reassign
1625            print ("\n\n*** Histogram "+str(hist)+" has repeated ranId. Fixing.\n")
1626            Histograms[hist]['ranId'] = ranId = ran.randint(0,sys.maxsize)
1627        hId = str(Histograms[hist]['hId'])
1628        HistIdLookup[hId] = (hist,ranId)
1629        HistRanIdLookup[ranId] = hId
1630        shortname = hist[:15]
1631        while shortname in ShortHistNames.values():
1632            shortname = hist[:11] + ' ('+ hId + ')'
1633        ShortHistNames[hId] = shortname
1634
1635def LookupAtomId(pId,ranId):
1636    '''Get the atom number from a phase and atom random Id
1637
1638    :param int/str pId: the sequential number of the phase
1639    :param int ranId: the random Id assigned to an atom
1640
1641    :returns: the index number of the atom (str)
1642    '''
1643    if not AtomRanIdLookup:
1644        raise Exception('Error: LookupAtomId called before IndexAllIds was run')
1645    if pId is None or pId == '':
1646        raise KeyError('Error: phase is invalid (None or blank)')
1647    pId = str(pId)
1648    if pId not in AtomRanIdLookup:
1649        raise KeyError('Error: LookupAtomId does not have phase '+pId)
1650    if ranId not in AtomRanIdLookup[pId]:
1651        raise KeyError('Error: LookupAtomId, ranId '+str(ranId)+' not in AtomRanIdLookup['+pId+']')
1652    return AtomRanIdLookup[pId][ranId]
1653
1654def LookupAtomLabel(pId,index):
1655    '''Get the atom label from a phase and atom index number
1656
1657    :param int/str pId: the sequential number of the phase
1658    :param int index: the index of the atom in the list of atoms
1659
1660    :returns: the label for the atom (str) and the random Id of the atom (int)
1661    '''
1662    if not AtomIdLookup:
1663        raise Exception('Error: LookupAtomLabel called before IndexAllIds was run')
1664    if pId is None or pId == '':
1665        raise KeyError('Error: phase is invalid (None or blank)')
1666    pId = str(pId)
1667    if pId not in AtomIdLookup:
1668        raise KeyError('Error: LookupAtomLabel does not have phase '+pId)
1669    if index not in AtomIdLookup[pId]:
1670        raise KeyError('Error: LookupAtomLabel, ranId '+str(index)+' not in AtomRanIdLookup['+pId+']')
1671    return AtomIdLookup[pId][index]
1672
1673def LookupPhaseId(ranId):
1674    '''Get the phase number and name from a phase random Id
1675
1676    :param int ranId: the random Id assigned to a phase
1677    :returns: the sequential Id (pId) number for the phase (str)
1678    '''
1679    if not PhaseRanIdLookup:
1680        raise Exception('Error: LookupPhaseId called before IndexAllIds was run')
1681    if ranId not in PhaseRanIdLookup:
1682        raise KeyError('Error: LookupPhaseId does not have ranId '+str(ranId))
1683    return PhaseRanIdLookup[ranId]
1684
1685def LookupPhaseName(pId):
1686    '''Get the phase number and name from a phase Id
1687
1688    :param int/str pId: the sequential assigned to a phase
1689    :returns:  (phase,ranId) where phase is the name of the phase (str)
1690      and ranId is the random # id for the phase (int)
1691    '''
1692    if not PhaseIdLookup:
1693        raise Exception('Error: LookupPhaseName called before IndexAllIds was run')
1694    if pId is None or pId == '':
1695        raise KeyError('Error: phase is invalid (None or blank)')
1696    pId = str(pId)
1697    if pId not in PhaseIdLookup:
1698        raise KeyError('Error: LookupPhaseName does not have index '+pId)
1699    return PhaseIdLookup[pId]
1700
1701def LookupHistId(ranId):
1702    '''Get the histogram number and name from a histogram random Id
1703
1704    :param int ranId: the random Id assigned to a histogram
1705    :returns: the sequential Id (hId) number for the histogram (str)
1706    '''
1707    if not HistRanIdLookup:
1708        raise Exception('Error: LookupHistId called before IndexAllIds was run')
1709    if ranId not in HistRanIdLookup:
1710        raise KeyError('Error: LookupHistId does not have ranId '+str(ranId))
1711    return HistRanIdLookup[ranId]
1712
1713def LookupHistName(hId):
1714    '''Get the histogram number and name from a histogram Id
1715
1716    :param int/str hId: the sequential assigned to a histogram
1717    :returns:  (hist,ranId) where hist is the name of the histogram (str)
1718      and ranId is the random # id for the histogram (int)
1719    '''
1720    if not HistIdLookup:
1721        raise Exception('Error: LookupHistName called before IndexAllIds was run')
1722    if hId is None or hId == '':
1723        raise KeyError('Error: histogram is invalid (None or blank)')
1724    hId = str(hId)
1725    if hId not in HistIdLookup:
1726        raise KeyError('Error: LookupHistName does not have index '+hId)
1727    return HistIdLookup[hId]
1728
1729def fmtVarDescr(varname):
1730    '''Return a string with a more complete description for a GSAS-II variable
1731
1732    :param str varname: A full G2 variable name with 2 or 3 or 4
1733       colons (<p>:<h>:name[:<a>] or <p>::RBname:<r>:<t>])
1734
1735    :returns: a string with the description
1736    '''
1737    s,l = VarDescr(varname)
1738    return s+": "+l
1739
1740def VarDescr(varname):
1741    '''Return two strings with a more complete description for a GSAS-II variable
1742
1743    :param str name: A full G2 variable name with 2 or 3 or 4
1744       colons (<p>:<h>:name[:<a>] or <p>::RBname:<r>:<t>])
1745
1746    :returns: (loc,meaning) where loc describes what item the variable is mapped
1747      (phase, histogram, etc.) and meaning describes what the variable does.
1748    '''
1749
1750    # special handling for parameter names without a colon
1751    # for now, assume self-defining
1752    if varname.find(':') == -1:
1753        return "Global",varname
1754
1755    l = getVarDescr(varname)
1756    if not l:
1757        return ("invalid variable name ("+str(varname)+")!"),""
1758#        return "invalid variable name!",""
1759
1760    if not l[-1]:
1761        l[-1] = "(variable needs a definition! Set it in CompileVarDesc)"
1762
1763    if len(l) == 3:         #SASD variable name!
1764        s = 'component:'+l[1]
1765        return s,l[-1]
1766    s = ""
1767    if l[0] is not None and l[1] is not None: # HAP: keep short
1768        if l[2] == "Scale": # fix up ambigous name
1769            l[5] = "Phase fraction"
1770        if l[0] == '*':
1771            lbl = 'Seq. ref.'
1772        else:
1773            lbl = ShortPhaseNames.get(l[0],'? #'+str(l[0]))
1774        if l[1] == '*':
1775            hlbl = 'Seq. ref.'
1776        else:
1777            hlbl = ShortHistNames.get(l[1],'? #'+str(l[1]))
1778        if hlbl[:4] == 'HKLF':
1779            hlbl = 'Xtl='+hlbl[5:]
1780        elif hlbl[:4] == 'PWDR':
1781            hlbl = 'Pwd='+hlbl[5:]
1782        else:
1783            hlbl = 'Hist='+hlbl
1784        s = "Ph="+str(lbl)+" * "+str(hlbl)
1785    else:
1786        if l[2] == "Scale": # fix up ambigous name: must be scale factor, since not HAP
1787            l[5] = "Scale factor"
1788        if l[2] == 'Back': # background parameters are "special", alas
1789            s = 'Hist='+ShortHistNames.get(l[1],'? #'+str(l[1]))
1790            l[-1] += ' #'+str(l[3])
1791        elif l[4] is not None: # rigid body parameter or modulation parm
1792            lbl = ShortPhaseNames.get(l[0],'phase?')
1793            if 'RB' in l[2]:    #rigid body parm
1794                s = "RB body #"+str(l[3])+" (type "+str(l[4])+") in "+str(lbl) + ','
1795            else: #modulation parm
1796                s = 'Atom %s wave %s in %s'%(LookupAtomLabel(l[0],l[3])[0],l[4],lbl)
1797        elif l[3] is not None: # atom parameter,
1798            lbl = ShortPhaseNames.get(l[0],'phase?')
1799            try:
1800                albl = LookupAtomLabel(l[0],l[3])[0]
1801            except KeyError:
1802                albl = 'Atom?'
1803            s = "Atom "+str(albl)+" in "+str(lbl)
1804        elif l[0] == '*':
1805            s = "All phases "
1806        elif l[0] is not None:
1807            lbl = ShortPhaseNames.get(l[0],'phase?')
1808            s = "Phase "+str(lbl)
1809        elif l[1] == '*':
1810            s = 'All hists'
1811        elif l[1] is not None:
1812            hlbl = ShortHistNames.get(l[1],'? #'+str(l[1]))
1813            if hlbl[:4] == 'HKLF':
1814                hlbl = 'Xtl='+hlbl[5:]
1815            elif hlbl[:4] == 'PWDR':
1816                hlbl = 'Pwd='+hlbl[5:]
1817            else:
1818                hlbl = 'Hist='+hlbl
1819            s = str(hlbl)
1820    if not s:
1821        s = 'Global'
1822    return s,l[-1]
1823
1824def getVarDescr(varname):
1825    '''Return a short description for a GSAS-II variable
1826
1827    :param str name: A full G2 variable name with 2 or 3 or 4
1828       colons (<p>:<h>:name[:<a1>][:<a2>])
1829
1830    :returns: a six element list as [`p`,`h`,`name`,`a1`,`a2`,`description`],
1831      where `p`, `h`, `a1`, `a2` are str values or `None`, for the phase number,
1832      the histogram number and the atom number; `name` will always be
1833      a str; and `description` is str or `None`.
1834      If the variable name is incorrectly formed (for example, wrong
1835      number of colons), `None` is returned instead of a list.
1836    '''
1837    l = varname.split(':')
1838    if len(l) == 2:     #SASD parameter name
1839        return varname,l[0],getDescr(l[1])
1840    if len(l) == 3:
1841        l += [None,None]
1842    elif len(l) == 4:
1843        l += [None]
1844    elif len(l) != 5:
1845        return None
1846    for i in (0,1,3,4):
1847        if l[i] == "":
1848            l[i] = None
1849    l += [getDescr(l[2])]
1850    return l
1851
1852def CompileVarDesc():
1853    '''Set the values in the variable lookup tables
1854    (:attr:`reVarDesc` and :attr:`reVarStep`).
1855    This is called in :func:`getDescr` and :func:`getVarStep` so this
1856    initialization is always done before use.
1857
1858    Note that keys may contain regular expressions, where '[xyz]'
1859    matches 'x' 'y' or 'z' (equivalently '[x-z]' describes this as range
1860    of values). '.*' matches any string. For example::
1861
1862    'AUiso':'Atomic isotropic displacement parameter',
1863
1864    will match variable ``'p::AUiso:a'``.
1865    If parentheses are used in the key, the contents of those parentheses can be
1866    used in the value, such as::
1867
1868    'AU([123][123])':'Atomic anisotropic displacement parameter U\\1',
1869
1870    will match ``AU11``, ``AU23``,.. and `U11`, `U23` etc will be displayed
1871    in the value when used.
1872
1873    '''
1874    if reVarDesc: return # already done
1875    for key,value in {
1876        # derived or other sequential vars
1877        '([abc])$' : 'Lattice parameter, \\1, from Ai and Djk', # N.B. '$' prevents match if any characters follow
1878        u'\u03B1' : u'Lattice parameter, \u03B1, from Ai and Djk',
1879        u'\u03B2' : u'Lattice parameter, \u03B2, from Ai and Djk',
1880        u'\u03B3' : u'Lattice parameter, \u03B3, from Ai and Djk',
1881        # ambiguous, alas:
1882        'Scale' : 'Phase or Histogram scale factor',
1883        # Phase vars (p::<var>)
1884        'A([0-5])' : ('Reciprocal metric tensor component \\1',1e-5),
1885        '[vV]ol' : 'Unit cell volume', # probably an error that both upper and lower case are used
1886        # Atom vars (p::<var>:a)
1887        'dA([xyz])$' : ('change to atomic coordinate, \\1',1e-6),
1888        'A([xyz])$' : '\\1 fractional atomic coordinate',
1889        'AUiso':('Atomic isotropic displacement parameter',1e-4),
1890        'AU([123][123])':('Atomic anisotropic displacement parameter U\\1',1e-4),
1891        'Afrac': ('Atomic site fraction parameter',1e-5),
1892        'Amul': 'Atomic site multiplicity value',
1893        'AM([xyz])$' : 'Atomic magnetic moment parameter, \\1',
1894        # Hist (:h:<var>) & Phase (HAP) vars (p:h:<var>)
1895        'Back': 'Background term',
1896        'BkPkint;(.*)':'Background peak #\\1 intensity',
1897        'BkPkpos;(.*)':'Background peak #\\1 position',
1898        'BkPksig;(.*)':'Background peak #\\1 Gaussian width',
1899        'BkPkgam;(.*)':'Background peak #\\1 Cauchy width',
1900#        'Back File' : 'Background file name',
1901        'BF mult' : 'Background file multiplier',
1902        'Bab([AU])': 'Babinet solvent scattering coef. \\1',
1903        'D([123][123])' : 'Anisotropic strain coef. \\1',
1904        'Extinction' : 'Extinction coef.',
1905        'MD' : 'March-Dollase coef.',
1906        'Mustrain;.*' : 'Microstrain coef.',
1907        'Size;.*' : 'Crystallite size value',
1908        'eA$' : 'Cubic mustrain value',
1909        'Ep$' : 'Primary extinction',
1910        'Es$' : 'Secondary type II extinction',
1911        'Eg$' : 'Secondary type I extinction',
1912        'Flack' : 'Flack parameter',
1913        'TwinFr' : 'Twin fraction',
1914        'Layer Disp'  : 'Layer displacement along beam',
1915        #Histogram vars (:h:<var>)
1916        'Absorption' : 'Absorption coef.',
1917        'LayerDisp'  : 'Bragg-Brentano Layer displacement',
1918        'Displace([XY])' : ('Debye-Scherrer sample displacement \\1',0.1),
1919        'Lam' : ('Wavelength',1e-6),
1920        'I\\(L2\\)\\/I\\(L1\\)' : ('Ka2/Ka1 intensity ratio',0.001),
1921        'Polariz\\.' : ('Polarization correction',1e-3),
1922        'SH/L' : ('FCJ peak asymmetry correction',1e-4),
1923        '([UVW])$' : ('Gaussian instrument broadening \\1',1e-5),
1924        '([XYZ])$' : ('Cauchy instrument broadening \\1',1e-5),
1925        'Zero' : 'Debye-Scherrer zero correction',
1926        'Shift' : 'Bragg-Brentano sample displ.',
1927        'SurfRoughA' : 'Bragg-Brenano surface roughness A',
1928        'SurfRoughB' : 'Bragg-Brenano surface roughness B',
1929        'Transparency' : 'Bragg-Brentano sample tranparency',
1930        'DebyeA' : 'Debye model amplitude',
1931        'DebyeR' : 'Debye model radius',
1932        'DebyeU' : 'Debye model Uiso',
1933        'RBV.*' : 'Vector rigid body parameter',
1934        'RBVO([aijk])' : 'Vector rigid body orientation parameter \\1',
1935        'RBVP([xyz])' : 'Vector rigid body \\1 position parameter',
1936        'RBVf' : 'Vector rigid body site fraction',
1937        'RBV([TLS])([123AB][123AB])' : 'Residue rigid body group disp. param.',
1938        'RBVU' : 'Residue rigid body group Uiso param.',
1939        'RBRO([aijk])' : 'Residue rigid body orientation parameter \\1',
1940        'RBRP([xyz])' : 'Residue rigid body \\1 position parameter',
1941        'RBRTr;.*' : 'Residue rigid body torsion parameter',
1942        'RBRf' : 'Residue rigid body site fraction',
1943        'RBR([TLS])([123AB][123AB])' : 'Residue rigid body group disp. param.',
1944        'RBRU' : 'Residue rigid body group Uiso param.',
1945        'constr([0-9]*)' : 'Parameter from constraint',
1946        'nv-([^_]+)_*' : 'New variable constraint parameter named \\1',
1947        # supersymmetry parameters  p::<var>:a:o 'Flen','Fcent'?
1948        'mV([0-2])$' : 'Modulation vector component \\1',
1949        'Fsin'  :   'Sin site fraction modulation',
1950        'Fcos'  :   'Cos site fraction modulation',
1951        'Fzero'  :   'Crenel function offset',      #may go away
1952        'Fwid'   :   'Crenel function width',
1953        'Tmin'   :   'ZigZag/Block min location',
1954        'Tmax'   :   'ZigZag/Block max location',
1955        '([XYZ])max': 'ZigZag/Block max value for \\1',
1956        '([XYZ])sin'  : 'Sin position wave for \\1',
1957        '([XYZ])cos'  : 'Cos position wave for \\1',
1958        'U([123][123])sin$' :  'Sin thermal wave for U\\1',
1959        'U([123][123])cos$' :  'Cos thermal wave for U\\1',
1960        'M([XYZ])sin$' :  'Sin mag. moment wave for \\1',
1961        'M([XYZ])cos$' :  'Cos mag. moment wave for \\1',
1962        # PDF peak parms (l:<var>;l = peak no.)
1963        'PDFpos'  : 'PDF peak position',
1964        'PDFmag'  : 'PDF peak magnitude',
1965        'PDFsig'  : 'PDF peak std. dev.',
1966        # SASD vars (l:<var>;l = component)
1967        'Aspect ratio' : 'Particle aspect ratio',
1968        'Length' : 'Cylinder length',
1969        'Diameter' : 'Cylinder/disk diameter',
1970        'Thickness' : 'Disk thickness',
1971        'Shell thickness' : 'Multiplier to get inner(<1) or outer(>1) sphere radius',
1972        'Dist' : 'Interparticle distance',
1973        'VolFr' : 'Dense scatterer volume fraction',
1974        'epis' : 'Sticky sphere epsilon',
1975        'Sticky' : 'Stickyness',
1976        'Depth' : 'Well depth',
1977        'Width' : 'Well width',
1978        'Volume' : 'Particle volume',
1979        'Radius' : 'Sphere/cylinder/disk radius',
1980        'Mean' : 'Particle mean radius',
1981        'StdDev' : 'Standard deviation in Mean',
1982        'G$': 'Guinier prefactor',
1983        'Rg$': 'Guinier radius of gyration',
1984        'B$': 'Porod prefactor',
1985        'P$': 'Porod power',
1986        'Cutoff': 'Porod cutoff',
1987        'PkInt': 'Bragg peak intensity',
1988        'PkPos': 'Bragg peak position',
1989        'PkSig': 'Bragg peak sigma',
1990        'PkGam': 'Bragg peak gamma',
1991        'e([12][12])' : 'strain tensor e\\1',   # strain vars e11, e22, e12
1992        'Dcalc': 'Calc. d-spacing',
1993        'Back$': 'background parameter',
1994        'pos$': 'peak position',
1995        'int$': 'peak intensity',
1996        'WgtFrac':'phase weight fraction',
1997        'alpha':'TOF profile term',
1998        'alpha-[01]':'Pink profile term',
1999        'beta-[01q]':'TOF/Pink profile term',
2000        'sig-[012q]':'TOF profile term',
2001        'dif[ABC]':'TOF to d-space calibration',
2002        'C\\([0-9]*,[0-9]*\\)' : 'spherical harmonics preferred orientation coef.',
2003        'Pressure': 'Pressure level for measurement in MPa',
2004        'Temperature': 'T value for measurement, K',
2005        'FreePrm([123])': 'User defined measurement parameter \\1',
2006        'Gonio. radius': 'Distance from sample to detector, mm',
2007        }.items():
2008        # Needs documentation: HAP: LeBail, newLeBail
2009        # hist: Azimuth, Chi, Omega, Phi, Bank, nDebye, nPeaks
2010       
2011        if len(value) == 2:
2012            #VarDesc[key] = value[0]
2013            reVarDesc[re.compile(key)] = value[0]
2014            reVarStep[re.compile(key)] = value[1]
2015        else:
2016            #VarDesc[key] = value
2017            reVarDesc[re.compile(key)] = value
2018
2019def removeNonRefined(parmList):
2020    '''Remove items from variable list that are not refined and should not
2021    appear as options for constraints
2022
2023    :param list parmList: a list of strings of form "p:h:VAR:a" where
2024      VAR is the variable name
2025
2026    :returns: a list after removing variables where VAR matches a
2027      entry in local variable NonRefinedList
2028    '''
2029    NonRefinedList = ['Omega','Type','Chi','Phi', 'Azimuth','Gonio. radius',
2030                          'Lam1','Lam2','Back','Temperature','Pressure',
2031                          'FreePrm1','FreePrm2','FreePrm3',
2032                          'Source','nPeaks','LeBail','newLeBail','Bank',
2033                          'nDebye', #'',
2034                    ]
2035    return [prm for prm in parmList if prm.split(':')[2] not in NonRefinedList]
2036       
2037def getDescr(name):
2038    '''Return a short description for a GSAS-II variable
2039
2040    :param str name: The descriptive part of the variable name without colons (:)
2041
2042    :returns: a short description or None if not found
2043    '''
2044
2045    CompileVarDesc() # compile the regular expressions, if needed
2046    for key in reVarDesc:
2047        m = key.match(name)
2048        if m:
2049            reVarDesc[key]
2050            return m.expand(reVarDesc[key])
2051    return None
2052
2053def getVarStep(name,parmDict=None):
2054    '''Return a step size for computing the derivative of a GSAS-II variable
2055
2056    :param str name: A complete variable name (with colons, :)
2057    :param dict parmDict: A dict with parameter values or None (default)
2058
2059    :returns: a float that should be an appropriate step size, either from
2060      the value supplied in :func:`CompileVarDesc` or based on the value for
2061      name in parmDict, if supplied. If not found or the value is zero,
2062      a default value of 1e-5 is used. If parmDict is None (default) and
2063      no value is provided in :func:`CompileVarDesc`, then None is returned.
2064    '''
2065    CompileVarDesc() # compile the regular expressions, if needed
2066    for key in reVarStep:
2067        m = key.match(name)
2068        if m:
2069            return reVarStep[key]
2070    if parmDict is None: return None
2071    val = parmDict.get(key,0.0)
2072    if abs(val) > 0.05:
2073        return abs(val)/1000.
2074    else:
2075        return 1e-5
2076
2077def GenWildCard(varlist):
2078    '''Generate wildcard versions of G2 variables. These introduce '*'
2079    for a phase, histogram or atom number (but only for one of these
2080    fields) but only when there is more than one matching variable in the
2081    input variable list. So if the input is this::
2082
2083      varlist = ['0::AUiso:0', '0::AUiso:1', '1::AUiso:0']
2084
2085    then the output will be this::
2086
2087       wildList = ['*::AUiso:0', '0::AUiso:*']
2088
2089    :param list varlist: an input list of GSAS-II variable names
2090      (such as 0::AUiso:0)
2091
2092    :returns: wildList, the generated list of wild card variable names.
2093    '''
2094    wild = []
2095    for i in (0,1,3):
2096        currentL = varlist[:]
2097        while currentL:
2098            item1 = currentL.pop(0)
2099            i1splt = item1.split(':')
2100            if i >= len(i1splt): continue
2101            if i1splt[i]:
2102                nextL = []
2103                i1splt[i] = '[0-9]+'
2104                rexp = re.compile(':'.join(i1splt))
2105                matchlist = [item1]
2106                for nxtitem in currentL:
2107                    if rexp.match(nxtitem):
2108                        matchlist += [nxtitem]
2109                    else:
2110                        nextL.append(nxtitem)
2111                if len(matchlist) > 1:
2112                    i1splt[i] = '*'
2113                    wild.append(':'.join(i1splt))
2114                currentL = nextL
2115    return wild
2116
2117def LookupWildCard(varname,varlist):
2118    '''returns a list of variable names from list varname
2119    that match wildcard name in varname
2120
2121    :param str varname: a G2 variable name containing a wildcard
2122      (such as \\*::var)
2123    :param list varlist: the list of all variable names used in
2124      the current project
2125    :returns: a list of matching GSAS-II variables (may be empty)
2126    '''
2127    rexp = re.compile(varname.replace('*','[0-9]+'))
2128    return sorted([var for var in varlist if rexp.match(var)])
2129
2130def prmLookup(name,prmDict):
2131    '''Looks for a parameter in a min/max dictionary, optionally
2132    considering a wild card for histogram or atom number (use of
2133    both will never occur at the same time).
2134
2135    :param name: a GSAS-II parameter name (str, see :func:`getVarDescr`
2136      and :func:`CompileVarDesc`) or a :class:`G2VarObj` object.
2137    :param dict prmDict: a min/max dictionary, (parmMinDict
2138      or parmMaxDict in Controls) where keys are :class:`G2VarObj`
2139      objects.
2140    :returns: Two values, (**matchname**, **value**), are returned where:
2141
2142       * **matchname** *(str)* is the :class:`G2VarObj` object
2143         corresponding to the actual matched name,
2144         which could contain a wildcard even if **name** does not; and
2145       * **value** *(float)* which contains the parameter limit.
2146    '''
2147    for key,value in prmDict.items():
2148        if str(key) == str(name): return key,value
2149        if key == name: return key,value
2150    return None,None
2151       
2152
2153def _lookup(dic,key):
2154    '''Lookup a key in a dictionary, where None returns an empty string
2155    but an unmatched key returns a question mark. Used in :class:`G2VarObj`
2156    '''
2157    if key is None:
2158        return ""
2159    elif key == "*":
2160        return "*"
2161    else:
2162        return dic.get(key,'?')
2163
2164def SortVariables(varlist):
2165    '''Sorts variable names in a sensible manner
2166    '''
2167    def cvnnums(var):
2168        v = []
2169        for i in var.split(':'):
2170            if i == '':
2171                v.append(-1)
2172                continue
2173            try:
2174                v.append(int(i))
2175            except:
2176                v.append(i)
2177        return v
2178    return sorted(varlist,key=cvnnums)
2179
2180class G2VarObj(object):
2181    '''Defines a GSAS-II variable either using the phase/atom/histogram
2182    unique Id numbers or using a character string that specifies
2183    variables by phase/atom/histogram number (which can change).
2184    Note that :func:`GSASIIstrIO.GetUsedHistogramsAndPhases`,
2185    which calls :func:`IndexAllIds` (or
2186    :func:`GSASIIscriptable.G2Project.index_ids`) should be used to
2187    (re)load the current Ids
2188    before creating or later using the G2VarObj object.
2189
2190    This can store rigid body variables, but does not translate the residue # and
2191    body # to/from random Ids
2192
2193    A :class:`G2VarObj` object can be created with a single parameter:
2194
2195    :param str/tuple varname: a single value can be used to create a :class:`G2VarObj`
2196      object. If a string, it must be of form "p:h:var" or "p:h:var:a", where
2197
2198     * p is the phase number (which may be left blank or may be '*' to indicate all phases);
2199     * h is the histogram number (which may be left blank or may be '*' to indicate all histograms);
2200     * a is the atom number (which may be left blank in which case the third colon is omitted).
2201       The atom number can be specified as '*' if a phase number is specified (not as '*').
2202       For rigid body variables, specify a will be a string of form "residue:body#"
2203
2204      Alternately a single tuple of form (Phase,Histogram,VarName,AtomID) can be used, where
2205      Phase, Histogram, and AtomID are None or are ranId values (or one can be '*')
2206      and VarName is a string. Note that if Phase is '*' then the AtomID is an atom number.
2207      For a rigid body variables, AtomID is a string of form "residue:body#".
2208
2209    If four positional arguments are supplied, they are:
2210
2211    :param str/int phasenum: The number for the phase (or None or '*')
2212    :param str/int histnum: The number for the histogram (or None or '*')
2213    :param str varname: a single value can be used to create a :class:`G2VarObj`
2214    :param str/int atomnum: The number for the atom (or None or '*')
2215
2216    '''
2217    IDdict = {}
2218    IDdict['phases'] = {}
2219    IDdict['hists'] = {}
2220    IDdict['atoms'] = {}
2221    def __init__(self,*args):
2222        self.phase = None
2223        self.histogram = None
2224        self.name = ''
2225        self.atom = None
2226        if len(args) == 1 and (type(args[0]) is list or type(args[0]) is tuple) and len(args[0]) == 4:
2227            # single arg with 4 values
2228            self.phase,self.histogram,self.name,self.atom = args[0]
2229        elif len(args) == 1 and ':' in args[0]:
2230            #parse a string
2231            lst = args[0].split(':')
2232            if lst[0] == '*':
2233                self.phase = '*'
2234                if len(lst) > 3:
2235                    self.atom = lst[3]
2236                self.histogram = HistIdLookup.get(lst[1],[None,None])[1]
2237            elif lst[1] == '*':
2238                self.histogram = '*'
2239                self.phase = PhaseIdLookup.get(lst[0],[None,None])[1]
2240            else:
2241                self.histogram = HistIdLookup.get(lst[1],[None,None])[1]
2242                self.phase = PhaseIdLookup.get(lst[0],[None,None])[1]
2243                if len(lst) == 4:
2244                    if lst[3] == '*':
2245                        self.atom = '*'
2246                    else:
2247                        self.atom = AtomIdLookup[lst[0]].get(lst[3],[None,None])[1]
2248                elif len(lst) == 5:
2249                    self.atom = lst[3]+":"+lst[4]
2250                elif len(lst) == 3:
2251                    pass
2252                else:
2253                    raise Exception("Incorrect number of colons in var name "+str(args[0]))
2254            self.name = lst[2]
2255        elif len(args) == 4:
2256            if args[0] == '*':
2257                self.phase = '*'
2258                self.atom = args[3]
2259            else:
2260                self.phase = PhaseIdLookup.get(str(args[0]),[None,None])[1]
2261                if args[3] == '*':
2262                    self.atom = '*'
2263                elif args[0] is not None:
2264                    self.atom = AtomIdLookup[args[0]].get(str(args[3]),[None,None])[1]
2265            if args[1] == '*':
2266                self.histogram = '*'
2267            else:
2268                self.histogram = HistIdLookup.get(str(args[1]),[None,None])[1]
2269            self.name = args[2]
2270        else:
2271            raise Exception("Incorrectly called GSAS-II parameter name")
2272
2273        #print "DEBUG: created ",self.phase,self.histogram,self.name,self.atom
2274
2275    def __str__(self):
2276        return self.varname()
2277
2278    def __hash__(self):
2279        'Allow G2VarObj to be a dict key by implementing hashing'
2280        return hash(self.varname())
2281
2282    def varname(self):
2283        '''Formats the GSAS-II variable name as a "traditional" GSAS-II variable
2284        string (p:h:<var>:a) or (p:h:<var>)
2285
2286        :returns: the variable name as a str
2287        '''
2288        a = ""
2289        if self.phase == "*":
2290            ph = "*"
2291            if self.atom:
2292                a = ":" + str(self.atom)
2293        else:
2294            ph = _lookup(PhaseRanIdLookup,self.phase)
2295            if self.atom == '*':
2296                a = ':*'
2297            elif self.atom:
2298                if ":" in str(self.atom):
2299                    a = ":" + str(self.atom)
2300                elif ph in AtomRanIdLookup:
2301                    a = ":" + AtomRanIdLookup[ph].get(self.atom,'?')
2302                else:
2303                    a = ":?"
2304        if self.histogram == "*":
2305            hist = "*"
2306        else:
2307            hist = _lookup(HistRanIdLookup,self.histogram)
2308        s = (ph + ":" + hist + ":" + str(self.name)) + a
2309        return s
2310
2311    def __repr__(self):
2312        '''Return the detailed contents of the object
2313        '''
2314        s = "<"
2315        if self.phase == '*':
2316            s += "Phases: all; "
2317            if self.atom is not None:
2318                if ":" in str(self.atom):
2319                    s += "Rigid body" + str(self.atom) + "; "
2320                else:
2321                    s += "Atom #" + str(self.atom) + "; "
2322        elif self.phase is not None:
2323            ph =  _lookup(PhaseRanIdLookup,self.phase)
2324            s += "Phase: rId=" + str(self.phase) + " (#"+ ph + "); "
2325            if self.atom == '*':
2326                s += "Atoms: all; "
2327            elif ":" in str(self.atom):
2328                s += "Rigid body" + str(self.atom) + "; "
2329            elif self.atom is not None:
2330                s += "Atom rId=" + str(self.atom)
2331                if ph in AtomRanIdLookup:
2332                    s += " (#" + AtomRanIdLookup[ph].get(self.atom,'?') + "); "
2333                else:
2334                    s += " (#? -- not found!); "
2335        if self.histogram == '*':
2336            s += "Histograms: all; "
2337        elif self.histogram is not None:
2338            hist = _lookup(HistRanIdLookup,self.histogram)
2339            s += "Histogram: rId=" + str(self.histogram) + " (#"+ hist + "); "
2340        s += 'Variable name="' + str(self.name) + '">'
2341        return s+" ("+self.varname()+")"
2342
2343    def __eq__(self, other):
2344        '''Allow comparison of G2VarObj to other G2VarObj objects or strings.
2345        If any field is a wildcard ('*') that field matches.
2346        '''
2347        if type(other) is str:
2348            other = G2VarObj(other)
2349        elif type(other) is not G2VarObj:
2350            raise Exception("Invalid type ({}) for G2VarObj comparison with {}"
2351                            .format(type(other),other))
2352        if self.phase != other.phase and self.phase != '*' and other.phase != '*':
2353            return False
2354        if self.histogram != other.histogram and self.histogram != '*' and other.histogram != '*':
2355            return False
2356        if self.atom != other.atom and self.atom != '*' and other.atom != '*':
2357            return False
2358        if self.name != other.name:
2359            return False
2360        return True
2361
2362    def _show(self):
2363        'For testing, shows the current lookup table'
2364        print ('phases'+ self.IDdict['phases'])
2365        print ('hists'+ self.IDdict['hists'])
2366        print ('atomDict'+ self.IDdict['atoms'])
2367
2368#==========================================================================
2369def SetDefaultSample():
2370    'Fills in default items for the Sample dictionary for Debye-Scherrer & SASD'
2371    return {
2372        'InstrName':'',
2373        'ranId':ran.randint(0,sys.maxsize),
2374        'Scale':[1.0,True],'Type':'Debye-Scherrer','Absorption':[0.0,False],
2375        'DisplaceX':[0.0,False],'DisplaceY':[0.0,False],
2376        'Temperature':300.,'Pressure':0.1,'Time':0.0,
2377        'FreePrm1':0.,'FreePrm2':0.,'FreePrm3':0.,
2378        'Gonio. radius':200.0,
2379        'Omega':0.0,'Chi':0.0,'Phi':0.0,'Azimuth':0.0,
2380#SASD items
2381        'Materials':[{'Name':'vacuum','VolFrac':1.0,},{'Name':'vacuum','VolFrac':0.0,}],
2382        'Thick':1.0,'Contrast':[0.0,0.0],       #contrast & anomalous contrast
2383        'Trans':1.0,                            #measured transmission
2384        'SlitLen':0.0,                          #Slit length - in Q(A-1)
2385        }
2386######################################################################
2387class ImportBaseclass(object):
2388    '''Defines a base class for the reading of input files (diffraction
2389    data, coordinates,...). See :ref:`Writing a Import Routine<import_routines>`
2390    for an explanation on how to use a subclass of this class.
2391    '''
2392    class ImportException(Exception):
2393        '''Defines an Exception that is used when an import routine hits an expected error,
2394        usually in .Reader.
2395
2396        Good practice is that the Reader should define a value in self.errors that
2397        tells the user some information about what is wrong with their file.
2398        '''
2399        pass
2400
2401    UseReader = True  # in __init__ set value of self.UseReader to False to skip use of current importer
2402    def __init__(self,formatName,longFormatName=None,
2403                 extensionlist=[],strictExtension=False,):
2404        self.formatName = formatName # short string naming file type
2405        if longFormatName: # longer string naming file type
2406            self.longFormatName = longFormatName
2407        else:
2408            self.longFormatName = formatName
2409        # define extensions that are allowed for the file type
2410        # for windows, remove any extensions that are duplicate, as case is ignored
2411        if sys.platform == 'windows' and extensionlist:
2412            extensionlist = list(set([s.lower() for s in extensionlist]))
2413        self.extensionlist = extensionlist
2414        # If strictExtension is True, the file will not be read, unless
2415        # the extension matches one in the extensionlist
2416        self.strictExtension = strictExtension
2417        self.errors = ''
2418        self.warnings = ''
2419        self.SciPy = False          #image reader needed scipy
2420        # used for readers that will use multiple passes to read
2421        # more than one data block
2422        self.repeat = False
2423        self.selections = []
2424        self.repeatcount = 0
2425        self.readfilename = '?'
2426        self.scriptable = False
2427        #print 'created',self.__class__
2428
2429    def ReInitialize(self):
2430        'Reinitialize the Reader to initial settings'
2431        self.errors = ''
2432        self.warnings = ''
2433        self.SciPy = False          #image reader needed scipy
2434        self.repeat = False
2435        self.repeatcount = 0
2436        self.readfilename = '?'
2437
2438
2439#    def Reader(self, filename, filepointer, ParentFrame=None, **unused):
2440#        '''This method must be supplied in the child class to read the file.
2441#        if the read fails either return False or raise an Exception
2442#        preferably of type ImportException.
2443#        '''
2444#        #start reading
2445#        raise ImportException("Error occurred while...")
2446#        self.errors += "Hint for user on why the error occur
2447#        return False # if an error occurs
2448#        return True # if read OK
2449
2450    def ExtensionValidator(self, filename):
2451        '''This methods checks if the file has the correct extension
2452       
2453        :returns:
2454       
2455          * False if this filename will not be supported by this reader (only
2456            when strictExtension is True)
2457          * True if the extension matches the list supplied by the reader
2458          * None if the reader allows un-registered extensions
2459         
2460        '''
2461        if filename:
2462            ext = ospath.splitext(filename)[1]
2463            if not ext and self.strictExtension: return False
2464            for ext in self.extensionlist:               
2465                if sys.platform == 'windows':
2466                    if filename.lower().endswith(ext): return True
2467                else:
2468                    if filename.endswith(ext): return True
2469        if self.strictExtension:
2470            return False
2471        else:
2472            return None
2473
2474    def ContentsValidator(self, filename):
2475        '''This routine will attempt to determine if the file can be read
2476        with the current format.
2477        This will typically be overridden with a method that
2478        takes a quick scan of [some of]
2479        the file contents to do a "sanity" check if the file
2480        appears to match the selected format.
2481        the file must be opened here with the correct format (binary/text)
2482        '''
2483        #filepointer.seek(0) # rewind the file pointer
2484        return True
2485
2486    def CIFValidator(self, filepointer):
2487        '''A :meth:`ContentsValidator` for use to validate CIF files.
2488        '''
2489        filepointer.seek(0)
2490        for i,l in enumerate(filepointer):
2491            if i >= 1000: return True
2492            '''Encountered only blank lines or comments in first 1000
2493            lines. This is unlikely, but assume it is CIF anyway, since we are
2494            even less likely to find a file with nothing but hashes and
2495            blank lines'''
2496            line = l.strip()
2497            if len(line) == 0: # ignore blank lines
2498                continue
2499            elif line.startswith('#'): # ignore comments
2500                continue
2501            elif line.startswith('data_'): # on the right track, accept this file
2502                return True
2503            else: # found something invalid
2504                self.errors = 'line '+str(i+1)+' contains unexpected data:\n'
2505                if all([ord(c) < 128 and ord(c) != 0 for c in str(l)]): # show only if ASCII
2506                    self.errors += '  '+str(l)
2507                else:
2508                    self.errors += '  (binary)'
2509                self.errors += '\n  Note: a CIF should only have blank lines or comments before'
2510                self.errors += '\n        a data_ statement begins a block.'
2511                return False
2512
2513######################################################################
2514class ImportPhase(ImportBaseclass):
2515    '''Defines a base class for the reading of files with coordinates
2516
2517    Objects constructed that subclass this (in import/G2phase_*.py etc.) will be used
2518    in :meth:`GSASIIdataGUI.GSASII.OnImportPhase` and in
2519    :func:`GSASIIscriptable.import_generic`.
2520    See :ref:`Writing a Import Routine<import_routines>`
2521    for an explanation on how to use this class.
2522
2523    '''
2524    def __init__(self,formatName,longFormatName=None,extensionlist=[],
2525        strictExtension=False,):
2526        # call parent __init__
2527        ImportBaseclass.__init__(self,formatName,longFormatName,
2528            extensionlist,strictExtension)
2529        self.Phase = None # a phase must be created with G2IO.SetNewPhase in the Reader
2530        self.SymOps = {} # specified when symmetry ops are in file (e.g. CIF)
2531        self.Constraints = None
2532
2533######################################################################
2534class ImportStructFactor(ImportBaseclass):
2535    '''Defines a base class for the reading of files with tables
2536    of structure factors.
2537
2538    Structure factors are read with a call to :meth:`GSASIIdataGUI.GSASII.OnImportSfact`
2539    which in turn calls :meth:`GSASIIdataGUI.GSASII.OnImportGeneric`, which calls
2540    methods :meth:`ExtensionValidator`, :meth:`ContentsValidator` and
2541    :meth:`Reader`.
2542
2543    See :ref:`Writing a Import Routine<import_routines>`
2544    for an explanation on how to use import classes in general. The specifics
2545    for reading a structure factor histogram require that
2546    the ``Reader()`` routine in the import
2547    class need to do only a few things: It
2548    should load :attr:`RefDict` item ``'RefList'`` with the reflection list,
2549    and set :attr:`Parameters` with the instrument parameters
2550    (initialized with :meth:`InitParameters` and set with :meth:`UpdateParameters`).
2551    '''
2552    def __init__(self,formatName,longFormatName=None,extensionlist=[],
2553        strictExtension=False,):
2554        ImportBaseclass.__init__(self,formatName,longFormatName,
2555            extensionlist,strictExtension)
2556
2557        # define contents of Structure Factor entry
2558        self.Parameters = []
2559        'self.Parameters is a list with two dicts for data parameter settings'
2560        self.InitParameters()
2561        self.RefDict = {'RefList':[],'FF':{},'Super':0}
2562        self.Banks = []             #for multi bank data (usually TOF)
2563        '''self.RefDict is a dict containing the reflection information, as read from the file.
2564        Item 'RefList' contains the reflection information. See the
2565        :ref:`Single Crystal Reflection Data Structure<XtalRefl_table>`
2566        for the contents of each row. Dict element 'FF'
2567        contains the form factor values for each element type; if this entry
2568        is left as initialized (an empty list) it will be initialized as needed later.
2569        '''
2570    def ReInitialize(self):
2571        'Reinitialize the Reader to initial settings'
2572        ImportBaseclass.ReInitialize(self)
2573        self.InitParameters()
2574        self.Banks = []             #for multi bank data (usually TOF)
2575        self.RefDict = {'RefList':[],'FF':{},'Super':0}
2576
2577    def InitParameters(self):
2578        'initialize the instrument parameters structure'
2579        Lambda = 0.70926
2580        HistType = 'SXC'
2581        self.Parameters = [{'Type':[HistType,HistType], # create the structure
2582                            'Lam':[Lambda,Lambda]
2583                            }, {}]
2584        'Parameters is a list with two dicts for data parameter settings'
2585
2586    def UpdateParameters(self,Type=None,Wave=None):
2587        'Revise the instrument parameters'
2588        if Type is not None:
2589            self.Parameters[0]['Type'] = [Type,Type]
2590        if Wave is not None:
2591            self.Parameters[0]['Lam'] = [Wave,Wave]
2592
2593######################################################################
2594class ImportPowderData(ImportBaseclass):
2595    '''Defines a base class for the reading of files with powder data.
2596
2597    Objects constructed that subclass this (in import/G2pwd_*.py etc.) will be used
2598    in :meth:`GSASIIdataGUI.GSASII.OnImportPowder` and in
2599    :func:`GSASIIscriptable.import_generic`.
2600    See :ref:`Writing a Import Routine<import_routines>`
2601    for an explanation on how to use this class.
2602    '''
2603    def __init__(self,formatName,longFormatName=None,
2604        extensionlist=[],strictExtension=False,):
2605        ImportBaseclass.__init__(self,formatName,longFormatName,
2606            extensionlist,strictExtension)
2607        self.clockWd = None  # used in TOF
2608        self.ReInitialize()
2609
2610    def ReInitialize(self):
2611        'Reinitialize the Reader to initial settings'
2612        ImportBaseclass.ReInitialize(self)
2613        self.powderentry = ['',None,None] #  (filename,Pos,Bank)
2614        self.powderdata = [] # Powder dataset
2615        '''A powder data set is a list with items [x,y,w,yc,yb,yd]:
2616                np.array(x), # x-axis values
2617                np.array(y), # powder pattern intensities
2618                np.array(w), # 1/sig(intensity)^2 values (weights)
2619                np.array(yc), # calc. intensities (zero)
2620                np.array(yb), # calc. background (zero)
2621                np.array(yd), # obs-calc profiles
2622        '''
2623        self.comments = []
2624        self.idstring = ''
2625        self.Sample = SetDefaultSample() # default sample parameters
2626        self.Controls = {}  # items to be placed in top-level Controls
2627        self.GSAS = None     # used in TOF
2628        self.repeat_instparm = True # Should a parm file be
2629        #                             used for multiple histograms?
2630        self.instparm = None # name hint from file of instparm to use
2631        self.instfile = '' # full path name to instrument parameter file
2632        self.instbank = '' # inst parm bank number
2633        self.instmsg = ''  # a label that gets printed to show
2634                           # where instrument parameters are from
2635        self.numbanks = 1
2636        self.instdict = {} # place items here that will be transferred to the instrument parameters
2637        self.pwdparms = {} # place parameters that are transferred directly to the tree
2638                           # here (typically from an existing GPX file)
2639######################################################################
2640class ImportSmallAngleData(ImportBaseclass):
2641    '''Defines a base class for the reading of files with small angle data.
2642    See :ref:`Writing a Import Routine<import_routines>`
2643    for an explanation on how to use this class.
2644    '''
2645    def __init__(self,formatName,longFormatName=None,extensionlist=[],
2646        strictExtension=False,):
2647
2648        ImportBaseclass.__init__(self,formatName,longFormatName,extensionlist,
2649            strictExtension)
2650        self.ReInitialize()
2651
2652    def ReInitialize(self):
2653        'Reinitialize the Reader to initial settings'
2654        ImportBaseclass.ReInitialize(self)
2655        self.smallangleentry = ['',None,None] #  (filename,Pos,Bank)
2656        self.smallangledata = [] # SASD dataset
2657        '''A small angle data set is a list with items [x,y,w,yc,yd]:
2658                np.array(x), # x-axis values
2659                np.array(y), # powder pattern intensities
2660                np.array(w), # 1/sig(intensity)^2 values (weights)
2661                np.array(yc), # calc. intensities (zero)
2662                np.array(yd), # obs-calc profiles
2663                np.array(yb), # preset bkg
2664        '''
2665        self.comments = []
2666        self.idstring = ''
2667        self.Sample = SetDefaultSample()
2668        self.GSAS = None     # used in TOF
2669        self.clockWd = None  # used in TOF
2670        self.numbanks = 1
2671        self.instdict = {} # place items here that will be transferred to the instrument parameters
2672
2673######################################################################
2674class ImportReflectometryData(ImportBaseclass):
2675    '''Defines a base class for the reading of files with reflectometry data.
2676    See :ref:`Writing a Import Routine<import_routines>`
2677    for an explanation on how to use this class.
2678    '''
2679    def __init__(self,formatName,longFormatName=None,extensionlist=[],
2680        strictExtension=False,):
2681
2682        ImportBaseclass.__init__(self,formatName,longFormatName,extensionlist,
2683            strictExtension)
2684        self.ReInitialize()
2685
2686    def ReInitialize(self):
2687        'Reinitialize the Reader to initial settings'
2688        ImportBaseclass.ReInitialize(self)
2689        self.reflectometryentry = ['',None,None] #  (filename,Pos,Bank)
2690        self.reflectometrydata = [] # SASD dataset
2691        '''A small angle data set is a list with items [x,y,w,yc,yd]:
2692                np.array(x), # x-axis values
2693                np.array(y), # powder pattern intensities
2694                np.array(w), # 1/sig(intensity)^2 values (weights)
2695                np.array(yc), # calc. intensities (zero)
2696                np.array(yd), # obs-calc profiles
2697                np.array(yb), # preset bkg
2698        '''
2699        self.comments = []
2700        self.idstring = ''
2701        self.Sample = SetDefaultSample()
2702        self.GSAS = None     # used in TOF
2703        self.clockWd = None  # used in TOF
2704        self.numbanks = 1
2705        self.instdict = {} # place items here that will be transferred to the instrument parameters
2706
2707######################################################################
2708class ImportPDFData(ImportBaseclass):
2709    '''Defines a base class for the reading of files with PDF G(R) data.
2710    See :ref:`Writing a Import Routine<import_routines>`
2711    for an explanation on how to use this class.
2712    '''
2713    def __init__(self,formatName,longFormatName=None,extensionlist=[],
2714        strictExtension=False,):
2715
2716        ImportBaseclass.__init__(self,formatName,longFormatName,extensionlist,
2717            strictExtension)
2718        self.ReInitialize()
2719
2720    def ReInitialize(self):
2721        'Reinitialize the Reader to initial settings'
2722        ImportBaseclass.ReInitialize(self)
2723        self.pdfentry = ['',None,None] #  (filename,Pos,Bank)
2724        self.pdfdata = [] # PDF G(R) dataset
2725        '''A pdf g(r) data set is a list with items [x,y]:
2726                np.array(x), # r-axis values
2727                np.array(y), # pdf g(r)
2728        '''
2729        self.comments = []
2730        self.idstring = ''
2731        self.numbanks = 1
2732
2733######################################################################
2734class ImportImage(ImportBaseclass):
2735    '''Defines a base class for the reading of images
2736
2737    Images are read in only these places:
2738
2739      * Initial reading is typically done from a menu item
2740        with a call to :meth:`GSASIIdataGUI.GSASII.OnImportImage`
2741        which in turn calls :meth:`GSASIIdataGUI.GSASII.OnImportGeneric`. That calls
2742        methods :meth:`ExtensionValidator`, :meth:`ContentsValidator` and
2743        :meth:`Reader`. This returns a list of reader objects for each read image.
2744        Also used in :func:`GSASIIscriptable.import_generic`.
2745
2746      * Images are read alternatively in :func:`GSASIIIO.ReadImages`, which puts image info
2747        directly into the data tree.
2748
2749      * Images are reloaded with :func:`GSASIIIO.GetImageData`.
2750
2751    When reading an image, the ``Reader()`` routine in the ImportImage class
2752    should set:
2753
2754      * :attr:`Comments`: a list of strings (str),
2755      * :attr:`Npix`: the number of pixels in the image (int),
2756      * :attr:`Image`: the actual image as a numpy array (np.array)
2757      * :attr:`Data`: a dict defining image parameters (dict). Within this dict the following
2758        data items are needed:
2759
2760         * 'pixelSize': size of each pixel in microns (such as ``[200.,200.]``.
2761         * 'wavelength': wavelength in :math:`\\AA`.
2762         * 'distance': distance of detector from sample in cm.
2763         * 'center': uncalibrated center of beam on detector (such as ``[204.8,204.8]``.
2764         * 'size': size of image (such as ``[2048,2048]``).
2765         * 'ImageTag': image number or other keyword used to retrieve image from
2766           a multi-image data file (defaults to ``1`` if not specified).
2767         * 'sumfile': holds sum image file name if a sum was produced from a multi image file
2768
2769    optional data items:
2770
2771      * :attr:`repeat`: set to True if there are additional images to
2772        read in the file, False otherwise
2773      * :attr:`repeatcount`: set to the number of the image.
2774
2775    Note that the above is initialized with :meth:`InitParameters`.
2776    (Also see :ref:`Writing a Import Routine<import_routines>`
2777    for an explanation on how to use import classes in general.)
2778    '''
2779    def __init__(self,formatName,longFormatName=None,extensionlist=[],
2780        strictExtension=False,):
2781        ImportBaseclass.__init__(self,formatName,longFormatName,
2782            extensionlist,strictExtension)
2783        self.InitParameters()
2784
2785    def ReInitialize(self):
2786        'Reinitialize the Reader to initial settings -- not used at present'
2787        ImportBaseclass.ReInitialize(self)
2788        self.InitParameters()
2789
2790    def InitParameters(self):
2791        'initialize the instrument parameters structure'
2792        self.Comments = ['No comments']
2793        self.Data = {'samplechangerpos':0.0,'det2theta':0.0,'Gain map':''}
2794        self.Npix = 0
2795        self.Image = None
2796        self.repeat = False
2797        self.repeatcount = 1
2798        self.sumfile = ''
2799
2800    def LoadImage(self,ParentFrame,imagefile,imagetag=None):
2801        '''Optionally, call this after reading in an image to load it into the tree.
2802        This saves time by preventing a reread of the same information.
2803        '''
2804        if ParentFrame:
2805            ParentFrame.ImageZ = self.Image   # store the image for plotting
2806            ParentFrame.oldImagefile = imagefile # save the name of the last image file read
2807            ParentFrame.oldImageTag = imagetag   # save the tag of the last image file read
2808
2809#################################################################################################
2810# shortcut routines
2811exp = np.exp
2812sind = sin = s = lambda x: np.sin(x*np.pi/180.)
2813cosd = cos = c = lambda x: np.cos(x*np.pi/180.)
2814tand = tan = t = lambda x: np.tan(x*np.pi/180.)
2815sqrt = sq = lambda x: np.sqrt(x)
2816pi = lambda: np.pi
2817
2818def FindFunction(f):
2819    '''Find the object corresponding to function f
2820
2821    :param str f: a function name such as 'numpy.exp'
2822    :returns: (pkgdict,pkgobj) where pkgdict contains a dict
2823      that defines the package location(s) and where pkgobj
2824      defines the object associated with the function.
2825      If the function is not found, pkgobj is None.
2826    '''
2827    df = f.split('.')
2828    pkgdict = {}
2829    # no listed module name, try in current namespace
2830    if len(df) == 1:
2831        try:
2832            fxnobj = eval(f)
2833            return pkgdict,fxnobj
2834        except (AttributeError, NameError):
2835            return None,None
2836
2837    # includes a package, see if package is already imported
2838    pkgnam = '.'.join(df[:-1])
2839    try:
2840        fxnobj = eval(f)
2841        pkgdict[pkgnam] = eval(pkgnam)
2842        return pkgdict,fxnobj
2843    except (AttributeError, NameError):
2844        pass
2845    # package not yet imported, so let's try
2846    if '.' not in sys.path: sys.path.append('.')
2847    pkgnam = '.'.join(df[:-1])
2848    #for pkg in f.split('.')[:-1]: # if needed, descend down the tree
2849    #    if pkgname:
2850    #        pkgname += '.' + pkg
2851    #    else:
2852    #        pkgname = pkg
2853    try:
2854        exec('import '+pkgnam)
2855        pkgdict[pkgnam] = eval(pkgnam)
2856        fxnobj = eval(f)
2857    except Exception as msg:
2858        print('load of '+pkgnam+' failed with error='+str(msg))
2859        return {},None
2860    # can we access the function? I am not exactly sure what
2861    #    I intended this to test originally (BHT)
2862    try:
2863        fxnobj = eval(f,globals(),pkgdict)
2864        return pkgdict,fxnobj
2865    except Exception as msg:
2866        print('call to',f,' failed with error=',str(msg))
2867        return None,None # not found
2868               
2869class ExpressionObj(object):
2870    '''Defines an object with a user-defined expression, to be used for
2871    secondary fits or restraints. Object is created null, but is changed
2872    using :meth:`LoadExpression`. This contains only the minimum
2873    information that needs to be stored to save and load the expression
2874    and how it is mapped to GSAS-II variables.
2875    '''
2876    def __init__(self):
2877        self.expression = ''
2878        'The expression as a text string'
2879        self.assgnVars = {}
2880        '''A dict where keys are label names in the expression mapping to a GSAS-II
2881        variable. The value a G2 variable name.
2882        Note that the G2 variable name may contain a wild-card and correspond to
2883        multiple values.
2884        '''
2885        self.freeVars = {}
2886        '''A dict where keys are label names in the expression mapping to a free
2887        parameter. The value is a list with:
2888
2889         * a name assigned to the parameter
2890         * a value for to the parameter and
2891         * a flag to determine if the variable is refined.
2892        '''
2893        self.depVar = None
2894
2895        self.lastError = ('','')
2896        '''Shows last encountered error in processing expression
2897        (list of 1-3 str values)'''
2898
2899        self.distance_dict  = None  # to be used for defining atom phase/symmetry info
2900        self.distance_atoms = None  # to be used for defining atom distances
2901
2902    def LoadExpression(self,expr,exprVarLst,varSelect,varName,varValue,varRefflag):
2903        '''Load the expression and associated settings into the object. Raises
2904        an exception if the expression is not parsed, if not all functions
2905        are defined or if not all needed parameter labels in the expression
2906        are defined.
2907
2908        This will not test if the variable referenced in these definitions
2909        are actually in the parameter dictionary. This is checked when the
2910        computation for the expression is done in :meth:`SetupCalc`.
2911
2912        :param str expr: the expression
2913        :param list exprVarLst: parameter labels found in the expression
2914        :param dict varSelect: this will be 0 for Free parameters
2915          and non-zero for expression labels linked to G2 variables.
2916        :param dict varName: Defines a name (str) associated with each free parameter
2917        :param dict varValue: Defines a value (float) associated with each free parameter
2918        :param dict varRefflag: Defines a refinement flag (bool)
2919          associated with each free parameter
2920        '''
2921        self.expression = expr
2922        self.compiledExpr = None
2923        self.freeVars = {}
2924        self.assgnVars = {}
2925        for v in exprVarLst:
2926            if varSelect[v] == 0:
2927                self.freeVars[v] = [
2928                    varName.get(v),
2929                    varValue.get(v),
2930                    varRefflag.get(v),
2931                    ]
2932            else:
2933                self.assgnVars[v] = varName[v]
2934        self.CheckVars()
2935
2936    def EditExpression(self,exprVarLst,varSelect,varName,varValue,varRefflag):
2937        '''Load the expression and associated settings from the object into
2938        arrays used for editing.
2939
2940        :param list exprVarLst: parameter labels found in the expression
2941        :param dict varSelect: this will be 0 for Free parameters
2942          and non-zero for expression labels linked to G2 variables.
2943        :param dict varName: Defines a name (str) associated with each free parameter
2944        :param dict varValue: Defines a value (float) associated with each free parameter
2945        :param dict varRefflag: Defines a refinement flag (bool)
2946          associated with each free parameter
2947
2948        :returns: the expression as a str
2949        '''
2950        for v in self.freeVars:
2951            varSelect[v] = 0
2952            varName[v] = self.freeVars[v][0]
2953            varValue[v] = self.freeVars[v][1]
2954            varRefflag[v] = self.freeVars[v][2]
2955        for v in self.assgnVars:
2956            varSelect[v] = 1
2957            varName[v] = self.assgnVars[v]
2958        return self.expression
2959
2960    def GetVaried(self):
2961        'Returns the names of the free parameters that will be refined'
2962        return ["::"+self.freeVars[v][0] for v in self.freeVars if self.freeVars[v][2]]
2963
2964    def GetVariedVarVal(self):
2965        'Returns the names and values of the free parameters that will be refined'
2966        return [("::"+self.freeVars[v][0],self.freeVars[v][1]) for v in self.freeVars if self.freeVars[v][2]]
2967
2968    def UpdateVariedVars(self,varyList,values):
2969        'Updates values for the free parameters (after a refinement); only updates refined vars'
2970        for v in self.freeVars:
2971            if not self.freeVars[v][2]: continue
2972            if "::"+self.freeVars[v][0] not in varyList: continue
2973            indx = list(varyList).index("::"+self.freeVars[v][0])
2974            self.freeVars[v][1] = values[indx]
2975
2976    def GetIndependentVars(self):
2977        'Returns the names of the required independent parameters used in expression'
2978        return [self.assgnVars[v] for v in self.assgnVars]
2979
2980    def CheckVars(self):
2981        '''Check that the expression can be parsed, all functions are
2982        defined and that input loaded into the object is internally
2983        consistent. If not an Exception is raised.
2984
2985        :returns: a dict with references to packages needed to
2986          find functions referenced in the expression.
2987        '''
2988        ret = self.ParseExpression(self.expression)
2989        if not ret:
2990            raise Exception("Expression parse error")
2991        exprLblList,fxnpkgdict = ret
2992        # check each var used in expression is defined
2993        defined = list(self.assgnVars.keys()) + list(self.freeVars.keys())
2994        notfound = []
2995        for var in exprLblList:
2996            if var not in defined:
2997                notfound.append(var)
2998        if notfound:
2999            msg = 'Not all variables defined'
3000            msg1 = 'The following variables were not defined: '
3001            msg2 = ''
3002            for var in notfound:
3003                if msg: msg += ', '
3004                msg += var
3005            self.lastError = (msg1,'  '+msg2)
3006            raise Exception(msg)
3007        return fxnpkgdict
3008
3009    def ParseExpression(self,expr):
3010        '''Parse an expression and return a dict of called functions and
3011        the variables used in the expression. Returns None in case an error
3012        is encountered. If packages are referenced in functions, they are loaded
3013        and the functions are looked up into the modules global
3014        workspace.
3015
3016        Note that no changes are made to the object other than
3017        saving an error message, so that this can be used for testing prior
3018        to the save.
3019
3020        :returns: a list of used variables
3021        '''
3022        self.lastError = ('','')
3023        import ast
3024        def ASTtransverse(node,fxn=False):
3025            '''Transverse a AST-parsed expresson, compiling a list of variables
3026            referenced in the expression. This routine is used recursively.
3027
3028            :returns: varlist,fxnlist where
3029              varlist is a list of referenced variable names and
3030              fxnlist is a list of used functions
3031            '''
3032            varlist = []
3033            fxnlist = []
3034            if isinstance(node, list):
3035                for b in node:
3036                    v,f = ASTtransverse(b,fxn)
3037                    varlist += v
3038                    fxnlist += f
3039            elif isinstance(node, ast.AST):
3040                for a, b in ast.iter_fields(node):
3041                    if isinstance(b, ast.AST):
3042                        if a == 'func':
3043                            fxnlist += ['.'.join(ASTtransverse(b,True)[0])]
3044                            continue
3045                        v,f = ASTtransverse(b,fxn)
3046                        varlist += v
3047                        fxnlist += f
3048                    elif isinstance(b, list):
3049                        v,f = ASTtransverse(b,fxn)
3050                        varlist += v
3051                        fxnlist += f
3052                    elif node.__class__.__name__ == "Name":
3053                        varlist += [b]
3054                    elif fxn and node.__class__.__name__ == "Attribute":
3055                        varlist += [b]
3056            return varlist,fxnlist
3057        try:
3058            exprast = ast.parse(expr)
3059        except SyntaxError:
3060            s = ''
3061            import traceback
3062            for i in traceback.format_exc().splitlines()[-3:-1]:
3063                if s: s += "\n"
3064                s += str(i)
3065            self.lastError = ("Error parsing expression:",s)
3066            return
3067        # find the variables & functions
3068        v,f = ASTtransverse(exprast)
3069        varlist = sorted(list(set(v)))
3070        fxnlist = list(set(f))
3071        pkgdict = {}
3072        # check the functions are defined
3073        for fxn in fxnlist:
3074            fxndict,fxnobj = FindFunction(fxn)
3075            if not fxnobj:
3076                self.lastError = ("Error: Invalid function",fxn,
3077                                  "is not defined")
3078                return
3079            if not hasattr(fxnobj,'__call__'):
3080                self.lastError = ("Error: Not a function.",fxn,
3081                                  "cannot be called as a function")
3082                return
3083            pkgdict.update(fxndict)
3084        return varlist,pkgdict
3085
3086    def GetDepVar(self):
3087        'return the dependent variable, or None'
3088        return self.depVar
3089
3090    def SetDepVar(self,var):
3091        'Set the dependent variable, if used'
3092        self.depVar = var
3093#==========================================================================
3094class ExpressionCalcObj(object):
3095    '''An object used to evaluate an expression from a :class:`ExpressionObj`
3096    object.
3097
3098    :param ExpressionObj exprObj: a :class:`~ExpressionObj` expression object with
3099      an expression string and mappings for the parameter labels in that object.
3100    '''
3101    def __init__(self,exprObj):
3102        self.eObj = exprObj
3103        'The expression and mappings; a :class:`ExpressionObj` object'
3104        self.compiledExpr = None
3105        'The expression as compiled byte-code'
3106        self.exprDict = {}
3107        '''dict that defines values for labels used in expression and packages
3108        referenced by functions
3109        '''
3110        self.lblLookup = {}
3111        '''Lookup table that specifies the expression label name that is
3112        tied to a particular GSAS-II parameters in the parmDict.
3113        '''
3114        self.fxnpkgdict = {}
3115        '''a dict with references to packages needed to
3116        find functions referenced in the expression.
3117        '''
3118        self.varLookup = {}
3119        '''Lookup table that specifies the GSAS-II variable(s)
3120        indexed by the expression label name. (Used for only for diagnostics
3121        not evaluation of expression.)
3122        '''
3123        self.su = None
3124        '''Standard error evaluation where supplied by the evaluator
3125        '''
3126        # Patch: for old-style expressions with a (now removed step size)
3127        if '2' in platform.python_version_tuple()[0]: 
3128            basestr = basestring
3129        else:
3130            basestr = str
3131        for v in self.eObj.assgnVars:
3132            if not isinstance(self.eObj.assgnVars[v], basestr):
3133                self.eObj.assgnVars[v] = self.eObj.assgnVars[v][0]
3134        self.parmDict = {}
3135        '''A copy of the parameter dictionary, for distance and angle computation
3136        '''
3137
3138    def SetupCalc(self,parmDict):
3139        '''Do all preparations to use the expression for computation.
3140        Adds the free parameter values to the parameter dict (parmDict).
3141        '''
3142        if self.eObj.expression.startswith('Dist') or self.eObj.expression.startswith('Angle'):
3143            return
3144        self.fxnpkgdict = self.eObj.CheckVars()
3145        # all is OK, compile the expression
3146        self.compiledExpr = compile(self.eObj.expression,'','eval')
3147
3148        # look at first value in parmDict to determine its type
3149        parmsInList = True
3150        if '2' in platform.python_version_tuple()[0]: 
3151            basestr = basestring
3152        else:
3153            basestr = str
3154        for key in parmDict:
3155            val = parmDict[key]
3156            if isinstance(val, basestr):
3157                parmsInList = False
3158                break
3159            try: # check if values are in lists
3160                val = parmDict[key][0]
3161            except (TypeError,IndexError):
3162                parmsInList = False
3163            break
3164
3165        # set up the dicts needed to speed computations
3166        self.exprDict = {}
3167        self.lblLookup = {}
3168        self.varLookup = {}
3169        for v in self.eObj.freeVars:
3170            varname = self.eObj.freeVars[v][0]
3171            varname = "::" + varname.lstrip(':').replace(' ','_').replace(':',';')
3172            self.lblLookup[varname] = v
3173            self.varLookup[v] = varname
3174            if parmsInList:
3175                parmDict[varname] = [self.eObj.freeVars[v][1],self.eObj.freeVars[v][2]]
3176            else:
3177                parmDict[varname] = self.eObj.freeVars[v][1]
3178            self.exprDict[v] = self.eObj.freeVars[v][1]
3179        for v in self.eObj.assgnVars:
3180            varname = self.eObj.assgnVars[v]
3181            if varname in parmDict:
3182                self.lblLookup[varname] = v
3183                self.varLookup[v] = varname
3184                if parmsInList:
3185                    self.exprDict[v] = parmDict[varname][0]
3186                else:
3187                    self.exprDict[v] = parmDict[varname]
3188            elif '*' in varname:
3189                varlist = LookupWildCard(varname,list(parmDict.keys()))
3190                if len(varlist) == 0:
3191                    raise Exception("No variables match "+str(v))
3192                for var in varlist:
3193                    self.lblLookup[var] = v
3194                if parmsInList:
3195                    self.exprDict[v] = np.array([parmDict[var][0] for var in varlist])
3196                else:
3197                    self.exprDict[v] = np.array([parmDict[var] for var in varlist])
3198                self.varLookup[v] = [var for var in varlist]
3199            else:
3200                self.exprDict[v] = None
3201#                raise Exception,"No value for variable "+str(v)
3202        self.exprDict.update(self.fxnpkgdict)
3203
3204    def UpdateVars(self,varList,valList):
3205        '''Update the dict for the expression with a set of values
3206        :param list varList: a list of variable names
3207        :param list valList: a list of corresponding values
3208        '''
3209        for var,val in zip(varList,valList):
3210            self.exprDict[self.lblLookup.get(var,'undefined: '+var)] = val
3211
3212    def UpdateDict(self,parmDict):
3213        '''Update the dict for the expression with values in a dict
3214        :param dict parmDict: a dict of values, items not in use are ignored
3215        '''
3216        if self.eObj.expression.startswith('Dist') or self.eObj.expression.startswith('Angle'):
3217            self.parmDict = parmDict
3218            return
3219        for var in parmDict:
3220            if var in self.lblLookup:
3221                self.exprDict[self.lblLookup[var]] = parmDict[var]
3222
3223    def EvalExpression(self):
3224        '''Evaluate an expression. Note that the expression
3225        and mapping are taken from the :class:`ExpressionObj` expression object
3226        and the parameter values were specified in :meth:`SetupCalc`.
3227        :returns: a single value for the expression. If parameter
3228        values are arrays (for example, from wild-carded variable names),
3229        the sum of the resulting expression is returned.
3230
3231        For example, if the expression is ``'A*B'``,
3232        where A is 2.0 and B maps to ``'1::Afrac:*'``, which evaluates to::
3233
3234        [0.5, 1, 0.5]
3235
3236        then the result will be ``4.0``.
3237        '''
3238        self.su = None
3239        if self.eObj.expression.startswith('Dist'):
3240#            GSASIIpath.IPyBreak()
3241            dist = G2mth.CalcDist(self.eObj.distance_dict, self.eObj.distance_atoms, self.parmDict)
3242            return dist
3243        elif self.eObj.expression.startswith('Angle'):
3244            angle = G2mth.CalcAngle(self.eObj.angle_dict, self.eObj.angle_atoms, self.parmDict)
3245            return angle
3246        if self.compiledExpr is None:
3247            raise Exception("EvalExpression called before SetupCalc")
3248        try:
3249            val = eval(self.compiledExpr,globals(),self.exprDict)
3250        except TypeError:
3251            val = None
3252        if not np.isscalar(val):
3253            val = np.sum(val)
3254        return val
3255
3256class G2Exception(Exception):
3257    'A generic GSAS-II exception class'
3258    def __init__(self,msg):
3259        self.msg = msg
3260    def __str__(self):
3261        return repr(self.msg)
3262
3263class G2RefineCancel(Exception):
3264    'Raised when Cancel is pressed in a refinement dialog'
3265    def __init__(self,msg):
3266        self.msg = msg
3267    def __str__(self):
3268        return repr(self.msg)
3269   
3270def HowDidIgetHere(wherecalledonly=False):
3271    '''Show a traceback with calls that brought us to the current location.
3272    Used for debugging.
3273    '''
3274    import traceback
3275    if wherecalledonly:
3276        i = traceback.format_list(traceback.extract_stack()[:-1])[-2]
3277        print(i.strip().rstrip())
3278    else:
3279        print (70*'*')
3280        for i in traceback.format_list(traceback.extract_stack()[:-1]): print(i.strip().rstrip())
3281        print (70*'*')
3282
3283# Note that this is GUI code and should be moved at somepoint
3284def CreatePDFitems(G2frame,PWDRtree,ElList,Qlimits,numAtm=1,FltBkg=0,PDFnames=[]):
3285    '''Create and initialize a new set of PDF tree entries
3286
3287    :param Frame G2frame: main GSAS-II tree frame object
3288    :param str PWDRtree: name of PWDR to be used to create PDF item
3289    :param dict ElList: data structure with composition
3290    :param list Qlimits: Q limits to be used for computing the PDF
3291    :param float numAtm: no. atom in chemical formula
3292    :param float FltBkg: flat background value
3293    :param list PDFnames: previously used PDF names
3294
3295    :returns: the Id of the newly created PDF entry
3296    '''
3297    PDFname = 'PDF '+PWDRtree[4:] # this places two spaces after PDF, which is needed is some places
3298    if PDFname in PDFnames:
3299        print('Skipping, entry already exists: '+PDFname)
3300        return None
3301    #PDFname = MakeUniqueLabel(PDFname,PDFnames)
3302    Id = G2frame.GPXtree.AppendItem(parent=G2frame.root,text=PDFname)
3303    Data = {
3304        'Sample':{'Name':PWDRtree,'Mult':1.0},
3305        'Sample Bkg.':{'Name':'','Mult':-1.0,'Refine':False},
3306        'Container':{'Name':'','Mult':-1.0,'Refine':False},
3307        'Container Bkg.':{'Name':'','Mult':-1.0},'ElList':ElList,
3308        'Geometry':'Cylinder','Diam':1.0,'Pack':0.50,'Form Vol':10.0*numAtm,'Flat Bkg':FltBkg,
3309        'DetType':'Area detector','ObliqCoeff':0.3,'Ruland':0.025,'QScaleLim':Qlimits,
3310        'Lorch':False,'BackRatio':0.0,'Rmax':100.,'noRing':False,'IofQmin':1.0,'Rmin':1.0,
3311        'I(Q)':[],'S(Q)':[],'F(Q)':[],'G(R)':[]}
3312    G2frame.GPXtree.SetItemPyData(G2frame.GPXtree.AppendItem(Id,text='PDF Controls'),Data)
3313    G2frame.GPXtree.SetItemPyData(G2frame.GPXtree.AppendItem(Id,text='PDF Peaks'),
3314        {'Limits':[1.,5.],'Background':[2,[0.,-0.2*np.pi],False],'Peaks':[]})
3315    return Id
3316
3317class ShowTiming(object):
3318    '''An object to use for timing repeated sections of code.
3319
3320    Create the object with::
3321       tim0 = ShowTiming()
3322
3323    Tag sections of code to be timed with::
3324       tim0.start('start')
3325       tim0.start('in section 1')
3326       tim0.start('in section 2')
3327       
3328    etc. (Note that each section should have a unique label.)
3329
3330    After the last section, end timing with::
3331       tim0.end()
3332
3333    Show timing results with::
3334       tim0.show()
3335       
3336    '''
3337    def __init__(self):
3338        self.timeSum =  []
3339        self.timeStart = []
3340        self.label = []
3341        self.prev = None
3342    def start(self,label):
3343        import time
3344        if label in self.label:
3345            i = self.label.index(label)
3346            self.timeStart[i] = time.time()
3347        else:
3348            i = len(self.label)
3349            self.timeSum.append(0.0)
3350            self.timeStart.append(time.time())
3351            self.label.append(label)
3352        if self.prev is not None:
3353            self.timeSum[self.prev] += self.timeStart[i] - self.timeStart[self.prev]
3354        self.prev = i
3355    def end(self):
3356        import time
3357        if self.prev is not None:
3358            self.timeSum[self.prev] += time.time() - self.timeStart[self.prev]
3359        self.prev = None
3360    def show(self):
3361        sumT = sum(self.timeSum)
3362        print('Timing results (total={:.2f} sec)'.format(sumT))
3363        for i,(lbl,val) in enumerate(zip(self.label,self.timeSum)):
3364            print('{} {:20} {:8.2f} ms {:5.2f}%'.format(i,lbl,1000.*val,100*val/sumT))
3365
3366def validateAtomDrawType(typ,generalData={}):
3367    '''Confirm that the selected Atom drawing type is valid for the current
3368    phase. If not, use 'vdW balls'. This is currently used only for setting a
3369    default when atoms are added to the atoms draw list.
3370    '''
3371    if typ in ('lines','vdW balls','sticks','balls & sticks','ellipsoids'):
3372        return typ
3373    # elif generalData.get('Type','') == 'macromolecular':
3374    #     if typ in ('backbone',):
3375    #         return typ
3376    return 'vdW balls'
3377
3378if __name__ == "__main__":
3379    # test variable descriptions
3380    for var in '0::Afrac:*',':1:Scale','1::dAx:0','::undefined':
3381        v = var.split(':')[2]
3382        print(var+':\t', getDescr(v),getVarStep(v))
3383    import sys; sys.exit()
3384    # test equation evaluation
3385    def showEQ(calcobj):
3386        print (50*'=')
3387        print (calcobj.eObj.expression+'='+calcobj.EvalExpression())
3388        for v in sorted(calcobj.varLookup):
3389            print ("  "+v+'='+calcobj.exprDict[v]+'='+calcobj.varLookup[v])
3390        # print '  Derivatives'
3391        # for v in calcobj.derivStep.keys():
3392        #     print '    d(Expr)/d('+v+') =',calcobj.EvalDeriv(v)
3393
3394    obj = ExpressionObj()
3395
3396    obj.expression = "A*np.exp(B)"
3397    obj.assgnVars =  {'B': '0::Afrac:1'}
3398    obj.freeVars =  {'A': [u'A', 0.5, True]}
3399    #obj.CheckVars()
3400    calcobj = ExpressionCalcObj(obj)
3401
3402    obj1 = ExpressionObj()
3403    obj1.expression = "A*np.exp(B)"
3404    obj1.assgnVars =  {'B': '0::Afrac:*'}
3405    obj1.freeVars =  {'A': [u'Free Prm A', 0.5, True]}
3406    #obj.CheckVars()
3407    calcobj1 = ExpressionCalcObj(obj1)
3408
3409    obj2 = ExpressionObj()
3410    obj2.distance_stuff = np.array([[0,1],[1,-1]])
3411    obj2.expression = "Dist(1,2)"
3412    GSASIIpath.InvokeDebugOpts()
3413    parmDict2 = {'0::Afrac:0':[0.0,True], '0::Afrac:1': [1.0,False]}
3414    calcobj2 = ExpressionCalcObj(obj2)
3415    calcobj2.SetupCalc(parmDict2)
3416    showEQ(calcobj2)
3417
3418    parmDict1 = {'0::Afrac:0':1.0, '0::Afrac:1': 1.0}
3419    print ('\nDict = '+parmDict1)
3420    calcobj.SetupCalc(parmDict1)
3421    showEQ(calcobj)
3422    calcobj1.SetupCalc(parmDict1)
3423    showEQ(calcobj1)
3424
3425    parmDict2 = {'0::Afrac:0':[0.0,True], '0::Afrac:1': [1.0,False]}
3426    print ('Dict = '+parmDict2)
3427    calcobj.SetupCalc(parmDict2)
3428    showEQ(calcobj)
3429    calcobj1.SetupCalc(parmDict2)
3430    showEQ(calcobj1)
3431    calcobj2.SetupCalc(parmDict2)
3432    showEQ(calcobj2)
Note: See TracBrowser for help on using the repository browser.