\(\renewcommand\AA{\text{Å}}\)

8. GSAS-II Structure Submodules

These modules are used to perform structure-related (mostly intensity) computations that are needed for refinement and related computations.

8.1. GSASIIstrMain: main structure routine

8.1.1. GSASIIstrMain Classes & Routines

GSASIIstrMain routines, used for refinement are found below.

GSASIIstrMain.AllPrmDerivs(Controls, Histograms, Phases, restraintDict, rigidbodyDict, parmDict, varyList, calcControls, pawleyLookup, symHold, dlg=None)[source]

Computes the derivative of the fitting function (total Chi**2) with respect to every parameter in the parameter dictionary (parmDict) by applying shift below the parameter value as well as above.

Returns:

a dict with the derivatives keyed by variable number. Derivatives are a list with three values: evaluated over v-d to v; v-d to v+d; v to v+d where v is the current value for the variable and d is a small delta value chosen for that variable type.

GSASIIstrMain.BestPlane(PlaneData)[source]

Needs a doc string

GSASIIstrMain.CheckLeBail(Phases)[source]

Check if there is a LeBail extraction in any histogram

Returns:

True if there is at least one LeBail flag turned on, False otherwise

GSASIIstrMain.DisAglTor(DATData)[source]

Needs a doc string

GSASIIstrMain.DoLeBail(GPXfile, dlg=None, cycles=10, refPlotUpdate=None, seqList=None)[source]

Fit LeBail intensities without changes to any other refined parameters. This is a stripped-down version of Refine() that does not perform any refinement cycles

Parameters:
  • GPXfile (str) – G2 .gpx file name

  • dlg (wx.ProgressDialog) – optional progress window to update. Default is None, which means no calls are made to this.

  • cycles (int) – Number of LeBail cycles to perform

  • refPlotUpdate (function) – Optional routine used to plot results. Default is None, which means no calls are made to this.

  • seqList (list) – List of histograms to be processed. Default is None which means that all used histograms in .gpx file are processed.

GSASIIstrMain.DoNoFit(GPXfile, key)[source]

Compute the diffraction pattern with no refinement of parameters.

TODO: At present, this will compute intensities all diffraction patterns in the project, but this likely can be made faster by dropping all the histograms except key from Histograms.

Parameters:
  • GPXfile (str) – G2 .gpx file name

  • key (str) – name of histogram to be computed

Returns:

the computed diffraction pattern for the selected histogram

GSASIIstrMain.PrintDistAngle(DisAglCtls, DisAglData, out=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>)[source]

Print distances and angles

Parameters:
  • DisAglCtls (dict) – contains distance/angle radii usually defined using GSASIIctrlGUI.DisAglDialog()

  • DisAglData (dict) – contains phase data: Items ‘OrigAtoms’ and ‘TargAtoms’ contain the atoms to be used for distance/angle origins and atoms to be used as targets. Item ‘SGData’ has the space group information (see Space Group object)

  • out (file) – file object for output. Defaults to sys.stdout.

GSASIIstrMain.Refine(GPXfile, dlg=None, makeBack=True, refPlotUpdate=None, newLeBail=False, allDerivs=False)[source]

Global refinement – refines to minimize against all histograms. This can be called in one of three ways, from GSASIIdataGUI.GSASII.OnRefine() in an interactive refinement, where dlg will be a wx.ProgressDialog, or non-interactively from GSASIIscriptable.G2Project.refine() or from do_refine(), where dlg will be None.

GSASIIstrMain.RefineCore(Controls, Histograms, Phases, restraintDict, rigidbodyDict, parmDict, varyList, calcControls, pawleyLookup, ifSeq, printFile, dlg, refPlotUpdate=None)[source]

Core optimization routines, shared between SeqRefine and Refine

Returns:

5-tuple of ifOk (bool), Rvals (dict), result, covMatrix, sig

GSASIIstrMain.ReportProblems(result, Rvals, varyList)[source]

Create a message based results from the refinement

GSASIIstrMain.RetDistAngle(DisAglCtls, DisAglData, dlg=None)[source]

Compute and return distances and angles

Parameters:
  • DisAglCtls (dict) – contains distance/angle radii usually defined using GSASIIctrlGUI.DisAglDialog()

  • DisAglData (dict) – contains phase data: Items ‘OrigAtoms’ and ‘TargAtoms’ contain the atoms to be used for distance/angle origins and atoms to be used as targets. Item ‘SGData’ has the space group information (see Space Group object)

Returns:

AtomLabels,DistArray,AngArray where:

AtomLabels is a dict of atom labels, keys are the atom number

DistArray is a dict keyed by the origin atom number where the value is a list of distance entries. The value for each distance is a list containing:

  1. the target atom number (int);

  2. the unit cell offsets added to x,y & z (tuple of int values)

  3. the symmetry operator number (which may be modified to indicate centering and center of symmetry)

  4. an interatomic distance in A (float)

  5. an uncertainty on the distance in A or 0.0 (float)

