Files
musrfit/doc/technical/html/classPFourier.html
Gitea Actions 8a8eac55b7 Deploy site
2025-12-27 12:57:44 +00:00

1467 lines
88 KiB
HTML
Raw Permalink Blame History

This file contains ambiguous Unicode characters
This file contains Unicode characters that might be confused with other characters. If you think that this is intentional, you can safely ignore this warning. Use the Escape button to reveal them.
<!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: PFourier 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">&#160;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&amp;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&amp;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&amp;dn=expat.txt MIT */
$(function(){initNavTree('classPFourier.html',''); initResizable(true); });
/* @license-end */
</script>
<div id="doc-content">
<div class="header">
<div class="summary">
<a href="#pub-methods">Public Member Functions</a> &#124;
<a href="#pub-static-methods">Static Public Member Functions</a> &#124;
<a href="#pri-methods">Private Member Functions</a> &#124;
<a href="#pri-attribs">Private Attributes</a> &#124;
<a href="classPFourier-members.html">List of all members</a> </div>
<div class="headertitle"><div class="title">PFourier Class Reference</div></div>
</div><!--header-->
<div class="contents">
<p><code>#include &lt;<a class="el" href="PFourier_8h_source.html">PFourier.h</a>&gt;</code></p>
<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:a027eb3e5a96b030add0e389e54d68333" id="r_a027eb3e5a96b030add0e389e54d68333"><td class="memItemLeft" align="right" valign="top">&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a027eb3e5a96b030add0e389e54d68333">PFourier</a> (TH1F *data, Int_t unitTag, Double_t startTime=0.0, Double_t endTime=0.0, Bool_t dcCorrected=false, UInt_t zeroPaddingPower=0)</td></tr>
<tr class="separator:a027eb3e5a96b030add0e389e54d68333"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:af890268d2e2f28a4fd6a02cd71be51fa" id="r_af890268d2e2f28a4fd6a02cd71be51fa"><td class="memItemLeft" align="right" valign="top">virtual&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#af890268d2e2f28a4fd6a02cd71be51fa">~PFourier</a> ()</td></tr>
<tr class="separator:af890268d2e2f28a4fd6a02cd71be51fa"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a5ff69ded125558e11c25570cf04e8c54" id="r_a5ff69ded125558e11c25570cf04e8c54"><td class="memItemLeft" align="right" valign="top">virtual void&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a5ff69ded125558e11c25570cf04e8c54">Transform</a> (UInt_t apodizationTag=0)</td></tr>
<tr class="separator:a5ff69ded125558e11c25570cf04e8c54"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a86496d50f5014fe731c723a156535902" id="r_a86496d50f5014fe731c723a156535902"><td class="memItemLeft" align="right" valign="top">virtual const char *&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a86496d50f5014fe731c723a156535902">GetDataTitle</a> ()</td></tr>
<tr class="separator:a86496d50f5014fe731c723a156535902"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ab09ed8afcb5a92cba298b41141b91723" id="r_ab09ed8afcb5a92cba298b41141b91723"><td class="memItemLeft" align="right" valign="top">virtual const Int_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#ab09ed8afcb5a92cba298b41141b91723">GetUnitTag</a> ()</td></tr>
<tr class="separator:ab09ed8afcb5a92cba298b41141b91723"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a4a650ee34155ce54f1f477e3c047a4f9" id="r_a4a650ee34155ce54f1f477e3c047a4f9"><td class="memItemLeft" align="right" valign="top">virtual Double_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a4a650ee34155ce54f1f477e3c047a4f9">GetResolution</a> ()</td></tr>
<tr class="separator:a4a650ee34155ce54f1f477e3c047a4f9"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a7ac730688ae85e4fefedfaec9ababada" id="r_a7ac730688ae85e4fefedfaec9ababada"><td class="memItemLeft" align="right" valign="top">virtual Double_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a7ac730688ae85e4fefedfaec9ababada">GetMaxFreq</a> ()</td></tr>
<tr class="separator:a7ac730688ae85e4fefedfaec9ababada"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a7060dcc81a0e9e0cb38fe3b4b4bbf44e" id="r_a7060dcc81a0e9e0cb38fe3b4b4bbf44e"><td class="memItemLeft" align="right" valign="top">virtual TH1F *&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a7060dcc81a0e9e0cb38fe3b4b4bbf44e">GetRealFourier</a> (const Double_t scale=1.0)</td></tr>
<tr class="separator:a7060dcc81a0e9e0cb38fe3b4b4bbf44e"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a58b268c1e9aa23a393cbb2f2a97d10e5" id="r_a58b268c1e9aa23a393cbb2f2a97d10e5"><td class="memItemLeft" align="right" valign="top">virtual TH1F *&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a58b268c1e9aa23a393cbb2f2a97d10e5">GetImaginaryFourier</a> (const Double_t scale=1.0)</td></tr>
<tr class="separator:a58b268c1e9aa23a393cbb2f2a97d10e5"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a7fc73c74743bf1e2631240417fa78e70" id="r_a7fc73c74743bf1e2631240417fa78e70"><td class="memItemLeft" align="right" valign="top">virtual TH1F *&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a7fc73c74743bf1e2631240417fa78e70">GetPowerFourier</a> (const Double_t scale=1.0)</td></tr>
<tr class="separator:a7fc73c74743bf1e2631240417fa78e70"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:abceffb7fff41008808e0c68a325b06d4" id="r_abceffb7fff41008808e0c68a325b06d4"><td class="memItemLeft" align="right" valign="top">virtual TH1F *&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#abceffb7fff41008808e0c68a325b06d4">GetPhaseFourier</a> (const Double_t scale=1.0)</td></tr>
<tr class="separator:abceffb7fff41008808e0c68a325b06d4"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ab0189b365f88b2dc9ea4db60b94fac97" id="r_ab0189b365f88b2dc9ea4db60b94fac97"><td class="memItemLeft" align="right" valign="top">virtual Bool_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#ab0189b365f88b2dc9ea4db60b94fac97">IsValid</a> ()</td></tr>
<tr class="separator:ab0189b365f88b2dc9ea4db60b94fac97"><td class="memSeparator" colspan="2">&#160;</td></tr>
</table><table class="memberdecls">
<tr class="heading"><td colspan="2"><h2 class="groupheader"><a id="pub-static-methods" name="pub-static-methods"></a>
Static Public Member Functions</h2></td></tr>
<tr class="memitem:a84b7f1f225ffd828ab04d016401a9fb9" id="r_a84b7f1f225ffd828ab04d016401a9fb9"><td class="memItemLeft" align="right" valign="top">static TH1F *&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a84b7f1f225ffd828ab04d016401a9fb9">GetPhaseOptRealFourier</a> (const TH1F *re, const TH1F *im, std::vector&lt; Double_t &gt; &amp;phase, const Double_t scale=1.0, const Double_t min=-1.0, const Double_t max=-1.0)</td></tr>
<tr class="separator:a84b7f1f225ffd828ab04d016401a9fb9"><td class="memSeparator" colspan="2">&#160;</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:abb91c0e97948c918e451349a5db7e2a1" id="r_abb91c0e97948c918e451349a5db7e2a1"><td class="memItemLeft" align="right" valign="top">virtual void&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#abb91c0e97948c918e451349a5db7e2a1">PrepareFFTwInputData</a> (UInt_t apodizationTag)</td></tr>
<tr class="separator:abb91c0e97948c918e451349a5db7e2a1"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a117de8385fde0e673e3a3ff8d06d3f76" id="r_a117de8385fde0e673e3a3ff8d06d3f76"><td class="memItemLeft" align="right" valign="top">virtual void&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a117de8385fde0e673e3a3ff8d06d3f76">ApodizeData</a> (Int_t apodizationTag)</td></tr>
<tr class="separator:a117de8385fde0e673e3a3ff8d06d3f76"><td class="memSeparator" colspan="2">&#160;</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:ab04950b3a00e1feffb3e207809ff35a0" id="r_ab04950b3a00e1feffb3e207809ff35a0"><td class="memItemLeft" align="right" valign="top">TH1F *&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#ab04950b3a00e1feffb3e207809ff35a0">fData</a></td></tr>
<tr class="memdesc:ab04950b3a00e1feffb3e207809ff35a0"><td class="mdescLeft">&#160;</td><td class="mdescRight">data histogram to be Fourier transformed. <br /></td></tr>
<tr class="separator:ab04950b3a00e1feffb3e207809ff35a0"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a2e9e489e46d77eb0646da26a1ebbec04" id="r_a2e9e489e46d77eb0646da26a1ebbec04"><td class="memItemLeft" align="right" valign="top">Bool_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a2e9e489e46d77eb0646da26a1ebbec04">fValid</a></td></tr>
<tr class="memdesc:a2e9e489e46d77eb0646da26a1ebbec04"><td class="mdescLeft">&#160;</td><td class="mdescRight">true = all boundary conditions fullfilled and hence a Fourier transform can be performed. <br /></td></tr>
<tr class="separator:a2e9e489e46d77eb0646da26a1ebbec04"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:abeef1583cb8be5d5f59598f2bee1d5e3" id="r_abeef1583cb8be5d5f59598f2bee1d5e3"><td class="memItemLeft" align="right" valign="top">Int_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#abeef1583cb8be5d5f59598f2bee1d5e3">fUnitTag</a></td></tr>
<tr class="memdesc:abeef1583cb8be5d5f59598f2bee1d5e3"><td class="mdescLeft">&#160;</td><td class="mdescRight">1=Field Units (G), 2=Field Units (T), 3=Frequency Units (MHz), 4=Angular Frequency Units (Mc/s) <br /></td></tr>
<tr class="separator:abeef1583cb8be5d5f59598f2bee1d5e3"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a2387d87514c271c59ecc6e08c92f008d" id="r_a2387d87514c271c59ecc6e08c92f008d"><td class="memItemLeft" align="right" valign="top">Int_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a2387d87514c271c59ecc6e08c92f008d">fApodization</a></td></tr>
<tr class="memdesc:a2387d87514c271c59ecc6e08c92f008d"><td class="mdescLeft">&#160;</td><td class="mdescRight">0=none, 1=weak, 2=medium, 3=strong <br /></td></tr>
<tr class="separator:a2387d87514c271c59ecc6e08c92f008d"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ad56792e88f211c12ed82bb696908cad1" id="r_ad56792e88f211c12ed82bb696908cad1"><td class="memItemLeft" align="right" valign="top">Double_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#ad56792e88f211c12ed82bb696908cad1">fTimeResolution</a></td></tr>
<tr class="memdesc:ad56792e88f211c12ed82bb696908cad1"><td class="mdescLeft">&#160;</td><td class="mdescRight">time resolution of the data histogram in (us) <br /></td></tr>
<tr class="separator:ad56792e88f211c12ed82bb696908cad1"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ad5354796ff7f532f8c27025c14528aa7" id="r_ad5354796ff7f532f8c27025c14528aa7"><td class="memItemLeft" align="right" valign="top">Double_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#ad5354796ff7f532f8c27025c14528aa7">fStartTime</a></td></tr>
<tr class="memdesc:ad5354796ff7f532f8c27025c14528aa7"><td class="mdescLeft">&#160;</td><td class="mdescRight">start time of the data histogram <br /></td></tr>
<tr class="separator:ad5354796ff7f532f8c27025c14528aa7"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ae55a05bb68c9205fd2795fc21a56d32a" id="r_ae55a05bb68c9205fd2795fc21a56d32a"><td class="memItemLeft" align="right" valign="top">Double_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#ae55a05bb68c9205fd2795fc21a56d32a">fEndTime</a></td></tr>
<tr class="memdesc:ae55a05bb68c9205fd2795fc21a56d32a"><td class="mdescLeft">&#160;</td><td class="mdescRight">end time of the data histogram <br /></td></tr>
<tr class="separator:ae55a05bb68c9205fd2795fc21a56d32a"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ae7c452d31e52add9d7d829c1fc34a420" id="r_ae7c452d31e52add9d7d829c1fc34a420"><td class="memItemLeft" align="right" valign="top">Bool_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#ae7c452d31e52add9d7d829c1fc34a420">fDCCorrected</a></td></tr>
<tr class="memdesc:ae7c452d31e52add9d7d829c1fc34a420"><td class="mdescLeft">&#160;</td><td class="mdescRight">if true, removed DC offset from signal before Fourier transformation, otherwise not <br /></td></tr>
<tr class="separator:ae7c452d31e52add9d7d829c1fc34a420"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a38f3445e133d2958124b9e0369391eb2" id="r_a38f3445e133d2958124b9e0369391eb2"><td class="memItemLeft" align="right" valign="top">UInt_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a38f3445e133d2958124b9e0369391eb2">fZeroPaddingPower</a></td></tr>
<tr class="memdesc:a38f3445e133d2958124b9e0369391eb2"><td class="mdescLeft">&#160;</td><td class="mdescRight">power for zero padding, if set &lt; 0 no zero padding will be done <br /></td></tr>
<tr class="separator:a38f3445e133d2958124b9e0369391eb2"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a5c136ca96c952bc48447d24c56450c0d" id="r_a5c136ca96c952bc48447d24c56450c0d"><td class="memItemLeft" align="right" valign="top">Double_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a5c136ca96c952bc48447d24c56450c0d">fResolution</a></td></tr>
<tr class="memdesc:a5c136ca96c952bc48447d24c56450c0d"><td class="mdescLeft">&#160;</td><td class="mdescRight">Fourier resolution (field, frequency, or angular frequency) <br /></td></tr>
<tr class="separator:a5c136ca96c952bc48447d24c56450c0d"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a59f32dfa39affb141ca0fc9330ef991a" id="r_a59f32dfa39affb141ca0fc9330ef991a"><td class="memItemLeft" align="right" valign="top">UInt_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a59f32dfa39affb141ca0fc9330ef991a">fNoOfData</a></td></tr>
<tr class="memdesc:a59f32dfa39affb141ca0fc9330ef991a"><td class="mdescLeft">&#160;</td><td class="mdescRight">number of bins in the time interval between fStartTime and fStopTime <br /></td></tr>
<tr class="separator:a59f32dfa39affb141ca0fc9330ef991a"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a557fd214ccff7b637e9fa47711f104a7" id="r_a557fd214ccff7b637e9fa47711f104a7"><td class="memItemLeft" align="right" valign="top">UInt_t&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a557fd214ccff7b637e9fa47711f104a7">fNoOfBins</a></td></tr>
<tr class="memdesc:a557fd214ccff7b637e9fa47711f104a7"><td class="mdescLeft">&#160;</td><td class="mdescRight">number of bins to be Fourier transformed. Might be different to fNoOfData due to zero padding <br /></td></tr>
<tr class="separator:a557fd214ccff7b637e9fa47711f104a7"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a11669439c5ab5710295b23afb2f41da7" id="r_a11669439c5ab5710295b23afb2f41da7"><td class="memItemLeft" align="right" valign="top">fftw_plan&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#a11669439c5ab5710295b23afb2f41da7">fFFTwPlan</a></td></tr>
<tr class="memdesc:a11669439c5ab5710295b23afb2f41da7"><td class="mdescLeft">&#160;</td><td class="mdescRight">fftw plan (see FFTW3 User Manual) <br /></td></tr>
<tr class="separator:a11669439c5ab5710295b23afb2f41da7"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ad91d86ae6b3b668bd0f8499b8c34ee67" id="r_ad91d86ae6b3b668bd0f8499b8c34ee67"><td class="memItemLeft" align="right" valign="top">fftw_complex *&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#ad91d86ae6b3b668bd0f8499b8c34ee67">fIn</a></td></tr>
<tr class="memdesc:ad91d86ae6b3b668bd0f8499b8c34ee67"><td class="mdescLeft">&#160;</td><td class="mdescRight">real part of the Fourier transform <br /></td></tr>
<tr class="separator:ad91d86ae6b3b668bd0f8499b8c34ee67"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ae2877b25fbca566be2f4f45404775025" id="r_ae2877b25fbca566be2f4f45404775025"><td class="memItemLeft" align="right" valign="top">fftw_complex *&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="#ae2877b25fbca566be2f4f45404775025">fOut</a></td></tr>
<tr class="memdesc:ae2877b25fbca566be2f4f45404775025"><td class="mdescLeft">&#160;</td><td class="mdescRight">imaginary part of the Fourier transform <br /></td></tr>
<tr class="separator:ae2877b25fbca566be2f4f45404775025"><td class="memSeparator" colspan="2">&#160;</td></tr>
</table>
<a name="details" id="details"></a><h2 class="groupheader">Detailed Description</h2>
<div class="textblock"><p>Fourier transform engine for μSR time-domain data.</p>
<p><a class="el" href="classPFourier.html">PFourier</a> converts time-domain μSR signals to frequency domain, revealing:</p><ul>
<li>Muon precession frequencies (field measurements)</li>
<li>Internal field distributions (superconductors, magnets)</li>
<li>Multiple muon stopping sites</li>
<li>Dynamic frequency fluctuations</li>
</ul>
<p><b>Key features:</b></p><ul>
<li>Uses FFTW3 library for efficient FFT computation</li>
<li>DC offset removal (for baseline correction)</li>
<li>Zero-padding (improves frequency interpolation)</li>
<li>Apodization/windowing (reduces spectral leakage)</li>
<li>Multiple output formats (real, imaginary, power, phase)</li>
<li>Unit conversion (field ↔ frequency)</li>
</ul>
<p><b>Workflow:</b></p><ol type="1">
<li>Create <a class="el" href="classPFourier.html">PFourier</a> with time histogram and settings</li>
<li>Call <a class="el" href="#a5ff69ded125558e11c25570cf04e8c54">Transform()</a> with desired apodization</li>
<li>Retrieve results: <a class="el" href="#a7060dcc81a0e9e0cb38fe3b4b4bbf44e">GetRealFourier()</a>, <a class="el" href="#a7fc73c74743bf1e2631240417fa78e70">GetPowerFourier()</a>, etc.</li>
</ol>
<p><b>Unit conversions:</b></p><ul>
<li>Gauss: ω(MHz) = γ_μ/(2π) × B(G) = 0.01355 × B(G)</li>
<li>Tesla: ω(MHz) = γ_μ/(2π) × B(T) = 135.54 × B(T)</li>
</ul>
<p><b>Example:</b> TF-μSR measurement at 100 G produces a peak at ~1.36 MHz in the Fourier spectrum. </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00201">201</a> of file <a class="el" href="PFourier_8h_source.html">PFourier.h</a>.</p>
</div><h2 class="groupheader">Constructor &amp; Destructor Documentation</h2>
<a id="a027eb3e5a96b030add0e389e54d68333" name="a027eb3e5a96b030add0e389e54d68333"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a027eb3e5a96b030add0e389e54d68333">&#9670;&#160;</a></span>PFourier()</h2>
<div class="memitem">
<div class="memproto">
<table class="memname">
<tr>
<td class="memname">PFourier::PFourier </td>
<td>(</td>
<td class="paramtype">TH1F *</td> <td class="paramname"><span class="paramname"><em>data</em></span>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">Int_t</td> <td class="paramname"><span class="paramname"><em>unitTag</em></span>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">Double_t</td> <td class="paramname"><span class="paramname"><em>startTime</em></span><span class="paramdefsep"> = </span><span class="paramdefval">0.0</span>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">Double_t</td> <td class="paramname"><span class="paramname"><em>endTime</em></span><span class="paramdefsep"> = </span><span class="paramdefval">0.0</span>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">Bool_t</td> <td class="paramname"><span class="paramname"><em>dcCorrected</em></span><span class="paramdefsep"> = </span><span class="paramdefval">false</span>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">UInt_t</td> <td class="paramname"><span class="paramname"><em>zeroPaddingPower</em></span><span class="paramdefsep"> = </span><span class="paramdefval">0</span>&#160;)</td>
</tr>
</table>
</div><div class="memdoc">
<p>Constructor for Fourier transformation.</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">data</td><td>Time histogram to transform </td></tr>
<tr><td class="paramname">unitTag</td><td>Output units (1=Gauss, 2=Tesla, 3=MHz, 4=Mc/s) </td></tr>
<tr><td class="paramname">startTime</td><td>Start time for transform in microseconds (0=from t0) </td></tr>
<tr><td class="paramname">endTime</td><td>End time for transform in microseconds (0=to end) </td></tr>
<tr><td class="paramname">dcCorrected</td><td>If true, remove DC offset before FFT </td></tr>
<tr><td class="paramname">zeroPaddingPower</td><td>Zero-pad to 2^N points (0=no padding)</td></tr>
</table>
</dd>
</dl>
<p>Constructor for Fourier transform engine.</p>
<p>Creates a <a class="el" href="classPFourier.html">PFourier</a> object configured to transform μSR time-domain data to frequency/field domain. The constructor validates inputs, allocates FFTW resources, and prepares for FFT computation.</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">data</td><td>TH1F histogram containing time-domain μSR signal. X-axis: time in microseconds, Y-axis: asymmetry or counts. </td></tr>
<tr><td class="paramname">unitTag</td><td>Output units for frequency axis:<ul>
<li>FOURIER_UNIT_GAUSS (1): Magnetic field in Gauss</li>
<li>FOURIER_UNIT_TESLA (2): Magnetic field in Tesla</li>
<li>FOURIER_UNIT_FREQ (3): Frequency in MHz</li>
<li>FOURIER_UNIT_CYCLES (4): Angular frequency in Mc/s (Mrad/s) </li>
</ul>
</td></tr>
<tr><td class="paramname">startTime</td><td>Start of transform window in μs (0.0 = from t₀). Default: 0.0 </td></tr>
<tr><td class="paramname">endTime</td><td>End of transform window in μs (0.0 = to end). Default: 0.0 </td></tr>
<tr><td class="paramname">dcCorrected</td><td>If true, remove DC offset before FFT to eliminate zero-frequency artifact. Default: false </td></tr>
<tr><td class="paramname">zeroPaddingPower</td><td>Zero-pad to 2^N points where N=zeroPaddingPower. Value 0 disables padding. Padding improves frequency interpolation but doesn't add information. Default: 0</td></tr>
</table>
</dd>
</dl>
<p><b>Unit conversions (using γ_μ/2π = 0.01355 MHz/G):</b></p><ul>
<li>Gauss: B(G) = ω(MHz) / 0.01355</li>
<li>Tesla: B(T) = ω(MHz) / 135.54</li>
</ul>
<p><b>Example usage:</b> </p><div class="fragment"><div class="line">TH1F *timeHisto = ...; <span class="comment">// μSR time spectrum</span></div>
<div class="line"><a class="code hl_function" href="#a027eb3e5a96b030add0e389e54d68333">PFourier</a> ft(timeHisto, <a class="code hl_define" href="PMusr_8h.html#a471cf3179039d0b28e02b9554eab3021">FOURIER_UNIT_GAUSS</a>, 0.1, 10.0, <span class="keyword">true</span>, 12);</div>
<div class="line">ft.Transform(<a class="code hl_define" href="PFourier_8h.html#a150af5e514df766e1b13a810242da5b9">F_APODIZATION_WEAK</a>);</div>
<div class="line">TH1F *powerSpec = ft.GetPowerFourier();</div>
<div class="ttc" id="aPFourier_8h_html_a150af5e514df766e1b13a810242da5b9"><div class="ttname"><a href="PFourier_8h.html#a150af5e514df766e1b13a810242da5b9">F_APODIZATION_WEAK</a></div><div class="ttdeci">#define F_APODIZATION_WEAK</div><div class="ttdoc">Weak apodization (gentle roll-off at edges)</div><div class="ttdef"><b>Definition</b> <a href="PFourier_8h_source.html#l00055">PFourier.h:55</a></div></div>
<div class="ttc" id="aPMusr_8h_html_a471cf3179039d0b28e02b9554eab3021"><div class="ttname"><a href="PMusr_8h.html#a471cf3179039d0b28e02b9554eab3021">FOURIER_UNIT_GAUSS</a></div><div class="ttdeci">#define FOURIER_UNIT_GAUSS</div><div class="ttdoc">Magnetic field in Gauss (G)</div><div class="ttdef"><b>Definition</b> <a href="PMusr_8h_source.html#l00272">PMusr.h:272</a></div></div>
<div class="ttc" id="aclassPFourier_html_a027eb3e5a96b030add0e389e54d68333"><div class="ttname"><a href="#a027eb3e5a96b030add0e389e54d68333">PFourier::PFourier</a></div><div class="ttdeci">PFourier(TH1F *data, Int_t unitTag, Double_t startTime=0.0, Double_t endTime=0.0, Bool_t dcCorrected=false, UInt_t zeroPaddingPower=0)</div><div class="ttdef"><b>Definition</b> <a href="PFourier_8cpp_source.html#l00482">PFourier.cpp:482</a></div></div>
</div><!-- fragment --><p>After construction, check <a class="el" href="#ab0189b365f88b2dc9ea4db60b94fac97">IsValid()</a> before calling <a class="el" href="#a5ff69ded125558e11c25570cf04e8c54">Transform()</a>.</p>
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#a5ff69ded125558e11c25570cf04e8c54">Transform()</a> </dd>
<dd>
<a class="el" href="#ab0189b365f88b2dc9ea4db60b94fac97">IsValid()</a> </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00482">482</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#l00053">F_APODIZATION_NONE</a>, <a class="el" href="PFourier_8h_source.html#l00311">fApodization</a>, <a class="el" href="PFourier_8h_source.html#l00306">fData</a>, <a class="el" href="PFourier_8h_source.html#l00316">fDCCorrected</a>, <a class="el" href="PFourier_8h_source.html#l00315">fEndTime</a>, <a class="el" href="PFourier_8h_source.html#l00322">fFFTwPlan</a>, <a class="el" href="PFourier_8h_source.html#l00323">fIn</a>, <a class="el" href="PFourier_8h_source.html#l00321">fNoOfBins</a>, <a class="el" href="PFourier_8h_source.html#l00320">fNoOfData</a>, <a class="el" href="PMusr_8h_source.html#l00278">FOURIER_UNIT_CYCLES</a>, <a class="el" href="PMusr_8h_source.html#l00276">FOURIER_UNIT_FREQ</a>, <a class="el" href="PMusr_8h_source.html#l00272">FOURIER_UNIT_GAUSS</a>, <a class="el" href="PMusr_8h_source.html#l00274">FOURIER_UNIT_TESLA</a>, <a class="el" href="PFourier_8h_source.html#l00324">fOut</a>, <a class="el" href="PFourier_8h_source.html#l00318">fResolution</a>, <a class="el" href="PFourier_8h_source.html#l00314">fStartTime</a>, <a class="el" href="PFourier_8h_source.html#l00313">fTimeResolution</a>, <a class="el" href="PFourier_8h_source.html#l00309">fUnitTag</a>, <a class="el" href="PFourier_8h_source.html#l00308">fValid</a>, <a class="el" href="PFourier_8h_source.html#l00317">fZeroPaddingPower</a>, <a class="el" href="PMusr_8h_source.html#l00138">GAMMA_BAR_MUON</a>, and <a class="el" href="PFourier_8cpp_source.html#l00045">PI</a>.</p>
</div>
</div>
<a id="af890268d2e2f28a4fd6a02cd71be51fa" name="af890268d2e2f28a4fd6a02cd71be51fa"></a>
<h2 class="memtitle"><span class="permalink"><a href="#af890268d2e2f28a4fd6a02cd71be51fa">&#9670;&#160;</a></span>~PFourier()</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">PFourier::~PFourier </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>Destructor - cleans up FFTW resources.</p>
<p>Releases all FFTW-allocated memory and destroys the FFT plan. Proper cleanup is essential to avoid memory leaks, as FFTW uses special aligned memory allocation.</p>
<p><b>Resources freed:</b></p><ul>
<li>fFFTwPlan: FFTW execution plan</li>
<li>fIn: Input array (time-domain data)</li>
<li>fOut: Output array (frequency-domain data) </li>
</ul>
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00584">584</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#l00322">fFFTwPlan</a>, <a class="el" href="PFourier_8h_source.html#l00323">fIn</a>, and <a class="el" href="PFourier_8h_source.html#l00324">fOut</a>.</p>
</div>
</div>
<h2 class="groupheader">Member Function Documentation</h2>
<a id="a117de8385fde0e673e3a3ff8d06d3f76" name="a117de8385fde0e673e3a3ff8d06d3f76"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a117de8385fde0e673e3a3ff8d06d3f76">&#9670;&#160;</a></span>ApodizeData()</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">void PFourier::ApodizeData </td>
<td>(</td>
<td class="paramtype">Int_t</td> <td class="paramname"><span class="paramname"><em>apodizationTag</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>Applies apodization (windowing) to time-domain data before FFT.</p>
<p><b>Purpose:</b> Apodization multiplies the time signal by a smooth window function that tapers to zero at the edges, reducing the Gibbs phenomenon (spectral leakage) caused by finite observation windows.</p>
<p><b>Trade-off:</b></p><ul>
<li><b>Advantage:</b> Sharper, cleaner frequency peaks with less sidelobes</li>
<li><b>Disadvantage:</b> Slightly broader main peaks, reduced amplitude accuracy</li>
</ul>
<p><b>Window functions:</b> Polynomial windows of the form: w(t) = Σ c_j·(t/T)^(2j) where t ∈ [0,T]</p>
<p>Three predefined window strengths with optimized coefficients:</p><ul>
<li><b>Weak:</b> Minimal smoothing, preserves amplitude</li>
<li><b>Medium:</b> Balanced smoothing for general use</li>
<li><b>Strong:</b> Maximum leakage suppression for crowded spectra</li>
</ul>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">apodizationTag</td><td>Window function type:<ul>
<li>F_APODIZATION_NONE: No windowing (rectangular)</li>
<li>F_APODIZATION_WEAK: Gentle taper</li>
<li>F_APODIZATION_MEDIUM: Moderate taper</li>
<li>F_APODIZATION_STRONG: Heavy taper</li>
</ul>
</td></tr>
</table>
</dd>
</dl>
<p>Window is applied in-place to fIn[] array.</p>
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#abb91c0e97948c918e451349a5db7e2a1">PrepareFFTwInputData()</a> </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l01226">1226</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#l00057">F_APODIZATION_MEDIUM</a>, <a class="el" href="PFourier_8h_source.html#l00053">F_APODIZATION_NONE</a>, <a class="el" href="PFourier_8h_source.html#l00059">F_APODIZATION_STRONG</a>, <a class="el" href="PFourier_8h_source.html#l00055">F_APODIZATION_WEAK</a>, <a class="el" href="PFourier_8h_source.html#l00323">fIn</a>, and <a class="el" href="PFourier_8h_source.html#l00320">fNoOfData</a>.</p>
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l01155">PrepareFFTwInputData()</a>.</p>
</div>
</div>
<a id="a86496d50f5014fe731c723a156535902" name="a86496d50f5014fe731c723a156535902"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a86496d50f5014fe731c723a156535902">&#9670;&#160;</a></span>GetDataTitle()</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">virtual const char * PFourier::GetDataTitle </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 original data histogram title </p><dl class="section return"><dt>Returns</dt><dd>Title string </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00232">232</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#l00306">fData</a>.</p>
</div>
</div>
<a id="a58b268c1e9aa23a393cbb2f2a97d10e5" name="a58b268c1e9aa23a393cbb2f2a97d10e5"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a58b268c1e9aa23a393cbb2f2a97d10e5">&#9670;&#160;</a></span>GetImaginaryFourier()</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">TH1F * PFourier::GetImaginaryFourier </td>
<td>(</td>
<td class="paramtype">const Double_t</td> <td class="paramname"><span class="paramname"><em>scale</em></span><span class="paramdefsep"> = </span><span class="paramdefval">1.0</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 imaginary part of Fourier transform as histogram.</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">scale</td><td>Scaling factor for amplitudes (default=1.0) </td></tr>
</table>
</dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>Pointer to TH1F histogram (caller must delete)</dd></dl>
<p>Returns imaginary part of Fourier transform as ROOT histogram.</p>
<p>The imaginary component represents the out-of-phase (sine) projection of the signal. Ideally zero for perfect phase correction, but contains signal information when phase is uncorrected.</p>
<p><b>Uses:</b></p><ul>
<li>Phase correction optimization (minimize imaginary component)</li>
<li>Dispersion-mode spectroscopy</li>
<li>Computing power spectrum: |F|² = Re² + Im²</li>
<li>Phase spectrum calculation: φ = atan(Im/Re)</li>
</ul>
<p>Histogram covers [0, f_Nyquist] with N/2 bins where N is the FFT size (including zero padding if applied).</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">scale</td><td>Multiplication factor applied to all amplitudes. Use for normalization or unit conversion. Default: 1.0</td></tr>
</table>
</dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>Pointer to newly allocated TH1F histogram. <b>Caller owns the histogram and must delete it.</b> Returns nullptr on error.</dd></dl>
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#a7060dcc81a0e9e0cb38fe3b4b4bbf44e">GetRealFourier()</a> </dd>
<dd>
<a class="el" href="#a84b7f1f225ffd828ab04d016401a9fb9">GetPhaseOptRealFourier()</a> </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00931">931</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#l00306">fData</a>, <a class="el" href="PFourier_8h_source.html#l00321">fNoOfBins</a>, <a class="el" href="PFourier_8h_source.html#l00324">fOut</a>, <a class="el" href="PFourier_8h_source.html#l00318">fResolution</a>, and <a class="el" href="PFourier_8h_source.html#l00308">fValid</a>.</p>
<p class="reference">Referenced by <a class="el" href="PMusrCanvas_8cpp_source.html#l03478">PMusrCanvas::HandleDifferenceFourier()</a>, and <a class="el" href="PMusrCanvas_8cpp_source.html#l03332">PMusrCanvas::HandleFourier()</a>.</p>
</div>
</div>
<a id="a7ac730688ae85e4fefedfaec9ababada" name="a7ac730688ae85e4fefedfaec9ababada"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a7ac730688ae85e4fefedfaec9ababada">&#9670;&#160;</a></span>GetMaxFreq()</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">Double_t PFourier::GetMaxFreq </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 maximum frequency (Nyquist frequency).</p>
<dl class="section return"><dt>Returns</dt><dd>Maximum frequency in output units</dd></dl>
<p>Returns the maximum frequency (Nyquist frequency) of the spectrum.</p>
<p>The Nyquist frequency is the highest frequency that can be unambiguously represented given the time resolution: f_Nyquist = 1 / (2·Δt)</p>
<p>Higher frequencies are aliased back into the [0, f_Nyquist] range, so the useful spectrum extends from 0 to f_Nyquist.</p>
<dl class="section return"><dt>Returns</dt><dd>Maximum frequency in units specified by fUnitTag:<ul>
<li>Gauss (if FOURIER_UNIT_GAUSS)</li>
<li>Tesla (if FOURIER_UNIT_TESLA)</li>
<li>MHz (if FOURIER_UNIT_FREQ)</li>
<li>Mc/s (if FOURIER_UNIT_CYCLES)</li>
</ul>
</dd></dl>
<p><b>Example:</b> For time resolution Δt = 0.01 μs: f_Nyquist = 1/(2×0.01) = 50 MHz ≈ 3690 Gauss</p>
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#a4a650ee34155ce54f1f477e3c047a4f9">GetResolution()</a> </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00685">685</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#l00321">fNoOfBins</a>, and <a class="el" href="PFourier_8h_source.html#l00318">fResolution</a>.</p>
</div>
</div>
<a id="abceffb7fff41008808e0c68a325b06d4" name="abceffb7fff41008808e0c68a325b06d4"></a>
<h2 class="memtitle"><span class="permalink"><a href="#abceffb7fff41008808e0c68a325b06d4">&#9670;&#160;</a></span>GetPhaseFourier()</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">TH1F * PFourier::GetPhaseFourier </td>
<td>(</td>
<td class="paramtype">const Double_t</td> <td class="paramname"><span class="paramname"><em>scale</em></span><span class="paramdefsep"> = </span><span class="paramdefval">1.0</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 spectrum arg(F(ω)) as histogram.</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">scale</td><td>Scaling factor (default=1.0) </td></tr>
</table>
</dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>Pointer to TH1F histogram (caller must delete)</dd></dl>
<p>Returns phase spectrum φ(ω) = arg(F) as ROOT histogram.</p>
<p>The phase spectrum shows the phase angle of the complex Fourier transform at each frequency: φ(ω) = atan2(Im, Re)</p>
<p><b>Range:</b> Phase values are in [-π, +π] radians, with proper quadrant handling using the two-argument arctangent.</p>
<p><b>Interpretation:</b></p><ul>
<li>Constant phase: All frequencies in phase (simple exponential decay)</li>
<li>Linear phase: Time-zero offset (shift theorem)</li>
<li>Quadratic phase: Dispersion or chirp</li>
<li>Random phase: Noise</li>
</ul>
<p><b>Applications:</b></p><ul>
<li>Diagnosing time-zero problems</li>
<li>Verifying phase correction quality</li>
<li>Identifying signal vs. noise regions</li>
</ul>
<p>Histogram covers [0, f_Nyquist] with N/2 bins where N is the FFT size (including zero padding if applied).</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">scale</td><td>Multiplication factor applied to phase values. Typically 1.0 (radians) or 180/π (degrees). Default: 1.0</td></tr>
</table>
</dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>Pointer to newly allocated TH1F histogram. <b>Caller owns the histogram and must delete it.</b> Returns nullptr on error.</dd></dl>
<p><b>Note:</b> Phase is undefined where amplitude is zero. In practice, phase values in noise regions are meaningless.</p>
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#a7060dcc81a0e9e0cb38fe3b4b4bbf44e">GetRealFourier()</a> </dd>
<dd>
<a class="el" href="#a58b268c1e9aa23a393cbb2f2a97d10e5">GetImaginaryFourier()</a> </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l01073">1073</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#l00306">fData</a>, <a class="el" href="PFourier_8h_source.html#l00321">fNoOfBins</a>, <a class="el" href="PFourier_8h_source.html#l00324">fOut</a>, <a class="el" href="PFourier_8h_source.html#l00318">fResolution</a>, <a class="el" href="PFourier_8h_source.html#l00308">fValid</a>, <a class="el" href="PFourier_8cpp_source.html#l00045">PI</a>, and <a class="el" href="PFourier_8cpp_source.html#l00046">PI_HALF</a>.</p>
<p class="reference">Referenced by <a class="el" href="PMusrCanvas_8cpp_source.html#l03478">PMusrCanvas::HandleDifferenceFourier()</a>, and <a class="el" href="PMusrCanvas_8cpp_source.html#l03332">PMusrCanvas::HandleFourier()</a>.</p>
</div>
</div>
<a id="a84b7f1f225ffd828ab04d016401a9fb9" name="a84b7f1f225ffd828ab04d016401a9fb9"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a84b7f1f225ffd828ab04d016401a9fb9">&#9670;&#160;</a></span>GetPhaseOptRealFourier()</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">TH1F * PFourier::GetPhaseOptRealFourier </td>
<td>(</td>
<td class="paramtype">const TH1F *</td> <td class="paramname"><span class="paramname"><em>re</em></span>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const TH1F *</td> <td class="paramname"><span class="paramname"><em>im</em></span>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">std::vector&lt; Double_t &gt; &amp;</td> <td class="paramname"><span class="paramname"><em>phase</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>scale</em></span><span class="paramdefsep"> = </span><span class="paramdefval">1.0</span>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Double_t</td> <td class="paramname"><span class="paramname"><em>min</em></span><span class="paramdefsep"> = </span><span class="paramdefval">-1.0</span>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Double_t</td> <td class="paramname"><span class="paramname"><em>max</em></span><span class="paramdefsep"> = </span><span class="paramdefval">-1.0</span>&#160;)</td>
</tr>
</table>
</td>
<td class="mlabels-right">
<span class="mlabels"><span class="mlabel static">static</span></span> </td>
</tr>
</table>
</div><div class="memdoc">
<p>Static method for phase-optimized real Fourier spectrum.</p>
<p>Applies phase correction to maximize real component using provided phase parameters.</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">re</td><td>Real part of Fourier transform </td></tr>
<tr><td class="paramname">im</td><td>Imaginary part of Fourier transform </td></tr>
<tr><td class="paramname">phase</td><td>Phase correction parameters [c₀, c₁] </td></tr>
<tr><td class="paramname">scale</td><td>Scaling factor (default=1.0) </td></tr>
<tr><td class="paramname">min</td><td>Minimum frequency for correction (-1=all) </td></tr>
<tr><td class="paramname">max</td><td>Maximum frequency for correction (-1=all) </td></tr>
</table>
</dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>Pointer to phase-corrected TH1F histogram</dd></dl>
<p>Returns phase-optimized real Fourier spectrum (static method).</p>
<p>This static method performs automatic phase correction on a complex Fourier spectrum to maximize the real component (absorption mode) while minimizing the imaginary component. The algorithm finds optimal phase parameters φ(ω) = c₀ + c₁·ω that produce the cleanest real spectrum.</p>
<p><b>Why phase correction?</b> Raw FFT produces mixed-phase spectra with arbitrary phase offsets due to:</p><ul>
<li>Uncertain time-zero (t₀) determination</li>
<li>Detector timing offsets</li>
<li>Signal propagation delays</li>
</ul>
<p><b>Algorithm:</b></p><ol type="1">
<li>Extract Re/Im data in specified frequency range [min, max]</li>
<li>Use <a class="el" href="classPFTPhaseCorrection.html">PFTPhaseCorrection</a> to minimize entropy+penalty functional</li>
<li>Apply optimal phase rotation to entire spectrum</li>
<li>Return phase-corrected real part</li>
</ol>
<p><b>Optimization range [min, max]:</b> Restrict to signal-containing regions for best results. Including noisy baseline degrades optimization.</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">re</td><td>Real part of Fourier transform (TH1F histogram) </td></tr>
<tr><td class="paramname">im</td><td>Imaginary part of Fourier transform (TH1F histogram) </td></tr>
<tr><td class="paramname">phase</td><td>Output parameter: phase[0]=c₀, phase[1]=c₁ (vector resized to 2) </td></tr>
<tr><td class="paramname">scale</td><td>Amplitude scaling factor. Default: 1.0 </td></tr>
<tr><td class="paramname">min</td><td>Minimum frequency/field for optimization window. Use -1.0 to include all frequencies. Default: -1.0 </td></tr>
<tr><td class="paramname">max</td><td>Maximum frequency/field for optimization window. Use -1.0 to include all frequencies. Default: -1.0</td></tr>
</table>
</dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>Pointer to newly allocated phase-corrected TH1F histogram. <b>Caller owns and must delete.</b> Returns nullptr on error.</dd></dl>
<p><b>Example:</b> </p><div class="fragment"><div class="line">TH1F *real = ft.GetRealFourier();</div>
<div class="line">TH1F *imag = ft.GetImaginaryFourier();</div>
<div class="line">std::vector&lt;Double_t&gt; phaseParams;</div>
<div class="line">TH1F *phCorrected = <a class="code hl_function" href="#a84b7f1f225ffd828ab04d016401a9fb9">PFourier::GetPhaseOptRealFourier</a>(real, imag, phaseParams, 1.0, 5.0, 50.0);</div>
<div class="line">std::cout &lt;&lt; <span class="stringliteral">&quot;Phase: &quot;</span> &lt;&lt; phaseParams[0] &lt;&lt; <span class="stringliteral">&quot; + &quot;</span> &lt;&lt; phaseParams[1] &lt;&lt; <span class="stringliteral">&quot;*w&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">phCorrected-&gt;Draw();</div>
<div class="line"><span class="keyword">delete</span> real; <span class="keyword">delete</span> imag; <span class="keyword">delete</span> phCorrected;</div>
<div class="ttc" id="aclassPFourier_html_a84b7f1f225ffd828ab04d016401a9fb9"><div class="ttname"><a href="#a84b7f1f225ffd828ab04d016401a9fb9">PFourier::GetPhaseOptRealFourier</a></div><div class="ttdeci">static TH1F * GetPhaseOptRealFourier(const TH1F *re, const TH1F *im, std::vector&lt; Double_t &gt; &amp;phase, const Double_t scale=1.0, const Double_t min=-1.0, const Double_t max=-1.0)</div><div class="ttdef"><b>Definition</b> <a href="PFourier_8cpp_source.html#l00817">PFourier.cpp:817</a></div></div>
</div><!-- fragment --><dl class="section see"><dt>See also</dt><dd><a class="el" href="classPFTPhaseCorrection.html">PFTPhaseCorrection</a> </dd>
<dd>
<a class="el" href="#a7060dcc81a0e9e0cb38fe3b4b4bbf44e">GetRealFourier()</a> </dd>
<dd>
<a class="el" href="#a58b268c1e9aa23a393cbb2f2a97d10e5">GetImaginaryFourier()</a> </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00817">817</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#l00195">PFTPhaseCorrection::GetPhaseCorrectionParam()</a>, <a class="el" href="PFourier_8h_source.html#l00109">PFTPhaseCorrection::IsValid()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00144">PFTPhaseCorrection::Minimize()</a>.</p>
<p class="reference">Referenced by <a class="el" href="PMusrCanvas_8cpp_source.html#l04259">PMusrCanvas::CalcPhaseOptReFT()</a>.</p>
</div>
</div>
<a id="a7fc73c74743bf1e2631240417fa78e70" name="a7fc73c74743bf1e2631240417fa78e70"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a7fc73c74743bf1e2631240417fa78e70">&#9670;&#160;</a></span>GetPowerFourier()</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">TH1F * PFourier::GetPowerFourier </td>
<td>(</td>
<td class="paramtype">const Double_t</td> <td class="paramname"><span class="paramname"><em>scale</em></span><span class="paramdefsep"> = </span><span class="paramdefval">1.0</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 power spectrum |F(ω)|² as histogram.</p>
<p>Power spectrum is always positive and shows signal strength at each frequency, useful for identifying dominant frequencies.</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">scale</td><td>Scaling factor for power (default=1.0) </td></tr>
</table>
</dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>Pointer to TH1F histogram (caller must delete)</dd></dl>
<p>Returns power spectrum |F(ω)| as ROOT histogram.</p>
<p>The power spectrum shows signal magnitude at each frequency, computed as: |F(ω)| = √(Re² + Im²)</p>
<p><b>Advantages:</b></p><ul>
<li>Always positive (easier interpretation)</li>
<li>Phase-independent (no phase correction needed)</li>
<li>Directly shows signal strength vs. frequency</li>
<li>Ideal for identifying all frequency components</li>
</ul>
<p><b>Applications:</b></p><ul>
<li>Quick survey of frequency content</li>
<li>Field distribution analysis in superconductors</li>
<li>Multiple muon site identification</li>
<li>TF-μSR field measurements</li>
</ul>
<p>Histogram covers [0, f_Nyquist] with N/2 bins where N is the FFT size (including zero padding if applied).</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">scale</td><td>Multiplication factor applied to all amplitudes. Use for normalization or unit conversion. Default: 1.0</td></tr>
</table>
</dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>Pointer to newly allocated TH1F histogram. <b>Caller owns the histogram and must delete it.</b> Returns nullptr on error.</dd></dl>
<p><b>Note:</b> This returns the amplitude |F|, not the power |F|². For power values, square the histogram contents.</p>
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#a7060dcc81a0e9e0cb38fe3b4b4bbf44e">GetRealFourier()</a> </dd>
<dd>
<a class="el" href="#abceffb7fff41008808e0c68a325b06d4">GetPhaseFourier()</a> </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l01001">1001</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#l00306">fData</a>, <a class="el" href="PFourier_8h_source.html#l00321">fNoOfBins</a>, <a class="el" href="PFourier_8h_source.html#l00324">fOut</a>, <a class="el" href="PFourier_8h_source.html#l00318">fResolution</a>, and <a class="el" href="PFourier_8h_source.html#l00308">fValid</a>.</p>
<p class="reference">Referenced by <a class="el" href="PMusrCanvas_8cpp_source.html#l03478">PMusrCanvas::HandleDifferenceFourier()</a>, and <a class="el" href="PMusrCanvas_8cpp_source.html#l03332">PMusrCanvas::HandleFourier()</a>.</p>
</div>
</div>
<a id="a7060dcc81a0e9e0cb38fe3b4b4bbf44e" name="a7060dcc81a0e9e0cb38fe3b4b4bbf44e"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a7060dcc81a0e9e0cb38fe3b4b4bbf44e">&#9670;&#160;</a></span>GetRealFourier()</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">TH1F * PFourier::GetRealFourier </td>
<td>(</td>
<td class="paramtype">const Double_t</td> <td class="paramname"><span class="paramname"><em>scale</em></span><span class="paramdefsep"> = </span><span class="paramdefval">1.0</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 real part of Fourier transform as histogram.</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">scale</td><td>Scaling factor for amplitudes (default=1.0) </td></tr>
</table>
</dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>Pointer to TH1F histogram (caller must delete)</dd></dl>
<p>Returns real part of Fourier transform as ROOT histogram.</p>
<p>The real component represents the in-phase (cosine) projection of the signal. Without phase correction, both positive and negative values typically appear due to arbitrary phase offsets.</p>
<p><b>Interpretation:</b></p><ul>
<li>Phase-corrected: Absorption-mode spectrum (predominantly positive)</li>
<li>Uncorrected: Mixed-phase spectrum (positive + negative regions)</li>
</ul>
<p>Histogram covers [0, f_Nyquist] with N/2 bins where N is the FFT size (including zero padding if applied).</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">scale</td><td>Multiplication factor applied to all amplitudes. Use for normalization or unit conversion. Default: 1.0</td></tr>
</table>
</dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>Pointer to newly allocated TH1F histogram. <b>Caller owns the histogram and must delete it.</b> Returns nullptr on error.</dd></dl>
<p><b>Example:</b> </p><div class="fragment"><div class="line">ft.Transform(<a class="code hl_define" href="PFourier_8h.html#a150af5e514df766e1b13a810242da5b9">F_APODIZATION_WEAK</a>);</div>
<div class="line">TH1F *realSpec = ft.GetRealFourier(1.0);</div>
<div class="line">realSpec-&gt;Draw();</div>
<div class="line"><span class="comment">// ... use histogram ...</span></div>
<div class="line"><span class="keyword">delete</span> realSpec;</div>
</div><!-- fragment --><dl class="section see"><dt>See also</dt><dd><a class="el" href="#a58b268c1e9aa23a393cbb2f2a97d10e5">GetImaginaryFourier()</a> </dd>
<dd>
<a class="el" href="#a84b7f1f225ffd828ab04d016401a9fb9">GetPhaseOptRealFourier()</a> </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00731">731</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#l00306">fData</a>, <a class="el" href="PFourier_8h_source.html#l00321">fNoOfBins</a>, <a class="el" href="PFourier_8h_source.html#l00324">fOut</a>, <a class="el" href="PFourier_8h_source.html#l00318">fResolution</a>, and <a class="el" href="PFourier_8h_source.html#l00308">fValid</a>.</p>
<p class="reference">Referenced by <a class="el" href="PMusrCanvas_8cpp_source.html#l03478">PMusrCanvas::HandleDifferenceFourier()</a>, and <a class="el" href="PMusrCanvas_8cpp_source.html#l03332">PMusrCanvas::HandleFourier()</a>.</p>
</div>
</div>
<a id="a4a650ee34155ce54f1f477e3c047a4f9" name="a4a650ee34155ce54f1f477e3c047a4f9"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a4a650ee34155ce54f1f477e3c047a4f9">&#9670;&#160;</a></span>GetResolution()</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 PFourier::GetResolution </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 frequency resolution (bin width in output units) </p><dl class="section return"><dt>Returns</dt><dd>Frequency resolution </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00240">240</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#l00318">fResolution</a>.</p>
</div>
</div>
<a id="ab09ed8afcb5a92cba298b41141b91723" name="ab09ed8afcb5a92cba298b41141b91723"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ab09ed8afcb5a92cba298b41141b91723">&#9670;&#160;</a></span>GetUnitTag()</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">virtual const Int_t PFourier::GetUnitTag </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 output unit tag (1=G, 2=T, 3=MHz, 4=Mc/s) </p><dl class="section return"><dt>Returns</dt><dd>Unit identifier </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00236">236</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#l00309">fUnitTag</a>.</p>
</div>
</div>
<a id="ab0189b365f88b2dc9ea4db60b94fac97" name="ab0189b365f88b2dc9ea4db60b94fac97"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ab0189b365f88b2dc9ea4db60b94fac97">&#9670;&#160;</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 PFourier::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 Fourier transform is ready </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#l00303">303</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#l00308">fValid</a>.</p>
<p class="reference">Referenced by <a class="el" href="PMusrCanvas_8cpp_source.html#l03478">PMusrCanvas::HandleDifferenceFourier()</a>, and <a class="el" href="PMusrCanvas_8cpp_source.html#l03332">PMusrCanvas::HandleFourier()</a>.</p>
</div>
</div>
<a id="abb91c0e97948c918e451349a5db7e2a1" name="abb91c0e97948c918e451349a5db7e2a1"></a>
<h2 class="memtitle"><span class="permalink"><a href="#abb91c0e97948c918e451349a5db7e2a1">&#9670;&#160;</a></span>PrepareFFTwInputData()</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">void PFourier::PrepareFFTwInputData </td>
<td>(</td>
<td class="paramtype">UInt_t</td> <td class="paramname"><span class="paramname"><em>apodizationTag</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>Prepares time-domain data for FFT with optional DC correction and apodization.</p>
<p>This method performs essential preprocessing steps:</p><ol type="1">
<li>Locates time-zero bin in the histogram</li>
<li>Extracts data from [fStartTime, fEndTime] window</li>
<li>Removes DC offset (if fDCCorrected=true) for baseline correction</li>
<li>Zero-pads to fNoOfBins (power of 2 for optimal FFT performance)</li>
<li>Applies apodization window to reduce spectral leakage</li>
</ol>
<p><b>DC Correction:</b> Subtracts mean value from time signal to remove constant offset, which would otherwise create a large artifact at zero frequency.</p>
<p><b>Zero Padding:</b> Increases frequency resolution by interpolation without adding information content. Often used for smoother spectra.</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">apodizationTag</td><td>Window function strength:<ul>
<li>F_APODIZATION_NONE (1): Rectangular window (no smoothing)</li>
<li>F_APODIZATION_WEAK (2): Gentle roll-off at edges</li>
<li>F_APODIZATION_MEDIUM (3): Moderate smoothing</li>
<li>F_APODIZATION_STRONG (4): Heavy smoothing for maximum leakage suppression</li>
</ul>
</td></tr>
</table>
</dd>
</dl>
<p>Result is stored in fIn[] array ready for FFTW execution.</p>
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#a117de8385fde0e673e3a3ff8d06d3f76">ApodizeData()</a> </dd>
<dd>
<a class="el" href="#a5ff69ded125558e11c25570cf04e8c54">Transform()</a> </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l01155">1155</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#l01226">ApodizeData()</a>, <a class="el" href="PFourier_8h_source.html#l00306">fData</a>, <a class="el" href="PFourier_8h_source.html#l00316">fDCCorrected</a>, <a class="el" href="PFourier_8h_source.html#l00323">fIn</a>, <a class="el" href="PFourier_8h_source.html#l00321">fNoOfBins</a>, <a class="el" href="PFourier_8h_source.html#l00320">fNoOfData</a>, <a class="el" href="PFourier_8h_source.html#l00314">fStartTime</a>, and <a class="el" href="PFourier_8h_source.html#l00313">fTimeResolution</a>.</p>
<p class="reference">Referenced by <a class="el" href="PFourier_8cpp_source.html#l00632">Transform()</a>.</p>
</div>
</div>
<a id="a5ff69ded125558e11c25570cf04e8c54" name="a5ff69ded125558e11c25570cf04e8c54"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a5ff69ded125558e11c25570cf04e8c54">&#9670;&#160;</a></span>Transform()</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">void PFourier::Transform </td>
<td>(</td>
<td class="paramtype">UInt_t</td> <td class="paramname"><span class="paramname"><em>apodizationTag</em></span><span class="paramdefsep"> = </span><span class="paramdefval">0</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 the Fourier transformation.</p>
<p>Applies optional apodization, computes FFT using FFTW3, and prepares output histograms in requested units.</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">apodizationTag</td><td>Apodization strength (0/1=none, 2=weak, 3=medium, 4=strong)</td></tr>
</table>
</dd>
</dl>
<p>Executes the Fourier transform with optional apodization.</p>
<p>This is the main workhorse method that performs the complete FFT pipeline:</p><ol type="1">
<li>Prepare input data (DC correction, zero padding)</li>
<li>Apply apodization window if requested</li>
<li>Execute FFTW plan (compute FFT)</li>
<li>Correct phase for non-zero start time</li>
</ol>
<p><b>Phase correction:</b> If fStartTime ≠ 0, applies phase shift e^(i·ω·t₀) to account for the time-zero offset. This ensures proper phase relationships in the frequency spectrum.</p>
<p><b>Apodization:</b> Window functions reduce spectral leakage (Gibbs phenomenon) at the cost of slight peak broadening:</p><ul>
<li>NONE: Best amplitude accuracy, worst leakage</li>
<li>WEAK: Minimal smoothing</li>
<li>MEDIUM: Balanced (recommended for most cases)</li>
<li>STRONG: Maximum leakage suppression</li>
</ul>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramname">apodizationTag</td><td>Window function strength:<ul>
<li>0 or F_APODIZATION_NONE: No windowing</li>
<li>1 or F_APODIZATION_WEAK: Gentle taper</li>
<li>2 or F_APODIZATION_MEDIUM: Moderate taper</li>
<li>3 or F_APODIZATION_STRONG: Heavy taper</li>
</ul>
</td></tr>
</table>
</dd>
</dl>
<p>After calling <a class="el" href="#a5ff69ded125558e11c25570cf04e8c54">Transform()</a>, retrieve results using: <a class="el" href="#a7060dcc81a0e9e0cb38fe3b4b4bbf44e">GetRealFourier()</a>, <a class="el" href="#a58b268c1e9aa23a393cbb2f2a97d10e5">GetImaginaryFourier()</a>, <a class="el" href="#a7fc73c74743bf1e2631240417fa78e70">GetPowerFourier()</a>, <a class="el" href="#abceffb7fff41008808e0c68a325b06d4">GetPhaseFourier()</a></p>
<dl class="section see"><dt>See also</dt><dd><a class="el" href="#abb91c0e97948c918e451349a5db7e2a1">PrepareFFTwInputData()</a> </dd>
<dd>
<a class="el" href="#a7fc73c74743bf1e2631240417fa78e70">GetPowerFourier()</a> </dd></dl>
<p class="definition">Definition at line <a class="el" href="PFourier_8cpp_source.html#l00632">632</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#l00306">fData</a>, <a class="el" href="PFourier_8h_source.html#l00322">fFFTwPlan</a>, <a class="el" href="PFourier_8h_source.html#l00321">fNoOfBins</a>, <a class="el" href="PFourier_8h_source.html#l00324">fOut</a>, <a class="el" href="PFourier_8h_source.html#l00314">fStartTime</a>, <a class="el" href="PFourier_8h_source.html#l00313">fTimeResolution</a>, <a class="el" href="PFourier_8h_source.html#l00308">fValid</a>, <a class="el" href="PFourier_8cpp_source.html#l00045">PI</a>, and <a class="el" href="PFourier_8cpp_source.html#l01155">PrepareFFTwInputData()</a>.</p>
<p class="reference">Referenced by <a class="el" href="PMusrCanvas_8cpp_source.html#l03478">PMusrCanvas::HandleDifferenceFourier()</a>, and <a class="el" href="PMusrCanvas_8cpp_source.html#l03332">PMusrCanvas::HandleFourier()</a>.</p>
</div>
</div>
<h2 class="groupheader">Member Data Documentation</h2>
<a id="a2387d87514c271c59ecc6e08c92f008d" name="a2387d87514c271c59ecc6e08c92f008d"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a2387d87514c271c59ecc6e08c92f008d">&#9670;&#160;</a></span>fApodization</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">Int_t PFourier::fApodization</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>0=none, 1=weak, 2=medium, 3=strong </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00311">311</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#l00482">PFourier()</a>.</p>
</div>
</div>
<a id="ab04950b3a00e1feffb3e207809ff35a0" name="ab04950b3a00e1feffb3e207809ff35a0"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ab04950b3a00e1feffb3e207809ff35a0">&#9670;&#160;</a></span>fData</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">TH1F* PFourier::fData</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>data histogram to be Fourier transformed. </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00306">306</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#l00232">GetDataTitle()</a>, <a class="el" href="PFourier_8cpp_source.html#l00931">GetImaginaryFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01073">GetPhaseFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01001">GetPowerFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l00731">GetRealFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l00482">PFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01155">PrepareFFTwInputData()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00632">Transform()</a>.</p>
</div>
</div>
<a id="ae7c452d31e52add9d7d829c1fc34a420" name="ae7c452d31e52add9d7d829c1fc34a420"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ae7c452d31e52add9d7d829c1fc34a420">&#9670;&#160;</a></span>fDCCorrected</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">Bool_t PFourier::fDCCorrected</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>if true, removed DC offset from signal before Fourier transformation, otherwise not </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00316">316</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#l00482">PFourier()</a>, and <a class="el" href="PFourier_8cpp_source.html#l01155">PrepareFFTwInputData()</a>.</p>
</div>
</div>
<a id="ae55a05bb68c9205fd2795fc21a56d32a" name="ae55a05bb68c9205fd2795fc21a56d32a"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ae55a05bb68c9205fd2795fc21a56d32a">&#9670;&#160;</a></span>fEndTime</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">Double_t PFourier::fEndTime</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>end time of the data histogram </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00315">315</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#l00482">PFourier()</a>.</p>
</div>
</div>
<a id="a11669439c5ab5710295b23afb2f41da7" name="a11669439c5ab5710295b23afb2f41da7"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a11669439c5ab5710295b23afb2f41da7">&#9670;&#160;</a></span>fFFTwPlan</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">fftw_plan PFourier::fFFTwPlan</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>fftw plan (see FFTW3 User Manual) </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00322">322</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#l00482">PFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l00632">Transform()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00584">~PFourier()</a>.</p>
</div>
</div>
<a id="ad91d86ae6b3b668bd0f8499b8c34ee67" name="ad91d86ae6b3b668bd0f8499b8c34ee67"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ad91d86ae6b3b668bd0f8499b8c34ee67">&#9670;&#160;</a></span>fIn</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">fftw_complex* PFourier::fIn</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>real part of the Fourier transform </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00323">323</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#l01226">ApodizeData()</a>, <a class="el" href="PFourier_8cpp_source.html#l00482">PFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01155">PrepareFFTwInputData()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00584">~PFourier()</a>.</p>
</div>
</div>
<a id="a557fd214ccff7b637e9fa47711f104a7" name="a557fd214ccff7b637e9fa47711f104a7"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a557fd214ccff7b637e9fa47711f104a7">&#9670;&#160;</a></span>fNoOfBins</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">UInt_t PFourier::fNoOfBins</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>number of bins to be Fourier transformed. Might be different to fNoOfData due to zero padding </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00321">321</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#l00931">GetImaginaryFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l00685">GetMaxFreq()</a>, <a class="el" href="PFourier_8cpp_source.html#l01073">GetPhaseFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01001">GetPowerFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l00731">GetRealFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l00482">PFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01155">PrepareFFTwInputData()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00632">Transform()</a>.</p>
</div>
</div>
<a id="a59f32dfa39affb141ca0fc9330ef991a" name="a59f32dfa39affb141ca0fc9330ef991a"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a59f32dfa39affb141ca0fc9330ef991a">&#9670;&#160;</a></span>fNoOfData</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">UInt_t PFourier::fNoOfData</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>number of bins in the time interval between fStartTime and fStopTime </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00320">320</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#l01226">ApodizeData()</a>, <a class="el" href="PFourier_8cpp_source.html#l00482">PFourier()</a>, and <a class="el" href="PFourier_8cpp_source.html#l01155">PrepareFFTwInputData()</a>.</p>
</div>
</div>
<a id="ae2877b25fbca566be2f4f45404775025" name="ae2877b25fbca566be2f4f45404775025"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ae2877b25fbca566be2f4f45404775025">&#9670;&#160;</a></span>fOut</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">fftw_complex* PFourier::fOut</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>imaginary part of the Fourier transform </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00324">324</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#l00931">GetImaginaryFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01073">GetPhaseFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01001">GetPowerFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l00731">GetRealFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l00482">PFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l00632">Transform()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00584">~PFourier()</a>.</p>
</div>
</div>
<a id="a5c136ca96c952bc48447d24c56450c0d" name="a5c136ca96c952bc48447d24c56450c0d"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a5c136ca96c952bc48447d24c56450c0d">&#9670;&#160;</a></span>fResolution</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">Double_t PFourier::fResolution</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>Fourier resolution (field, frequency, or angular frequency) </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00318">318</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#l00931">GetImaginaryFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l00685">GetMaxFreq()</a>, <a class="el" href="PFourier_8cpp_source.html#l01073">GetPhaseFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01001">GetPowerFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l00731">GetRealFourier()</a>, <a class="el" href="PFourier_8h_source.html#l00240">GetResolution()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00482">PFourier()</a>.</p>
</div>
</div>
<a id="ad5354796ff7f532f8c27025c14528aa7" name="ad5354796ff7f532f8c27025c14528aa7"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ad5354796ff7f532f8c27025c14528aa7">&#9670;&#160;</a></span>fStartTime</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">Double_t PFourier::fStartTime</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>start time of the data histogram </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00314">314</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#l00482">PFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01155">PrepareFFTwInputData()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00632">Transform()</a>.</p>
</div>
</div>
<a id="ad56792e88f211c12ed82bb696908cad1" name="ad56792e88f211c12ed82bb696908cad1"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ad56792e88f211c12ed82bb696908cad1">&#9670;&#160;</a></span>fTimeResolution</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">Double_t PFourier::fTimeResolution</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>time resolution of the data histogram in (us) </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00313">313</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#l00482">PFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01155">PrepareFFTwInputData()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00632">Transform()</a>.</p>
</div>
</div>
<a id="abeef1583cb8be5d5f59598f2bee1d5e3" name="abeef1583cb8be5d5f59598f2bee1d5e3"></a>
<h2 class="memtitle"><span class="permalink"><a href="#abeef1583cb8be5d5f59598f2bee1d5e3">&#9670;&#160;</a></span>fUnitTag</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">Int_t PFourier::fUnitTag</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>1=Field Units (G), 2=Field Units (T), 3=Frequency Units (MHz), 4=Angular Frequency Units (Mc/s) </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00309">309</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#l00236">GetUnitTag()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00482">PFourier()</a>.</p>
</div>
</div>
<a id="a2e9e489e46d77eb0646da26a1ebbec04" name="a2e9e489e46d77eb0646da26a1ebbec04"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a2e9e489e46d77eb0646da26a1ebbec04">&#9670;&#160;</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 PFourier::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>true = all boundary conditions fullfilled and hence a Fourier transform can be performed. </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00308">308</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#l00931">GetImaginaryFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01073">GetPhaseFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l01001">GetPowerFourier()</a>, <a class="el" href="PFourier_8cpp_source.html#l00731">GetRealFourier()</a>, <a class="el" href="PFourier_8h_source.html#l00303">IsValid()</a>, <a class="el" href="PFourier_8cpp_source.html#l00482">PFourier()</a>, and <a class="el" href="PFourier_8cpp_source.html#l00632">Transform()</a>.</p>
</div>
</div>
<a id="a38f3445e133d2958124b9e0369391eb2" name="a38f3445e133d2958124b9e0369391eb2"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a38f3445e133d2958124b9e0369391eb2">&#9670;&#160;</a></span>fZeroPaddingPower</h2>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
<tr>
<td class="mlabels-left">
<table class="memname">
<tr>
<td class="memname">UInt_t PFourier::fZeroPaddingPower</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>power for zero padding, if set &lt; 0 no zero padding will be done </p>
<p class="definition">Definition at line <a class="el" href="PFourier_8h_source.html#l00317">317</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#l00482">PFourier()</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="classPFourier.html">PFourier</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>