source: sphinxdocs/build/html/_modules/GSASIIstrMain.html @ 1709

Last change on this file since 1709 was 1709, checked in by toby, 8 years ago

change sphinx docs links; rebuild & fixing minor formatting

  • Property svn:mime-type set to text/html
File size: 161.2 KB
Line 
1<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN"
2  "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
3
4
5<html xmlns="http://www.w3.org/1999/xhtml">
6  <head>
7    <meta http-equiv="Content-Type" content="text/html; charset=utf-8" />
8   
9    <title>GSASIIstrMain &mdash; GSAS-II 0.2.0 documentation</title>
10   
11    <link rel="stylesheet" href="../_static/default.css" type="text/css" />
12    <link rel="stylesheet" href="../_static/pygments.css" type="text/css" />
13   
14    <script type="text/javascript">
15      var DOCUMENTATION_OPTIONS = {
16        URL_ROOT:    '../',
17        VERSION:     '0.2.0',
18        COLLAPSE_INDEX: false,
19        FILE_SUFFIX: '.html',
20        HAS_SOURCE:  true
21      };
22    </script>
23    <script type="text/javascript" src="../_static/jquery.js"></script>
24    <script type="text/javascript" src="../_static/underscore.js"></script>
25    <script type="text/javascript" src="../_static/doctools.js"></script>
26    <script type="text/javascript" src="http://cdn.mathjax.org/mathjax/latest/MathJax.js?config=TeX-AMS-MML_HTMLorMML"></script>
27    <link rel="top" title="GSAS-II 0.2.0 documentation" href="../index.html" />
28    <link rel="up" title="Module code" href="index.html" /> 
29  </head>
30  <body>
31    <div class="related">
32      <h3>Navigation</h3>
33      <ul>
34        <li class="right" style="margin-right: 10px">
35          <a href="../genindex.html" title="General Index"
36             accesskey="I">index</a></li>
37        <li class="right" >
38          <a href="../py-modindex.html" title="Python Module Index"
39             >modules</a> |</li>
40        <li><a href="../index.html">GSAS-II 0.2.0 documentation</a> &raquo;</li>
41          <li><a href="index.html" accesskey="U">Module code</a> &raquo;</li> 
42      </ul>
43    </div> 
44
45    <div class="document">
46      <div class="documentwrapper">
47        <div class="bodywrapper">
48          <div class="body">
49           
50  <h1>Source code for GSASIIstrMain</h1><div class="highlight"><pre>
51<span class="c"># -*- coding: utf-8 -*-</span>
52<span class="sd">&#39;&#39;&#39;</span>
53<span class="sd">*GSASIIstrMain: main structure routine*</span>
54<span class="sd">---------------------------------------</span>
55
56<span class="sd">&#39;&#39;&#39;</span>
57<span class="c">########### SVN repository information ###################</span>
58<span class="c"># $Date: 2015-03-03 15:01:18 -0600 (Tue, 03 Mar 2015) $</span>
59<span class="c"># $Author: toby $</span>
60<span class="c"># $Revision: 1688 $</span>
61<span class="c"># $URL: https://subversion.xray.aps.anl.gov/pyGSAS/trunk/GSASIIstrMain.py $</span>
62<span class="c"># $Id: GSASIIstrMain.py 1688 2015-03-03 21:01:18Z toby $</span>
63<span class="c">########### SVN repository information ###################</span>
64<span class="kn">import</span> <span class="nn">sys</span>
65<span class="kn">import</span> <span class="nn">os</span>
66<span class="kn">import</span> <span class="nn">os.path</span> <span class="kn">as</span> <span class="nn">ospath</span>
67<span class="kn">import</span> <span class="nn">time</span>
68<span class="kn">import</span> <span class="nn">math</span>
69<span class="kn">import</span> <span class="nn">copy</span>
70<span class="kn">import</span> <span class="nn">random</span>
71<span class="kn">import</span> <span class="nn">cPickle</span>
72<span class="kn">import</span> <span class="nn">numpy</span> <span class="kn">as</span> <span class="nn">np</span>
73<span class="kn">import</span> <span class="nn">numpy.ma</span> <span class="kn">as</span> <span class="nn">ma</span>
74<span class="kn">import</span> <span class="nn">numpy.linalg</span> <span class="kn">as</span> <span class="nn">nl</span>
75<span class="kn">import</span> <span class="nn">scipy.optimize</span> <span class="kn">as</span> <span class="nn">so</span>
76<span class="kn">import</span> <span class="nn">GSASIIpath</span>
77<span class="n">GSASIIpath</span><span class="o">.</span><span class="n">SetVersionNumber</span><span class="p">(</span><span class="s">&quot;$Revision: 1688 $&quot;</span><span class="p">)</span>
78<span class="kn">import</span> <span class="nn">GSASIIlattice</span> <span class="kn">as</span> <span class="nn">G2lat</span>
79<span class="kn">import</span> <span class="nn">GSASIIspc</span> <span class="kn">as</span> <span class="nn">G2spc</span>
80<span class="kn">import</span> <span class="nn">GSASIImapvars</span> <span class="kn">as</span> <span class="nn">G2mv</span>
81<span class="kn">import</span> <span class="nn">GSASIImath</span> <span class="kn">as</span> <span class="nn">G2mth</span>
82<span class="kn">import</span> <span class="nn">GSASIIstrIO</span> <span class="kn">as</span> <span class="nn">G2stIO</span>
83<span class="kn">import</span> <span class="nn">GSASIIstrMath</span> <span class="kn">as</span> <span class="nn">G2stMth</span>
84<span class="kn">import</span> <span class="nn">GSASIIobj</span> <span class="kn">as</span> <span class="nn">G2obj</span>
85
86<span class="n">sind</span> <span class="o">=</span> <span class="k">lambda</span> <span class="n">x</span><span class="p">:</span> <span class="n">np</span><span class="o">.</span><span class="n">sin</span><span class="p">(</span><span class="n">x</span><span class="o">*</span><span class="n">np</span><span class="o">.</span><span class="n">pi</span><span class="o">/</span><span class="mf">180.</span><span class="p">)</span>
87<span class="n">cosd</span> <span class="o">=</span> <span class="k">lambda</span> <span class="n">x</span><span class="p">:</span> <span class="n">np</span><span class="o">.</span><span class="n">cos</span><span class="p">(</span><span class="n">x</span><span class="o">*</span><span class="n">np</span><span class="o">.</span><span class="n">pi</span><span class="o">/</span><span class="mf">180.</span><span class="p">)</span>
88<span class="n">tand</span> <span class="o">=</span> <span class="k">lambda</span> <span class="n">x</span><span class="p">:</span> <span class="n">np</span><span class="o">.</span><span class="n">tan</span><span class="p">(</span><span class="n">x</span><span class="o">*</span><span class="n">np</span><span class="o">.</span><span class="n">pi</span><span class="o">/</span><span class="mf">180.</span><span class="p">)</span>
89<span class="n">asind</span> <span class="o">=</span> <span class="k">lambda</span> <span class="n">x</span><span class="p">:</span> <span class="mf">180.</span><span class="o">*</span><span class="n">np</span><span class="o">.</span><span class="n">arcsin</span><span class="p">(</span><span class="n">x</span><span class="p">)</span><span class="o">/</span><span class="n">np</span><span class="o">.</span><span class="n">pi</span>
90<span class="n">acosd</span> <span class="o">=</span> <span class="k">lambda</span> <span class="n">x</span><span class="p">:</span> <span class="mf">180.</span><span class="o">*</span><span class="n">np</span><span class="o">.</span><span class="n">arccos</span><span class="p">(</span><span class="n">x</span><span class="p">)</span><span class="o">/</span><span class="n">np</span><span class="o">.</span><span class="n">pi</span>
91<span class="n">atan2d</span> <span class="o">=</span> <span class="k">lambda</span> <span class="n">y</span><span class="p">,</span><span class="n">x</span><span class="p">:</span> <span class="mf">180.</span><span class="o">*</span><span class="n">np</span><span class="o">.</span><span class="n">arctan2</span><span class="p">(</span><span class="n">y</span><span class="p">,</span><span class="n">x</span><span class="p">)</span><span class="o">/</span><span class="n">np</span><span class="o">.</span><span class="n">pi</span>
92   
93<span class="n">ateln2</span> <span class="o">=</span> <span class="mf">8.0</span><span class="o">*</span><span class="n">math</span><span class="o">.</span><span class="n">log</span><span class="p">(</span><span class="mf">2.0</span><span class="p">)</span>
94<span class="n">DEBUG</span> <span class="o">=</span> <span class="bp">False</span>
95
96<div class="viewcode-block" id="RefineCore"><a class="viewcode-back" href="../GSASIIstruc.html#GSASIIstrMain.RefineCore">[docs]</a><span class="k">def</span> <span class="nf">RefineCore</span><span class="p">(</span><span class="n">Controls</span><span class="p">,</span><span class="n">Histograms</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">restraintDict</span><span class="p">,</span><span class="n">rigidbodyDict</span><span class="p">,</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">,</span>
97    <span class="n">calcControls</span><span class="p">,</span><span class="n">pawleyLookup</span><span class="p">,</span><span class="n">ifPrint</span><span class="p">,</span><span class="n">printFile</span><span class="p">,</span><span class="n">dlg</span><span class="p">):</span>
98    <span class="s">&#39;Core optimization routines, shared between SeqRefine and Refine&#39;</span>
99<span class="c">#    print &#39;current&#39;,varyList</span>
100<span class="c">#    for item in parmDict: print item,parmDict[item] ######### show dict just before refinement</span>
101    <span class="n">G2mv</span><span class="o">.</span><span class="n">Map2Dict</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">)</span>
102    <span class="n">Rvals</span> <span class="o">=</span> <span class="p">{}</span>
103    <span class="k">while</span> <span class="bp">True</span><span class="p">:</span>
104        <span class="n">begin</span> <span class="o">=</span> <span class="n">time</span><span class="o">.</span><span class="n">time</span><span class="p">()</span>
105        <span class="n">values</span> <span class="o">=</span>  <span class="n">np</span><span class="o">.</span><span class="n">array</span><span class="p">(</span><span class="n">G2stMth</span><span class="o">.</span><span class="n">Dict2Values</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span> <span class="n">varyList</span><span class="p">))</span>
106        <span class="c"># test code to compute GOF and save for external repeat</span>
107        <span class="c">#args = ([Histograms,Phases,restraintDict,rigidbodyDict],parmDict,varyList,calcControls,pawleyLookup,dlg)</span>
108        <span class="c">#print &#39;*** before fit chi**2&#39;,np.sum(G2stMth.errRefine(values,*args)**2)            </span>
109        <span class="c">#fl = open(&#39;beforeFit.cpickle&#39;,&#39;wb&#39;)</span>
110        <span class="c">#import cPickle</span>
111        <span class="c">#cPickle.dump(values,fl,1)</span>
112        <span class="c">#cPickle.dump(args[:-1],fl,1)</span>
113        <span class="c">#fl.close()</span>
114        <span class="n">Ftol</span> <span class="o">=</span> <span class="n">Controls</span><span class="p">[</span><span class="s">&#39;min dM/M&#39;</span><span class="p">]</span>
115        <span class="n">Factor</span> <span class="o">=</span> <span class="n">Controls</span><span class="p">[</span><span class="s">&#39;shift factor&#39;</span><span class="p">]</span>
116        <span class="k">if</span> <span class="s">&#39;Jacobian&#39;</span> <span class="ow">in</span> <span class="n">Controls</span><span class="p">[</span><span class="s">&#39;deriv type&#39;</span><span class="p">]:</span>           
117            <span class="n">result</span> <span class="o">=</span> <span class="n">so</span><span class="o">.</span><span class="n">leastsq</span><span class="p">(</span><span class="n">G2stMth</span><span class="o">.</span><span class="n">errRefine</span><span class="p">,</span><span class="n">values</span><span class="p">,</span><span class="n">Dfun</span><span class="o">=</span><span class="n">G2stMth</span><span class="o">.</span><span class="n">dervRefine</span><span class="p">,</span><span class="n">full_output</span><span class="o">=</span><span class="bp">True</span><span class="p">,</span>
118                <span class="n">ftol</span><span class="o">=</span><span class="n">Ftol</span><span class="p">,</span><span class="n">col_deriv</span><span class="o">=</span><span class="bp">True</span><span class="p">,</span><span class="n">factor</span><span class="o">=</span><span class="n">Factor</span><span class="p">,</span>
119                <span class="n">args</span><span class="o">=</span><span class="p">([</span><span class="n">Histograms</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">restraintDict</span><span class="p">,</span><span class="n">rigidbodyDict</span><span class="p">],</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">,</span><span class="n">calcControls</span><span class="p">,</span><span class="n">pawleyLookup</span><span class="p">,</span><span class="n">dlg</span><span class="p">))</span>
120            <span class="n">ncyc</span> <span class="o">=</span> <span class="nb">int</span><span class="p">(</span><span class="n">result</span><span class="p">[</span><span class="mi">2</span><span class="p">][</span><span class="s">&#39;nfev&#39;</span><span class="p">]</span><span class="o">/</span><span class="mi">2</span><span class="p">)</span>
121        <span class="k">elif</span> <span class="s">&#39;Hessian&#39;</span> <span class="ow">in</span> <span class="n">Controls</span><span class="p">[</span><span class="s">&#39;deriv type&#39;</span><span class="p">]:</span>
122            <span class="n">maxCyc</span> <span class="o">=</span> <span class="n">Controls</span><span class="p">[</span><span class="s">&#39;max cyc&#39;</span><span class="p">]</span>
123            <span class="n">result</span> <span class="o">=</span> <span class="n">G2mth</span><span class="o">.</span><span class="n">HessianLSQ</span><span class="p">(</span><span class="n">G2stMth</span><span class="o">.</span><span class="n">errRefine</span><span class="p">,</span><span class="n">values</span><span class="p">,</span><span class="n">Hess</span><span class="o">=</span><span class="n">G2stMth</span><span class="o">.</span><span class="n">HessRefine</span><span class="p">,</span><span class="n">ftol</span><span class="o">=</span><span class="n">Ftol</span><span class="p">,</span><span class="n">maxcyc</span><span class="o">=</span><span class="n">maxCyc</span><span class="p">,</span><span class="n">Print</span><span class="o">=</span><span class="n">ifPrint</span><span class="p">,</span>
124                <span class="n">args</span><span class="o">=</span><span class="p">([</span><span class="n">Histograms</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">restraintDict</span><span class="p">,</span><span class="n">rigidbodyDict</span><span class="p">],</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">,</span><span class="n">calcControls</span><span class="p">,</span><span class="n">pawleyLookup</span><span class="p">,</span><span class="n">dlg</span><span class="p">))</span>
125            <span class="n">ncyc</span> <span class="o">=</span> <span class="n">result</span><span class="p">[</span><span class="mi">2</span><span class="p">][</span><span class="s">&#39;num cyc&#39;</span><span class="p">]</span><span class="o">+</span><span class="mi">1</span>
126            <span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;lamMax&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">result</span><span class="p">[</span><span class="mi">2</span><span class="p">][</span><span class="s">&#39;lamMax&#39;</span><span class="p">]</span>
127        <span class="k">else</span><span class="p">:</span>           <span class="c">#&#39;numeric&#39;</span>
128            <span class="n">result</span> <span class="o">=</span> <span class="n">so</span><span class="o">.</span><span class="n">leastsq</span><span class="p">(</span><span class="n">G2stMth</span><span class="o">.</span><span class="n">errRefine</span><span class="p">,</span><span class="n">values</span><span class="p">,</span><span class="n">full_output</span><span class="o">=</span><span class="bp">True</span><span class="p">,</span><span class="n">ftol</span><span class="o">=</span><span class="n">Ftol</span><span class="p">,</span><span class="n">epsfcn</span><span class="o">=</span><span class="mf">1.e-8</span><span class="p">,</span><span class="n">factor</span><span class="o">=</span><span class="n">Factor</span><span class="p">,</span>
129                <span class="n">args</span><span class="o">=</span><span class="p">([</span><span class="n">Histograms</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">restraintDict</span><span class="p">,</span><span class="n">rigidbodyDict</span><span class="p">],</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">,</span><span class="n">calcControls</span><span class="p">,</span><span class="n">pawleyLookup</span><span class="p">,</span><span class="n">dlg</span><span class="p">))</span>
130            <span class="n">ncyc</span> <span class="o">=</span> <span class="mi">1</span>
131            <span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">varyList</span><span class="p">):</span>
132                <span class="n">ncyc</span> <span class="o">=</span> <span class="nb">int</span><span class="p">(</span><span class="n">result</span><span class="p">[</span><span class="mi">2</span><span class="p">][</span><span class="s">&#39;nfev&#39;</span><span class="p">]</span><span class="o">/</span><span class="nb">len</span><span class="p">(</span><span class="n">varyList</span><span class="p">))</span>
133<span class="c">#        table = dict(zip(varyList,zip(values,result[0],(result[0]-values))))</span>
134<span class="c">#        for item in table: print item,table[item]               #useful debug - are things shifting?</span>
135        <span class="n">runtime</span> <span class="o">=</span> <span class="n">time</span><span class="o">.</span><span class="n">time</span><span class="p">()</span><span class="o">-</span><span class="n">begin</span>
136        <span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;converged&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">result</span><span class="p">[</span><span class="mi">2</span><span class="p">]</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s">&#39;Converged&#39;</span><span class="p">)</span>
137        <span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;DelChi2&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">result</span><span class="p">[</span><span class="mi">2</span><span class="p">]</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s">&#39;DelChi2&#39;</span><span class="p">,</span><span class="o">-</span><span class="mf">1.</span><span class="p">)</span>
138        <span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;chisq&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">sum</span><span class="p">(</span><span class="n">result</span><span class="p">[</span><span class="mi">2</span><span class="p">][</span><span class="s">&#39;fvec&#39;</span><span class="p">]</span><span class="o">**</span><span class="mi">2</span><span class="p">)</span>
139        <span class="n">G2stMth</span><span class="o">.</span><span class="n">Values2Dict</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span> <span class="n">varyList</span><span class="p">,</span> <span class="n">result</span><span class="p">[</span><span class="mi">0</span><span class="p">])</span>
140        <span class="n">G2mv</span><span class="o">.</span><span class="n">Dict2Map</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">)</span>
141        <span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;Nobs&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">Histograms</span><span class="p">[</span><span class="s">&#39;Nobs&#39;</span><span class="p">]</span>
142        <span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;Rwp&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">sqrt</span><span class="p">(</span><span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;chisq&#39;</span><span class="p">]</span><span class="o">/</span><span class="n">Histograms</span><span class="p">[</span><span class="s">&#39;sumwYo&#39;</span><span class="p">])</span><span class="o">*</span><span class="mf">100.</span>      <span class="c">#to %</span>
143        <span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;GOF&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">sqrt</span><span class="p">(</span><span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;chisq&#39;</span><span class="p">]</span><span class="o">/</span><span class="p">(</span><span class="n">Histograms</span><span class="p">[</span><span class="s">&#39;Nobs&#39;</span><span class="p">]</span><span class="o">-</span><span class="nb">len</span><span class="p">(</span><span class="n">varyList</span><span class="p">)))</span>
144        <span class="k">print</span> <span class="o">&gt;&gt;</span><span class="n">printFile</span><span class="p">,</span><span class="s">&#39; Number of function calls:&#39;</span><span class="p">,</span><span class="n">result</span><span class="p">[</span><span class="mi">2</span><span class="p">][</span><span class="s">&#39;nfev&#39;</span><span class="p">],</span><span class="s">&#39; Number of observations: &#39;</span><span class="p">,</span><span class="n">Histograms</span><span class="p">[</span><span class="s">&#39;Nobs&#39;</span><span class="p">],</span><span class="s">&#39; Number of parameters: &#39;</span><span class="p">,</span><span class="nb">len</span><span class="p">(</span><span class="n">varyList</span><span class="p">)</span>
145        <span class="k">print</span> <span class="o">&gt;&gt;</span><span class="n">printFile</span><span class="p">,</span><span class="s">&#39; Refinement time = </span><span class="si">%8.3f</span><span class="s">s, </span><span class="si">%8.3f</span><span class="s">s/cycle, for </span><span class="si">%d</span><span class="s"> cycles&#39;</span><span class="o">%</span><span class="p">(</span><span class="n">runtime</span><span class="p">,</span><span class="n">runtime</span><span class="o">/</span><span class="n">ncyc</span><span class="p">,</span><span class="n">ncyc</span><span class="p">)</span>
146        <span class="k">print</span> <span class="o">&gt;&gt;</span><span class="n">printFile</span><span class="p">,</span><span class="s">&#39; wR = </span><span class="si">%7.2f%%</span><span class="s">, chi**2 = </span><span class="si">%12.6g</span><span class="s">, reduced chi**2 = </span><span class="si">%6.2f</span><span class="s">&#39;</span><span class="o">%</span><span class="p">(</span><span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;Rwp&#39;</span><span class="p">],</span><span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;chisq&#39;</span><span class="p">],</span><span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;GOF&#39;</span><span class="p">]</span><span class="o">**</span><span class="mi">2</span><span class="p">)</span>
147        <span class="n">IfOK</span> <span class="o">=</span> <span class="bp">True</span>
148        <span class="k">try</span><span class="p">:</span>
149            <span class="n">covMatrix</span> <span class="o">=</span> <span class="n">result</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span><span class="o">*</span><span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;GOF&#39;</span><span class="p">]</span><span class="o">**</span><span class="mi">2</span>
150            <span class="n">sig</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">sqrt</span><span class="p">(</span><span class="n">np</span><span class="o">.</span><span class="n">diag</span><span class="p">(</span><span class="n">covMatrix</span><span class="p">))</span>
151            <span class="k">if</span> <span class="n">np</span><span class="o">.</span><span class="n">any</span><span class="p">(</span><span class="n">np</span><span class="o">.</span><span class="n">isnan</span><span class="p">(</span><span class="n">sig</span><span class="p">)):</span>
152                <span class="k">print</span> <span class="s">&#39;*** Least squares aborted - some invalid esds possible ***&#39;</span>
153<span class="c">#            table = dict(zip(varyList,zip(values,result[0],(result[0]-values)/sig)))</span>
154<span class="c">#            for item in table: print item,table[item]               #useful debug - are things shifting?</span>
155            <span class="k">break</span>                   <span class="c">#refinement succeeded - finish up!</span>
156        <span class="k">except</span> <span class="ne">TypeError</span><span class="p">,</span><span class="ne">FloatingPointError</span><span class="p">:</span>          <span class="c">#result[1] is None on singular matrix</span>
157            <span class="n">IfOK</span> <span class="o">=</span> <span class="bp">False</span>
158            <span class="k">if</span> <span class="ow">not</span> <span class="nb">len</span><span class="p">(</span><span class="n">varyList</span><span class="p">):</span>
159                <span class="n">covMatrix</span> <span class="o">=</span> <span class="p">[]</span>
160                <span class="n">sig</span> <span class="o">=</span> <span class="p">[]</span>
161                <span class="k">break</span>
162            <span class="k">print</span> <span class="s">&#39;**** Refinement failed - singular matrix ****&#39;</span>
163            <span class="k">if</span> <span class="s">&#39;Hessian&#39;</span> <span class="ow">in</span> <span class="n">Controls</span><span class="p">[</span><span class="s">&#39;deriv type&#39;</span><span class="p">]:</span>
164                <span class="n">num</span> <span class="o">=</span> <span class="nb">len</span><span class="p">(</span><span class="n">varyList</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span>
165                <span class="k">for</span> <span class="n">i</span><span class="p">,</span><span class="n">val</span> <span class="ow">in</span> <span class="nb">enumerate</span><span class="p">(</span><span class="n">np</span><span class="o">.</span><span class="n">flipud</span><span class="p">(</span><span class="n">result</span><span class="p">[</span><span class="mi">2</span><span class="p">][</span><span class="s">&#39;psing&#39;</span><span class="p">])):</span>
166                    <span class="k">if</span> <span class="n">val</span><span class="p">:</span>
167                        <span class="k">print</span> <span class="s">&#39;Removing parameter: &#39;</span><span class="p">,</span><span class="n">varyList</span><span class="p">[</span><span class="n">num</span><span class="o">-</span><span class="n">i</span><span class="p">]</span>
168                        <span class="k">del</span><span class="p">(</span><span class="n">varyList</span><span class="p">[</span><span class="n">num</span><span class="o">-</span><span class="n">i</span><span class="p">])</span>                   
169            <span class="k">else</span><span class="p">:</span>
170                <span class="n">Ipvt</span> <span class="o">=</span> <span class="n">result</span><span class="p">[</span><span class="mi">2</span><span class="p">][</span><span class="s">&#39;ipvt&#39;</span><span class="p">]</span>
171                <span class="k">for</span> <span class="n">i</span><span class="p">,</span><span class="n">ipvt</span> <span class="ow">in</span> <span class="nb">enumerate</span><span class="p">(</span><span class="n">Ipvt</span><span class="p">):</span>
172                    <span class="k">if</span> <span class="ow">not</span> <span class="n">np</span><span class="o">.</span><span class="n">sum</span><span class="p">(</span><span class="n">result</span><span class="p">[</span><span class="mi">2</span><span class="p">][</span><span class="s">&#39;fjac&#39;</span><span class="p">],</span><span class="n">axis</span><span class="o">=</span><span class="mi">1</span><span class="p">)[</span><span class="n">i</span><span class="p">]:</span>
173                        <span class="k">print</span> <span class="s">&#39;Removing parameter: &#39;</span><span class="p">,</span><span class="n">varyList</span><span class="p">[</span><span class="n">ipvt</span><span class="o">-</span><span class="mi">1</span><span class="p">]</span>
174                        <span class="k">del</span><span class="p">(</span><span class="n">varyList</span><span class="p">[</span><span class="n">ipvt</span><span class="o">-</span><span class="mi">1</span><span class="p">])</span>
175                        <span class="k">break</span>
176    <span class="n">G2stMth</span><span class="o">.</span><span class="n">GetFobsSq</span><span class="p">(</span><span class="n">Histograms</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">parmDict</span><span class="p">,</span><span class="n">calcControls</span><span class="p">)</span>
177    <span class="k">return</span> <span class="n">IfOK</span><span class="p">,</span><span class="n">Rvals</span><span class="p">,</span><span class="n">result</span><span class="p">,</span><span class="n">covMatrix</span><span class="p">,</span><span class="n">sig</span>
178</div>
179<div class="viewcode-block" id="Refine"><a class="viewcode-back" href="../GSASIIstruc.html#GSASIIstrMain.Refine">[docs]</a><span class="k">def</span> <span class="nf">Refine</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">,</span><span class="n">dlg</span><span class="p">):</span>
180    <span class="s">&#39;Global refinement -- refines to minimize against all histograms&#39;</span>
181    <span class="kn">import</span> <span class="nn">pytexture</span> <span class="kn">as</span> <span class="nn">ptx</span>
182    <span class="n">ptx</span><span class="o">.</span><span class="n">pyqlmninit</span><span class="p">()</span>            <span class="c">#initialize fortran arrays for spherical harmonics</span>
183   
184    <span class="n">printFile</span> <span class="o">=</span> <span class="nb">open</span><span class="p">(</span><span class="n">ospath</span><span class="o">.</span><span class="n">splitext</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)[</span><span class="mi">0</span><span class="p">]</span><span class="o">+</span><span class="s">&#39;.lst&#39;</span><span class="p">,</span><span class="s">&#39;w&#39;</span><span class="p">)</span>
185    <span class="n">G2stIO</span><span class="o">.</span><span class="n">ShowBanner</span><span class="p">(</span><span class="n">printFile</span><span class="p">)</span>
186    <span class="n">varyList</span> <span class="o">=</span> <span class="p">[]</span>
187    <span class="n">parmDict</span> <span class="o">=</span> <span class="p">{}</span>
188    <span class="n">G2mv</span><span class="o">.</span><span class="n">InitVars</span><span class="p">()</span>   
189    <span class="n">Controls</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetControls</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)</span>
190    <span class="n">G2stIO</span><span class="o">.</span><span class="n">ShowControls</span><span class="p">(</span><span class="n">Controls</span><span class="p">,</span><span class="n">printFile</span><span class="p">)</span>
191    <span class="n">calcControls</span> <span class="o">=</span> <span class="p">{}</span>
192    <span class="n">calcControls</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">Controls</span><span class="p">)</span>           
193    <span class="n">constrDict</span><span class="p">,</span><span class="n">fixedList</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetConstraints</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)</span>
194    <span class="n">restraintDict</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetRestraints</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)</span>
195    <span class="n">Histograms</span><span class="p">,</span><span class="n">Phases</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetUsedHistogramsAndPhases</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)</span>
196    <span class="k">if</span> <span class="ow">not</span> <span class="n">Phases</span><span class="p">:</span>
197        <span class="k">print</span> <span class="s">&#39; *** ERROR - you have no phases! ***&#39;</span>
198        <span class="k">print</span> <span class="s">&#39; *** Refine aborted ***&#39;</span>
199        <span class="k">raise</span> <span class="ne">Exception</span>
200    <span class="k">if</span> <span class="ow">not</span> <span class="n">Histograms</span><span class="p">:</span>
201        <span class="k">print</span> <span class="s">&#39; *** ERROR - you have no data to refine with! ***&#39;</span>
202        <span class="k">print</span> <span class="s">&#39; *** Refine aborted ***&#39;</span>
203        <span class="k">raise</span> <span class="ne">Exception</span>       
204    <span class="n">rigidbodyDict</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetRigidBodies</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)</span>
205    <span class="n">rbIds</span> <span class="o">=</span> <span class="n">rigidbodyDict</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s">&#39;RBIds&#39;</span><span class="p">,{</span><span class="s">&#39;Vector&#39;</span><span class="p">:[],</span><span class="s">&#39;Residue&#39;</span><span class="p">:[]})</span>
206    <span class="n">rbVary</span><span class="p">,</span><span class="n">rbDict</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetRigidBodyModels</span><span class="p">(</span><span class="n">rigidbodyDict</span><span class="p">,</span><span class="n">pFile</span><span class="o">=</span><span class="n">printFile</span><span class="p">)</span>
207    <span class="n">Natoms</span><span class="p">,</span><span class="n">atomIndx</span><span class="p">,</span><span class="n">phaseVary</span><span class="p">,</span><span class="n">phaseDict</span><span class="p">,</span><span class="n">pawleyLookup</span><span class="p">,</span><span class="n">FFtables</span><span class="p">,</span><span class="n">BLtables</span><span class="p">,</span><span class="n">maxSSwave</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetPhaseData</span><span class="p">(</span><span class="n">Phases</span><span class="p">,</span><span class="n">restraintDict</span><span class="p">,</span><span class="n">rbIds</span><span class="p">,</span><span class="n">pFile</span><span class="o">=</span><span class="n">printFile</span><span class="p">)</span>
208    <span class="n">calcControls</span><span class="p">[</span><span class="s">&#39;atomIndx&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">atomIndx</span>
209    <span class="n">calcControls</span><span class="p">[</span><span class="s">&#39;Natoms&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">Natoms</span>
210    <span class="n">calcControls</span><span class="p">[</span><span class="s">&#39;FFtables&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">FFtables</span>
211    <span class="n">calcControls</span><span class="p">[</span><span class="s">&#39;BLtables&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">BLtables</span>
212    <span class="n">calcControls</span><span class="p">[</span><span class="s">&#39;maxSSwave&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">maxSSwave</span>
213    <span class="n">hapVary</span><span class="p">,</span><span class="n">hapDict</span><span class="p">,</span><span class="n">controlDict</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetHistogramPhaseData</span><span class="p">(</span><span class="n">Phases</span><span class="p">,</span><span class="n">Histograms</span><span class="p">,</span><span class="n">pFile</span><span class="o">=</span><span class="n">printFile</span><span class="p">)</span>
214    <span class="n">calcControls</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">controlDict</span><span class="p">)</span>
215    <span class="n">histVary</span><span class="p">,</span><span class="n">histDict</span><span class="p">,</span><span class="n">controlDict</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetHistogramData</span><span class="p">(</span><span class="n">Histograms</span><span class="p">,</span><span class="n">pFile</span><span class="o">=</span><span class="n">printFile</span><span class="p">)</span>
216    <span class="n">calcControls</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">controlDict</span><span class="p">)</span>
217    <span class="n">varyList</span> <span class="o">=</span> <span class="n">rbVary</span><span class="o">+</span><span class="n">phaseVary</span><span class="o">+</span><span class="n">hapVary</span><span class="o">+</span><span class="n">histVary</span>
218    <span class="n">parmDict</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">rbDict</span><span class="p">)</span>
219    <span class="n">parmDict</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">phaseDict</span><span class="p">)</span>
220    <span class="n">parmDict</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">hapDict</span><span class="p">)</span>
221    <span class="n">parmDict</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">histDict</span><span class="p">)</span>
222    <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetFprime</span><span class="p">(</span><span class="n">calcControls</span><span class="p">,</span><span class="n">Histograms</span><span class="p">)</span>
223    <span class="c"># do constraint processing</span>
224    <span class="n">varyListStart</span> <span class="o">=</span> <span class="nb">tuple</span><span class="p">(</span><span class="n">varyList</span><span class="p">)</span> <span class="c"># save the original varyList before dependent vars are removed</span>
225    <span class="k">try</span><span class="p">:</span>
226        <span class="n">groups</span><span class="p">,</span><span class="n">parmlist</span> <span class="o">=</span> <span class="n">G2mv</span><span class="o">.</span><span class="n">GroupConstraints</span><span class="p">(</span><span class="n">constrDict</span><span class="p">)</span>
227        <span class="n">G2mv</span><span class="o">.</span><span class="n">GenerateConstraints</span><span class="p">(</span><span class="n">groups</span><span class="p">,</span><span class="n">parmlist</span><span class="p">,</span><span class="n">varyList</span><span class="p">,</span><span class="n">constrDict</span><span class="p">,</span><span class="n">fixedList</span><span class="p">,</span><span class="n">parmDict</span><span class="p">)</span>
228    <span class="k">except</span><span class="p">:</span>
229        <span class="k">print</span> <span class="s">&#39; *** ERROR - your constraints are internally inconsistent ***&#39;</span>
230        <span class="c">#errmsg, warnmsg = G2mv.CheckConstraints(varyList,constrDict,fixedList)</span>
231        <span class="c">#print &#39;Errors&#39;,errmsg</span>
232        <span class="c">#if warnmsg: print &#39;Warnings&#39;,warnmsg</span>
233        <span class="k">raise</span> <span class="ne">Exception</span><span class="p">(</span><span class="s">&#39; *** Refine aborted ***&#39;</span><span class="p">)</span>
234<span class="c">#    print G2mv.VarRemapShow(varyList)</span>
235   
236    <span class="n">ifPrint</span> <span class="o">=</span> <span class="bp">True</span>
237    <span class="k">print</span> <span class="o">&gt;&gt;</span><span class="n">printFile</span><span class="p">,</span><span class="s">&#39;</span><span class="se">\n</span><span class="s"> Refinement results:&#39;</span>
238    <span class="k">print</span> <span class="o">&gt;&gt;</span><span class="n">printFile</span><span class="p">,</span><span class="mi">135</span><span class="o">*</span><span class="s">&#39;-&#39;</span>
239    <span class="n">IfOK</span><span class="p">,</span><span class="n">Rvals</span><span class="p">,</span><span class="n">result</span><span class="p">,</span><span class="n">covMatrix</span><span class="p">,</span><span class="n">sig</span> <span class="o">=</span> <span class="n">RefineCore</span><span class="p">(</span><span class="n">Controls</span><span class="p">,</span><span class="n">Histograms</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">restraintDict</span><span class="p">,</span>
240        <span class="n">rigidbodyDict</span><span class="p">,</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">,</span><span class="n">calcControls</span><span class="p">,</span><span class="n">pawleyLookup</span><span class="p">,</span><span class="n">ifPrint</span><span class="p">,</span><span class="n">printFile</span><span class="p">,</span><span class="n">dlg</span><span class="p">)</span>
241    <span class="n">sigDict</span> <span class="o">=</span> <span class="nb">dict</span><span class="p">(</span><span class="nb">zip</span><span class="p">(</span><span class="n">varyList</span><span class="p">,</span><span class="n">sig</span><span class="p">))</span>
242    <span class="n">newCellDict</span> <span class="o">=</span> <span class="n">G2stMth</span><span class="o">.</span><span class="n">GetNewCellParms</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">)</span>
243    <span class="n">newAtomDict</span> <span class="o">=</span> <span class="n">G2stMth</span><span class="o">.</span><span class="n">ApplyXYZshifts</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">)</span>
244    <span class="n">covData</span> <span class="o">=</span> <span class="p">{</span><span class="s">&#39;variables&#39;</span><span class="p">:</span><span class="n">result</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span><span class="s">&#39;varyList&#39;</span><span class="p">:</span><span class="n">varyList</span><span class="p">,</span><span class="s">&#39;sig&#39;</span><span class="p">:</span><span class="n">sig</span><span class="p">,</span><span class="s">&#39;Rvals&#39;</span><span class="p">:</span><span class="n">Rvals</span><span class="p">,</span>
245               <span class="s">&#39;varyListStart&#39;</span><span class="p">:</span><span class="n">varyListStart</span><span class="p">,</span>
246               <span class="s">&#39;covMatrix&#39;</span><span class="p">:</span><span class="n">covMatrix</span><span class="p">,</span><span class="s">&#39;title&#39;</span><span class="p">:</span><span class="n">GPXfile</span><span class="p">,</span><span class="s">&#39;newAtomDict&#39;</span><span class="p">:</span><span class="n">newAtomDict</span><span class="p">,</span>
247               <span class="s">&#39;newCellDict&#39;</span><span class="p">:</span><span class="n">newCellDict</span><span class="p">,</span><span class="s">&#39;freshCOV&#39;</span><span class="p">:</span><span class="bp">True</span><span class="p">}</span>
248    <span class="c"># add the uncertainties into the esd dictionary (sigDict)</span>
249    <span class="n">sigDict</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">G2mv</span><span class="o">.</span><span class="n">ComputeDepESD</span><span class="p">(</span><span class="n">covMatrix</span><span class="p">,</span><span class="n">varyList</span><span class="p">,</span><span class="n">parmDict</span><span class="p">))</span>
250    <span class="n">G2mv</span><span class="o">.</span><span class="n">PrintIndependentVars</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">,</span><span class="n">sigDict</span><span class="p">,</span><span class="n">pFile</span><span class="o">=</span><span class="n">printFile</span><span class="p">)</span>
251    <span class="n">G2stMth</span><span class="o">.</span><span class="n">ApplyRBModels</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">rigidbodyDict</span><span class="p">,</span><span class="bp">True</span><span class="p">)</span>
252    <span class="n">G2stIO</span><span class="o">.</span><span class="n">SetRigidBodyModels</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">sigDict</span><span class="p">,</span><span class="n">rigidbodyDict</span><span class="p">,</span><span class="n">printFile</span><span class="p">)</span>
253    <span class="n">G2stIO</span><span class="o">.</span><span class="n">SetPhaseData</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">sigDict</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">rbIds</span><span class="p">,</span><span class="n">covData</span><span class="p">,</span><span class="n">restraintDict</span><span class="p">,</span><span class="n">printFile</span><span class="p">)</span>
254    <span class="n">G2stIO</span><span class="o">.</span><span class="n">SetHistogramPhaseData</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">sigDict</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">Histograms</span><span class="p">,</span><span class="n">pFile</span><span class="o">=</span><span class="n">printFile</span><span class="p">)</span>
255    <span class="n">G2stIO</span><span class="o">.</span><span class="n">SetHistogramData</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">sigDict</span><span class="p">,</span><span class="n">Histograms</span><span class="p">,</span><span class="n">pFile</span><span class="o">=</span><span class="n">printFile</span><span class="p">)</span>
256    <span class="n">G2stIO</span><span class="o">.</span><span class="n">SetUsedHistogramsAndPhases</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">,</span><span class="n">Histograms</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">rigidbodyDict</span><span class="p">,</span><span class="n">covData</span><span class="p">)</span>
257    <span class="n">printFile</span><span class="o">.</span><span class="n">close</span><span class="p">()</span>
258    <span class="k">print</span> <span class="s">&#39; Refinement results are in file: &#39;</span><span class="o">+</span><span class="n">ospath</span><span class="o">.</span><span class="n">splitext</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)[</span><span class="mi">0</span><span class="p">]</span><span class="o">+</span><span class="s">&#39;.lst&#39;</span>
259    <span class="k">print</span> <span class="s">&#39; ***** Refinement successful *****&#39;</span>
260   
261<span class="c">#for testing purposes!!!</span>
262    <span class="k">if</span> <span class="n">DEBUG</span><span class="p">:</span>
263<span class="c">#needs: values,HistoPhases,parmDict,varylist,calcControls,pawleyLookup</span>
264        <span class="kn">import</span> <span class="nn">cPickle</span>
265        <span class="n">fl</span> <span class="o">=</span> <span class="nb">open</span><span class="p">(</span><span class="s">&#39;testDeriv.dat&#39;</span><span class="p">,</span><span class="s">&#39;wb&#39;</span><span class="p">)</span>
266        <span class="n">cPickle</span><span class="o">.</span><span class="n">dump</span><span class="p">(</span><span class="n">result</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span><span class="n">fl</span><span class="p">,</span><span class="mi">1</span><span class="p">)</span>
267        <span class="n">cPickle</span><span class="o">.</span><span class="n">dump</span><span class="p">([</span><span class="n">Histograms</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">restraintDict</span><span class="p">,</span><span class="n">rigidbodyDict</span><span class="p">],</span><span class="n">fl</span><span class="p">,</span><span class="mi">1</span><span class="p">)</span>
268        <span class="n">cPickle</span><span class="o">.</span><span class="n">dump</span><span class="p">([</span><span class="n">G2mv</span><span class="o">.</span><span class="n">dependentParmList</span><span class="p">,</span><span class="n">G2mv</span><span class="o">.</span><span class="n">arrayList</span><span class="p">,</span><span class="n">G2mv</span><span class="o">.</span><span class="n">invarrayList</span><span class="p">,</span>
269            <span class="n">G2mv</span><span class="o">.</span><span class="n">indParmList</span><span class="p">,</span><span class="n">G2mv</span><span class="o">.</span><span class="n">invarrayList</span><span class="p">],</span><span class="n">fl</span><span class="p">,</span><span class="mi">1</span><span class="p">)</span>
270        <span class="n">cPickle</span><span class="o">.</span><span class="n">dump</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">fl</span><span class="p">,</span><span class="mi">1</span><span class="p">)</span>
271        <span class="n">cPickle</span><span class="o">.</span><span class="n">dump</span><span class="p">(</span><span class="n">varyList</span><span class="p">,</span><span class="n">fl</span><span class="p">,</span><span class="mi">1</span><span class="p">)</span>
272        <span class="n">cPickle</span><span class="o">.</span><span class="n">dump</span><span class="p">(</span><span class="n">calcControls</span><span class="p">,</span><span class="n">fl</span><span class="p">,</span><span class="mi">1</span><span class="p">)</span>
273        <span class="n">cPickle</span><span class="o">.</span><span class="n">dump</span><span class="p">(</span><span class="n">pawleyLookup</span><span class="p">,</span><span class="n">fl</span><span class="p">,</span><span class="mi">1</span><span class="p">)</span>
274        <span class="n">fl</span><span class="o">.</span><span class="n">close</span><span class="p">()</span>
275
276    <span class="k">if</span> <span class="n">dlg</span><span class="p">:</span>
277        <span class="k">return</span> <span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;Rwp&#39;</span><span class="p">]</span>
278</div>
279<div class="viewcode-block" id="SeqRefine"><a class="viewcode-back" href="../GSASIIstruc.html#GSASIIstrMain.SeqRefine">[docs]</a><span class="k">def</span> <span class="nf">SeqRefine</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">,</span><span class="n">dlg</span><span class="p">):</span>
280    <span class="sd">&#39;&#39;&#39;Perform a sequential refinement -- cycles through all selected histgrams,</span>
281<span class="sd">    one at a time</span>
282<span class="sd">    &#39;&#39;&#39;</span>
283    <span class="kn">import</span> <span class="nn">pytexture</span> <span class="kn">as</span> <span class="nn">ptx</span>
284    <span class="n">ptx</span><span class="o">.</span><span class="n">pyqlmninit</span><span class="p">()</span>            <span class="c">#initialize fortran arrays for spherical harmonics</span>
285   
286    <span class="n">printFile</span> <span class="o">=</span> <span class="nb">open</span><span class="p">(</span><span class="n">ospath</span><span class="o">.</span><span class="n">splitext</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)[</span><span class="mi">0</span><span class="p">]</span><span class="o">+</span><span class="s">&#39;.lst&#39;</span><span class="p">,</span><span class="s">&#39;w&#39;</span><span class="p">)</span>
287    <span class="k">print</span> <span class="s">&#39;Starting Sequential Refinement&#39;</span>
288    <span class="n">G2stIO</span><span class="o">.</span><span class="n">ShowBanner</span><span class="p">(</span><span class="n">printFile</span><span class="p">)</span>
289    <span class="n">Controls</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetControls</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)</span>
290    <span class="n">G2stIO</span><span class="o">.</span><span class="n">ShowControls</span><span class="p">(</span><span class="n">Controls</span><span class="p">,</span><span class="n">printFile</span><span class="p">,</span><span class="n">SeqRef</span><span class="o">=</span><span class="bp">True</span><span class="p">)</span>           
291    <span class="n">restraintDict</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetRestraints</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)</span>
292    <span class="n">Histograms</span><span class="p">,</span><span class="n">Phases</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetUsedHistogramsAndPhases</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)</span>
293    <span class="k">if</span> <span class="ow">not</span> <span class="n">Phases</span><span class="p">:</span>
294        <span class="k">print</span> <span class="s">&#39; *** ERROR - you have no phases to refine! ***&#39;</span>
295        <span class="k">print</span> <span class="s">&#39; *** Refine aborted ***&#39;</span>
296        <span class="k">raise</span> <span class="ne">Exception</span>
297    <span class="k">if</span> <span class="ow">not</span> <span class="n">Histograms</span><span class="p">:</span>
298        <span class="k">print</span> <span class="s">&#39; *** ERROR - you have no data to refine with! ***&#39;</span>
299        <span class="k">print</span> <span class="s">&#39; *** Refine aborted ***&#39;</span>
300        <span class="k">raise</span> <span class="ne">Exception</span>
301    <span class="n">rigidbodyDict</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetRigidBodies</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)</span>
302    <span class="n">rbIds</span> <span class="o">=</span> <span class="n">rigidbodyDict</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s">&#39;RBIds&#39;</span><span class="p">,{</span><span class="s">&#39;Vector&#39;</span><span class="p">:[],</span><span class="s">&#39;Residue&#39;</span><span class="p">:[]})</span>
303    <span class="n">rbVary</span><span class="p">,</span><span class="n">rbDict</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetRigidBodyModels</span><span class="p">(</span><span class="n">rigidbodyDict</span><span class="p">,</span><span class="n">pFile</span><span class="o">=</span><span class="n">printFile</span><span class="p">)</span>
304    <span class="n">Natoms</span><span class="p">,</span><span class="n">atomIndx</span><span class="p">,</span><span class="n">phaseVary</span><span class="p">,</span><span class="n">phaseDict</span><span class="p">,</span><span class="n">pawleyLookup</span><span class="p">,</span><span class="n">FFtables</span><span class="p">,</span><span class="n">BLtables</span><span class="p">,</span><span class="n">maxSSwave</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetPhaseData</span><span class="p">(</span><span class="n">Phases</span><span class="p">,</span><span class="n">restraintDict</span><span class="p">,</span><span class="n">rbIds</span><span class="p">,</span><span class="bp">False</span><span class="p">,</span><span class="n">printFile</span><span class="p">)</span>
305    <span class="k">for</span> <span class="n">item</span> <span class="ow">in</span> <span class="n">phaseVary</span><span class="p">:</span>
306        <span class="k">if</span> <span class="s">&#39;::A0&#39;</span> <span class="ow">in</span> <span class="n">item</span><span class="p">:</span>
307            <span class="k">print</span> <span class="s">&#39;**** WARNING - lattice parameters should not be refined in a sequential refinement ****&#39;</span>
308            <span class="k">print</span> <span class="s">&#39;****           instead use the Dij parameters for each powder histogram            ****&#39;</span>
309    <span class="k">if</span> <span class="s">&#39;Seq Data&#39;</span> <span class="ow">in</span> <span class="n">Controls</span><span class="p">:</span>
310        <span class="n">histNames</span> <span class="o">=</span> <span class="n">Controls</span><span class="p">[</span><span class="s">&#39;Seq Data&#39;</span><span class="p">]</span>
311    <span class="k">else</span><span class="p">:</span>
312        <span class="n">histNames</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetHistogramNames</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">,[</span><span class="s">&#39;PWDR&#39;</span><span class="p">,])</span>
313    <span class="k">if</span> <span class="s">&#39;Reverse Seq&#39;</span> <span class="ow">in</span> <span class="n">Controls</span><span class="p">:</span>
314        <span class="k">if</span> <span class="n">Controls</span><span class="p">[</span><span class="s">&#39;Reverse Seq&#39;</span><span class="p">]:</span>
315            <span class="n">histNames</span><span class="o">.</span><span class="n">reverse</span><span class="p">()</span>
316    <span class="n">SeqResult</span> <span class="o">=</span> <span class="p">{</span><span class="s">&#39;histNames&#39;</span><span class="p">:</span><span class="n">histNames</span><span class="p">}</span>
317    <span class="n">makeBack</span> <span class="o">=</span> <span class="bp">True</span>
318    <span class="n">Histo</span> <span class="o">=</span> <span class="p">{}</span>
319    <span class="n">NewparmDict</span> <span class="o">=</span> <span class="p">{}</span>
320    <span class="k">for</span> <span class="n">ihst</span><span class="p">,</span><span class="n">histogram</span> <span class="ow">in</span> <span class="nb">enumerate</span><span class="p">(</span><span class="n">histNames</span><span class="p">):</span>
321        <span class="k">print</span><span class="p">(</span><span class="s">&#39;Refining with &#39;</span><span class="o">+</span><span class="nb">str</span><span class="p">(</span><span class="n">histogram</span><span class="p">))</span>
322        <span class="n">ifPrint</span> <span class="o">=</span> <span class="bp">False</span>
323        <span class="k">if</span> <span class="n">dlg</span><span class="p">:</span>
324            <span class="n">dlg</span><span class="o">.</span><span class="n">SetTitle</span><span class="p">(</span><span class="s">&#39;Residual for histogram &#39;</span><span class="o">+</span><span class="nb">str</span><span class="p">(</span><span class="n">ihst</span><span class="p">))</span>
325        <span class="n">calcControls</span> <span class="o">=</span> <span class="p">{}</span>
326        <span class="n">calcControls</span><span class="p">[</span><span class="s">&#39;atomIndx&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">atomIndx</span>
327        <span class="n">calcControls</span><span class="p">[</span><span class="s">&#39;Natoms&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">Natoms</span>
328        <span class="n">calcControls</span><span class="p">[</span><span class="s">&#39;FFtables&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">FFtables</span>
329        <span class="n">calcControls</span><span class="p">[</span><span class="s">&#39;BLtables&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">BLtables</span>
330        <span class="n">calcControls</span><span class="p">[</span><span class="s">&#39;maxSSwave&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">maxSSwave</span>
331        <span class="n">Histo</span> <span class="o">=</span> <span class="p">{</span><span class="n">histogram</span><span class="p">:</span><span class="n">Histograms</span><span class="p">[</span><span class="n">histogram</span><span class="p">],}</span>
332        <span class="n">hapVary</span><span class="p">,</span><span class="n">hapDict</span><span class="p">,</span><span class="n">controlDict</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetHistogramPhaseData</span><span class="p">(</span><span class="n">Phases</span><span class="p">,</span><span class="n">Histo</span><span class="p">,</span><span class="n">Print</span><span class="o">=</span><span class="bp">False</span><span class="p">)</span>
333        <span class="n">calcControls</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">controlDict</span><span class="p">)</span>
334        <span class="n">histVary</span><span class="p">,</span><span class="n">histDict</span><span class="p">,</span><span class="n">controlDict</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetHistogramData</span><span class="p">(</span><span class="n">Histo</span><span class="p">,</span><span class="bp">False</span><span class="p">)</span>
335        <span class="n">calcControls</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">controlDict</span><span class="p">)</span>
336        <span class="n">varyList</span> <span class="o">=</span> <span class="n">rbVary</span><span class="o">+</span><span class="n">phaseVary</span><span class="o">+</span><span class="n">hapVary</span><span class="o">+</span><span class="n">histVary</span>
337        <span class="k">if</span> <span class="ow">not</span> <span class="n">ihst</span><span class="p">:</span>
338            <span class="c"># save the initial vary list, but without histogram numbers on parameters</span>
339            <span class="n">saveVaryList</span> <span class="o">=</span> <span class="n">varyList</span><span class="p">[:]</span>
340            <span class="k">for</span> <span class="n">i</span><span class="p">,</span><span class="n">item</span> <span class="ow">in</span> <span class="nb">enumerate</span><span class="p">(</span><span class="n">saveVaryList</span><span class="p">):</span>
341                <span class="n">items</span> <span class="o">=</span> <span class="n">item</span><span class="o">.</span><span class="n">split</span><span class="p">(</span><span class="s">&#39;:&#39;</span><span class="p">)</span>
342                <span class="k">if</span> <span class="n">items</span><span class="p">[</span><span class="mi">1</span><span class="p">]:</span>
343                    <span class="n">items</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">=</span> <span class="s">&#39;&#39;</span>
344                <span class="n">item</span> <span class="o">=</span> <span class="s">&#39;:&#39;</span><span class="o">.</span><span class="n">join</span><span class="p">(</span><span class="n">items</span><span class="p">)</span>
345                <span class="n">saveVaryList</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">=</span> <span class="n">item</span>
346            <span class="n">SeqResult</span><span class="p">[</span><span class="s">&#39;varyList&#39;</span><span class="p">]</span> <span class="o">=</span> <span class="n">saveVaryList</span>
347        <span class="n">origvaryList</span> <span class="o">=</span> <span class="n">varyList</span><span class="p">[:]</span>
348        <span class="n">parmDict</span> <span class="o">=</span> <span class="p">{}</span>
349        <span class="n">parmDict</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">phaseDict</span><span class="p">)</span>
350        <span class="n">parmDict</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">hapDict</span><span class="p">)</span>
351        <span class="n">parmDict</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">histDict</span><span class="p">)</span>
352        <span class="k">if</span> <span class="n">Controls</span><span class="p">[</span><span class="s">&#39;Copy2Next&#39;</span><span class="p">]:</span>
353            <span class="n">parmDict</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">NewparmDict</span><span class="p">)</span>
354        <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetFprime</span><span class="p">(</span><span class="n">calcControls</span><span class="p">,</span><span class="n">Histo</span><span class="p">)</span>
355        <span class="c"># do constraint processing</span>
356        <span class="c">#reload(G2mv) # debug</span>
357        <span class="n">G2mv</span><span class="o">.</span><span class="n">InitVars</span><span class="p">()</span>   
358        <span class="n">constrDict</span><span class="p">,</span><span class="n">fixedList</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">GetConstraints</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)</span>
359        <span class="n">varyListStart</span> <span class="o">=</span> <span class="nb">tuple</span><span class="p">(</span><span class="n">varyList</span><span class="p">)</span> <span class="c"># save the original varyList before dependent vars are removed</span>
360        <span class="k">try</span><span class="p">:</span>
361            <span class="n">groups</span><span class="p">,</span><span class="n">parmlist</span> <span class="o">=</span> <span class="n">G2mv</span><span class="o">.</span><span class="n">GroupConstraints</span><span class="p">(</span><span class="n">constrDict</span><span class="p">)</span>
362            <span class="n">G2mv</span><span class="o">.</span><span class="n">GenerateConstraints</span><span class="p">(</span><span class="n">groups</span><span class="p">,</span><span class="n">parmlist</span><span class="p">,</span><span class="n">varyList</span><span class="p">,</span><span class="n">constrDict</span><span class="p">,</span><span class="n">fixedList</span><span class="p">,</span><span class="n">parmDict</span><span class="p">,</span><span class="n">SeqHist</span><span class="o">=</span><span class="n">ihst</span><span class="p">)</span>
363            <span class="n">constraintInfo</span> <span class="o">=</span> <span class="p">(</span><span class="n">groups</span><span class="p">,</span><span class="n">parmlist</span><span class="p">,</span><span class="n">constrDict</span><span class="p">,</span><span class="n">fixedList</span><span class="p">,</span><span class="n">ihst</span><span class="p">)</span>
364        <span class="k">except</span><span class="p">:</span>
365            <span class="k">print</span> <span class="s">&#39; *** ERROR - your constraints are internally inconsistent ***&#39;</span>
366            <span class="c">#errmsg, warnmsg = G2mv.CheckConstraints(varyList,constrDict,fixedList)</span>
367            <span class="c">#print &#39;Errors&#39;,errmsg</span>
368            <span class="c">#if warnmsg: print &#39;Warnings&#39;,warnmsg</span>
369            <span class="k">raise</span> <span class="ne">Exception</span><span class="p">(</span><span class="s">&#39; *** Refine aborted ***&#39;</span><span class="p">)</span>
370        <span class="c">#print G2mv.VarRemapShow(varyList)</span>
371        <span class="k">if</span> <span class="ow">not</span> <span class="n">ihst</span><span class="p">:</span>
372            <span class="c"># first histogram to refine against</span>
373            <span class="n">firstVaryList</span> <span class="o">=</span> <span class="p">[]</span>
374            <span class="k">for</span> <span class="n">item</span> <span class="ow">in</span> <span class="n">varyList</span><span class="p">:</span>
375                <span class="n">items</span> <span class="o">=</span> <span class="n">item</span><span class="o">.</span><span class="n">split</span><span class="p">(</span><span class="s">&#39;:&#39;</span><span class="p">)</span>
376                <span class="k">if</span> <span class="n">items</span><span class="p">[</span><span class="mi">1</span><span class="p">]:</span>
377                    <span class="n">items</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">=</span> <span class="s">&#39;&#39;</span>
378                <span class="n">item</span> <span class="o">=</span> <span class="s">&#39;:&#39;</span><span class="o">.</span><span class="n">join</span><span class="p">(</span><span class="n">items</span><span class="p">)</span>
379                <span class="n">firstVaryList</span><span class="o">.</span><span class="n">append</span><span class="p">(</span><span class="n">item</span><span class="p">)</span>
380            <span class="n">newVaryList</span> <span class="o">=</span> <span class="n">firstVaryList</span>
381        <span class="k">else</span><span class="p">:</span>
382            <span class="n">newVaryList</span> <span class="o">=</span> <span class="p">[]</span>
383            <span class="k">for</span> <span class="n">item</span> <span class="ow">in</span> <span class="n">varyList</span><span class="p">:</span>
384                <span class="n">items</span> <span class="o">=</span> <span class="n">item</span><span class="o">.</span><span class="n">split</span><span class="p">(</span><span class="s">&#39;:&#39;</span><span class="p">)</span>
385                <span class="k">if</span> <span class="n">items</span><span class="p">[</span><span class="mi">1</span><span class="p">]:</span>
386                    <span class="n">items</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">=</span> <span class="s">&#39;&#39;</span>
387                <span class="n">item</span> <span class="o">=</span> <span class="s">&#39;:&#39;</span><span class="o">.</span><span class="n">join</span><span class="p">(</span><span class="n">items</span><span class="p">)</span>
388                <span class="n">newVaryList</span><span class="o">.</span><span class="n">append</span><span class="p">(</span><span class="n">item</span><span class="p">)</span>
389        <span class="k">if</span> <span class="n">newVaryList</span> <span class="o">!=</span> <span class="n">firstVaryList</span> <span class="ow">and</span> <span class="n">Controls</span><span class="p">[</span><span class="s">&#39;Copy2Next&#39;</span><span class="p">]:</span>
390            <span class="c"># variable lists are expected to match between sequential refinements when Copy2Next is on</span>
391            <span class="k">print</span> <span class="s">&#39;**** ERROR - variable list for this histogram does not match previous&#39;</span>
392            <span class="k">print</span> <span class="s">&#39;     Copy of variables is not possible&#39;</span>
393            <span class="k">print</span> <span class="s">&#39;</span><span class="se">\n</span><span class="s">current histogram&#39;</span><span class="p">,</span><span class="n">histogram</span><span class="p">,</span><span class="s">&#39;has&#39;</span><span class="p">,</span><span class="nb">len</span><span class="p">(</span><span class="n">newVaryList</span><span class="p">),</span><span class="s">&#39;variables&#39;</span>
394            <span class="n">combined</span> <span class="o">=</span> <span class="nb">list</span><span class="p">(</span><span class="nb">set</span><span class="p">(</span><span class="n">firstVaryList</span><span class="o">+</span><span class="n">newVaryList</span><span class="p">))</span>
395            <span class="n">c</span> <span class="o">=</span> <span class="p">[</span><span class="n">var</span> <span class="k">for</span> <span class="n">var</span> <span class="ow">in</span> <span class="n">combined</span> <span class="k">if</span> <span class="n">var</span> <span class="ow">not</span> <span class="ow">in</span> <span class="n">newVaryList</span><span class="p">]</span>
396            <span class="n">p</span> <span class="o">=</span> <span class="p">[</span><span class="n">var</span> <span class="k">for</span> <span class="n">var</span> <span class="ow">in</span> <span class="n">combined</span> <span class="k">if</span> <span class="n">var</span> <span class="ow">not</span> <span class="ow">in</span> <span class="n">firstVaryList</span><span class="p">]</span>
397            <span class="n">line</span> <span class="o">=</span> <span class="s">&#39;Variables in previous but not in current: &#39;</span>
398            <span class="k">if</span> <span class="n">c</span><span class="p">:</span>
399                <span class="k">for</span> <span class="n">var</span> <span class="ow">in</span> <span class="n">c</span><span class="p">:</span>
400                    <span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">line</span><span class="p">)</span> <span class="o">&gt;</span> <span class="mi">100</span><span class="p">:</span>
401                        <span class="k">print</span> <span class="n">line</span>
402                        <span class="n">line</span> <span class="o">=</span> <span class="s">&#39;    &#39;</span>
403                    <span class="n">line</span> <span class="o">+=</span> <span class="n">var</span> <span class="o">+</span> <span class="s">&#39;, &#39;</span>
404            <span class="k">else</span><span class="p">:</span>
405                <span class="n">line</span> <span class="o">+=</span> <span class="s">&#39;none&#39;</span>
406            <span class="k">print</span> <span class="n">line</span>
407            <span class="k">print</span> <span class="s">&#39;</span><span class="se">\n</span><span class="s">Previous refinement has&#39;</span><span class="p">,</span><span class="nb">len</span><span class="p">(</span><span class="n">firstVaryList</span><span class="p">),</span><span class="s">&#39;variables&#39;</span>
408            <span class="n">line</span> <span class="o">=</span> <span class="s">&#39;Variables in current but not in previous: &#39;</span>
409            <span class="k">if</span> <span class="n">p</span><span class="p">:</span>
410                <span class="k">for</span> <span class="n">var</span> <span class="ow">in</span> <span class="n">p</span><span class="p">:</span>
411                    <span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">line</span><span class="p">)</span> <span class="o">&gt;</span> <span class="mi">100</span><span class="p">:</span>
412                        <span class="k">print</span> <span class="n">line</span>
413                        <span class="n">line</span> <span class="o">=</span> <span class="s">&#39;    &#39;</span>
414                    <span class="n">line</span> <span class="o">+=</span> <span class="n">var</span> <span class="o">+</span> <span class="s">&#39;, &#39;</span>
415            <span class="k">else</span><span class="p">:</span>
416                <span class="n">line</span> <span class="o">+=</span> <span class="s">&#39;none&#39;</span>
417            <span class="k">print</span> <span class="n">line</span>
418            <span class="k">raise</span> <span class="ne">Exception</span>
419       
420        <span class="n">ifPrint</span> <span class="o">=</span> <span class="bp">False</span>
421        <span class="k">print</span> <span class="o">&gt;&gt;</span><span class="n">printFile</span><span class="p">,</span><span class="s">&#39;</span><span class="se">\n</span><span class="s"> Refinement results for histogram: v&#39;</span><span class="o">+</span><span class="n">histogram</span>
422        <span class="k">print</span> <span class="o">&gt;&gt;</span><span class="n">printFile</span><span class="p">,</span><span class="mi">135</span><span class="o">*</span><span class="s">&#39;-&#39;</span>
423        <span class="n">IfOK</span><span class="p">,</span><span class="n">Rvals</span><span class="p">,</span><span class="n">result</span><span class="p">,</span><span class="n">covMatrix</span><span class="p">,</span><span class="n">sig</span> <span class="o">=</span> <span class="n">RefineCore</span><span class="p">(</span><span class="n">Controls</span><span class="p">,</span><span class="n">Histo</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">restraintDict</span><span class="p">,</span>
424            <span class="n">rigidbodyDict</span><span class="p">,</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">,</span><span class="n">calcControls</span><span class="p">,</span><span class="n">pawleyLookup</span><span class="p">,</span><span class="n">ifPrint</span><span class="p">,</span><span class="n">printFile</span><span class="p">,</span><span class="n">dlg</span><span class="p">)</span>
425
426        <span class="k">print</span> <span class="s">&#39;  wR = </span><span class="si">%7.2f%%</span><span class="s">, chi**2 = </span><span class="si">%12.6g</span><span class="s">, reduced chi**2 = </span><span class="si">%6.2f</span><span class="s">, last delta chi = </span><span class="si">%.4f</span><span class="s">&#39;</span><span class="o">%</span><span class="p">(</span>
427            <span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;Rwp&#39;</span><span class="p">],</span><span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;chisq&#39;</span><span class="p">],</span><span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;GOF&#39;</span><span class="p">]</span><span class="o">**</span><span class="mi">2</span><span class="p">,</span><span class="n">Rvals</span><span class="p">[</span><span class="s">&#39;DelChi2&#39;</span><span class="p">])</span>
428        <span class="c"># add the uncertainties into the esd dictionary (sigDict)</span>
429        <span class="n">sigDict</span> <span class="o">=</span> <span class="nb">dict</span><span class="p">(</span><span class="nb">zip</span><span class="p">(</span><span class="n">varyList</span><span class="p">,</span><span class="n">sig</span><span class="p">))</span>
430        <span class="c"># the uncertainties for dependent constrained parms into the esd dict</span>
431        <span class="n">sigDict</span><span class="o">.</span><span class="n">update</span><span class="p">(</span><span class="n">G2mv</span><span class="o">.</span><span class="n">ComputeDepESD</span><span class="p">(</span><span class="n">covMatrix</span><span class="p">,</span><span class="n">varyList</span><span class="p">,</span><span class="n">parmDict</span><span class="p">))</span>
432
433        <span class="c"># a dict with values &amp; esds for dependent (constrained) parameters</span>
434        <span class="n">depParmDict</span> <span class="o">=</span> <span class="p">{</span><span class="n">i</span><span class="p">:(</span><span class="n">parmDict</span><span class="p">[</span><span class="n">i</span><span class="p">],</span><span class="n">sigDict</span><span class="p">[</span><span class="n">i</span><span class="p">])</span> <span class="k">for</span> <span class="n">i</span> <span class="ow">in</span> <span class="n">varyListStart</span>
435                       <span class="k">if</span> <span class="n">i</span> <span class="ow">not</span> <span class="ow">in</span> <span class="n">varyList</span><span class="p">}</span>
436        <span class="n">newCellDict</span> <span class="o">=</span> <span class="n">copy</span><span class="o">.</span><span class="n">deepcopy</span><span class="p">(</span><span class="n">G2stMth</span><span class="o">.</span><span class="n">GetNewCellParms</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">))</span>
437        <span class="n">newAtomDict</span> <span class="o">=</span> <span class="n">copy</span><span class="o">.</span><span class="n">deepcopy</span><span class="p">(</span><span class="n">G2stMth</span><span class="o">.</span><span class="n">ApplyXYZshifts</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">varyList</span><span class="p">))</span>
438        <span class="n">histRefData</span> <span class="o">=</span> <span class="p">{</span>
439            <span class="s">&#39;variables&#39;</span><span class="p">:</span><span class="n">result</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span><span class="s">&#39;varyList&#39;</span><span class="p">:</span><span class="n">varyList</span><span class="p">,</span><span class="s">&#39;sig&#39;</span><span class="p">:</span><span class="n">sig</span><span class="p">,</span><span class="s">&#39;Rvals&#39;</span><span class="p">:</span><span class="n">Rvals</span><span class="p">,</span>
440            <span class="s">&#39;varyListStart&#39;</span><span class="p">:</span><span class="n">varyListStart</span><span class="p">,</span>
441            <span class="s">&#39;covMatrix&#39;</span><span class="p">:</span><span class="n">covMatrix</span><span class="p">,</span><span class="s">&#39;title&#39;</span><span class="p">:</span><span class="n">histogram</span><span class="p">,</span><span class="s">&#39;newAtomDict&#39;</span><span class="p">:</span><span class="n">newAtomDict</span><span class="p">,</span>
442            <span class="s">&#39;newCellDict&#39;</span><span class="p">:</span><span class="n">newCellDict</span><span class="p">,</span><span class="s">&#39;depParmDict&#39;</span><span class="p">:</span><span class="n">depParmDict</span><span class="p">,</span>
443            <span class="s">&#39;constraintInfo&#39;</span><span class="p">:</span><span class="n">constraintInfo</span><span class="p">,</span>
444            <span class="s">&#39;parmDict&#39;</span><span class="p">:</span><span class="n">parmDict</span><span class="p">}</span>
445        <span class="n">SeqResult</span><span class="p">[</span><span class="n">histogram</span><span class="p">]</span> <span class="o">=</span> <span class="n">histRefData</span>
446        <span class="n">G2stMth</span><span class="o">.</span><span class="n">ApplyRBModels</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">rigidbodyDict</span><span class="p">,</span><span class="bp">True</span><span class="p">)</span>
447<span class="c">#        G2stIO.SetRigidBodyModels(parmDict,sigDict,rigidbodyDict,printFile)</span>
448        <span class="n">G2stIO</span><span class="o">.</span><span class="n">SetHistogramPhaseData</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">sigDict</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">Histo</span><span class="p">,</span><span class="n">ifPrint</span><span class="p">,</span><span class="n">printFile</span><span class="p">)</span>
449        <span class="n">G2stIO</span><span class="o">.</span><span class="n">SetHistogramData</span><span class="p">(</span><span class="n">parmDict</span><span class="p">,</span><span class="n">sigDict</span><span class="p">,</span><span class="n">Histo</span><span class="p">,</span><span class="n">ifPrint</span><span class="p">,</span><span class="n">printFile</span><span class="p">)</span>
450        <span class="n">G2stIO</span><span class="o">.</span><span class="n">SetUsedHistogramsAndPhases</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">,</span><span class="n">Histo</span><span class="p">,</span><span class="n">Phases</span><span class="p">,</span><span class="n">rigidbodyDict</span><span class="p">,</span><span class="n">histRefData</span><span class="p">,</span><span class="n">makeBack</span><span class="p">)</span>
451        <span class="n">makeBack</span> <span class="o">=</span> <span class="bp">False</span>
452        <span class="n">NewparmDict</span> <span class="o">=</span> <span class="p">{}</span>
453        <span class="c"># make dict of varied parameters in current histogram, renamed to</span>
454        <span class="c"># next histogram, for use in next refinement. </span>
455        <span class="k">if</span> <span class="n">Controls</span><span class="p">[</span><span class="s">&#39;Copy2Next&#39;</span><span class="p">]</span> <span class="ow">and</span> <span class="n">ihst</span> <span class="o">&lt;</span> <span class="nb">len</span><span class="p">(</span><span class="n">histNames</span><span class="p">)</span><span class="o">-</span><span class="mi">1</span><span class="p">:</span>
456            <span class="n">hId</span> <span class="o">=</span> <span class="n">Histo</span><span class="p">[</span><span class="n">histogram</span><span class="p">][</span><span class="s">&#39;hId&#39;</span><span class="p">]</span> <span class="c"># current histogram</span>
457            <span class="n">nexthId</span> <span class="o">=</span> <span class="n">Histograms</span><span class="p">[</span><span class="n">histNames</span><span class="p">[</span><span class="n">ihst</span><span class="o">+</span><span class="mi">1</span><span class="p">]][</span><span class="s">&#39;hId&#39;</span><span class="p">]</span>
458            <span class="k">for</span> <span class="n">parm</span> <span class="ow">in</span> <span class="nb">set</span><span class="p">(</span><span class="nb">list</span><span class="p">(</span><span class="n">varyList</span><span class="p">)</span><span class="o">+</span><span class="nb">list</span><span class="p">(</span><span class="n">varyListStart</span><span class="p">)):</span>
459                <span class="n">items</span> <span class="o">=</span> <span class="n">parm</span><span class="o">.</span><span class="n">split</span><span class="p">(</span><span class="s">&#39;:&#39;</span><span class="p">)</span>
460                <span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">items</span><span class="p">)</span> <span class="o">&lt;</span> <span class="mi">3</span><span class="p">:</span> <span class="k">continue</span>
461                <span class="k">if</span> <span class="nb">str</span><span class="p">(</span><span class="n">hId</span><span class="p">)</span> <span class="ow">in</span> <span class="n">items</span><span class="p">[</span><span class="mi">1</span><span class="p">]:</span>
462                    <span class="n">items</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">=</span> <span class="nb">str</span><span class="p">(</span><span class="n">nexthId</span><span class="p">)</span>
463                    <span class="n">newparm</span> <span class="o">=</span> <span class="s">&#39;:&#39;</span><span class="o">.</span><span class="n">join</span><span class="p">(</span><span class="n">items</span><span class="p">)</span>
464                    <span class="n">NewparmDict</span><span class="p">[</span><span class="n">newparm</span><span class="p">]</span> <span class="o">=</span> <span class="n">parmDict</span><span class="p">[</span><span class="n">parm</span><span class="p">]</span>
465    <span class="n">G2stIO</span><span class="o">.</span><span class="n">SetSeqResult</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">,</span><span class="n">Histograms</span><span class="p">,</span><span class="n">SeqResult</span><span class="p">)</span>
466    <span class="n">printFile</span><span class="o">.</span><span class="n">close</span><span class="p">()</span>
467    <span class="k">print</span> <span class="s">&#39; Sequential refinement results are in file: &#39;</span><span class="o">+</span><span class="n">ospath</span><span class="o">.</span><span class="n">splitext</span><span class="p">(</span><span class="n">GPXfile</span><span class="p">)[</span><span class="mi">0</span><span class="p">]</span><span class="o">+</span><span class="s">&#39;.lst&#39;</span>
468    <span class="k">print</span> <span class="s">&#39; ***** Sequential refinement successful *****&#39;</span>
469</div>
470<div class="viewcode-block" id="RetDistAngle"><a class="viewcode-back" href="../GSASIIstruc.html#GSASIIstrMain.RetDistAngle">[docs]</a><span class="k">def</span> <span class="nf">RetDistAngle</span><span class="p">(</span><span class="n">DisAglCtls</span><span class="p">,</span><span class="n">DisAglData</span><span class="p">):</span>
471    <span class="sd">&#39;&#39;&#39;Compute and return distances and angles</span>
472
473<span class="sd">    :param dict DisAglCtls: contains distance/angle radii usually defined using</span>
474<span class="sd">       :func:`GSASIIgrid.DisAglDialog`</span>
475<span class="sd">    :param dict DisAglData: contains phase data: </span>
476<span class="sd">       Items &#39;OrigAtoms&#39; and &#39;TargAtoms&#39; contain the atoms to be used</span>
477<span class="sd">       for distance/angle origins and atoms to be used as targets.</span>
478<span class="sd">       Item &#39;SGData&#39; has the space group information (see :ref:`Space Group object&lt;SGData_table&gt;`)</span>
479
480<span class="sd">    :returns: AtomLabels,DistArray,AngArray where:</span>
481
482<span class="sd">      **AtomLabels** is a dict of atom labels, keys are the atom number</span>
483
484<span class="sd">      **DistArray** is a dict keyed by the origin atom number where the value is a list</span>
485<span class="sd">      of distance entries. The value for each distance is a list containing:</span>
486<span class="sd">      </span>
487<span class="sd">        0) the target atom number (int);</span>
488<span class="sd">        1) the unit cell offsets added to x,y &amp; z (tuple of int values)</span>
489<span class="sd">        2) the symmetry operator number (which may be modified to indicate centering and center of symmetry)</span>
490<span class="sd">        3) an interatomic distance in A (float)</span>
491<span class="sd">        4) an uncertainty on the distance in A or 0.0 (float)</span>
492
493<span class="sd">      **AngArray** is a dict keyed by the origin (central) atom number where</span>
494<span class="sd">      the value is a list of</span>
495<span class="sd">      angle entries. The value for each angle entry consists of three values:</span>
496
497<span class="sd">        0) a distance item reference for one neighbor (int)</span>
498<span class="sd">        1) a distance item reference for a second neighbor (int)</span>
499<span class="sd">        2) a angle, uncertainty pair; the s.u. may be zero (tuple of two floats)</span>
500
501<span class="sd">      The AngArray distance reference items refer directly to the index of the items in the</span>
502<span class="sd">      DistArray item for the list of distances for the central atom. </span>
503<span class="sd">    &#39;&#39;&#39;</span>
504    <span class="kn">import</span> <span class="nn">numpy.ma</span> <span class="kn">as</span> <span class="nn">ma</span>
505   
506    <span class="n">SGData</span> <span class="o">=</span> <span class="n">DisAglData</span><span class="p">[</span><span class="s">&#39;SGData&#39;</span><span class="p">]</span>
507    <span class="n">Cell</span> <span class="o">=</span> <span class="n">DisAglData</span><span class="p">[</span><span class="s">&#39;Cell&#39;</span><span class="p">]</span>
508   
509    <span class="n">Amat</span><span class="p">,</span><span class="n">Bmat</span> <span class="o">=</span> <span class="n">G2lat</span><span class="o">.</span><span class="n">cell2AB</span><span class="p">(</span><span class="n">Cell</span><span class="p">[:</span><span class="mi">6</span><span class="p">])</span>
510    <span class="n">covData</span> <span class="o">=</span> <span class="p">{}</span>
511    <span class="k">if</span> <span class="s">&#39;covData&#39;</span> <span class="ow">in</span> <span class="n">DisAglData</span><span class="p">:</span>   
512        <span class="n">covData</span> <span class="o">=</span> <span class="n">DisAglData</span><span class="p">[</span><span class="s">&#39;covData&#39;</span><span class="p">]</span>
513        <span class="n">covMatrix</span> <span class="o">=</span> <span class="n">covData</span><span class="p">[</span><span class="s">&#39;covMatrix&#39;</span><span class="p">]</span>
514        <span class="n">varyList</span> <span class="o">=</span> <span class="n">covData</span><span class="p">[</span><span class="s">&#39;varyList&#39;</span><span class="p">]</span>
515        <span class="n">pfx</span> <span class="o">=</span> <span class="nb">str</span><span class="p">(</span><span class="n">DisAglData</span><span class="p">[</span><span class="s">&#39;pId&#39;</span><span class="p">])</span><span class="o">+</span><span class="s">&#39;::&#39;</span>
516        <span class="n">A</span> <span class="o">=</span> <span class="n">G2lat</span><span class="o">.</span><span class="n">cell2A</span><span class="p">(</span><span class="n">Cell</span><span class="p">[:</span><span class="mi">6</span><span class="p">])</span>
517        <span class="n">cellSig</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">getCellEsd</span><span class="p">(</span><span class="n">pfx</span><span class="p">,</span><span class="n">SGData</span><span class="p">,</span><span class="n">A</span><span class="p">,</span><span class="n">covData</span><span class="p">)</span>
518        <span class="n">names</span> <span class="o">=</span> <span class="p">[</span><span class="s">&#39; a = &#39;</span><span class="p">,</span><span class="s">&#39; b = &#39;</span><span class="p">,</span><span class="s">&#39; c = &#39;</span><span class="p">,</span><span class="s">&#39; alpha = &#39;</span><span class="p">,</span><span class="s">&#39; beta = &#39;</span><span class="p">,</span><span class="s">&#39; gamma = &#39;</span><span class="p">,</span><span class="s">&#39; Volume = &#39;</span><span class="p">]</span>
519        <span class="n">valEsd</span> <span class="o">=</span> <span class="p">[</span><span class="n">G2mth</span><span class="o">.</span><span class="n">ValEsd</span><span class="p">(</span><span class="n">Cell</span><span class="p">[</span><span class="n">i</span><span class="p">],</span><span class="n">cellSig</span><span class="p">[</span><span class="n">i</span><span class="p">],</span><span class="bp">True</span><span class="p">)</span> <span class="k">for</span> <span class="n">i</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">7</span><span class="p">)]</span>
520
521    <span class="n">Factor</span> <span class="o">=</span> <span class="n">DisAglCtls</span><span class="p">[</span><span class="s">&#39;Factors&#39;</span><span class="p">]</span>
522    <span class="n">Radii</span> <span class="o">=</span> <span class="nb">dict</span><span class="p">(</span><span class="nb">zip</span><span class="p">(</span><span class="n">DisAglCtls</span><span class="p">[</span><span class="s">&#39;AtomTypes&#39;</span><span class="p">],</span><span class="nb">zip</span><span class="p">(</span><span class="n">DisAglCtls</span><span class="p">[</span><span class="s">&#39;BondRadii&#39;</span><span class="p">],</span><span class="n">DisAglCtls</span><span class="p">[</span><span class="s">&#39;AngleRadii&#39;</span><span class="p">])))</span>
523    <span class="n">indices</span> <span class="o">=</span> <span class="p">(</span><span class="o">-</span><span class="mi">1</span><span class="p">,</span><span class="mi">0</span><span class="p">,</span><span class="mi">1</span><span class="p">)</span>
524    <span class="n">Units</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">array</span><span class="p">([[</span><span class="n">h</span><span class="p">,</span><span class="n">k</span><span class="p">,</span><span class="n">l</span><span class="p">]</span> <span class="k">for</span> <span class="n">h</span> <span class="ow">in</span> <span class="n">indices</span> <span class="k">for</span> <span class="n">k</span> <span class="ow">in</span> <span class="n">indices</span> <span class="k">for</span> <span class="n">l</span> <span class="ow">in</span> <span class="n">indices</span><span class="p">])</span>
525    <span class="n">origAtoms</span> <span class="o">=</span> <span class="n">DisAglData</span><span class="p">[</span><span class="s">&#39;OrigAtoms&#39;</span><span class="p">]</span>
526    <span class="n">targAtoms</span> <span class="o">=</span> <span class="n">DisAglData</span><span class="p">[</span><span class="s">&#39;TargAtoms&#39;</span><span class="p">]</span>
527    <span class="n">AtomLabels</span> <span class="o">=</span> <span class="p">{}</span>
528    <span class="k">for</span> <span class="n">Oatom</span> <span class="ow">in</span> <span class="n">origAtoms</span><span class="p">:</span>
529        <span class="n">AtomLabels</span><span class="p">[</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">0</span><span class="p">]]</span> <span class="o">=</span> <span class="n">Oatom</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span>
530    <span class="k">for</span> <span class="n">Oatom</span> <span class="ow">in</span> <span class="n">targAtoms</span><span class="p">:</span>
531        <span class="n">AtomLabels</span><span class="p">[</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">0</span><span class="p">]]</span> <span class="o">=</span> <span class="n">Oatom</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span>
532    <span class="n">DistArray</span> <span class="o">=</span> <span class="p">{}</span>
533    <span class="n">AngArray</span> <span class="o">=</span> <span class="p">{}</span>
534    <span class="k">for</span> <span class="n">Oatom</span> <span class="ow">in</span> <span class="n">origAtoms</span><span class="p">:</span>
535        <span class="n">DistArray</span><span class="p">[</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">0</span><span class="p">]]</span> <span class="o">=</span> <span class="p">[]</span>
536        <span class="n">AngArray</span><span class="p">[</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">0</span><span class="p">]]</span> <span class="o">=</span> <span class="p">[]</span>
537        <span class="n">OxyzNames</span> <span class="o">=</span> <span class="s">&#39;&#39;</span>
538        <span class="n">IndBlist</span> <span class="o">=</span> <span class="p">[]</span>
539        <span class="n">Dist</span> <span class="o">=</span> <span class="p">[]</span>
540        <span class="n">Vect</span> <span class="o">=</span> <span class="p">[]</span>
541        <span class="n">VectA</span> <span class="o">=</span> <span class="p">[]</span>
542        <span class="n">angles</span> <span class="o">=</span> <span class="p">[]</span>
543        <span class="k">for</span> <span class="n">Tatom</span> <span class="ow">in</span> <span class="n">targAtoms</span><span class="p">:</span>
544            <span class="n">Xvcov</span> <span class="o">=</span> <span class="p">[]</span>
545            <span class="n">TxyzNames</span> <span class="o">=</span> <span class="s">&#39;&#39;</span>
546            <span class="k">if</span> <span class="s">&#39;covData&#39;</span> <span class="ow">in</span> <span class="n">DisAglData</span><span class="p">:</span>
547                <span class="n">OxyzNames</span> <span class="o">=</span> <span class="p">[</span><span class="n">pfx</span><span class="o">+</span><span class="s">&#39;dAx:</span><span class="si">%d</span><span class="s">&#39;</span><span class="o">%</span><span class="p">(</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">0</span><span class="p">]),</span><span class="n">pfx</span><span class="o">+</span><span class="s">&#39;dAy:</span><span class="si">%d</span><span class="s">&#39;</span><span class="o">%</span><span class="p">(</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">0</span><span class="p">]),</span><span class="n">pfx</span><span class="o">+</span><span class="s">&#39;dAz:</span><span class="si">%d</span><span class="s">&#39;</span><span class="o">%</span><span class="p">(</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">0</span><span class="p">])]</span>
548                <span class="n">TxyzNames</span> <span class="o">=</span> <span class="p">[</span><span class="n">pfx</span><span class="o">+</span><span class="s">&#39;dAx:</span><span class="si">%d</span><span class="s">&#39;</span><span class="o">%</span><span class="p">(</span><span class="n">Tatom</span><span class="p">[</span><span class="mi">0</span><span class="p">]),</span><span class="n">pfx</span><span class="o">+</span><span class="s">&#39;dAy:</span><span class="si">%d</span><span class="s">&#39;</span><span class="o">%</span><span class="p">(</span><span class="n">Tatom</span><span class="p">[</span><span class="mi">0</span><span class="p">]),</span><span class="n">pfx</span><span class="o">+</span><span class="s">&#39;dAz:</span><span class="si">%d</span><span class="s">&#39;</span><span class="o">%</span><span class="p">(</span><span class="n">Tatom</span><span class="p">[</span><span class="mi">0</span><span class="p">])]</span>
549                <span class="n">Xvcov</span> <span class="o">=</span> <span class="n">G2mth</span><span class="o">.</span><span class="n">getVCov</span><span class="p">(</span><span class="n">OxyzNames</span><span class="o">+</span><span class="n">TxyzNames</span><span class="p">,</span><span class="n">varyList</span><span class="p">,</span><span class="n">covMatrix</span><span class="p">)</span>
550            <span class="n">result</span> <span class="o">=</span> <span class="n">G2spc</span><span class="o">.</span><span class="n">GenAtom</span><span class="p">(</span><span class="n">Tatom</span><span class="p">[</span><span class="mi">3</span><span class="p">:</span><span class="mi">6</span><span class="p">],</span><span class="n">SGData</span><span class="p">,</span><span class="bp">False</span><span class="p">,</span><span class="n">Move</span><span class="o">=</span><span class="bp">False</span><span class="p">)</span>
551            <span class="n">BsumR</span> <span class="o">=</span> <span class="p">(</span><span class="n">Radii</span><span class="p">[</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">2</span><span class="p">]][</span><span class="mi">0</span><span class="p">]</span><span class="o">+</span><span class="n">Radii</span><span class="p">[</span><span class="n">Tatom</span><span class="p">[</span><span class="mi">2</span><span class="p">]][</span><span class="mi">0</span><span class="p">])</span><span class="o">*</span><span class="n">Factor</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span>
552            <span class="n">AsumR</span> <span class="o">=</span> <span class="p">(</span><span class="n">Radii</span><span class="p">[</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">2</span><span class="p">]][</span><span class="mi">1</span><span class="p">]</span><span class="o">+</span><span class="n">Radii</span><span class="p">[</span><span class="n">Tatom</span><span class="p">[</span><span class="mi">2</span><span class="p">]][</span><span class="mi">1</span><span class="p">])</span><span class="o">*</span><span class="n">Factor</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span>
553            <span class="k">for</span> <span class="n">Txyz</span><span class="p">,</span><span class="n">Top</span><span class="p">,</span><span class="n">Tunit</span> <span class="ow">in</span> <span class="n">result</span><span class="p">:</span>
554                <span class="n">Dx</span> <span class="o">=</span> <span class="p">(</span><span class="n">Txyz</span><span class="o">-</span><span class="n">np</span><span class="o">.</span><span class="n">array</span><span class="p">(</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">3</span><span class="p">:</span><span class="mi">6</span><span class="p">]))</span><span class="o">+</span><span class="n">Units</span>
555                <span class="n">dx</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">inner</span><span class="p">(</span><span class="n">Amat</span><span class="p">,</span><span class="n">Dx</span><span class="p">)</span>
556                <span class="n">dist</span> <span class="o">=</span> <span class="n">ma</span><span class="o">.</span><span class="n">masked_less</span><span class="p">(</span><span class="n">np</span><span class="o">.</span><span class="n">sqrt</span><span class="p">(</span><span class="n">np</span><span class="o">.</span><span class="n">sum</span><span class="p">(</span><span class="n">dx</span><span class="o">**</span><span class="mi">2</span><span class="p">,</span><span class="n">axis</span><span class="o">=</span><span class="mi">0</span><span class="p">)),</span><span class="mf">0.5</span><span class="p">)</span>
557                <span class="n">IndB</span> <span class="o">=</span> <span class="n">ma</span><span class="o">.</span><span class="n">nonzero</span><span class="p">(</span><span class="n">ma</span><span class="o">.</span><span class="n">masked_greater</span><span class="p">(</span><span class="n">dist</span><span class="o">-</span><span class="n">BsumR</span><span class="p">,</span><span class="mf">0.</span><span class="p">))</span>
558                <span class="k">if</span> <span class="n">np</span><span class="o">.</span><span class="n">any</span><span class="p">(</span><span class="n">IndB</span><span class="p">):</span>
559                    <span class="k">for</span> <span class="n">indb</span> <span class="ow">in</span> <span class="n">IndB</span><span class="p">:</span>
560                        <span class="k">for</span> <span class="n">i</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="nb">len</span><span class="p">(</span><span class="n">indb</span><span class="p">)):</span>
561                            <span class="k">if</span> <span class="nb">str</span><span class="p">(</span><span class="n">dx</span><span class="o">.</span><span class="n">T</span><span class="p">[</span><span class="n">indb</span><span class="p">][</span><span class="n">i</span><span class="p">])</span> <span class="ow">not</span> <span class="ow">in</span> <span class="n">IndBlist</span><span class="p">:</span>
562                                <span class="n">IndBlist</span><span class="o">.</span><span class="n">append</span><span class="p">(</span><span class="nb">str</span><span class="p">(</span><span class="n">dx</span><span class="o">.</span><span class="n">T</span><span class="p">[</span><span class="n">indb</span><span class="p">][</span><span class="n">i</span><span class="p">]))</span>
563                                <span class="n">unit</span> <span class="o">=</span> <span class="n">Units</span><span class="p">[</span><span class="n">indb</span><span class="p">][</span><span class="n">i</span><span class="p">]</span>
564                                <span class="n">tunit</span> <span class="o">=</span> <span class="p">(</span><span class="n">unit</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="o">+</span><span class="n">Tunit</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span><span class="n">unit</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span><span class="o">+</span><span class="n">Tunit</span><span class="p">[</span><span class="mi">1</span><span class="p">],</span><span class="n">unit</span><span class="p">[</span><span class="mi">2</span><span class="p">]</span><span class="o">+</span><span class="n">Tunit</span><span class="p">[</span><span class="mi">2</span><span class="p">])</span>
565                                <span class="n">pdpx</span> <span class="o">=</span> <span class="n">G2mth</span><span class="o">.</span><span class="n">getDistDerv</span><span class="p">(</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">3</span><span class="p">:</span><span class="mi">6</span><span class="p">],</span><span class="n">Tatom</span><span class="p">[</span><span class="mi">3</span><span class="p">:</span><span class="mi">6</span><span class="p">],</span><span class="n">Amat</span><span class="p">,</span><span class="n">unit</span><span class="p">,</span><span class="n">Top</span><span class="p">,</span><span class="n">SGData</span><span class="p">)</span>
566                                <span class="n">sig</span> <span class="o">=</span> <span class="mf">0.0</span>
567                                <span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">Xvcov</span><span class="p">):</span>
568                                    <span class="n">sig</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">sqrt</span><span class="p">(</span><span class="n">np</span><span class="o">.</span><span class="n">inner</span><span class="p">(</span><span class="n">pdpx</span><span class="p">,</span><span class="n">np</span><span class="o">.</span><span class="n">inner</span><span class="p">(</span><span class="n">Xvcov</span><span class="p">,</span><span class="n">pdpx</span><span class="p">)))</span>
569                                <span class="n">Dist</span><span class="o">.</span><span class="n">append</span><span class="p">([</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span><span class="n">Tatom</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span><span class="n">tunit</span><span class="p">,</span><span class="n">Top</span><span class="p">,</span><span class="n">ma</span><span class="o">.</span><span class="n">getdata</span><span class="p">(</span><span class="n">dist</span><span class="p">[</span><span class="n">indb</span><span class="p">])[</span><span class="n">i</span><span class="p">],</span><span class="n">sig</span><span class="p">])</span>
570                                <span class="k">if</span> <span class="p">(</span><span class="n">Dist</span><span class="p">[</span><span class="o">-</span><span class="mi">1</span><span class="p">][</span><span class="o">-</span><span class="mi">2</span><span class="p">]</span><span class="o">-</span><span class="n">AsumR</span><span class="p">)</span> <span class="o">&lt;=</span> <span class="mf">0.</span><span class="p">:</span>
571                                    <span class="n">Vect</span><span class="o">.</span><span class="n">append</span><span class="p">(</span><span class="n">dx</span><span class="o">.</span><span class="n">T</span><span class="p">[</span><span class="n">indb</span><span class="p">][</span><span class="n">i</span><span class="p">]</span><span class="o">/</span><span class="n">Dist</span><span class="p">[</span><span class="o">-</span><span class="mi">1</span><span class="p">][</span><span class="o">-</span><span class="mi">2</span><span class="p">])</span>
572                                    <span class="n">VectA</span><span class="o">.</span><span class="n">append</span><span class="p">([</span><span class="n">OxyzNames</span><span class="p">,</span><span class="n">np</span><span class="o">.</span><span class="n">array</span><span class="p">(</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">3</span><span class="p">:</span><span class="mi">6</span><span class="p">]),</span><span class="n">TxyzNames</span><span class="p">,</span><span class="n">np</span><span class="o">.</span><span class="n">array</span><span class="p">(</span><span class="n">Tatom</span><span class="p">[</span><span class="mi">3</span><span class="p">:</span><span class="mi">6</span><span class="p">]),</span><span class="n">unit</span><span class="p">,</span><span class="n">Top</span><span class="p">])</span>
573                                <span class="k">else</span><span class="p">:</span>
574                                    <span class="n">Vect</span><span class="o">.</span><span class="n">append</span><span class="p">([</span><span class="mf">0.</span><span class="p">,</span><span class="mf">0.</span><span class="p">,</span><span class="mf">0.</span><span class="p">])</span>
575                                    <span class="n">VectA</span><span class="o">.</span><span class="n">append</span><span class="p">([])</span>
576        <span class="k">for</span> <span class="n">D</span> <span class="ow">in</span> <span class="n">Dist</span><span class="p">:</span>
577            <span class="n">DistArray</span><span class="p">[</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">0</span><span class="p">]]</span><span class="o">.</span><span class="n">append</span><span class="p">(</span><span class="n">D</span><span class="p">[</span><span class="mi">1</span><span class="p">:])</span>
578        <span class="n">Vect</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">array</span><span class="p">(</span><span class="n">Vect</span><span class="p">)</span>
579        <span class="n">angles</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">zeros</span><span class="p">((</span><span class="nb">len</span><span class="p">(</span><span class="n">Vect</span><span class="p">),</span><span class="nb">len</span><span class="p">(</span><span class="n">Vect</span><span class="p">)))</span>
580        <span class="n">angsig</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">zeros</span><span class="p">((</span><span class="nb">len</span><span class="p">(</span><span class="n">Vect</span><span class="p">),</span><span class="nb">len</span><span class="p">(</span><span class="n">Vect</span><span class="p">)))</span>
581        <span class="k">for</span> <span class="n">i</span><span class="p">,</span><span class="n">veca</span> <span class="ow">in</span> <span class="nb">enumerate</span><span class="p">(</span><span class="n">Vect</span><span class="p">):</span>
582            <span class="k">if</span> <span class="n">np</span><span class="o">.</span><span class="n">any</span><span class="p">(</span><span class="n">veca</span><span class="p">):</span>
583                <span class="k">for</span> <span class="n">j</span><span class="p">,</span><span class="n">vecb</span> <span class="ow">in</span> <span class="nb">enumerate</span><span class="p">(</span><span class="n">Vect</span><span class="p">):</span>
584                    <span class="k">if</span> <span class="n">np</span><span class="o">.</span><span class="n">any</span><span class="p">(</span><span class="n">vecb</span><span class="p">):</span>
585                        <span class="n">angles</span><span class="p">[</span><span class="n">i</span><span class="p">][</span><span class="n">j</span><span class="p">],</span><span class="n">angsig</span><span class="p">[</span><span class="n">i</span><span class="p">][</span><span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">G2mth</span><span class="o">.</span><span class="n">getAngSig</span><span class="p">(</span><span class="n">VectA</span><span class="p">[</span><span class="n">i</span><span class="p">],</span><span class="n">VectA</span><span class="p">[</span><span class="n">j</span><span class="p">],</span><span class="n">Amat</span><span class="p">,</span><span class="n">SGData</span><span class="p">,</span><span class="n">covData</span><span class="p">)</span>
586                        <span class="k">if</span> <span class="n">i</span> <span class="o">&lt;=</span> <span class="n">j</span><span class="p">:</span> <span class="k">continue</span>
587                        <span class="n">AngArray</span><span class="p">[</span><span class="n">Oatom</span><span class="p">[</span><span class="mi">0</span><span class="p">]]</span><span class="o">.</span><span class="n">append</span><span class="p">((</span><span class="n">i</span><span class="p">,</span><span class="n">j</span><span class="p">,</span>
588                            <span class="n">G2mth</span><span class="o">.</span><span class="n">getAngSig</span><span class="p">(</span><span class="n">VectA</span><span class="p">[</span><span class="n">i</span><span class="p">],</span><span class="n">VectA</span><span class="p">[</span><span class="n">j</span><span class="p">],</span><span class="n">Amat</span><span class="p">,</span><span class="n">SGData</span><span class="p">,</span><span class="n">covData</span><span class="p">)))</span>
589    <span class="k">return</span> <span class="n">AtomLabels</span><span class="p">,</span><span class="n">DistArray</span><span class="p">,</span><span class="n">AngArray</span>
590</div>
591<div class="viewcode-block" id="PrintDistAngle"><a class="viewcode-back" href="../GSASIIstruc.html#GSASIIstrMain.PrintDistAngle">[docs]</a><span class="k">def</span> <span class="nf">PrintDistAngle</span><span class="p">(</span><span class="n">DisAglCtls</span><span class="p">,</span><span class="n">DisAglData</span><span class="p">,</span><span class="n">out</span><span class="o">=</span><span class="n">sys</span><span class="o">.</span><span class="n">stdout</span><span class="p">):</span>
592    <span class="sd">&#39;&#39;&#39;Print distances and angles</span>
593
594<span class="sd">    :param dict DisAglCtls: contains distance/angle radii usually defined using</span>
595<span class="sd">       :func:`GSASIIgrid.DisAglDialog`</span>
596<span class="sd">    :param dict DisAglData: contains phase data: </span>
597<span class="sd">       Items &#39;OrigAtoms&#39; and &#39;TargAtoms&#39; contain the atoms to be used</span>
598<span class="sd">       for distance/angle origins and atoms to be used as targets.</span>
599<span class="sd">       Item &#39;SGData&#39; has the space group information (see :ref:`Space Group object&lt;SGData_table&gt;`)</span>
600<span class="sd">    :param file out: file object for output. Defaults to sys.stdout.    </span>
601<span class="sd">    &#39;&#39;&#39;</span>
602    <span class="kn">import</span> <span class="nn">numpy.ma</span> <span class="kn">as</span> <span class="nn">ma</span>
603    <span class="k">def</span> <span class="nf">MyPrint</span><span class="p">(</span><span class="n">s</span><span class="p">):</span>
604        <span class="n">out</span><span class="o">.</span><span class="n">write</span><span class="p">(</span><span class="n">s</span><span class="o">+</span><span class="s">&#39;</span><span class="se">\n</span><span class="s">&#39;</span><span class="p">)</span>
605        <span class="c"># print(s,file=out) # use in Python 3</span>
606   
607    <span class="k">def</span> <span class="nf">ShowBanner</span><span class="p">(</span><span class="n">name</span><span class="p">):</span>
608        <span class="n">MyPrint</span><span class="p">(</span><span class="mi">80</span><span class="o">*</span><span class="s">&#39;*&#39;</span><span class="p">)</span>
609        <span class="n">MyPrint</span><span class="p">(</span><span class="s">&#39;   Interatomic Distances and Angles for phase &#39;</span><span class="o">+</span><span class="n">name</span><span class="p">)</span>
610        <span class="n">MyPrint</span><span class="p">((</span><span class="mi">80</span><span class="o">*</span><span class="s">&#39;*&#39;</span><span class="p">)</span><span class="o">+</span><span class="s">&#39;</span><span class="se">\n</span><span class="s">&#39;</span><span class="p">)</span>
611
612    <span class="n">ShowBanner</span><span class="p">(</span><span class="n">DisAglCtls</span><span class="p">[</span><span class="s">&#39;Name&#39;</span><span class="p">])</span>
613    <span class="n">SGData</span> <span class="o">=</span> <span class="n">DisAglData</span><span class="p">[</span><span class="s">&#39;SGData&#39;</span><span class="p">]</span>
614    <span class="n">SGtext</span><span class="p">,</span><span class="n">SGtable</span> <span class="o">=</span> <span class="n">G2spc</span><span class="o">.</span><span class="n">SGPrint</span><span class="p">(</span><span class="n">SGData</span><span class="p">)</span>
615    <span class="k">for</span> <span class="n">line</span> <span class="ow">in</span> <span class="n">SGtext</span><span class="p">:</span> <span class="n">MyPrint</span><span class="p">(</span><span class="n">line</span><span class="p">)</span>
616    <span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">SGtable</span><span class="p">):</span>
617        <span class="k">for</span> <span class="n">i</span><span class="p">,</span><span class="n">item</span> <span class="ow">in</span> <span class="nb">enumerate</span><span class="p">(</span><span class="n">SGtable</span><span class="p">[::</span><span class="mi">2</span><span class="p">]):</span>
618            <span class="n">line</span> <span class="o">=</span> <span class="s">&#39; </span><span class="si">%s</span><span class="s"> </span><span class="si">%s</span><span class="s">&#39;</span><span class="o">%</span><span class="p">(</span><span class="n">item</span><span class="o">.</span><span class="n">ljust</span><span class="p">(</span><span class="mi">30</span><span class="p">),</span><span class="n">SGtable</span><span class="p">[</span><span class="mi">2</span><span class="o">*</span><span class="n">i</span><span class="o">+</span><span class="mi">1</span><span class="p">]</span><span class="o">.</span><span class="n">ljust</span><span class="p">(</span><span class="mi">30</span><span class="p">))</span>
619            <span class="n">MyPrint</span><span class="p">(</span><span class="n">line</span><span class="p">)</span>   
620    <span class="k">else</span><span class="p">:</span>
621        <span class="n">MyPrint</span><span class="p">(</span><span class="s">&#39; ( 1)    </span><span class="si">%s</span><span class="s">&#39;</span><span class="o">%</span><span class="p">(</span><span class="n">SGtable</span><span class="p">[</span><span class="mi">0</span><span class="p">]))</span> 
622    <span class="n">Cell</span> <span class="o">=</span> <span class="n">DisAglData</span><span class="p">[</span><span class="s">&#39;Cell&#39;</span><span class="p">]</span>
623   
624    <span class="n">Amat</span><span class="p">,</span><span class="n">Bmat</span> <span class="o">=</span> <span class="n">G2lat</span><span class="o">.</span><span class="n">cell2AB</span><span class="p">(</span><span class="n">Cell</span><span class="p">[:</span><span class="mi">6</span><span class="p">])</span>
625    <span class="n">covData</span> <span class="o">=</span> <span class="p">{}</span>
626    <span class="k">if</span> <span class="s">&#39;covData&#39;</span> <span class="ow">in</span> <span class="n">DisAglData</span><span class="p">:</span>   
627        <span class="n">covData</span> <span class="o">=</span> <span class="n">DisAglData</span><span class="p">[</span><span class="s">&#39;covData&#39;</span><span class="p">]</span>
628        <span class="n">covMatrix</span> <span class="o">=</span> <span class="n">covData</span><span class="p">[</span><span class="s">&#39;covMatrix&#39;</span><span class="p">]</span>
629        <span class="n">varyList</span> <span class="o">=</span> <span class="n">covData</span><span class="p">[</span><span class="s">&#39;varyList&#39;</span><span class="p">]</span>
630        <span class="n">pfx</span> <span class="o">=</span> <span class="nb">str</span><span class="p">(</span><span class="n">DisAglData</span><span class="p">[</span><span class="s">&#39;pId&#39;</span><span class="p">])</span><span class="o">+</span><span class="s">&#39;::&#39;</span>
631        <span class="n">A</span> <span class="o">=</span> <span class="n">G2lat</span><span class="o">.</span><span class="n">cell2A</span><span class="p">(</span><span class="n">Cell</span><span class="p">[:</span><span class="mi">6</span><span class="p">])</span>
632        <span class="n">cellSig</span> <span class="o">=</span> <span class="n">G2stIO</span><span class="o">.</span><span class="n">getCellEsd</span><span class="p">(</span><span class="n">pfx</span><span class="p">,</span><span class="n">SGData</span><span class="p">,</span><span class="n">A</span><span class="p">,</span><span class="n">covData</span><span class="p">)</span>
633        <span class="n">names</span