AngArray is a dict keyed by the origin (central) atom number where the value is a list of angle entries. The value for each angle entry consists of three values:

  1. a distance item reference for one neighbor (int)

  2. a distance item reference for a second neighbor (int)

  3. a angle, uncertainty pair; the s.u. may be zero (tuple of two floats)

The AngArray distance reference items refer directly to the index of the items in the DistArray item for the list of distances for the central atom.

GSASIIstrMain.SeqRefine(GPXfile, dlg, refPlotUpdate=None)[source]

Perform a sequential refinement – cycles through all selected histgrams, one at a time

GSASIIstrMain.do_refine(*args)[source]

Called to run a refinement when this module is executed

GSASIIstrMain.dropOOBvars(varyList, parmDict, sigDict, Controls, parmFrozenList)[source]

Find variables in the parameters dict that are outside the ranges (in parmMinDict and parmMaxDict) and set them to the limits values. Add any such variables into the list of frozen variable (parmFrozenList). Returns a list of newly frozen variables, if any.

GSASIIstrMain.phaseCheck(phaseVary, Phases, histogram)[source]

Removes unused parameters from phase varylist if phase not in histogram for seq refinement removes vars in “Fix FXU” and “FixedSeqVars” here

8.2. GSASIIstrMath - structure math routines

8.2.1. GSASIIstrMath Classes & Routines

GSASIIstrMath routines, used for refinement computations are found below.

GSASIIstrMath.ApplyRBModelDervs(dFdvDict, parmDict, rigidbodyDict, Phase)[source]

Computes rigid body derivatives; there are none for Spin RBs

GSASIIstrMath.ApplyRBModels(parmDict, Phases, rigidbodyDict, Update=False)[source]

Takes RB info from RBModels in Phase and RB data in rigidbodyDict along with current RB values in parmDict & modifies atom contents (fxyz & Uij) of parmDict

GSASIIstrMath.ApplyXYZshifts(parmDict, varyList)[source]

takes atom x,y,z shift and applies it to corresponding atom x,y,z value

Parameters:
  • parmDict (dict) – parameter dictionary

  • varyList (list) – list of variables (not used!)

Returns:

newAtomDict - dictionary of new atomic coordinate names & values; key is parameter shift name

GSASIIstrMath.Dict2Values(parmdict, varylist)[source]

Use before call to leastsq to setup list of values for the parameters in parmdict, as selected by key in varylist

GSASIIstrMath.GetAbsorb(refl, im, hfx, calcControls, parmDict)[source]

Needs a doc string

GSASIIstrMath.GetAbsorbDerv(refl, im, hfx, calcControls, parmDict)[source]

Needs a doc string

GSASIIstrMath.GetAtomFXU(pfx, calcControls, parmDict)[source]

Needs a doc string

GSASIIstrMath.GetAtomSSFXU(pfx, calcControls, parmDict)[source]

Needs a doc string

GSASIIstrMath.GetFobsSq(Histograms, Phases, parmDict, calcControls)[source]

Compute the observed structure factors for Powder histograms and store in reflection array Multiprocessing support added

GSASIIstrMath.GetHStrainShift(refl, im, SGData, phfx, hfx, calcControls, parmDict)[source]

Computes the shifts in peak position due to the Hydrostatic strain (HStrain, Dij terms). This routine is not used anywhere

GSASIIstrMath.GetHStrainShiftDerv(refl, im, SGData, phfx, hfx, calcControls, parmDict)[source]

Computes the derivatives due to the shifts in peak position from Hydrostatic strain (HStrain, Dij terms).

GSASIIstrMath.GetIntensityCorr(refl, im, uniq, G, g, pfx, phfx, hfx, SGData, calcControls, parmDict)[source]

Needs a doc string

GSASIIstrMath.GetIntensityDerv(refl, im, wave, uniq, G, g, pfx, phfx, hfx, SGData, calcControls, parmDict)[source]

Needs a doc string

GSASIIstrMath.GetNewCellParms(parmDict, varyList)[source]

Compute unit cell tensor terms from varied Aij and Dij values. Terms are included in the dict only if Aij or Dij is varied.

GSASIIstrMath.GetPrefOri(uniq, G, g, phfx, hfx, SGData, calcControls, parmDict)[source]

March-Dollase preferred orientation correction

GSASIIstrMath.GetPrefOriDerv(refl, im, uniq, G, g, phfx, hfx, SGData, calcControls, parmDict)[source]

Needs a doc string

GSASIIstrMath.GetPwdrExt(refl, im, pfx, phfx, hfx, calcControls, parmDict)[source]

Needs a doc string

GSASIIstrMath.GetPwdrExtDerv(refl, im, pfx, phfx, hfx, calcControls, parmDict)[source]

Needs a doc string

GSASIIstrMath.GetReflPos(refl, im, wave, A, pfx, hfx, phfx, calcControls, parmDict)[source]

