1242 lines
70 KiB
HTML
1242 lines
70 KiB
HTML
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "https://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
|
||
<html xmlns="http://www.w3.org/1999/xhtml" lang="en-US">
|
||
<head>
|
||
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
|
||
<meta http-equiv="X-UA-Compatible" content="IE=11"/>
|
||
<meta name="generator" content="Doxygen 1.13.2"/>
|
||
<meta name="viewport" content="width=device-width, initial-scale=1"/>
|
||
<title>musrfit: PFTPhaseCorrection Class Reference</title>
|
||
<link href="tabs.css" rel="stylesheet" type="text/css"/>
|
||
<script type="text/javascript" src="jquery.js"></script>
|
||
<script type="text/javascript" src="dynsections.js"></script>
|
||
<script type="text/javascript" src="clipboard.js"></script>
|
||
<link href="navtree.css" rel="stylesheet" type="text/css"/>
|
||
<script type="text/javascript" src="navtreedata.js"></script>
|
||
<script type="text/javascript" src="navtree.js"></script>
|
||
<script type="text/javascript" src="resize.js"></script>
|
||
<script type="text/javascript" src="cookie.js"></script>
|
||
<link href="doxygen.css" rel="stylesheet" type="text/css" />
|
||
</head>
|
||
<body>
|
||
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
|
||
<div id="titlearea">
|
||
<table cellspacing="0" cellpadding="0">
|
||
<tbody>
|
||
<tr id="projectrow">
|
||
<td id="projectalign">
|
||
<div id="projectname">musrfit<span id="projectnumber"> 1.9.9</span>
|
||
</div>
|
||
</td>
|
||
</tr>
|
||
</tbody>
|
||
</table>
|
||
</div>
|
||
<!-- end header part -->
|
||
<!-- Generated by Doxygen 1.13.2 -->
|
||
<script type="text/javascript">
|
||
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
|
||
$(function() { codefold.init(0); });
|
||
/* @license-end */
|
||
</script>
|
||
<script type="text/javascript" src="menudata.js"></script>
|
||
<script type="text/javascript" src="menu.js"></script>
|
||
<script type="text/javascript">
|
||
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
|
||
$(function() {
|
||
initMenu('',false,false,'search.php','Search',true);
|
||
});
|
||
/* @license-end */
|
||
</script>
|
||
<div id="main-nav"></div>
|
||
</div><!-- top -->
|
||
<div id="side-nav" class="ui-resizable side-nav-resizable">
|
||
<div id="nav-tree">
|
||
<div id="nav-tree-contents">
|
||
<div id="nav-sync" class="sync"></div>
|
||
</div>
|
||
</div>
|
||
<div id="splitbar" style="-moz-user-select:none;"
|
||
class="ui-resizable-handle">
|
||
</div>
|
||
</div>
|
||
<script type="text/javascript">
|
||
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
|
||
$(function(){initNavTree('classPFTPhaseCorrection.html',''); initResizable(true); });
|
||
/* @license-end */
|
||
</script>
|
||
<div id="doc-content">
|
||
<div class="header">
|
||
<div class="summary">
|
||
<a href="#pub-methods">Public Member Functions</a> |
|
||
<a href="#pri-methods">Private Member Functions</a> |
|
||
<a href="#pri-attribs">Private Attributes</a> |
|
||
<a href="classPFTPhaseCorrection-members.html">List of all members</a> </div>
|
||
<div class="headertitle"><div class="title">PFTPhaseCorrection Class Reference</div></div>
|
||
</div><!--header-->
|
||
<div class="contents">
|
||
|
||
<p><code>#include <<a class="el" href="PFourier_8h_source.html">PFourier.h</a>></code></p>
|
||
<div class="dynheader">
|
||
Inheritance diagram for PFTPhaseCorrection:</div>
|
||
<div class="dyncontent">
|
||
<div class="center"><img src="classPFTPhaseCorrection__inherit__graph.png" border="0" usemap="#aPFTPhaseCorrection_inherit__map" alt="Inheritance graph"/></div>
|
||
<map name="aPFTPhaseCorrection_inherit__map" id="aPFTPhaseCorrection_inherit__map">
|
||
<area shape="rect" title=" " alt="" coords="17,81,162,109"/>
|
||
<area shape="rect" title=" " alt="" coords="5,5,174,33"/>
|
||
<area shape="poly" title=" " alt="" coords="92,49,92,81,87,81,87,49"/>
|
||
</map>
|
||
<center><span class="legend">[<a target="top" href="graph_legend.html">legend</a>]</span></center></div>
|
||
<div class="dynheader">
|
||
Collaboration diagram for PFTPhaseCorrection:</div>
|
||
<div class="dyncontent">
|
||
<div class="center"><img src="classPFTPhaseCorrection__coll__graph.png" border="0" usemap="#aPFTPhaseCorrection_coll__map" alt="Collaboration graph"/></div>
|
||
<map name="aPFTPhaseCorrection_coll__map" id="aPFTPhaseCorrection_coll__map">
|
||
<area shape="rect" title=" " alt="" coords="17,81,162,109"/>
|
||
<area shape="rect" title=" " alt="" coords="5,5,174,33"/>
|
||
<area shape="poly" title=" " alt="" coords="92,49,92,81,87,81,87,49"/>
|
||
</map>
|
||
<center><span class="legend">[<a target="top" href="graph_legend.html">legend</a>]</span></center></div>
|
||
<table class="memberdecls">
|
||
<tr class="heading"><td colspan="2"><h2 class="groupheader"><a id="pub-methods" name="pub-methods"></a>
|
||
Public Member Functions</h2></td></tr>
|
||
<tr class="memitem:a073fbb2695d7132a157862c6faead778" id="r_a073fbb2695d7132a157862c6faead778"><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="#a073fbb2695d7132a157862c6faead778">PFTPhaseCorrection</a> (const Int_t minBin=-1, const Int_t maxBin=-1)</td></tr>
|
||
<tr class="separator:a073fbb2695d7132a157862c6faead778"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:ada87329c4c8775cc000eeb4268781ffc" id="r_ada87329c4c8775cc000eeb4268781ffc"><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="#ada87329c4c8775cc000eeb4268781ffc">PFTPhaseCorrection</a> (std::vector< Double_t > &reFT, std::vector< Double_t > &imFT, const Int_t minBin=-1, const Int_t maxBin=-1)</td></tr>
|
||
<tr class="separator:ada87329c4c8775cc000eeb4268781ffc"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a46a54ad30dbde0ce9df26eadc17860c1" id="r_a46a54ad30dbde0ce9df26eadc17860c1"><td class="memItemLeft" align="right" valign="top">virtual </td><td class="memItemRight" valign="bottom"><a class="el" href="#a46a54ad30dbde0ce9df26eadc17860c1">~PFTPhaseCorrection</a> ()</td></tr>
|
||
<tr class="separator:a46a54ad30dbde0ce9df26eadc17860c1"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a4a16bfd3cf495e2ee83026541ce14621" id="r_a4a16bfd3cf495e2ee83026541ce14621"><td class="memItemLeft" align="right" valign="top">virtual Bool_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#a4a16bfd3cf495e2ee83026541ce14621">IsValid</a> ()</td></tr>
|
||
<tr class="separator:a4a16bfd3cf495e2ee83026541ce14621"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a05863bb58d111803ab36b922700bb667" id="r_a05863bb58d111803ab36b922700bb667"><td class="memItemLeft" align="right" valign="top">virtual void </td><td class="memItemRight" valign="bottom"><a class="el" href="#a05863bb58d111803ab36b922700bb667">Minimize</a> ()</td></tr>
|
||
<tr class="separator:a05863bb58d111803ab36b922700bb667"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a3c4e1295338f8110e91bef28abe430c0" id="r_a3c4e1295338f8110e91bef28abe430c0"><td class="memItemLeft" align="right" valign="top">virtual void </td><td class="memItemRight" valign="bottom"><a class="el" href="#a3c4e1295338f8110e91bef28abe430c0">SetGamma</a> (const Double_t gamma)</td></tr>
|
||
<tr class="separator:a3c4e1295338f8110e91bef28abe430c0"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a70069a47f75e2ec309d278f9de2def7c" id="r_a70069a47f75e2ec309d278f9de2def7c"><td class="memItemLeft" align="right" valign="top">virtual void </td><td class="memItemRight" valign="bottom"><a class="el" href="#a70069a47f75e2ec309d278f9de2def7c">SetPh</a> (const Double_t c0, const Double_t c1)</td></tr>
|
||
<tr class="separator:a70069a47f75e2ec309d278f9de2def7c"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:ae84c2e25aa0d68e214c57f2e3066a38e" id="r_ae84c2e25aa0d68e214c57f2e3066a38e"><td class="memItemLeft" align="right" valign="top">virtual Double_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#ae84c2e25aa0d68e214c57f2e3066a38e">GetGamma</a> ()</td></tr>
|
||
<tr class="separator:ae84c2e25aa0d68e214c57f2e3066a38e"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:abdac531b17bda660e3648eee992ed4ca" id="r_abdac531b17bda660e3648eee992ed4ca"><td class="memItemLeft" align="right" valign="top">virtual Double_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#abdac531b17bda660e3648eee992ed4ca">GetPhaseCorrectionParam</a> (UInt_t idx)</td></tr>
|
||
<tr class="separator:abdac531b17bda660e3648eee992ed4ca"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:ab40b23635a454d645fb4861b8417319f" id="r_ab40b23635a454d645fb4861b8417319f"><td class="memItemLeft" align="right" valign="top">virtual Double_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#ab40b23635a454d645fb4861b8417319f">GetMinimum</a> ()</td></tr>
|
||
<tr class="separator:ab40b23635a454d645fb4861b8417319f"><td class="memSeparator" colspan="2"> </td></tr>
|
||
</table><table class="memberdecls">
|
||
<tr class="heading"><td colspan="2"><h2 class="groupheader"><a id="pri-methods" name="pri-methods"></a>
|
||
Private Member Functions</h2></td></tr>
|
||
<tr class="memitem:a2551fcabb21bc94ec1dd999ad4d81313" id="r_a2551fcabb21bc94ec1dd999ad4d81313"><td class="memItemLeft" align="right" valign="top">virtual void </td><td class="memItemRight" valign="bottom"><a class="el" href="#a2551fcabb21bc94ec1dd999ad4d81313">Init</a> ()</td></tr>
|
||
<tr class="memdesc:a2551fcabb21bc94ec1dd999ad4d81313"><td class="mdescLeft"> </td><td class="mdescRight">keeps the minimum of the entropy/penalty minimization <br /></td></tr>
|
||
<tr class="separator:a2551fcabb21bc94ec1dd999ad4d81313"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a66fe9679aac2280c7adcdaf5e2c36b54" id="r_a66fe9679aac2280c7adcdaf5e2c36b54"><td class="memItemLeft" align="right" valign="top">virtual void </td><td class="memItemRight" valign="bottom"><a class="el" href="#a66fe9679aac2280c7adcdaf5e2c36b54">CalcPhasedFT</a> () const</td></tr>
|
||
<tr class="separator:a66fe9679aac2280c7adcdaf5e2c36b54"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a10d666f39e1ca8b602fdb295ef2ba7a9" id="r_a10d666f39e1ca8b602fdb295ef2ba7a9"><td class="memItemLeft" align="right" valign="top">virtual void </td><td class="memItemRight" valign="bottom"><a class="el" href="#a10d666f39e1ca8b602fdb295ef2ba7a9">CalcRealPhFTDerivative</a> () const</td></tr>
|
||
<tr class="separator:a10d666f39e1ca8b602fdb295ef2ba7a9"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:aa65f5dd4a52f0cb2e3fd8bba7b4a2945" id="r_aa65f5dd4a52f0cb2e3fd8bba7b4a2945"><td class="memItemLeft" align="right" valign="top">virtual Double_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#aa65f5dd4a52f0cb2e3fd8bba7b4a2945">Penalty</a> () const</td></tr>
|
||
<tr class="separator:aa65f5dd4a52f0cb2e3fd8bba7b4a2945"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a4e8ddb15191632e59f4f3fcd86932199" id="r_a4e8ddb15191632e59f4f3fcd86932199"><td class="memItemLeft" align="right" valign="top">virtual Double_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#a4e8ddb15191632e59f4f3fcd86932199">Entropy</a> () const</td></tr>
|
||
<tr class="separator:a4e8ddb15191632e59f4f3fcd86932199"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a1e2255b2b4d7644ff6ee489d5ce3082f" id="r_a1e2255b2b4d7644ff6ee489d5ce3082f"><td class="memItemLeft" align="right" valign="top">virtual Double_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#a1e2255b2b4d7644ff6ee489d5ce3082f">Up</a> () const</td></tr>
|
||
<tr class="separator:a1e2255b2b4d7644ff6ee489d5ce3082f"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a6e4cc9da6c28cb1816b14a174e70694a" id="r_a6e4cc9da6c28cb1816b14a174e70694a"><td class="memItemLeft" align="right" valign="top">virtual Double_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#a6e4cc9da6c28cb1816b14a174e70694a">operator()</a> (const std::vector< Double_t > &) const</td></tr>
|
||
<tr class="separator:a6e4cc9da6c28cb1816b14a174e70694a"><td class="memSeparator" colspan="2"> </td></tr>
|
||
</table><table class="memberdecls">
|
||
<tr class="heading"><td colspan="2"><h2 class="groupheader"><a id="pri-attribs" name="pri-attribs"></a>
|
||
Private Attributes</h2></td></tr>
|
||
<tr class="memitem:a1d20df49ab2b5bba162eab9e0088396a" id="r_a1d20df49ab2b5bba162eab9e0088396a"><td class="memItemLeft" align="right" valign="top">Bool_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#a1d20df49ab2b5bba162eab9e0088396a">fValid</a></td></tr>
|
||
<tr class="separator:a1d20df49ab2b5bba162eab9e0088396a"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a4f19633e8a329c6f54e12683d1b357e2" id="r_a4f19633e8a329c6f54e12683d1b357e2"><td class="memItemLeft" align="right" valign="top">std::vector< Double_t > </td><td class="memItemRight" valign="bottom"><a class="el" href="#a4f19633e8a329c6f54e12683d1b357e2">fReal</a></td></tr>
|
||
<tr class="separator:a4f19633e8a329c6f54e12683d1b357e2"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:af424a8bd0dc0ef16a2537d0ac2ec0abb" id="r_af424a8bd0dc0ef16a2537d0ac2ec0abb"><td class="memItemLeft" align="right" valign="top">std::vector< Double_t > </td><td class="memItemRight" valign="bottom"><a class="el" href="#af424a8bd0dc0ef16a2537d0ac2ec0abb">fImag</a></td></tr>
|
||
<tr class="memdesc:af424a8bd0dc0ef16a2537d0ac2ec0abb"><td class="mdescLeft"> </td><td class="mdescRight">original real Fourier data set <br /></td></tr>
|
||
<tr class="separator:af424a8bd0dc0ef16a2537d0ac2ec0abb"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a1fa993c00f55291ed2b58c17ecb96d4b" id="r_a1fa993c00f55291ed2b58c17ecb96d4b"><td class="memItemLeft" align="right" valign="top">std::vector< Double_t > </td><td class="memItemRight" valign="bottom"><a class="el" href="#a1fa993c00f55291ed2b58c17ecb96d4b">fRealPh</a></td></tr>
|
||
<tr class="memdesc:a1fa993c00f55291ed2b58c17ecb96d4b"><td class="mdescLeft"> </td><td class="mdescRight">original imag Fourier data set <br /></td></tr>
|
||
<tr class="separator:a1fa993c00f55291ed2b58c17ecb96d4b"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a4eb80579fb232cdd5dc76cc99101192b" id="r_a4eb80579fb232cdd5dc76cc99101192b"><td class="memItemLeft" align="right" valign="top">std::vector< Double_t > </td><td class="memItemRight" valign="bottom"><a class="el" href="#a4eb80579fb232cdd5dc76cc99101192b">fImagPh</a></td></tr>
|
||
<tr class="memdesc:a4eb80579fb232cdd5dc76cc99101192b"><td class="mdescLeft"> </td><td class="mdescRight">phased real Fourier data set <br /></td></tr>
|
||
<tr class="separator:a4eb80579fb232cdd5dc76cc99101192b"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a48c629e27e8f5f8832c921d14d9a1b62" id="r_a48c629e27e8f5f8832c921d14d9a1b62"><td class="memItemLeft" align="right" valign="top">std::vector< Double_t > </td><td class="memItemRight" valign="bottom"><a class="el" href="#a48c629e27e8f5f8832c921d14d9a1b62">fRealPhD</a></td></tr>
|
||
<tr class="memdesc:a48c629e27e8f5f8832c921d14d9a1b62"><td class="mdescLeft"> </td><td class="mdescRight">phased imag Fourier data set <br /></td></tr>
|
||
<tr class="separator:a48c629e27e8f5f8832c921d14d9a1b62"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a1f1c220a61b47d89df38a57a400867cc" id="r_a1f1c220a61b47d89df38a57a400867cc"><td class="memItemLeft" align="right" valign="top">Int_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#a1f1c220a61b47d89df38a57a400867cc">fMinBin</a></td></tr>
|
||
<tr class="memdesc:a1f1c220a61b47d89df38a57a400867cc"><td class="mdescLeft"> </td><td class="mdescRight">1st derivative of fRealPh <br /></td></tr>
|
||
<tr class="separator:a1f1c220a61b47d89df38a57a400867cc"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:abcd9de24b2045dba8af6c534976d4417" id="r_abcd9de24b2045dba8af6c534976d4417"><td class="memItemLeft" align="right" valign="top">Int_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#abcd9de24b2045dba8af6c534976d4417">fMaxBin</a></td></tr>
|
||
<tr class="memdesc:abcd9de24b2045dba8af6c534976d4417"><td class="mdescLeft"> </td><td class="mdescRight">minimum bin from Fourier range to be used for the phase correction estimate <br /></td></tr>
|
||
<tr class="separator:abcd9de24b2045dba8af6c534976d4417"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a9e911366ef26bbd1b510ab85e2496e7e" id="r_a9e911366ef26bbd1b510ab85e2496e7e"><td class="memItemLeft" align="right" valign="top">Double_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#a9e911366ef26bbd1b510ab85e2496e7e">fPh_c0</a></td></tr>
|
||
<tr class="memdesc:a9e911366ef26bbd1b510ab85e2496e7e"><td class="mdescLeft"> </td><td class="mdescRight">maximum bin from Fourier range to be used for the phase correction estimate <br /></td></tr>
|
||
<tr class="separator:a9e911366ef26bbd1b510ab85e2496e7e"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a6f8a487092c8d8e5b2e10c797ae6368f" id="r_a6f8a487092c8d8e5b2e10c797ae6368f"><td class="memItemLeft" align="right" valign="top">Double_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#a6f8a487092c8d8e5b2e10c797ae6368f">fPh_c1</a></td></tr>
|
||
<tr class="memdesc:a6f8a487092c8d8e5b2e10c797ae6368f"><td class="mdescLeft"> </td><td class="mdescRight">constant part of the phase dispersion used for the phase correction <br /></td></tr>
|
||
<tr class="separator:a6f8a487092c8d8e5b2e10c797ae6368f"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a8e5d6c8081d596f892034ce74081c34c" id="r_a8e5d6c8081d596f892034ce74081c34c"><td class="memItemLeft" align="right" valign="top">Double_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#a8e5d6c8081d596f892034ce74081c34c">fGamma</a></td></tr>
|
||
<tr class="memdesc:a8e5d6c8081d596f892034ce74081c34c"><td class="mdescLeft"> </td><td class="mdescRight">linear part of the phase dispersion used for the phase correction <br /></td></tr>
|
||
<tr class="separator:a8e5d6c8081d596f892034ce74081c34c"><td class="memSeparator" colspan="2"> </td></tr>
|
||
<tr class="memitem:a7d2ee1776fd5cee6cfc61cfd48409750" id="r_a7d2ee1776fd5cee6cfc61cfd48409750"><td class="memItemLeft" align="right" valign="top">Double_t </td><td class="memItemRight" valign="bottom"><a class="el" href="#a7d2ee1776fd5cee6cfc61cfd48409750">fMin</a></td></tr>
|
||
<tr class="memdesc:a7d2ee1776fd5cee6cfc61cfd48409750"><td class="mdescLeft"> </td><td class="mdescRight">gamma parameter to balance between entropy and penalty <br /></td></tr>
|
||
<tr class="separator:a7d2ee1776fd5cee6cfc61cfd48409750"><td class="memSeparator" colspan="2"> </td></tr>
|
||
</table>
|
||
<a name="details" id="details"></a><h2 class="groupheader">Detailed Description</h2>
|
||
<div class="textblock"><p>Phase correction optimizer for Fourier transforms.</p>
|
||
<p>This class performs automatic phase correction on complex Fourier spectra to maximize the real component and minimize the imaginary component. Phase errors arise from:</p><ul>
|
||
<li>Uncertain time-zero determination</li>
|
||
<li>Detector time offsets</li>
|
||
<li>Signal dispersion</li>
|
||
</ul>
|
||
<p><b>Algorithm:</b> Minimizes a combined entropy-penalty functional using Minuit2, finding optimal phase parameters (constant + linear dispersion): φ(ω) = c₀ + c₁·ω</p>
|
||
<p><b>Applications:</b></p><ul>
|
||
<li>Improving signal clarity in real Fourier spectra</li>
|
||
<li>Identifying field distributions in vortex lattices</li>
|
||
<li>Resolving closely-spaced frequency components</li>
|
||
</ul>
|
||
<p><b>Usage:</b> Specify frequency range for optimization to focus on signal peaks while avoiding noise regions. </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00084">84</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
</div><h2 class="groupheader">Constructor & Destructor Documentation</h2>
|
||
<a id="a073fbb2695d7132a157862c6faead778" name="a073fbb2695d7132a157862c6faead778"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a073fbb2695d7132a157862c6faead778">◆ </a></span>PFTPhaseCorrection() <span class="overload">[1/2]</span></h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">PFTPhaseCorrection::PFTPhaseCorrection </td>
|
||
<td>(</td>
|
||
<td class="paramtype">const Int_t</td> <td class="paramname"><span class="paramname"><em>minBin</em></span><span class="paramdefsep"> = </span><span class="paramdefval">-1</span>, </td>
|
||
</tr>
|
||
<tr>
|
||
<td class="paramkey"></td>
|
||
<td></td>
|
||
<td class="paramtype">const Int_t</td> <td class="paramname"><span class="paramname"><em>maxBin</em></span><span class="paramdefsep"> = </span><span class="paramdefval">-1</span> )</td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Constructor for phase correction with default Fourier data.</p>
|
||
<dl class="params"><dt>Parameters</dt><dd>
|
||
<table class="params">
|
||
<tr><td class="paramname">minBin</td><td>Minimum frequency bin for optimization (-1 = use all) </td></tr>
|
||
<tr><td class="paramname">maxBin</td><td>Maximum frequency bin for optimization (-1 = use all)</td></tr>
|
||
</table>
|
||
</dd>
|
||
</dl>
|
||
<p>Default constructor for phase correction optimizer.</p>
|
||
<p>Initializes a phase correction object with default Fourier data. The actual Fourier data must be set externally before minimization.</p>
|
||
<dl class="params"><dt>Parameters</dt><dd>
|
||
<table class="params">
|
||
<tr><td class="paramname">minBin</td><td>Minimum frequency bin index for optimization region. Use -1 to optimize over all bins. Default: -1. </td></tr>
|
||
<tr><td class="paramname">maxBin</td><td>Maximum frequency bin index for optimization region. Use -1 to optimize over all bins. Default: -1.</td></tr>
|
||
</table>
|
||
</dd>
|
||
</dl>
|
||
<p><b>Note:</b> Restricting the optimization range to signal-containing regions improves results by excluding noisy baseline regions. </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00067">67</a> of file <a class="el" href="PFourier_8cpp_source.html">PFourier.cpp</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00154">fMaxBin</a>, <a class="el" href="PFourier_8h_source.html#l00153">fMinBin</a>, and <a class="el" href="PFourier_8cpp_source.html#l00248">Init()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="ada87329c4c8775cc000eeb4268781ffc" name="ada87329c4c8775cc000eeb4268781ffc"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#ada87329c4c8775cc000eeb4268781ffc">◆ </a></span>PFTPhaseCorrection() <span class="overload">[2/2]</span></h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">PFTPhaseCorrection::PFTPhaseCorrection </td>
|
||
<td>(</td>
|
||
<td class="paramtype">std::vector< Double_t > &</td> <td class="paramname"><span class="paramname"><em>reFT</em></span>, </td>
|
||
</tr>
|
||
<tr>
|
||
<td class="paramkey"></td>
|
||
<td></td>
|
||
<td class="paramtype">std::vector< Double_t > &</td> <td class="paramname"><span class="paramname"><em>imFT</em></span>, </td>
|
||
</tr>
|
||
<tr>
|
||
<td class="paramkey"></td>
|
||
<td></td>
|
||
<td class="paramtype">const Int_t</td> <td class="paramname"><span class="paramname"><em>minBin</em></span><span class="paramdefsep"> = </span><span class="paramdefval">-1</span>, </td>
|
||
</tr>
|
||
<tr>
|
||
<td class="paramkey"></td>
|
||
<td></td>
|
||
<td class="paramtype">const Int_t</td> <td class="paramname"><span class="paramname"><em>maxBin</em></span><span class="paramdefsep"> = </span><span class="paramdefval">-1</span> )</td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Constructor with explicit Fourier data.</p>
|
||
<dl class="params"><dt>Parameters</dt><dd>
|
||
<table class="params">
|
||
<tr><td class="paramname">reFT</td><td>Real part of Fourier transform </td></tr>
|
||
<tr><td class="paramname">imFT</td><td>Imaginary part of Fourier transform </td></tr>
|
||
<tr><td class="paramname">minBin</td><td>Minimum frequency bin for optimization </td></tr>
|
||
<tr><td class="paramname">maxBin</td><td>Maximum frequency bin for optimization</td></tr>
|
||
</table>
|
||
</dd>
|
||
</dl>
|
||
<p>Constructor with explicit Fourier transform data.</p>
|
||
<p>Creates a phase correction optimizer with pre-computed Fourier data. The optimizer will find phase parameters φ(ω) = c₀ + c₁·ω that maximize the real component while minimizing the imaginary component.</p>
|
||
<dl class="params"><dt>Parameters</dt><dd>
|
||
<table class="params">
|
||
<tr><td class="paramname">reFT</td><td>Real part of the Fourier transform (vector of amplitudes). </td></tr>
|
||
<tr><td class="paramname">imFT</td><td>Imaginary part of the Fourier transform (vector of amplitudes). </td></tr>
|
||
<tr><td class="paramname">minBin</td><td>Minimum frequency bin index for optimization region. Use -1 to optimize over all bins. Default: -1. </td></tr>
|
||
<tr><td class="paramname">maxBin</td><td>Maximum frequency bin index for optimization region. Use -1 to optimize over all bins. Default: -1.</td></tr>
|
||
</table>
|
||
</dd>
|
||
</dl>
|
||
<p><b>Usage example:</b> </p><div class="fragment"><div class="line">std::vector<Double_t> real, imag;</div>
|
||
<div class="line"><span class="comment">// ... fill real and imag with FFT results ...</span></div>
|
||
<div class="line"><a class="code hl_function" href="#a073fbb2695d7132a157862c6faead778">PFTPhaseCorrection</a> phCorr(real, imag, 10, 500);</div>
|
||
<div class="line">phCorr.Minimize();</div>
|
||
<div class="line">Double_t c0 = phCorr.GetPhaseCorrectionParam(0);</div>
|
||
<div class="line">Double_t c1 = phCorr.GetPhaseCorrectionParam(1);</div>
|
||
<div class="ttc" id="aclassPFTPhaseCorrection_html_a073fbb2695d7132a157862c6faead778"><div class="ttname"><a href="#a073fbb2695d7132a157862c6faead778">PFTPhaseCorrection::PFTPhaseCorrection</a></div><div class="ttdeci">PFTPhaseCorrection(const Int_t minBin=-1, const Int_t maxBin=-1)</div><div class="ttdef"><b>Definition</b> <a href="PFourier_8cpp_source.html#l00067">PFourier.cpp:67</a></div></div>
|
||
</div><!-- fragment -->
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00100">100</a> of file <a class="el" href="PFourier_8cpp_source.html">PFourier.cpp</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00148">fImag</a>, <a class="el" href="PFourier_8h_source.html#l00150">fImagPh</a>, <a class="el" href="PFourier_8h_source.html#l00154">fMaxBin</a>, <a class="el" href="PFourier_8h_source.html#l00153">fMinBin</a>, <a class="el" href="PFourier_8h_source.html#l00147">fReal</a>, <a class="el" href="PFourier_8h_source.html#l00149">fRealPh</a>, and <a class="el" href="PFourier_8cpp_source.html#l00248">Init()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a46a54ad30dbde0ce9df26eadc17860c1" name="a46a54ad30dbde0ce9df26eadc17860c1"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a46a54ad30dbde0ce9df26eadc17860c1">◆ </a></span>~PFTPhaseCorrection()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">virtual PFTPhaseCorrection::~PFTPhaseCorrection </td>
|
||
<td>(</td>
|
||
<td class="paramname"><span class="paramname"><em></em></span></td><td>)</td>
|
||
<td></td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel inline">inline</span><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00105">105</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<h2 class="groupheader">Member Function Documentation</h2>
|
||
<a id="a66fe9679aac2280c7adcdaf5e2c36b54" name="a66fe9679aac2280c7adcdaf5e2c36b54"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a66fe9679aac2280c7adcdaf5e2c36b54">◆ </a></span>CalcPhasedFT()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">void PFTPhaseCorrection::CalcPhasedFT </td>
|
||
<td>(</td>
|
||
<td class="paramname"><span class="paramname"><em></em></span></td><td>)</td>
|
||
<td> const</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Calculates phase-corrected Fourier transform.</p>
|
||
<p>Applies the phase correction φ(ω) = c₀ + c₁·(i/N) to the complex Fourier spectrum using rotation in the complex plane:</p><ul>
|
||
<li>F'_re(ω) = F_re(ω)·cos(φ) - F_im(ω)·sin(φ)</li>
|
||
<li>F'_im(ω) = F_re(ω)·sin(φ) + F_im(ω)·cos(φ)</li>
|
||
</ul>
|
||
<p>The corrected spectra are stored in fRealPh and fImagPh.</p>
|
||
<p><b>Physics:</b> This rotation compensates for:</p><ul>
|
||
<li>Time-zero offset (constant phase c₀)</li>
|
||
<li>Timing dispersion (linear phase c₁)</li>
|
||
</ul>
|
||
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#a10d666f39e1ca8b602fdb295ef2ba7a9">CalcRealPhFTDerivative()</a> </dd></dl>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00276">276</a> of file <a class="el" href="PFourier_8cpp_source.html">PFourier.cpp</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00148">fImag</a>, <a class="el" href="PFourier_8h_source.html#l00150">fImagPh</a>, <a class="el" href="PFourier_8h_source.html#l00155">fPh_c0</a>, <a class="el" href="PFourier_8h_source.html#l00156">fPh_c1</a>, <a class="el" href="PFourier_8h_source.html#l00147">fReal</a>, and <a class="el" href="PFourier_8h_source.html#l00149">fRealPh</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00425">operator()()</a>, and <a class="el" href="PFourier_8h_source.html#l00126">SetPh()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a10d666f39e1ca8b602fdb295ef2ba7a9" name="a10d666f39e1ca8b602fdb295ef2ba7a9"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a10d666f39e1ca8b602fdb295ef2ba7a9">◆ </a></span>CalcRealPhFTDerivative()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">void PFTPhaseCorrection::CalcRealPhFTDerivative </td>
|
||
<td>(</td>
|
||
<td class="paramname"><span class="paramname"><em></em></span></td><td>)</td>
|
||
<td> const</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Calculates first derivative of phase-corrected real Fourier spectrum.</p>
|
||
<p>Computes the finite difference derivative: dF/dω ≈ F(i+1) - F(i) for all interior points. Boundary points are set to 1.0.</p>
|
||
<p>The derivative is used in the entropy calculation, where smooth spectra (small derivatives) have lower entropy and are preferred by the optimization algorithm.</p>
|
||
<p><b>Purpose:</b> Entropy minimization favors phase corrections that produce smooth, well-resolved real spectra with minimal oscillations.</p>
|
||
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#a4e8ddb15191632e59f4f3fcd86932199">Entropy()</a> </dd></dl>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00307">307</a> of file <a class="el" href="PFourier_8cpp_source.html">PFourier.cpp</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00149">fRealPh</a>, and <a class="el" href="PFourier_8h_source.html#l00151">fRealPhD</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00425">operator()()</a>, and <a class="el" href="PFourier_8h_source.html#l00126">SetPh()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a4e8ddb15191632e59f4f3fcd86932199" name="a4e8ddb15191632e59f4f3fcd86932199"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a4e8ddb15191632e59f4f3fcd86932199">◆ </a></span>Entropy()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">Double_t PFTPhaseCorrection::Entropy </td>
|
||
<td>(</td>
|
||
<td class="paramname"><span class="paramname"><em></em></span></td><td>)</td>
|
||
<td> const</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Calculates Shannon entropy of the real spectrum derivative.</p>
|
||
<p>Computes: S = -Σ[p_i·ln(p_i)] where p_i = |dF/dω|_i / Σ|dF/dω|</p>
|
||
<p>This entropy measure quantifies the "smoothness" of the real Fourier spectrum:</p><ul>
|
||
<li><b>Low entropy:</b> Concentrated, smooth spectrum (desired)</li>
|
||
<li><b>High entropy:</b> Dispersed, noisy spectrum (undesired)</li>
|
||
</ul>
|
||
<p><b>Physical interpretation:</b> Correct phase alignment produces sharp, well-defined peaks with smooth baselines, resulting in concentrated derivatives and low entropy.</p>
|
||
<dl class="section return"><dt>Returns</dt><dd>Shannon entropy of normalized derivative distribution</dd></dl>
|
||
<p><b>Note:</b> Small values (< 10⁻¹⁵) are skipped to avoid numerical issues with log(0).</p>
|
||
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#a10d666f39e1ca8b602fdb295ef2ba7a9">CalcRealPhFTDerivative()</a> </dd>
|
||
<dd>
|
||
<a class="el" href="#aa65f5dd4a52f0cb2e3fd8bba7b4a2945">Penalty()</a> </dd></dl>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00376">376</a> of file <a class="el" href="PFourier_8cpp_source.html">PFourier.cpp</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00154">fMaxBin</a>, <a class="el" href="PFourier_8h_source.html#l00153">fMinBin</a>, and <a class="el" href="PFourier_8h_source.html#l00151">fRealPhD</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00425">operator()()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="ae84c2e25aa0d68e214c57f2e3066a38e" name="ae84c2e25aa0d68e214c57f2e3066a38e"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#ae84c2e25aa0d68e214c57f2e3066a38e">◆ </a></span>GetGamma()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">virtual Double_t PFTPhaseCorrection::GetGamma </td>
|
||
<td>(</td>
|
||
<td class="paramname"><span class="paramname"><em></em></span></td><td>)</td>
|
||
<td></td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel inline">inline</span><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Returns the gamma parameter </p><dl class="section return"><dt>Returns</dt><dd>Balancing factor between entropy and penalty </dd></dl>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00130">130</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00157">fGamma</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="ab40b23635a454d645fb4861b8417319f" name="ab40b23635a454d645fb4861b8417319f"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#ab40b23635a454d645fb4861b8417319f">◆ </a></span>GetMinimum()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">Double_t PFTPhaseCorrection::GetMinimum </td>
|
||
<td>(</td>
|
||
<td class="paramname"><span class="paramname"><em></em></span></td><td>)</td>
|
||
<td></td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Returns the minimum value of the optimization functional </p><dl class="section return"><dt>Returns</dt><dd>Minimum value achieved</dd></dl>
|
||
<p>Returns the minimum value of the optimization functional.</p>
|
||
<p>This value represents the combined entropy-penalty functional at the optimal phase parameters. Lower values indicate better phase correction. A value of -1.0 indicates failed minimization.</p>
|
||
<dl class="section return"><dt>Returns</dt><dd>Minimum functional value, or -1.0 if minimization failed.</dd></dl>
|
||
<p><b>Note:</b> Always check <a class="el" href="#a4a16bfd3cf495e2ee83026541ce14621">IsValid()</a> before using this value. </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00223">223</a> of file <a class="el" href="PFourier_8cpp_source.html">PFourier.cpp</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00158">fMin</a>, and <a class="el" href="PFourier_8h_source.html#l00145">fValid</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="abdac531b17bda660e3648eee992ed4ca" name="abdac531b17bda660e3648eee992ed4ca"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#abdac531b17bda660e3648eee992ed4ca">◆ </a></span>GetPhaseCorrectionParam()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">Double_t PFTPhaseCorrection::GetPhaseCorrectionParam </td>
|
||
<td>(</td>
|
||
<td class="paramtype">UInt_t</td> <td class="paramname"><span class="paramname"><em>idx</em></span></td><td>)</td>
|
||
<td></td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Gets phase correction parameter.</p>
|
||
<dl class="params"><dt>Parameters</dt><dd>
|
||
<table class="params">
|
||
<tr><td class="paramname">idx</td><td>Parameter index (0=c₀, 1=c₁) </td></tr>
|
||
</table>
|
||
</dd>
|
||
</dl>
|
||
<dl class="section return"><dt>Returns</dt><dd>Phase parameter value</dd></dl>
|
||
<p>Retrieves optimized phase correction parameter by index.</p>
|
||
<p>Returns the phase parameter determined by <a class="el" href="#a05863bb58d111803ab36b922700bb667">Minimize()</a>. The phase correction formula is: φ(ω) = c₀ + c₁·(ω/ω_max)</p>
|
||
<dl class="params"><dt>Parameters</dt><dd>
|
||
<table class="params">
|
||
<tr><td class="paramname">idx</td><td>Parameter index:<ul>
|
||
<li>0: Returns c₀ (constant phase offset in radians)</li>
|
||
<li>1: Returns c₁ (linear phase dispersion coefficient)</li>
|
||
</ul>
|
||
</td></tr>
|
||
</table>
|
||
</dd>
|
||
</dl>
|
||
<dl class="section return"><dt>Returns</dt><dd>Phase correction parameter value, or 0.0 if idx is invalid.</dd></dl>
|
||
<p><b>Example:</b> </p><div class="fragment"><div class="line">phCorr.Minimize();</div>
|
||
<div class="line"><span class="keywordflow">if</span> (phCorr.IsValid()) {</div>
|
||
<div class="line"> Double_t c0 = phCorr.GetPhaseCorrectionParam(0);</div>
|
||
<div class="line"> Double_t c1 = phCorr.GetPhaseCorrectionParam(1);</div>
|
||
<div class="line"> std::cout << <span class="stringliteral">"Phase: "</span> << c0 << <span class="stringliteral">" + "</span> << c1 << <span class="stringliteral">"*w"</span> << std::endl;</div>
|
||
<div class="line">}</div>
|
||
</div><!-- fragment -->
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00195">195</a> of file <a class="el" href="PFourier_8cpp_source.html">PFourier.cpp</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00155">fPh_c0</a>, and <a class="el" href="PFourier_8h_source.html#l00156">fPh_c1</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00817">PFourier::GetPhaseOptRealFourier()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a2551fcabb21bc94ec1dd999ad4d81313" name="a2551fcabb21bc94ec1dd999ad4d81313"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a2551fcabb21bc94ec1dd999ad4d81313">◆ </a></span>Init()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">void PFTPhaseCorrection::Init </td>
|
||
<td>(</td>
|
||
<td class="paramname"><span class="paramname"><em></em></span></td><td>)</td>
|
||
<td></td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p>keeps the minimum of the entropy/penalty minimization </p>
|
||
<p>Initializes phase correction object with default values.</p>
|
||
<p>Sets initial state:</p><ul>
|
||
<li>fValid = true (object ready for use)</li>
|
||
<li>fPh_c0 = 0.0 (no constant phase offset)</li>
|
||
<li>fPh_c1 = 0.0 (no linear phase dispersion)</li>
|
||
<li>fGamma = 1.0 (equal weighting of entropy and penalty)</li>
|
||
<li>fMin = -1.0 (no minimization performed yet)</li>
|
||
</ul>
|
||
<p>Called by constructors to establish consistent initial state. </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00248">248</a> of file <a class="el" href="PFourier_8cpp_source.html">PFourier.cpp</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00157">fGamma</a>, <a class="el" href="PFourier_8h_source.html#l00158">fMin</a>, <a class="el" href="PFourier_8h_source.html#l00155">fPh_c0</a>, <a class="el" href="PFourier_8h_source.html#l00156">fPh_c1</a>, and <a class="el" href="PFourier_8h_source.html#l00145">fValid</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00067">PFTPhaseCorrection()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00100">PFTPhaseCorrection()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a4a16bfd3cf495e2ee83026541ce14621" name="a4a16bfd3cf495e2ee83026541ce14621"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a4a16bfd3cf495e2ee83026541ce14621">◆ </a></span>IsValid()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">virtual Bool_t PFTPhaseCorrection::IsValid </td>
|
||
<td>(</td>
|
||
<td class="paramname"><span class="paramname"><em></em></span></td><td>)</td>
|
||
<td></td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel inline">inline</span><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Returns true if phase correction initialized successfully </p><dl class="section return"><dt>Returns</dt><dd>Validity status </dd></dl>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00109">109</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00145">fValid</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00817">PFourier::GetPhaseOptRealFourier()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a05863bb58d111803ab36b922700bb667" name="a05863bb58d111803ab36b922700bb667"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a05863bb58d111803ab36b922700bb667">◆ </a></span>Minimize()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">void PFTPhaseCorrection::Minimize </td>
|
||
<td>(</td>
|
||
<td class="paramname"><span class="paramname"><em></em></span></td><td>)</td>
|
||
<td></td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Performs phase correction minimization.</p>
|
||
<p>Uses Minuit2 to find optimal phase parameters that maximize the real spectrum while minimizing imaginary components.</p>
|
||
<p>Performs phase correction minimization using Minuit2.</p>
|
||
<p>This method finds optimal phase parameters (c₀, c₁) by minimizing a combined entropy-penalty functional:</p><ul>
|
||
<li><b>Entropy term:</b> Promotes smooth, concentrated real spectra</li>
|
||
<li><b>Penalty term:</b> Penalizes negative values in real spectrum</li>
|
||
</ul>
|
||
<p>The phase correction is: φ(ω) = c₀ + c₁·ω where:</p><ul>
|
||
<li>c₀ corrects constant phase offset (time-zero uncertainty)</li>
|
||
<li>c₁ corrects linear phase dispersion (detector timing effects)</li>
|
||
</ul>
|
||
<p><b>Algorithm:</b></p><ol type="1">
|
||
<li>Initialize Minuit2 parameters with starting values (c₀=0, c₁=0)</li>
|
||
<li>Minimize entropy + γ×penalty using MnMinimize</li>
|
||
<li>Store optimal parameters in fPh_c0, fPh_c1</li>
|
||
</ol>
|
||
<p>After calling this method, use <a class="el" href="#abdac531b17bda660e3648eee992ed4ca">GetPhaseCorrectionParam()</a> to retrieve the optimal phase parameters and <a class="el" href="#a4a16bfd3cf495e2ee83026541ce14621">IsValid()</a> to check success.</p>
|
||
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#abdac531b17bda660e3648eee992ed4ca">GetPhaseCorrectionParam()</a> </dd>
|
||
<dd>
|
||
<a class="el" href="#a3c4e1295338f8110e91bef28abe430c0">SetGamma()</a> </dd>
|
||
<dd>
|
||
<a class="el" href="#a4a16bfd3cf495e2ee83026541ce14621">IsValid()</a> </dd></dl>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00144">144</a> of file <a class="el" href="PFourier_8cpp_source.html">PFourier.cpp</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00158">fMin</a>, <a class="el" href="PFourier_8h_source.html#l00155">fPh_c0</a>, <a class="el" href="PFourier_8h_source.html#l00156">fPh_c1</a>, and <a class="el" href="PFourier_8h_source.html#l00145">fValid</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00817">PFourier::GetPhaseOptRealFourier()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a6e4cc9da6c28cb1816b14a174e70694a" name="a6e4cc9da6c28cb1816b14a174e70694a"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a6e4cc9da6c28cb1816b14a174e70694a">◆ </a></span>operator()()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">double PFTPhaseCorrection::operator() </td>
|
||
<td>(</td>
|
||
<td class="paramtype">const std::vector< Double_t > &</td> <td class="paramname"><span class="paramname"><em></em></span></td><td>)</td>
|
||
<td> const</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Objective function for Minuit2 minimization (FCNBase interface).</p>
|
||
<p>Evaluates the combined entropy-penalty functional: f(c₀, c₁) = S + γ·P where:</p><ul>
|
||
<li>S = entropy of phase-corrected real spectrum derivative</li>
|
||
<li>P = penalty for negative values in real spectrum</li>
|
||
<li>γ = balancing parameter (fGamma)</li>
|
||
</ul>
|
||
<dl class="params"><dt>Parameters</dt><dd>
|
||
<table class="params">
|
||
<tr><td class="paramname">par</td><td>Parameter vector: [0]=c₀ (constant phase), [1]=c₁ (linear phase) </td></tr>
|
||
</table>
|
||
</dd>
|
||
</dl>
|
||
<dl class="section return"><dt>Returns</dt><dd>Combined functional value to be minimized</dd></dl>
|
||
<p><b>Algorithm flow:</b></p><ol type="1">
|
||
<li>Apply phase correction with parameters par</li>
|
||
<li>Calculate phase-corrected Fourier transform</li>
|
||
<li>Compute spectrum derivative</li>
|
||
<li>Return entropy + penalty</li>
|
||
</ol>
|
||
<p>Called repeatedly by Minuit2 during optimization.</p>
|
||
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#a4e8ddb15191632e59f4f3fcd86932199">Entropy()</a> </dd>
|
||
<dd>
|
||
<a class="el" href="#aa65f5dd4a52f0cb2e3fd8bba7b4a2945">Penalty()</a> </dd>
|
||
<dd>
|
||
<a class="el" href="#a05863bb58d111803ab36b922700bb667">Minimize()</a> </dd></dl>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00425">425</a> of file <a class="el" href="PFourier_8cpp_source.html">PFourier.cpp</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8cpp_source.html#l00276">CalcPhasedFT()</a>, <a class="el" href="PFourier_8cpp_source.html#l00307">CalcRealPhFTDerivative()</a>, <a class="el" href="PFourier_8cpp_source.html#l00376">Entropy()</a>, <a class="el" href="PFourier_8h_source.html#l00155">fPh_c0</a>, <a class="el" href="PFourier_8h_source.html#l00156">fPh_c1</a>, and <a class="el" href="PFourier_8cpp_source.html#l00339">Penalty()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="aa65f5dd4a52f0cb2e3fd8bba7b4a2945" name="aa65f5dd4a52f0cb2e3fd8bba7b4a2945"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#aa65f5dd4a52f0cb2e3fd8bba7b4a2945">◆ </a></span>Penalty()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">Double_t PFTPhaseCorrection::Penalty </td>
|
||
<td>(</td>
|
||
<td class="paramname"><span class="paramname"><em></em></span></td><td>)</td>
|
||
<td> const</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Calculates penalty term for negative real spectrum values.</p>
|
||
<p>Computes: P = γ·Σ[F_re(ω)² for all F_re(ω) < 0]</p>
|
||
<p>The penalty term penalizes negative values in the real Fourier spectrum, since physically meaningful absorption-mode spectra should be predominantly positive. The γ parameter controls the relative importance of this constraint.</p>
|
||
<dl class="section return"><dt>Returns</dt><dd>Weighted penalty value (γ×sum of squared negative amplitudes)</dd></dl>
|
||
<p><b>Note:</b> Only bins in range [fMinBin, fMaxBin) are considered, allowing focus on signal-containing regions.</p>
|
||
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#a3c4e1295338f8110e91bef28abe430c0">SetGamma()</a> </dd>
|
||
<dd>
|
||
<a class="el" href="#a4e8ddb15191632e59f4f3fcd86932199">Entropy()</a> </dd></dl>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00339">339</a> of file <a class="el" href="PFourier_8cpp_source.html">PFourier.cpp</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00157">fGamma</a>, <a class="el" href="PFourier_8h_source.html#l00154">fMaxBin</a>, <a class="el" href="PFourier_8h_source.html#l00153">fMinBin</a>, and <a class="el" href="PFourier_8h_source.html#l00149">fRealPh</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00425">operator()()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a3c4e1295338f8110e91bef28abe430c0" name="a3c4e1295338f8110e91bef28abe430c0"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a3c4e1295338f8110e91bef28abe430c0">◆ </a></span>SetGamma()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">virtual void PFTPhaseCorrection::SetGamma </td>
|
||
<td>(</td>
|
||
<td class="paramtype">const Double_t</td> <td class="paramname"><span class="paramname"><em>gamma</em></span></td><td>)</td>
|
||
<td></td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel inline">inline</span><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Sets the gamma balancing parameter between entropy and penalty </p><dl class="params"><dt>Parameters</dt><dd>
|
||
<table class="params">
|
||
<tr><td class="paramname">gamma</td><td>Balancing factor (typical range: 0.1 to 10) </td></tr>
|
||
</table>
|
||
</dd>
|
||
</dl>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00121">121</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8h_source.html#l00157">fGamma</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a70069a47f75e2ec309d278f9de2def7c" name="a70069a47f75e2ec309d278f9de2def7c"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a70069a47f75e2ec309d278f9de2def7c">◆ </a></span>SetPh()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">virtual void PFTPhaseCorrection::SetPh </td>
|
||
<td>(</td>
|
||
<td class="paramtype">const Double_t</td> <td class="paramname"><span class="paramname"><em>c0</em></span>, </td>
|
||
</tr>
|
||
<tr>
|
||
<td class="paramkey"></td>
|
||
<td></td>
|
||
<td class="paramtype">const Double_t</td> <td class="paramname"><span class="paramname"><em>c1</em></span> )</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel inline">inline</span><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
<p>Sets phase correction parameters manually </p><dl class="params"><dt>Parameters</dt><dd>
|
||
<table class="params">
|
||
<tr><td class="paramname">c0</td><td>Constant phase offset in degrees </td></tr>
|
||
<tr><td class="paramname">c1</td><td>Linear phase dispersion coefficient </td></tr>
|
||
</table>
|
||
</dd>
|
||
</dl>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00126">126</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">References <a class="el" href="PFourier_8cpp_source.html#l00276">CalcPhasedFT()</a>, <a class="el" href="PFourier_8cpp_source.html#l00307">CalcRealPhFTDerivative()</a>, <a class="el" href="PFourier_8h_source.html#l00155">fPh_c0</a>, and <a class="el" href="PFourier_8h_source.html#l00156">fPh_c1</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a1e2255b2b4d7644ff6ee489d5ce3082f" name="a1e2255b2b4d7644ff6ee489d5ce3082f"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a1e2255b2b4d7644ff6ee489d5ce3082f">◆ </a></span>Up()</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">virtual Double_t PFTPhaseCorrection::Up </td>
|
||
<td>(</td>
|
||
<td class="paramname"><span class="paramname"><em></em></span></td><td>)</td>
|
||
<td> const</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel inline">inline</span><span class="mlabel private">private</span><span class="mlabel virtual">virtual</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00166">166</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<h2 class="groupheader">Member Data Documentation</h2>
|
||
<a id="a8e5d6c8081d596f892034ce74081c34c" name="a8e5d6c8081d596f892034ce74081c34c"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a8e5d6c8081d596f892034ce74081c34c">◆ </a></span>fGamma</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">Double_t PFTPhaseCorrection::fGamma</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p>linear part of the phase dispersion used for the phase correction </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00157">157</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8h_source.html#l00130">GetGamma()</a>, <a class="el" href="PFourier_8cpp_source.html#l00248">Init()</a>, <a class="el" href="PFourier_8cpp_source.html#l00339">Penalty()</a>, and <a class="el" href="PFourier_8h_source.html#l00121">SetGamma()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="af424a8bd0dc0ef16a2537d0ac2ec0abb" name="af424a8bd0dc0ef16a2537d0ac2ec0abb"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#af424a8bd0dc0ef16a2537d0ac2ec0abb">◆ </a></span>fImag</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">std::vector<Double_t> PFTPhaseCorrection::fImag</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p>original real Fourier data set </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00148">148</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00276">CalcPhasedFT()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00100">PFTPhaseCorrection()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a4eb80579fb232cdd5dc76cc99101192b" name="a4eb80579fb232cdd5dc76cc99101192b"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a4eb80579fb232cdd5dc76cc99101192b">◆ </a></span>fImagPh</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">std::vector<Double_t> PFTPhaseCorrection::fImagPh</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel mutable">mutable</span><span class="mlabel private">private</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p>phased real Fourier data set </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00150">150</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00276">CalcPhasedFT()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00100">PFTPhaseCorrection()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="abcd9de24b2045dba8af6c534976d4417" name="abcd9de24b2045dba8af6c534976d4417"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#abcd9de24b2045dba8af6c534976d4417">◆ </a></span>fMaxBin</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">Int_t PFTPhaseCorrection::fMaxBin</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p>minimum bin from Fourier range to be used for the phase correction estimate </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00154">154</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00376">Entropy()</a>, <a class="el" href="PFourier_8cpp_source.html#l00339">Penalty()</a>, <a class="el" href="PFourier_8cpp_source.html#l00067">PFTPhaseCorrection()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00100">PFTPhaseCorrection()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a7d2ee1776fd5cee6cfc61cfd48409750" name="a7d2ee1776fd5cee6cfc61cfd48409750"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a7d2ee1776fd5cee6cfc61cfd48409750">◆ </a></span>fMin</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">Double_t PFTPhaseCorrection::fMin</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p>gamma parameter to balance between entropy and penalty </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00158">158</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00223">GetMinimum()</a>, <a class="el" href="PFourier_8cpp_source.html#l00248">Init()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00144">Minimize()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a1f1c220a61b47d89df38a57a400867cc" name="a1f1c220a61b47d89df38a57a400867cc"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a1f1c220a61b47d89df38a57a400867cc">◆ </a></span>fMinBin</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">Int_t PFTPhaseCorrection::fMinBin</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p>1st derivative of fRealPh </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00153">153</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00376">Entropy()</a>, <a class="el" href="PFourier_8cpp_source.html#l00339">Penalty()</a>, <a class="el" href="PFourier_8cpp_source.html#l00067">PFTPhaseCorrection()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00100">PFTPhaseCorrection()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a9e911366ef26bbd1b510ab85e2496e7e" name="a9e911366ef26bbd1b510ab85e2496e7e"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a9e911366ef26bbd1b510ab85e2496e7e">◆ </a></span>fPh_c0</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">Double_t PFTPhaseCorrection::fPh_c0</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel mutable">mutable</span><span class="mlabel private">private</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p>maximum bin from Fourier range to be used for the phase correction estimate </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00155">155</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00276">CalcPhasedFT()</a>, <a class="el" href="PFourier_8cpp_source.html#l00195">GetPhaseCorrectionParam()</a>, <a class="el" href="PFourier_8cpp_source.html#l00248">Init()</a>, <a class="el" href="PFourier_8cpp_source.html#l00144">Minimize()</a>, <a class="el" href="PFourier_8cpp_source.html#l00425">operator()()</a>, and <a class="el" href="PFourier_8h_source.html#l00126">SetPh()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a6f8a487092c8d8e5b2e10c797ae6368f" name="a6f8a487092c8d8e5b2e10c797ae6368f"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a6f8a487092c8d8e5b2e10c797ae6368f">◆ </a></span>fPh_c1</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">Double_t PFTPhaseCorrection::fPh_c1</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel mutable">mutable</span><span class="mlabel private">private</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p>constant part of the phase dispersion used for the phase correction </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00156">156</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00276">CalcPhasedFT()</a>, <a class="el" href="PFourier_8cpp_source.html#l00195">GetPhaseCorrectionParam()</a>, <a class="el" href="PFourier_8cpp_source.html#l00248">Init()</a>, <a class="el" href="PFourier_8cpp_source.html#l00144">Minimize()</a>, <a class="el" href="PFourier_8cpp_source.html#l00425">operator()()</a>, and <a class="el" href="PFourier_8h_source.html#l00126">SetPh()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a4f19633e8a329c6f54e12683d1b357e2" name="a4f19633e8a329c6f54e12683d1b357e2"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a4f19633e8a329c6f54e12683d1b357e2">◆ </a></span>fReal</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">std::vector<Double_t> PFTPhaseCorrection::fReal</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00147">147</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00276">CalcPhasedFT()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00100">PFTPhaseCorrection()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a1fa993c00f55291ed2b58c17ecb96d4b" name="a1fa993c00f55291ed2b58c17ecb96d4b"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a1fa993c00f55291ed2b58c17ecb96d4b">◆ </a></span>fRealPh</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">std::vector<Double_t> PFTPhaseCorrection::fRealPh</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel mutable">mutable</span><span class="mlabel private">private</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p>original imag Fourier data set </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00149">149</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00276">CalcPhasedFT()</a>, <a class="el" href="PFourier_8cpp_source.html#l00307">CalcRealPhFTDerivative()</a>, <a class="el" href="PFourier_8cpp_source.html#l00339">Penalty()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00100">PFTPhaseCorrection()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a48c629e27e8f5f8832c921d14d9a1b62" name="a48c629e27e8f5f8832c921d14d9a1b62"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a48c629e27e8f5f8832c921d14d9a1b62">◆ </a></span>fRealPhD</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">std::vector<Double_t> PFTPhaseCorrection::fRealPhD</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel mutable">mutable</span><span class="mlabel private">private</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p>phased imag Fourier data set </p>
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00151">151</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00307">CalcRealPhFTDerivative()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00376">Entropy()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<a id="a1d20df49ab2b5bba162eab9e0088396a" name="a1d20df49ab2b5bba162eab9e0088396a"></a>
|
||
<h2 class="memtitle"><span class="permalink"><a href="#a1d20df49ab2b5bba162eab9e0088396a">◆ </a></span>fValid</h2>
|
||
|
||
<div class="memitem">
|
||
<div class="memproto">
|
||
<table class="mlabels">
|
||
<tr>
|
||
<td class="mlabels-left">
|
||
<table class="memname">
|
||
<tr>
|
||
<td class="memname">Bool_t PFTPhaseCorrection::fValid</td>
|
||
</tr>
|
||
</table>
|
||
</td>
|
||
<td class="mlabels-right">
|
||
<span class="mlabels"><span class="mlabel private">private</span></span> </td>
|
||
</tr>
|
||
</table>
|
||
</div><div class="memdoc">
|
||
|
||
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00145">145</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
|
||
|
||
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00223">GetMinimum()</a>, <a class="el" href="PFourier_8cpp_source.html#l00248">Init()</a>, <a class="el" href="PFourier_8h_source.html#l00109">IsValid()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00144">Minimize()</a>.</p>
|
||
|
||
</div>
|
||
</div>
|
||
<hr/>The documentation for this class was generated from the following files:<ul>
|
||
<li>/workspace/LMU/musrfit/src/include/<a class="el" href="PFourier_8h_source.html">PFourier.h</a></li>
|
||
<li>/workspace/LMU/musrfit/src/classes/<a class="el" href="PFourier_8cpp_source.html">PFourier.cpp</a></li>
|
||
</ul>
|
||
</div><!-- contents -->
|
||
</div><!-- doc-content -->
|
||
<!-- start footer part -->
|
||
<div id="nav-path" class="navpath"><!-- id is needed for treeview function! -->
|
||
<ul>
|
||
<li class="navelem"><a class="el" href="classPFTPhaseCorrection.html">PFTPhaseCorrection</a></li>
|
||
<li class="footer">Generated by <a href="https://www.doxygen.org/index.html"><img class="footer" src="doxygen.svg" width="104" height="31" alt="doxygen"/></a> 1.13.2 </li>
|
||
</ul>
|
||
</div>
|
||
</body>
|
||
</html>
|