\(\renewcommand\AA{\text{Å}}\)
4. GSASIIobj: Data objects & Docs
This module defines many data structures used in GSAS-II, as well as provides misc. support routines for accessing them.
4.1. GSASIIobj Classes and routines
Classes and routines defined in GSASIIobj
follow.
- GSASIIobj.AddPhase2Index(rdObj, filename)[source]
Add a phase to the index during reading Used where constraints are generated during import (ISODISTORT CIFs)
- GSASIIobj.AtomIdLookup = {}
dict listing for each phase index as a str, the atom label and atom random Id, keyed by atom sequential index as a str; best to access this using
LookupAtomLabel()
- GSASIIobj.AtomRanIdLookup = {}
dict listing for each phase the atom sequential index keyed by atom random Id; best to access this using
LookupAtomId()
- GSASIIobj.CompileVarDesc()[source]
Set the values in the variable lookup tables (
reVarDesc
andreVarStep
). This is called ingetDescr()
andgetVarStep()
so this initialization is always done before use. These variables are also used in script makeVarTbl.py which creates the table in section 3.2 of the Sphinx docs (Parameter names in GSAS-II).Note that keys may contain regular expressions, where ‘[xyz]’ matches ‘x’ ‘y’ or ‘z’ (equivalently ‘[x-z]’ describes this as range of values). ‘.*’ matches any string. For example:
'AUiso':'Atomic isotropic displacement parameter',
will match variable
'p::AUiso:a'
. If parentheses are used in the key, the contents of those parentheses can be used in the value, such as:'AU([123][123])':'Atomic anisotropic displacement parameter U\1',
will match
AU11
,AU23
,… and U11, U23 etc will be displayed in the value when used.
- GSASIIobj.CreatePDFitems(G2frame, PWDRtree, ElList, Qlimits, numAtm=1, FltBkg=0, PDFnames=[])[source]
Create and initialize a new set of PDF tree entries
- Parameters:
G2frame (Frame) – main GSAS-II tree frame object
PWDRtree (str) – name of PWDR to be used to create PDF item
ElList (dict) – data structure with composition
Qlimits (list) – Q limits to be used for computing the PDF
numAtm (float) – no. atom in chemical formula
FltBkg (float) – flat background value
PDFnames (list) – previously used PDF names
- Returns:
the Id of the newly created PDF entry
- GSASIIobj.DefaultControls = {'Author': 'no name', 'Copy2Next': False, 'F**2': False, 'FreePrm1': 'Sample humidity (%)', 'FreePrm2': 'Sample voltage (V)', 'FreePrm3': 'Applied load (MN)', 'HatomFix': False, 'Reverse Seq': False, 'SVDtol': 1e-06, 'ShowCell': False, 'UsrReject': {'MaxD': 500.0, 'MaxDF/F': 100.0, 'MinD': 0.05, 'MinExt': 0.01, 'minF/sig': 0.0}, 'deriv type': 'analytic Hessian', 'max cyc': 3, 'min dM/M': 0.001, 'newLeBail': False, 'shift factor': 1.0}
Values to be used as defaults for the initial contents of the
Controls
data tree item.
- class GSASIIobj.ExpressionCalcObj(exprObj)[source]
An object used to evaluate an expression from a
ExpressionObj
object.- Parameters:
exprObj (ExpressionObj) – a
ExpressionObj
expression object with an expression string and mappings for the parameter labels in that object.
- EvalExpression()[source]
Evaluate an expression. Note that the expression and mapping are taken from the
ExpressionObj
expression object and the parameter values were specified inSetupCalc()
. :returns: a single value for the expression. If parameter values are arrays (for example, from wild-carded variable names), the sum of the resulting expression is returned.For example, if the expression is
'A*B'
, where A is 2.0 and B maps to'1::Afrac:*'
, which evaluates to:[0.5, 1, 0.5]
then the result will be
4.0
.
- SetupCalc(parmDict)[source]
Do all preparations to use the expression for computation. Adds the free parameter values to the parameter dict (parmDict).
- UpdateDict(parmDict)[source]
Update the dict for the expression with values in a dict :param dict parmDict: a dict of values, items not in use are ignored
- UpdateVars(varList, valList)[source]
Update the dict for the expression with a set of values :param list varList: a list of variable names :param list valList: a list of corresponding values
- __weakref__
list of weak references to the object
- compiledExpr
The expression as compiled byte-code
- eObj
The expression and mappings; a
ExpressionObj
object
- exprDict
dict that defines values for labels used in expression and packages referenced by functions
- fxnpkgdict
a dict with references to packages needed to find functions referenced in the expression.
- lblLookup
Lookup table that specifies the expression label name that is tied to a particular GSAS-II parameters in the parmDict.
- parmDict
A copy of the parameter dictionary, for distance and angle computation
- su
Standard error evaluation where supplied by the evaluator
- varLookup
Lookup table that specifies the GSAS-II variable(s) indexed by the expression label name. (Used for only for diagnostics not evaluation of expression.)
- class GSASIIobj.ExpressionObj[source]
Defines an object with a user-defined expression, to be used for secondary fits or restraints. Object is created null, but is changed using
LoadExpression()
. This contains only the minimum information that needs to be stored to save and load the expression and how it is mapped to GSAS-II variables.- CheckVars()[source]
Check that the expression can be parsed, all functions are defined and that input loaded into the object is internally consistent. If not an Exception is raised.
- Returns:
a dict with references to packages needed to find functions referenced in the expression.
- EditExpression(exprVarLst, varSelect, varName, varValue, varRefflag)[source]
Load the expression and associated settings from the object into arrays used for editing.
- Parameters:
exprVarLst (list) – parameter labels found in the expression
varSelect (dict) – this will be 0 for Free parameters and non-zero for expression labels linked to G2 variables.
varName (dict) – Defines a name (str) associated with each free parameter
varValue (dict) – Defines a value (float) associated with each free parameter
varRefflag (dict) – Defines a refinement flag (bool) associated with each free parameter
- Returns:
the expression as a str
- GetIndependentVars()[source]
Returns the names of the required independent parameters used in expression
- LoadExpression(expr, exprVarLst, varSelect, varName, varValue, varRefflag)[source]
Load the expression and associated settings into the object. Raises an exception if the expression is not parsed, if not all functions are defined or if not all needed parameter labels in the expression are defined.
This will not test if the variable referenced in these definitions are actually in the parameter dictionary. This is checked when the computation for the expression is done in
SetupCalc()
.- Parameters:
expr (str) – the expression
exprVarLst (list) – parameter labels found in the expression
varSelect (dict) – this will be 0 for Free parameters and non-zero for expression labels linked to G2 variables.
varName (dict) – Defines a name (str) associated with each free parameter
varValue (dict) – Defines a value (float) associated with each free parameter
varRefflag (dict) – Defines a refinement flag (bool) associated with each free parameter
- ParseExpression(expr)[source]
Parse an expression and return a dict of called functions and the variables used in the expression. Returns None in case an error is encountered. If packages are referenced in functions, they are loaded and the functions are looked up into the modules global workspace.
Note that no changes are made to the object other than saving an error message, so that this can be used for testing prior to the save.
- Returns:
a list of used variables
- UpdateVariedVars(varyList, values)[source]
Updates values for the free parameters (after a refinement); only updates refined vars
- __weakref__
list of weak references to the object
- assgnVars
A dict where keys are label names in the expression mapping to a GSAS-II variable. The value a G2 variable name. Note that the G2 variable name may contain a wild-card and correspond to multiple values.
- expression
The expression as a text string
- freeVars
A dict where keys are label names in the expression mapping to a free parameter. The value is a list with:
a name assigned to the parameter
a value for to the parameter and
a flag to determine if the variable is refined.
- lastError
Shows last encountered error in processing expression (list of 1-3 str values)
- GSASIIobj.FindFunction(f)[source]
Find the object corresponding to function f
- Parameters:
f (str) – a function name such as ‘numpy.exp’
- Returns:
(pkgdict,pkgobj) where pkgdict contains a dict that defines the package location(s) and where pkgobj defines the object associated with the function. If the function is not found, pkgobj is None.
- exception GSASIIobj.G2Exception(msg)[source]
A generic GSAS-II exception class
- __weakref__
list of weak references to the object
- exception GSASIIobj.G2RefineCancel(msg)[source]
Raised when Cancel is pressed in a refinement dialog
- __weakref__
list of weak references to the object
- class GSASIIobj.G2VarObj(*args)[source]
Defines a GSAS-II variable either using the phase/atom/histogram unique Id numbers or using a character string that specifies variables by phase/atom/histogram number (which can change). Note that
GSASIIstrIO.GetUsedHistogramsAndPhases()
, which callsIndexAllIds()
(orGSASIIscriptable.G2Project.index_ids()
) should be used to (re)load the current Ids before creating or later using the G2VarObj object.This can store rigid body variables, but does not translate the residue # and body # to/from random Ids
A
G2VarObj
object can be created with a single parameter:- Parameters:
varname (str/tuple) –
- a single value can be used to create a
G2VarObj
object. If a string, it must be of form “p:h:var” or “p:h:var:a”, where
p is the phase number (which may be left blank or may be ‘*’ to indicate all phases);
h is the histogram number (which may be left blank or may be ‘*’ to indicate all histograms);
a is the atom number (which may be left blank in which case the third colon is omitted). The atom number can be specified as ‘*’ if a phase number is specified (not as ‘*’). For rigid body variables, specify a will be a string of form “residue:body#”
Alternately a single tuple of form (Phase,Histogram,VarName,AtomID) can be used, where Phase, Histogram, and AtomID are None or are ranId values (or one can be ‘*’) and VarName is a string. Note that if Phase is ‘*’ then the AtomID is an atom number. For a rigid body variables, AtomID is a string of form “residue:body#”.
- a single value can be used to create a
If four positional arguments are supplied, they are:
- Parameters:
phasenum (str/int) – The number for the phase (or None or ‘*’)
histnum (str/int) – The number for the histogram (or None or ‘*’)
varname (str) – a single value can be used to create a
G2VarObj
atomnum (str/int) – The number for the atom (or None or ‘*’)
- __eq__(other)[source]
Allow comparison of G2VarObj to other G2VarObj objects or strings. If any field is a wildcard (‘*’) that field matches.
- __weakref__
list of weak references to the object
- fmtVarByMode(seqmode, note, warnmsg)[source]
Format a parameter object for display. Note that these changes are only temporary and are only shown only when the Constraints data tree is selected.
In a non-sequential refinement or where the mode is ‘use-all’, the name is converted unchanged to a str
In a sequential refinement when the mode is ‘wildcards-only’ the name is converted unchanged to a str but a warning is added for non-wildcarded HAP or Histogram parameters
In a sequential refinement or where the mode is ‘auto-wildcard’, a histogram number is converted to a wildcard (*) and then converted to str
- Parameters:
mode (str) – the sequential mode (see above)
note (str) – value displayed on the line of the constraint/equiv.
warnmsg (str) – a message saying the constraint is not used
- Returns:
varname, explain, note, warnmsg (all str values) where:
varname is the parameter expressed as a string,
explain is blank unless there is a warning explanation about the parameter or blank
note is the previous value unless overridden
warnmsg is the previous value unless overridden
- GSASIIobj.GenWildCard(varlist)[source]
Generate wildcard versions of G2 variables. These introduce ‘*’ for a phase, histogram or atom number (but only for one of these fields) but only when there is more than one matching variable in the input variable list. So if the input is this:
varlist = ['0::AUiso:0', '0::AUiso:1', '1::AUiso:0']
then the output will be this:
wildList = ['*::AUiso:0', '0::AUiso:*']
- Parameters:
varlist (list) – an input list of GSAS-II variable names (such as 0::AUiso:0)
- Returns:
wildList, the generated list of wild card variable names.
- GSASIIobj.GetPhaseNames(fl)[source]
Returns a list of phase names found under ‘Phases’ in GSASII gpx file NB: there is another one of these in GSASIIstrIO.py that uses the gpx filename
- Parameters:
fl (file) – opened .gpx file
- Returns:
list of phase names
- GSASIIobj.HistIdLookup = {}
dict listing histogram name and random Id, keyed by sequential histogram index as a str; best to access this using
LookupHistName()
- GSASIIobj.HistRanIdLookup = {}
dict listing histogram sequential index keyed by histogram random Id; best to access this using
LookupHistId()
- GSASIIobj.HowDidIgetHere(wherecalledonly=False)[source]
Show a traceback with calls that brought us to the current location. Used for debugging.
- Parameters:
wherecalledonly (bool) – When True, the entire calling stack is shown. When False (default), only the 2nd to last stack entry (the routine that called the calling routine is shown.
- class GSASIIobj.ImportBaseclass(formatName, longFormatName=None, extensionlist=[], strictExtension=False)[source]
Defines a base class for the reading of input files (diffraction data, coordinates,…). See Writing a Import Routine for an explanation on how to use a subclass of this class.
- CIFValidator(filepointer)[source]
A
ContentsValidator()
for use to validate CIF files.
- ContentsValidator(filename)[source]
This routine will attempt to determine if the file can be read with the current format. This will typically be overridden with a method that takes a quick scan of [some of] the file contents to do a “sanity” check if the file appears to match the selected format. the file must be opened here with the correct format (binary/text)
- ExtensionValidator(filename)[source]
This methods checks if the file has the correct extension
- Returns:
False if this filename will not be supported by this reader (only when strictExtension is True)
True if the extension matches the list supplied by the reader
None if the reader allows un-registered extensions
- exception ImportException[source]
Defines an Exception that is used when an import routine hits an expected error, usually in .Reader.
Good practice is that the Reader should define a value in self.errors that tells the user some information about what is wrong with their file.
- __weakref__
list of weak references to the object
- __weakref__
list of weak references to the object
- class GSASIIobj.ImportImage(formatName, longFormatName=None, extensionlist=[], strictExtension=False)[source]
Defines a base class for the reading of images
Images are read in only these places:
Initial reading is typically done from a menu item with a call to
GSASIIdataGUI.GSASII.OnImportImage()
which in turn callsGSASIIdataGUI.GSASII.OnImportGeneric()
. That calls methodsExtensionValidator()
,ContentsValidator()
andReader()
. This returns a list of reader objects for each read image. Also used inGSASIIscriptable.import_generic()
.Images are read alternatively in
GSASIImiscGUI.ReadImages()
, which puts image info directly into the data tree.Unlike all other data types read by GSAS-II, images are only kept in memory as they are used and function
GSASIIfiles.GetImageData()
orGSASIIfiles.RereadImageData()
is used to reread images if they are reloaded. For quick retrieval of previously read images, it may be useful to save sums of images or save a keyword (seeImageTag
, below
When reading an image, the
Reader()
routine in the ImportImage class should set:Comments
: a list of strings (str),Npix
: the number of pixels in the image (int),Image
: the actual image as a numpy array (np.array)Data
: a dict defining image parameters (dict). Within this dict the following data items are used:pixelSize
: size of each pixel (x,y) in microns (such as[200.,200.]
.wavelength
: wavelength in \(\AA\).distance
: distance of detector from sample in cm.center
: uncalibrated center of beam on detector (such as[204.8,204.8]
, in mm measured from top left corner of the detectorsize
: size of image in pixels (x,y) (such as[2048,2048]
).ImageTag
: image number or other keyword used to retrieve image from a multi-image data file (defaults to1
if not specified).sumfile
: holds sum image file name if a sum was produced from a multi image filePolaVal
: has two values, the polarization fraction (typically 0.95-0.99 for synchrotrons, 0.5 for lab instruments) and a refinement flag (such as[0.99, False]
).setdist
: nominal distance from sample to detector. Note thatdistance
may be changed during calibration, butsetdist
will not be, so that calibration may be repeated.
optional data items:
repeat
: set to True if there are additional images to read in the file, False otherwiserepeatcount
: set to the number of the image.
Note that the above is initialized with
InitParameters()
. (Also see Writing a Import Routine for an explanation on how to use import classes in general.)
- class GSASIIobj.ImportPDFData(formatName, longFormatName=None, extensionlist=[], strictExtension=False)[source]
Defines a base class for the reading of files with PDF G(R) data. See Writing a Import Routine for an explanation on how to use this class.
- class GSASIIobj.ImportPhase(formatName, longFormatName=None, extensionlist=[], strictExtension=False)[source]
Defines a base class for the reading of files with coordinates
Objects constructed that subclass this (in import/G2phase_*.py etc.) will be used in
GSASIIdataGUI.GSASII.OnImportPhase()
and inGSASIIscriptable.import_generic()
. See Writing a Import Routine for an explanation on how to use this class.
- class GSASIIobj.ImportPowderData(formatName, longFormatName=None, extensionlist=[], strictExtension=False)[source]
Defines a base class for the reading of files with powder data.
Objects constructed that subclass this (in import/G2pwd_*.py etc.) will be used in
GSASIIdataGUI.GSASII.OnImportPowder()
and inGSASIIscriptable.import_generic()
. See Writing a Import Routine for an explanation on how to use this class.
- class GSASIIobj.ImportReflectometryData(formatName, longFormatName=None, extensionlist=[], strictExtension=False)[source]
Defines a base class for the reading of files with reflectometry data. See Writing a Import Routine for an explanation on how to use this class.
- class GSASIIobj.ImportSmallAngleData(formatName, longFormatName=None, extensionlist=[], strictExtension=False)[source]
Defines a base class for the reading of files with small angle data. See Writing a Import Routine for an explanation on how to use this class.
- class GSASIIobj.ImportStructFactor(formatName, longFormatName=None, extensionlist=[], strictExtension=False)[source]
Defines a base class for the reading of files with tables of structure factors.
Structure factors are read with a call to
GSASIIdataGUI.GSASII.OnImportSfact()
which in turn callsGSASIIdataGUI.GSASII.OnImportGeneric()
, which calls methodsExtensionValidator()
,ContentsValidator()
andReader()
.See Writing a Import Routine for an explanation on how to use import classes in general. The specifics for reading a structure factor histogram require that the
Reader()
routine in the import class need to do only a few things: It should loadRefDict
item'RefList'
with the reflection list, and setParameters
with the instrument parameters (initialized withInitParameters()
and set withUpdateParameters()
).- Banks
self.RefDict is a dict containing the reflection information, as read from the file. Item ‘RefList’ contains the reflection information. See the Single Crystal Reflection Data Structure for the contents of each row. Dict element ‘FF’ contains the form factor values for each element type; if this entry is left as initialized (an empty list) it will be initialized as needed later.
- Parameters
self.Parameters is a list with two dicts for data parameter settings
- GSASIIobj.IndexAllIds(Histograms, Phases)[source]
Scan through the used phases & histograms and create an index to the random numbers of phases, histograms and atoms. While doing this, confirm that assigned random numbers are unique – just in case lightning strikes twice in the same place.
Note: this code assumes that the atom random Id (ranId) is the last element each atom record.
This is called when phases & histograms are looked up in these places (only):
GSASIIstrIO.GetUsedHistogramsAndPhases()
(which loads the histograms and phases from a GPX file),GetUsedHistogramsAndPhasesfromTree()
(which does the same thing but from the data tree.)OnFileClose()
(clears out an old project)
Note that globals
PhaseIdLookup
andPhaseRanIdLookup
are also set inAddPhase2Index()
to temporarily assign a phase number as a phase is being imported.TODO: do we need a lookup for rigid body variables?
- GSASIIobj.LookupAtomId(pId, ranId)[source]
Get the atom number from a phase and atom random Id
- Parameters:
pId (int/str) – the sequential number of the phase
ranId (int) – the random Id assigned to an atom
- Returns:
the index number of the atom (str)
- GSASIIobj.LookupAtomLabel(pId, index)[source]
Get the atom label from a phase and atom index number
- Parameters:
pId (int/str) – the sequential number of the phase
index (int) – the index of the atom in the list of atoms
- Returns:
the label for the atom (str) and the random Id of the atom (int)
- GSASIIobj.LookupHistId(ranId)[source]
Get the histogram number and name from a histogram random Id
- Parameters:
ranId (int) – the random Id assigned to a histogram
- Returns:
the sequential Id (hId) number for the histogram (str)
- GSASIIobj.LookupHistName(hId)[source]
Get the histogram number and name from a histogram Id
- Parameters:
hId (int/str) – the sequential assigned to a histogram
- Returns:
(hist,ranId) where hist is the name of the histogram (str) and ranId is the random # id for the histogram (int)
- GSASIIobj.LookupPhaseId(ranId)[source]
Get the phase number and name from a phase random Id
- Parameters:
ranId (int) – the random Id assigned to a phase
- Returns:
the sequential Id (pId) number for the phase (str)
- GSASIIobj.LookupPhaseName(pId)[source]
Get the phase number and name from a phase Id
- Parameters:
pId (int/str) – the sequential assigned to a phase
- Returns:
(phase,ranId) where phase is the name of the phase (str) and ranId is the random # id for the phase (int)
- GSASIIobj.LookupWildCard(varname, varlist)[source]
returns a list of variable names from list varname that match wildcard name in varname
- Parameters:
varname (str) – a G2 variable name containing a wildcard (such as *::var)
varlist (list) – the list of all variable names used in the current project
- Returns:
a list of matching GSAS-II variables (may be empty)
- GSASIIobj.MakeUniqueLabel(lbl, labellist)[source]
Make sure that every a label is unique against a list by adding digits at the end until it is not found in list.
- Parameters:
lbl (str) – the input label
labellist (list) – the labels that have already been encountered
- Returns:
lbl if not found in labellist or lbl with
_1-9
(or_10-99
, etc.) appended at the end
- GSASIIobj.PhaseIdLookup = {}
dict listing phase name and random Id keyed by sequential phase index as a str; best to access this using
LookupPhaseName()
- GSASIIobj.PhaseRanIdLookup = {}
dict listing phase sequential index keyed by phase random Id; best to access this using
LookupPhaseId()
- GSASIIobj.ReadCIF(URLorFile)[source]
Open a CIF, which may be specified as a file name or as a URL using PyCifRW (from James Hester). The open routine gets confused with DOS names that begin with a letter and colon “C:dir” so this routine will try to open the passed name as a file and if that fails, try it as a URL
- Parameters:
URLorFile (str) – string containing a URL or a file name. Code will try first to open it as a file and then as a URL.
- Returns:
a PyCifRW CIF object.
- GSASIIobj.SetDefaultSample()[source]
Fills in default items for the Sample dictionary for Debye-Scherrer & SASD
- GSASIIobj.SetNewPhase(Name='New Phase', SGData=None, cell=None, Super=None)[source]
Create a new phase dict with default values for various parameters
- Parameters:
Name (str) – Name for new Phase
SGData (dict) – space group data from
GSASIIspc:SpcGroup()
; defaults to data for P 1cell (list) – unit cell parameter list; defaults to [1.0,1.0,1.0,90.,90,90.,1.]
- GSASIIobj.ShortHistNames = {}
a dict containing a possibly shortened and when non-unique numbered version of the histogram name. Keyed by the histogram sequential index.
- GSASIIobj.ShortPhaseNames = {}
a dict containing a possibly shortened and when non-unique numbered version of the phase name. Keyed by the phase sequential index.
- class GSASIIobj.ShowTiming[source]
An object to use for timing repeated sections of code.
- Create the object with::
tim0 = ShowTiming()
- Tag sections of code to be timed with::
tim0.start(‘start’) tim0.start(‘in section 1’) tim0.start(‘in section 2’)
etc. (Note that each section should have a unique label.)
- After the last section, end timing with::
tim0.end()
- Show timing results with::
tim0.show()
- __weakref__
list of weak references to the object
- GSASIIobj.StripUnicode(string, subs='.')[source]
Strip non-ASCII characters from strings
- Parameters:
string (str) – string to strip Unicode characters from
subs (str) – character(s) to place into string in place of each Unicode character. Defaults to ‘.’
- Returns:
a new string with only ASCII characters
- GSASIIobj.TestIndexAll()[source]
Test if
IndexAllIds()
has been called to index all phases and histograms (this is needed beforeG2VarObj()
can be used.- Returns:
Returns True if indexing is needed.
- GSASIIobj.VarDescr(varname)[source]
Return two strings with a more complete description for a GSAS-II variable
- Parameters:
name (str) – A full G2 variable name with 2 or 3 or 4 colons (<p>:<h>:name[:<a>] or <p>::RBname:<r>:<t>])
- Returns:
(loc,meaning) where loc describes what item the variable is mapped (phase, histogram, etc.) and meaning describes what the variable does.
- GSASIIobj._lookup(dic, key)[source]
Lookup a key in a dictionary, where None returns an empty string but an unmatched key returns a question mark. Used in
G2VarObj
- GSASIIobj.fmtVarDescr(varname)[source]
Return a string with a more complete description for a GSAS-II variable
- Parameters:
varname (str) – A full G2 variable name with 2 or 3 or 4 colons (<p>:<h>:name[:<a>] or <p>::RBname:<r>:<t>])
- Returns:
a string with the description
- GSASIIobj.getDescr(name)[source]
Return a short description for a GSAS-II variable
- Parameters:
name (str) – The descriptive part of the variable name without colons (:)
- Returns:
a short description or None if not found
- GSASIIobj.getVarDescr(varname)[source]
Return a short description for a GSAS-II variable
- Parameters:
name (str) – A full G2 variable name with 2 or 3 or 4 colons (<p>:<h>:name[:<a1>][:<a2>])
- Returns:
a six element list as [p,`h`,`name`,`a1`,`a2`,`description`], where p, h, a1, a2 are str values or None, for the phase number, the histogram number and the atom number; name will always be a str; and description is str or None. If the variable name is incorrectly formed (for example, wrong number of colons), None is returned instead of a list.
- GSASIIobj.getVarStep(name, parmDict=None)[source]
Return a step size for computing the derivative of a GSAS-II variable
- Parameters:
name (str) – A complete variable name (with colons, :)
parmDict (dict) – A dict with parameter values or None (default)
- Returns:
a float that should be an appropriate step size, either from the value supplied in
CompileVarDesc()
or based on the value for name in parmDict, if supplied. If not found or the value is zero, a default value of 1e-5 is used. If parmDict is None (default) and no value is provided inCompileVarDesc()
, then None is returned.
- GSASIIobj.patchControls(Controls)[source]
patch routine to convert variable names used in parameter limits to G2VarObj objects (See Parameter Limits description.)
- GSASIIobj.prmLookup(name, prmDict)[source]
Looks for a parameter in a min/max dictionary, optionally considering a wild card for histogram or atom number (use of both will never occur at the same time).
- Parameters:
name – a GSAS-II parameter name (str, see
getVarDescr()
andCompileVarDesc()
) or aG2VarObj
object.prmDict (dict) – a min/max dictionary, (parmMinDict or parmMaxDict in Controls) where keys are
G2VarObj
objects.
- Returns:
Two values, (matchname, value), are returned where:
matchname (str) is the
G2VarObj
object corresponding to the actual matched name, which could contain a wildcard even if name does not; andvalue (float) which contains the parameter limit.
- GSASIIobj.reVarDesc = {re.compile('([UVW])$'): 'Gaussian instrument broadening \\1', re.compile('([XYZ])$'): 'Cauchy instrument broadening \\1', re.compile('([XYZ])cos'): 'Cos position wave for \\1', re.compile('([XYZ])max'): 'ZigZag/Block max value for \\1', re.compile('([XYZ])sin'): 'Sin position wave for \\1', re.compile('([abc])$'): 'Lattice parameter, \\1, from Ai and Djk', re.compile('([vV]ol)'): 'Unit cell volume', re.compile('A([0-5])'): 'Reciprocal metric tensor component \\1', re.compile('A([xyz])$'): 'Fractional atomic coordinate, \\1', re.compile('AD\\([0-6],-[0-6]\\)([0-6])'): ' Atomic sp. harm. coeff for orbital, \\1', re.compile('AD\\([0-6],[0-6]\\)([0-6])'): ' Atomic sp. harm. coeff for orbital, \\1', re.compile('AM([xyz])$'): 'Atomic magnetic moment parameter, \\1', re.compile('ANe([01])'): ' Atomic <j0> orbital population for orbital, \\1', re.compile('AU([123][123])'): 'Atomic anisotropic displacement parameter U\\1', re.compile('AUiso'): 'Atomic isotropic displacement parameter', re.compile('Absorption'): 'Absorption coef.', re.compile('Afrac'): 'Atomic site fraction parameter', re.compile('Akappa([0-6])'): ' Atomic orbital softness for orbital, \\1', re.compile('Amul'): 'Atomic site multiplicity value', re.compile('Aspect ratio'): 'Particle aspect ratio', re.compile('B$'): 'Porod prefactor', re.compile('BF mult'): 'Background file multiplier', re.compile('Bab([AU])'): 'Babinet solvent scattering coef. \\1', re.compile('Back$'): 'background parameter', re.compile('Back(.*)'): 'Background term #\\1', re.compile('BkPkgam;(.*)'): 'Background peak #\\1 Cauchy width', re.compile('BkPkint;(.*)'): 'Background peak #\\1 intensity', re.compile('BkPkpos;(.*)'): 'Background peak #\\1 position', re.compile('BkPksig;(.*)'): 'Background peak #\\1 Gaussian width', re.compile('C\\([0-9]*,[0-9]*\\)'): 'spherical harmonics preferred orientation coef.', re.compile('Cutoff'): 'Porod cutoff', re.compile('D([123][123])'): 'Anisotropic strain coef. \\1', re.compile('Dcalc'): 'Calc. d-spacing', re.compile('DebyeA'): 'Debye model amplitude', re.compile('DebyeR'): 'Debye model radius', re.compile('DebyeU'): 'Debye model Uiso', re.compile('Depth'): 'Well depth', re.compile('Diameter'): 'Cylinder/disk diameter', re.compile('Displace([XY])'): 'Debye-Scherrer sample displacement \\1', re.compile('Dist'): 'Interparticle distance', re.compile('Eg$'): 'Secondary type I extinction', re.compile('Ep$'): 'Primary extinction', re.compile('Es$'): 'Secondary type II extinction', re.compile('Extinction'): 'Extinction coef.', re.compile('Fcos'): 'Cos site fraction modulation', re.compile('Flack'): 'Flack parameter', re.compile('FreePrm([123])'): 'User defined measurement parameter \\1', re.compile('Fsin'): 'Sin site fraction modulation', re.compile('Fwid'): 'Crenel function width', re.compile('Fzero'): 'Crenel function offset', re.compile('G$'): 'Guinier prefactor', re.compile('Gonio. radius'): 'Distance from sample to detector, mm', re.compile('I\\(L2\\)\\/I\\(L1\\)'): 'Ka2/Ka1 intensity ratio', re.compile('Lam'): 'Wavelength', re.compile('Layer Disp'): 'Layer displacement along beam', re.compile('LayerDisp'): 'Bragg-Brentano Layer displacement', re.compile('Length'): 'Cylinder length', re.compile('M([XYZ])cos$'): 'Cos mag. moment wave for \\1', re.compile('M([XYZ])sin$'): 'Sin mag. moment wave for \\1', re.compile('MD'): 'March-Dollase coef.', re.compile('Mean'): 'Particle mean radius', re.compile('Mustrain;.*'): 'Microstrain coefficient (delta Q/Q x 10**6)', re.compile('P$'): 'Porod power', re.compile('PDFmag'): 'PDF peak magnitude', re.compile('PDFpos'): 'PDF peak position', re.compile('PDFsig'): 'PDF peak std. dev.', re.compile('PkGam'): 'Bragg peak gamma', re.compile('PkInt'): 'Bragg peak intensity', re.compile('PkPos'): 'Bragg peak position', re.compile('PkSig'): 'Bragg peak sigma', re.compile('Polariz.'): 'Polarization correction', re.compile('Pressure'): 'Pressure level for measurement in MPa', re.compile('RBR([TLS])([123AB][123AB])'): 'Residue rigid body group disp. param.', re.compile('RBRO([aijk])'): 'Residue rigid body orientation parameter \\1', re.compile('RBRP([xyz])'): 'Residue rigid body \\1 position parameter', re.compile('RBRTr;.*'): 'Residue rigid body torsion parameter', re.compile('RBRU'): 'Residue rigid body group Uiso param.', re.compile('RBRf'): 'Residue rigid body site fraction', re.compile('RBSAtNo'): 'Atom number for spinning rigid body', re.compile('RBSO([aijk])'): 'Spinning rigid body orientation parameter \\1', re.compile('RBSP([xyz])'): 'Spinning rigid body \\1 position parameter', re.compile('RBSShC([1-20,1-20])'): 'Spinning rigid body sph. harmonics term', re.compile('RBSShRadius'): 'Spinning rigid body shell radius', re.compile('RBV([TLS])([123AB][123AB])'): 'Residue rigid body group disp. param.', re.compile('RBV.*'): 'Vector rigid body parameter', re.compile('RBVO([aijk])'): 'Vector rigid body orientation parameter \\1', re.compile('RBVP([xyz])'): 'Vector rigid body \\1 position parameter', re.compile('RBVU'): 'Residue rigid body group Uiso param.', re.compile('RBVf'): 'Vector rigid body site fraction', re.compile('Radius'): 'Sphere/cylinder/disk radius', re.compile('Rg$'): 'Guinier radius of gyration', re.compile('SH/L'): 'FCJ peak asymmetry correction', re.compile('Scale'): 'Phase fraction (as p:h:Scale) or Histogram scale factor (as :h:Scale)', re.compile('Shell thickness'): 'Multiplier to get inner(<1) or outer(>1) sphere radius', re.compile('Shift'): 'Bragg-Brentano sample displ.', re.compile('Size;.*'): 'Crystallite size value (in microns)', re.compile('StdDev'): 'Standard deviation in Mean', re.compile('Sticky'): 'Stickyness', re.compile('SurfRoughA'): 'Bragg-Brenano surface roughness A', re.compile('SurfRoughB'): 'Bragg-Brenano surface roughness B', re.compile('Temperature'): 'T value for measurement, K', re.compile('Thickness'): 'Disk thickness', re.compile('Tmax'): 'ZigZag/Block max location', re.compile('Tmin'): 'ZigZag/Block min location', re.compile('Transparency'): 'Bragg-Brentano sample tranparency', re.compile('TwinFr'): 'Twin fraction', re.compile('U([123][123])cos$'): 'Cos thermal wave for U\\1', re.compile('U([123][123])sin$'): 'Sin thermal wave for U\\1', re.compile('VolFr'): 'Dense scatterer volume fraction', re.compile('Volume'): 'Particle volume', re.compile('WgtFrac'): 'phase weight fraction', re.compile('Width'): 'Well width', re.compile('Zero'): 'Debye-Scherrer zero correction', re.compile('alpha'): 'TOF profile term', re.compile('alpha-([01])'): 'Pink profile term', re.compile('beta-([01q])'): 'TOF/Pink profile term', re.compile('constr([0-9]*)'): 'Generated degree of freedom from constraint', re.compile('dA([xyz])$'): 'Refined change to atomic coordinate, \\1', re.compile('dif([ABC])'): 'TOF to d-space calibration', re.compile('e([12][12])'): 'strain tensor e\\1', re.compile('eA$'): 'Cubic mustrain value', re.compile('epis'): 'Sticky sphere epsilon', re.compile('int$'): 'peak intensity', re.compile('mV([0-2])$'): 'Modulation vector component \\1', re.compile('nv-(.+)'): 'New variable assignment with name \\1', re.compile('pos$'): 'peak position', re.compile('sig-([012q])'): 'TOF profile term', re.compile('α'): 'Lattice parameter, α, computed from both Ai and Djk', re.compile('β'): 'Lattice parameter, β, computed from both Ai and Djk', re.compile('γ'): 'Lattice parameter, γ, computed from both Ai and Djk'}
This dictionary lists descriptions for GSAS-II variables where keys are compiled regular expressions that will match the name portion of a parameter name. Initialized in
CompileVarDesc()
.
- GSASIIobj.reVarStep = {re.compile('([UVW])$'): 1e-05, re.compile('([XYZ])$'): 1e-05, re.compile('A([0-5])'): 1e-05, re.compile('AU([123][123])'): 0.0001, re.compile('AUiso'): 0.0001, re.compile('Afrac'): 1e-05, re.compile('Displace([XY])'): 0.1, re.compile('I\\(L2\\)\\/I\\(L1\\)'): 0.001, re.compile('Lam'): 1e-06, re.compile('Polariz.'): 0.001, re.compile('SH/L'): 0.0001, re.compile('dA([xyz])$'): 1e-06}
This dictionary lists the preferred step size for numerical derivative computation w/r to a GSAS-II variable. Keys are compiled regular expressions and values are the step size for that parameter. Initialized in
CompileVarDesc()
.
- GSASIIobj.removeNonRefined(parmList)[source]
Remove items from variable list that are not refined and should not appear as options for constraints
- Parameters:
parmList (list) – a list of strings of form “p:h:VAR:a” where VAR is the variable name
- Returns:
a list after removing variables where VAR matches a entry in local variable NonRefinedList
- GSASIIobj.restraintNames = [['Bond', 'Bonds'], ['Angle', 'Angles'], ['Plane', 'Planes'], ['Chiral', 'Volumes'], ['Torsion', 'Torsions'], ['Rama', 'Ramas'], ['ChemComp', 'Sites'], ['Texture', 'HKLs'], ['Moments', 'Moments'], ['General', 'General']]
Names of restraint keys for the restraint dict and the location of the restraints in each dict