Needs a doc string

GSASIIstrMath.GetReflPosDerv(refl, im, wave, A, pfx, hfx, phfx, calcControls, parmDict)[source]

Needs a doc string

GSASIIstrMath.GetSampleSigGam(refl, im, wave, G, GB, SGData, hfx, phfx, calcControls, parmDict)[source]

Computes the sample-dependent Lorentzian & Gaussian peak width contributions from size & microstrain parameters :param float wave: wavelength for CW data, 2-theta for EDX data

GSASIIstrMath.GetSampleSigGamDerv(refl, im, wave, G, GB, SGData, hfx, phfx, calcControls, parmDict)[source]

Computes the derivatives on sample-dependent Lorentzian & Gaussian peak widths contributions from size & microstrain parameters :param float wave: wavelength for CW data, 2-theta for EDX data

GSASIIstrMath.HessRefine(values, HistoPhases, parmDict, varylist, calcControls, pawleyLookup, dlg)[source]

Loop over histograms and compute derivatives of the fitting model (M) with respect to all parameters. For each histogram, the Jacobian matrix, dMdv, with dimensions (n by m) where n is the number of parameters and m is the number of data points in the histogram. The (n by n) Hessian is computed from each Jacobian and it is returned. This routine is used when refinement derivatives are selected as “analtytic Hessian” in Controls.

Returns:

Vec,Hess where Vec is the least-squares vector and Hess is the Hessian

GSASIIstrMath.MagStructureFactor2(refDict, G, hfx, pfx, SGData, calcControls, parmDict)[source]

Compute neutron magnetic structure factors for all h,k,l for phase puts the result, F^2, in each ref[8] in refList operates on blocks of 100 reflections for speed input:

Parameters:
  • refDict (dict) – where ‘RefList’ list where each ref = h,k,l,it,d,… ‘FF’ dict of form factors - filed in below

  • G (np.array) – reciprocal metric tensor

  • pfx (str) – phase id string

  • SGData (dict) – space group info. dictionary output from SpcGroup

  • calcControls (dict)

  • ParmDict (dict)

Returns:

copy of new refList - used in calculating numerical derivatives

GSASIIstrMath.MagStructureFactorDerv(refDict, G, hfx, pfx, SGData, calcControls, parmDict)[source]

Compute nonmagnetic structure factor derivatives on blocks of reflections in magnetic structures - for powders/nontwins only input:

Parameters:
  • refDict (dict) – where ‘RefList’ list where each ref = h,k,l,it,d,… ‘FF’ dict of form factors - filled in below

  • G (np.array) – reciprocal metric tensor

  • hfx (str) – histogram id string

  • pfx (str) – phase id string

  • SGData (dict) – space group info. dictionary output from SpcGroup

  • calcControls (dict)

  • parmDict (dict)

Returns:

dict dFdvDict: dictionary of derivatives

GSASIIstrMath.MagStructureFactorDerv2(refDict, G, hfx, pfx, SGData, calcControls, parmDict)[source]

Compute magnetic structure factor derivatives numerically - for powders/nontwins only input:

Parameters:
  • refDict (dict) – where ‘RefList’ list where each ref = h,k,l,it,d,… ‘FF’ dict of form factors - filled in below

  • G (np.array) – reciprocal metric tensor

  • hfx (str) – histogram id string

  • pfx (str) – phase id string

  • SGData (dict) – space group info. dictionary output from SpcGroup

  • calcControls (dict)

  • parmDict (dict)

Returns:

dict dFdvDict: dictionary of magnetic derivatives

GSASIIstrMath.MakeSpHarmFF(HKL, Bmat, SHCdict, Tdata, hType, FFtables, ORBtables, BLtables, FF, SQ, ifDeriv=False)[source]

Computes hkl dependent form factors & derivatives from spinning rigid bodies :param array HKL: reflection hkl set to be considered :param array Bmat: inv crystal to Cartesian transfomation matrix :param dict SHCdict: RB spin/deformation data :param array Tdata: atom type info :param str hType: histogram type :param dict FFtables: x-ray form factor tables :param dict ORBtables: x-ray orbital form factor tables :param dict BLtables: neutron scattering lenghts :param array FF: form factors - will be modified by adding the spin/deformation RB spherical harmonics terms :param array SQ: 1/4d^2 for the HKL set :param bool ifDeriv: True if dFF/dcoff to be returned

Returns:

dict dFFdS of derivatives if ifDeriv = True

GSASIIstrMath.SCExtinction(ref, im, phfx, hfx, pfx, calcControls, parmDict, varyList)[source]

Single crystal extinction function; returns extinction & derivative

GSASIIstrMath.SHPOcal(refl, im, g, phfx, hfx, SGData, calcControls, parmDict)[source]

spherical harmonics preferred orientation (cylindrical symmetry only)

GSASIIstrMath.SHPOcalDerv(refl, im, g, phfx, hfx, SGData, calcControls, parmDict)[source]

spherical harmonics preferred orientation derivatives (cylindrical symmetry only)

GSASIIstrMath.SHTXcal(refl, im, g, pfx, hfx, SGData, calcControls, parmDict)[source]

Spherical harmonics texture

GSASIIstrMath.SHTXcalDerv(refl, im, g, pfx, hfx, SGData, calcControls, parmDict)[source]

Spherical harmonics texture derivatives

GSASIIstrMath.SStructureFactor(refDict, G, hfx, pfx, SGData, SSGData, calcControls, parmDict)[source]

Compute super structure factors for all h,k,l,m for phase - no twins puts the result, F^2, in each ref[9] in refList works on blocks of 32 reflections for speed input:

Parameters:
  • refDict (dict) – where ‘RefList’ list where each ref = h,k,l,m,it,d,… ‘FF’ dict of form factors - filed in below

  • G (np.array) – reciprocal metric tensor

  • pfx (str) – phase id string

  • SGData (dict) – space group info. dictionary output from SpcGroup

  • calcControls (dict)

  • ParmDict (dict)

GSASIIstrMath.SStructureFactorDerv(refDict, im, G, hfx, pfx, SGData, SSGData, calcControls, parmDict)[source]

Compute super structure factor derivatives for all h,k,l,m for phase - no twins Only Fourier component are done analytically here input:

Parameters:
  • refDict (dict) – where ‘RefList’ list where each ref = h,k,l,m,it,d,… ‘FF’ dict of form factors - filled in below

  • im (int) – = 1 (could be eliminated)

  • G (np.array) – reciprocal metric tensor

  • hfx (str) – histogram id string

  • pfx (str) – phase id string

  • SGData (dict) – space group info. dictionary output from SpcGroup

  • SSGData (dict) – super space group info.

  • calcControls (dict)

  • ParmDict (dict)

Returns:

dict dFdvDict: dictionary of derivatives

GSASIIstrMath.SStructureFactorDerv2(refDict, im, G, hfx, pfx, SGData, SSGData, calcControls, parmDict)[source]

Compute super structure factor derivatives for all h,k,l,m for phase - no twins input:

Parameters:
  • refDict (dict) – where ‘RefList’ list where each ref = h,k,l,m,it,d,… ‘FF’ dict of form factors - filled in below

  • im (int) – = 1 (could be eliminated)

  • G (np.array) – reciprocal metric tensor

  • hfx (str) – histogram id string

  • pfx (str) – phase id string

  • SGData (dict) – space group info. dictionary output from SpcGroup

  • SSGData (dict) – super space group info.

  • calcControls (dict)

  • ParmDict (dict)

Returns:

dict dFdvDict: dictionary of derivatives

GSASIIstrMath.SStructureFactorDervTw(refDict, im, G, hfx, pfx, SGData, SSGData, calcControls, parmDict)[source]

Needs a doc string

GSASIIstrMath.SStructureFactorTw(refDict, G, hfx, pfx, SGData, SSGData, calcControls, parmDict)[source]

Compute super structure factors for all h,k,l,m for phase - twins only puts the result, F^2, in each ref[8+im] in refList works on blocks of 32 reflections for speed input:

Parameters:
  • refDict (dict) – where ‘RefList’ list where each ref = h,k,l,m,it,d,… ‘FF’ dict of form factors - filed in below

  • G (np.array) – reciprocal metric tensor

  • pfx (str) – phase id string

  • SGData (dict) – space group info. dictionary output from SpcGroup

  • calcControls (dict)

  • ParmDict (dict)

GSASIIstrMath.StructureFactor2(refDict, G, hfx, pfx, SGData, calcControls, parmDict)[source]

Compute structure factors for all h,k,l for phase puts the result, F^2, in each ref[8] in refList operates on blocks of 100 reflections for speed input:

Parameters:
  • refDict (dict) – where ‘RefList’ list where each ref = h,k,l,it,d,… ‘FF’ dict of form factors - filed in below

  • G (np.array) – reciprocal metric tensor

  • pfx (str) – phase id string

  • SGData (dict) – space group info. dictionary output from SpcGroup

  • calcControls (dict)

  • ParmDict (dict)

GSASIIstrMath.StructureFactorDerv2(refDict, G, hfx, pfx, SGData, calcControls, parmDict)[source]

Compute structure factor derivatives on blocks of reflections - for powders/nontwins only faster than StructureFactorDerv - correct for powders/nontwins!! input:

Parameters:
  • refDict (dict) – where ‘RefList’ list where each ref = h,k,l,it,d,… ‘FF’ dict of form factors - filled in below

  • G (np.array) – reciprocal metric tensor

  • hfx (str) – histogram id string

  • pfx (str) – phase id string

  • SGData (dict) – space group info. dictionary output from SpcGroup

  • calcControls (dict)

  • parmDict (dict)

Returns:

dict dFdvDict: dictionary of derivatives

GSASIIstrMath.StructureFactorDervTw2(refDict, G, hfx, pfx, SGData, calcControls, parmDict)[source]

Compute structure factor derivatives on blocks of reflections - for twins only faster than StructureFactorDervTw input:

Parameters:
  • refDict (dict) – where ‘RefList’ list where each ref = h,k,l,it,d,… ‘FF’ dict of form factors - filled in below

  • G (np.array) – reciprocal metric tensor

  • hfx (str) – histogram id string

  • pfx (str) – phase id string

  • SGData (dict) – space group info. dictionary output from SpcGroup

  • calcControls (dict)

  • parmDict (dict)

Returns:

dict dFdvDict: dictionary of derivatives

GSASIIstrMath.Values2Dict(parmdict, varylist, values)[source]

Use after call to leastsq to update the parameter dictionary with values corresponding to keys in varylist

GSASIIstrMath.calcMassFracs(varyList, covMatrix, Phases, hist, hId)[source]

Compute mass fractions and their uncertainties for all phases in a histogram. Computed using the covariance matrix, along with the derivatives for the mass fraction equations.

Parameters:
  • varyList (list) – list of varied parameters

  • covMatrix (np.array) – covariance matrix, order of rows and columns must match varyList

  • Phases (dict) – data structure (from tree or .gpx) with all phase information

  • hist (str) – name of selected histogram

  • hId (int) – number of current histogram

Returns:

valDict,sigDict which contain the mass fraction values and sigmas, keyed by “ph:h:WgtFrac”

GSASIIstrMath.dervHKLF(Histogram, Phase, calcControls, varylist, parmDict, rigidbodyDict)[source]

Loop over reflections in a HKLF histogram and compute derivatives of the fitting model (M) with respect to all parameters. Independent and dependant dM/dp arrays are returned to either dervRefine or HessRefine.

Returns:

GSASIIstrMath.dervRefine(values, HistoPhases, parmDict, varylist, calcControls, pawleyLookup, dlg)[source]

Loop over histograms and compute derivatives of the fitting model (M) with respect to all parameters. Results are returned in a Jacobian matrix (aka design matrix) of dimensions (n by m) where n is the number of parameters and m is the number of data points. This can exceed memory when m gets large. This routine is used when refinement derivatives are selected as “analtytic Jacobian” in Controls.

Returns:

Jacobian numpy.array dMdv for all histograms concatinated

GSASIIstrMath.errRefine(values, HistoPhases, parmDict, varylist, calcControls, pawleyLookup, dlg=None)[source]

Computes the point-by-point discrepancies between every data point in every histogram and the observed value. Used in the Jacobian, Hessian & numeric least-squares to compute function

Returns:

an np array of differences between observed and computed diffraction values.

GSASIIstrMath.getPowderProfile(parmDict, x, varylist, Histogram, Phases, calcControls, pawleyLookup, histogram=None)[source]

Computes the powder pattern for a histogram based on contributions from all used phases

GSASIIstrMath.getPowderProfileDerv(args)[source]

Computes the derivatives of the computed powder pattern with respect to all refined parameters. Used for single processor & Multiprocessor versions

GSASIIstrMath.penaltyDeriv(pNames, pVal, HistoPhases, calcControls, parmDict, varyList)[source]

Compute derivatives on user-supplied and built-in restraint (penalty) functions

where pNames is list of restraint labels

Returns:

array pDerv: partial derivatives by variable# in varList and restraint# in pNames (pDerv[variable#][restraint#])

GSASIIstrMath.penaltyFxn(HistoPhases, calcControls, parmDict, varyList)[source]

Compute user-supplied and built-in restraint functions

8.3. GSASIIstrIO: structure I/O routines

8.3.1. GSASIIstrIO Classes & Routines

GSASIIstrIO routines, used for refinement to read from GPX files and print to the .LST file. Used for refinements and in G2scriptable.

This file should not contain any wxpython references as this must be used in non-GUI settings.

GSASIIstrIO.GPXBackup(GPXfile, makeBack=True)[source]

makes a backup of the specified .gpx file

Parameters:
  • GPXfile (str) – full .gpx file name

  • makeBack (bool) – if True (default), the backup is written to a new file; if False, the last backup is overwritten

Returns:

the name of the backup file that was written

GSASIIstrIO.GetAllPhaseData(GPXfile, PhaseName)[source]

Returns the entire dictionary for PhaseName from GSASII gpx file

Parameters:
  • GPXfile (str) – full .gpx file name

  • PhaseName (str) – phase name

Returns:

phase dictionary or None if PhaseName is not present

GSASIIstrIO.GetControls(GPXfile)[source]

Returns dictionary of control items found in GSASII gpx file

Parameters:

GPXfile (str) – full .gpx file name

Returns:

dictionary of control items

GSASIIstrIO.GetFprime(controlDict, Histograms)[source]

Needs a doc string

GSASIIstrIO.GetFullGPX(GPXfile)[source]

Returns complete contents of GSASII gpx file. Used in GSASIIscriptable.LoadDictFromProjFile().

Parameters:

GPXfile (str) – full .gpx file name

Returns:

Project,nameList, where

  • Project (dict) is a representation of gpx file following the GSAS-II tree structure for each item: key = tree name (e.g. ‘Controls’, ‘Restraints’, etc.), data is dict

  • nameList (list) has names of main tree entries & subentries used to reconstruct project file

GSASIIstrIO.GetHistogramData(Histograms, Print=True, pFile=None)[source]

needs a doc string

GSASIIstrIO.GetHistogramNames(GPXfile, hTypes)[source]

Returns a list of histogram names found in a GSAS-II .gpx file that match specifed histogram types. Names are returned in the order they appear in the file.

Parameters:
  • GPXfile (str) – full .gpx file name

  • hTypes (str) – list of histogram types

Returns:

list of histogram names (types = PWDR & HKLF)

GSASIIstrIO.GetHistogramPhaseData(Phases, Histograms, Controls={}, Print=True, pFile=None, resetRefList=True)[source]

Loads the HAP histogram/phase information into dicts

Parameters:
  • Phases (dict) – phase information

  • Histograms (dict) – Histogram information

  • Print (bool) – prints information as it is read

  • pFile (file) – file object to print to (the default, None causes printing to the console)

  • resetRefList (bool) – Should the contents of the Reflection List be initialized on loading. The default, True, initializes the Reflection List as it is loaded.

Returns:

(hapVary,hapDict,controlDict) * hapVary: list of refined variables * hapDict: dict with refined variables and their values * controlDict: dict with fixed parameters

GSASIIstrIO.GetHistograms(GPXfile, hNames)[source]

Returns a dictionary of histograms found in GSASII gpx file

Parameters:
  • GPXfile (str) – full .gpx file name

  • hNames (str) – list of histogram names

Returns:

dictionary of histograms (types = PWDR & HKLF)

GSASIIstrIO.GetPawleyConstr(SGLaue, PawleyRef, im, pawleyVary)[source]

needs a doc string

GSASIIstrIO.GetPhaseData(PhaseData, RestraintDict={}, rbIds={}, Print=True, pFile=None, seqHistName=None, symHold=None)[source]

Setup the phase information for a structural refinement, used for regular and sequential refinements, optionally printing information to the .lst file (if Print is True). Used as part of refinements but also to generate information on refinement settings. Can be used with dicts from data tree or read from the GPX file. Note that this routine shares a name with routine G2frame.GetPhaseData() (GSASIIdata.GSASII.GetPhaseData()) that instead returns the phase dict(s) from the tree.

Parameters:
  • PhaseData (dict) – the contents of the Phase tree item (may be read from .gpx file) with information on all phases

  • RestraintDict (dict) – an optional dict with restraint information

  • rbIds (dict) – an optional dict with rigid body information

  • Print (bool) – a flag that determines if information will be formatted and printed to the .lst file

  • pFile (file) – a file object (created by open) where print information is sent when Print is True

  • seqHistName (str) – will be None, except for sequential fits. For sequential fits, this can be the name of the current histogram or ‘All’. If a histogram name is supplied, only the phases used in the current histogram are loaded. If ‘All’ is specified, all phases are loaded (used for error checking).

  • symHold (list) – if not None (None is the default) the names of parameters held due to symmetry are placed in this list even if not varied. (Used in G2constrGUI and for parameter impact estimates in AllPrmDerivs).

Returns:

lots of stuff: Natoms,atomIndx,phaseVary,phaseDict,pawleyLookup, FFtables,EFtables,ORBtables,BLtables,MFtables,maxSSwave (see code for details).

GSASIIstrIO.GetPhaseNames(GPXfile)[source]

Returns a list of phase names found under ‘Phases’ in GSASII gpx file

Parameters:

GPXfile (str) – full .gpx file name

Returns:

list of phase names

GSASIIstrIO.GetRestraints(GPXfile)[source]

Read the restraints from the GPX file. Throws an exception if not found in the .GPX file

GSASIIstrIO.GetRigidBodies(GPXfile)[source]

Read the rigid body models from the GPX file

GSASIIstrIO.GetRigidBodyModels(rigidbodyDict, Print=True, pFile=None)[source]

Get Rigid body info from tree entry and print it to .LST file Adds variables and dict items for vector RBs, but for Residue bodies this is done in GetPhaseData().

GSASIIstrIO.GetSeqResult(GPXfile)[source]

Returns the sequential results table information from a GPX file. Called at the beginning of GSASIIstrMain.SeqRefine()

Parameters:

GPXfile (str) – full .gpx file name

Returns:

a dict containing the sequential results table

GSASIIstrIO.GetUsedHistogramsAndPhases(GPXfile)[source]

Returns all histograms that are found in any phase and any phase that uses a histogram. This also assigns numbers to used phases and histograms by the order they appear in the file.

Parameters:

GPXfile (str) – full .gpx file name

Returns:

(Histograms,Phases)

  • Histograms = dictionary of histograms as {name:data,…}

  • Phases = dictionary of phases that use histograms

GSASIIstrIO.IndexGPX(GPXfile, read=False)[source]

Create an index to a GPX file, optionally the file into memory. The byte size of the GPX file is saved. If this routine is called again, and if this size does not change, indexing is not repeated since it is assumed the file has not changed (this can be overriden by setting read=True).

Parameters:

GPXfile (str) – full .gpx file name

Returns:

Project,nameList if read=, where

  • Project (dict) is a representation of gpx file following the GSAS-II tree structure for each item: key = tree name (e.g. ‘Controls’, ‘Restraints’, etc.), data is dict

  • nameList (list) has names of main tree entries & subentries used to reconstruct project file

GSASIIstrIO.PrintISOmodes(pFile, Phases, parmDict, sigDict)[source]

Prints the values for the ISODISTORT modes into the project’s .lst file after a refinement.

GSASIIstrIO.PrintIndependentVars(parmDict, varyList, sigDict, PrintAll=False, pFile=None)[source]

Print the values and uncertainties on the independent parameters

GSASIIstrIO.PrintRestraints(cell, SGData, AtPtrs, Atoms, AtLookup, textureData, phaseRest, pFile)[source]

needs a doc string

GSASIIstrIO.ReadCheckConstraints(GPXfile, seqHist=None, Histograms=None, Phases=None)[source]

Load constraints and related info and return any error or warning messages This is done from the GPX file rather than the tree.

Parameters:
  • GPXfile (str) – specifies the path to a .gpx file.

  • seqHist (str) – specifies a histogram to be loaded for a sequential refinement. If None (default) all are loaded.

  • Histograms (dict) – output from GetUsedHistogramsAndPhases(), can optionally be supplied to save time for sequential refinements

  • Phases (dict) – output from GetUsedHistogramsAndPhases(), can optionally be supplied to save time for sequential refinements

GSASIIstrIO.ReadConstraints(GPXfile, seqHist=None)[source]

Read the constraints from the GPX file and interpret them

called in ReadCheckConstraints(), GSASIIstrMain.Refine() and GSASIIstrMain.SeqRefine().

GSASIIstrIO.SaveUpdatedHistogramsAndPhases(GPXfile, Histograms, Phases, RigidBodies, CovData, parmFrozen)[source]

Save phase and histogram information into “pseudo-gpx” files. The phase information is overwritten each time this is called, but histogram information is appended after each sequential step.

Parameters:
  • GPXfile (str) – full .gpx file name

  • Histograms (dict) – dictionary of histograms as {name:data,…}

  • Phases (dict) – dictionary of phases that use histograms

  • RigidBodies (dict) – dictionary of rigid bodies

  • CovData (dict) – dictionary of refined variables, varyList, & covariance matrix

  • parmFrozen (dict) – dict with frozen parameters for all phases and histograms (specified as str values)

GSASIIstrIO.SetHistogramData(parmDict, sigDict, Histograms, calcControls, Print=True, pFile=None, seq=False)[source]

Shows histogram data after a refinement

GSASIIstrIO.SetHistogramPhaseData(parmDict, sigDict, Phases, Histograms, calcControls, Print=True, pFile=None, covMatrix=[], varyList=[])[source]

Updates parmDict with HAP results from refinement and prints a summary if Print is True

GSASIIstrIO.SetISOmodes(parmDict, sigDict, Phases, pFile=None)[source]

After a refinement, sets the values for the ISODISTORT modes into the parameter and s.u. dicts. Also, in the case of a non-sequential refinement, prints them into the project’s .lst file.

Parameters:
  • parmDict (dict) – parameter dict

  • sigDict (dict) – s.u. (uncertainty) dict

  • Phases (dict) – Phase info from tree/.gpx

  • pFile (file) – file for .lst info or None (None for sequential fits)

GSASIIstrIO.SetPhaseData(parmDict, sigDict, Phases, RBIds, covData, RestraintDict=None, pFile=None)[source]

Called after a refinement to transfer parameters from the parameter dict to the phase(s) information read from a GPX file. Also prints values to the .lst file

GSASIIstrIO.SetRigidBodyModels(parmDict, sigDict, rigidbodyDict, pFile=None)[source]

needs a doc string

GSASIIstrIO.SetSeqResult(GPXfile, Histograms, SeqResult)[source]

Places the sequential results information into a GPX file after a refinement has been completed. Called at the end of GSASIIstrMain.SeqRefine()

Parameters:

GPXfile (str) – full .gpx file name

GSASIIstrIO.SetUsedHistogramsAndPhases(GPXfile, Histograms, Phases, RigidBodies, CovData, parmFrozenList, makeBack=True)[source]

Updates gpxfile from all histograms that are found in any phase and any phase that used a histogram. Also updates rigid body definitions. This is used for non-sequential fits, but not for sequential fitting.

Parameters:
  • GPXfile (str) – full .gpx file name

  • Histograms (dict) – dictionary of histograms as {name:data,…}

  • Phases (dict) – dictionary of phases that use histograms

  • RigidBodies (dict) – dictionary of rigid bodies

  • CovData (dict) – dictionary of refined variables, varyList, & covariance matrix

  • parmFrozenList (list) – list of parameters (as str) that are frozen due to limits; converted to GSASIIobj.G2VarObj objects.

  • makeBack (bool) – True if new backup of .gpx file is to be made; else use the last one made

GSASIIstrIO.SetupSeqSavePhases(GPXfile)[source]

Initialize the files used to save intermediate results from sequential fits.

GSASIIstrIO.ShowBanner(pFile=None)[source]

Print authorship, copyright and citation notice

GSASIIstrIO.ShowControls(Controls, pFile=None, SeqRef=False, preFrozenCount=0)[source]

Print controls information

GSASIIstrIO.WriteRBObjPOAndSig(pfx, rbfx, rbsx, parmDict, sigDict)[source]

Cribbed version of PrintRBObjPOAndSig but returns lists of strings. Moved so it can be used in ExportCIF

GSASIIstrIO.WriteRBObjSHCAndSig(pfx, rbfx, rbsx, parmDict, sigDict, SHC)[source]

Cribbed version of PrintRBObjTorAndSig but returns lists of strings. Moved so it can be used in ExportCIF

GSASIIstrIO.WriteRBObjTLSAndSig(pfx, rbfx, rbsx, TLS, parmDict, sigDict)[source]

Cribbed version of PrintRBObjTLSAndSig but returns lists of strings. Moved so it can be used in ExportCIF

GSASIIstrIO.WriteRBObjTorAndSig(pfx, rbsx, parmDict, sigDict, nTors)[source]

Cribbed version of PrintRBObjTorAndSig but returns lists of strings. Moved so it can be used in ExportCIF

GSASIIstrIO.WriteResRBModel(RBModel)[source]

Write description of a residue rigid body. Code shifted from PrintResRBModel to make usable from G2export_CIF

GSASIIstrIO.WriteVecRBModel(RBModel, sigDict={}, irb=None)[source]

Write description of a vector rigid body. Code shifted from PrintVecRBModel to make usable from G2export_CIF

GSASIIstrIO.cellFill(pfx, SGData, parmDict, sigDict)[source]

Returns the filled-out reciprocal cell (A) terms and their uncertainties from the parameter and sig dictionaries.

Parameters:
  • pfx (str) – parameter prefix (“n::”, where n is a phase number)

  • SGdata (dict) – a symmetry object

  • parmDict (dict) – a dictionary of parameters

  • sigDict (dict) – a dictionary of uncertainties on parameters

Returns:

A,sigA where each is a list of six terms with the A terms

GSASIIstrIO.cellVary(pfx, SGData)[source]

Creates equivalences for a phase based on the Laue class. Returns a list of A tensor terms that are non-zero.

GSASIIstrIO.fmtESD(varname, SigDict, fmtcode, ndig=None, ndec=None)[source]

Format an uncertainty value as requested, but surround the number by () if the parameter is set by an equivalence or by [] if the parameter is set by an constraint

Parameters:

fmtcode (str) – can be a single letter such as ‘g’ or ‘f’, or a format string, such as ‘%10.5f’. For the latter, leave ndig and ndec as None.

GSASIIstrIO.getBackupName(GPXfile, makeBack)[source]

Get the name for the backup .gpx file name

Parameters:
  • GPXfile (str) – full .gpx file name

  • makeBack (bool) – if True the name of a new file is returned, if False the name of the last file that exists is returned

Returns:

the name of a backup file

GSASIIstrIO.getCellEsd(pfx, SGData, A, covData, unique=False)[source]

Compute the standard uncertainty on cell parameters

Parameters:
  • pfx (str) – prefix of form p::

  • SGdata – space group information

  • A (list) – Reciprocal cell Ai terms

  • covData (dict) – covariance tree item

  • unique (bool) – when True, only directly refined parameters (a in cubic, a & alpha in rhombohedral cells) are assigned positive s.u. values. Used for CIF generation.

GSASIIstrIO.getCellSU(pId, hId, SGData, parmDict, covData)[source]

Compute the unit cell parameters and standard uncertainties where lattice parameters and Hstrain (Dij) may be refined. This is called only for generation of CIFs.

Parameters:
  • pId – phase index

  • hId – histogram index

  • SGdata – space group information

  • parmDict (dict) – parameter dict, must have all non-zero Dij and Ai terms

  • covData (dict) – covariance tree item

GSASIIstrIO.gpxSize = -1

Global variables used in IndexGPX() to see if file has changed (gpxSize) and to index where to find each 1st-level tree item in the file.